/[MITgcm]/MITgcm/model/src/thermodynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/thermodynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.21 - (show annotations) (download)
Tue Apr 30 22:53:17 2002 UTC (22 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint45a_post
Changes since 1.20: +11 -9 lines
extend the computation domain of Rho to cover the whole domain.
 (necessary when using GM-Advect form + stagger time step with Olx=Oly=3)

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.20 2002/03/24 02:36:39 heimbach Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #ifdef ALLOW_AUTODIFF_TAMC
6 # ifdef ALLOW_GMREDI
7 # include "GMREDI_OPTIONS.h"
8 # endif
9 # ifdef ALLOW_KPP
10 # include "KPP_OPTIONS.h"
11 # endif
12 #endif /* ALLOW_AUTODIFF_TAMC */
13
14 CBOP
15 C !ROUTINE: THERMODYNAMICS
16 C !INTERFACE:
17 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
18 C !DESCRIPTION: \bv
19 C *==========================================================*
20 C | SUBROUTINE THERMODYNAMICS
21 C | o Controlling routine for the prognostic part of the
22 C | thermo-dynamics.
23 C *===========================================================
24 C | The algorithm...
25 C |
26 C | "Correction Step"
27 C | =================
28 C | Here we update the horizontal velocities with the surface
29 C | pressure such that the resulting flow is either consistent
30 C | with the free-surface evolution or the rigid-lid:
31 C | U[n] = U* + dt x d/dx P
32 C | V[n] = V* + dt x d/dy P
33 C |
34 C | "Calculation of Gs"
35 C | ===================
36 C | This is where all the accelerations and tendencies (ie.
37 C | physics, parameterizations etc...) are calculated
38 C | rho = rho ( theta[n], salt[n] )
39 C | b = b(rho, theta)
40 C | K31 = K31 ( rho )
41 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
42 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
43 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
44 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
45 C |
46 C | "Time-stepping" or "Prediction"
47 C | ================================
48 C | The models variables are stepped forward with the appropriate
49 C | time-stepping scheme (currently we use Adams-Bashforth II)
50 C | - For momentum, the result is always *only* a "prediction"
51 C | in that the flow may be divergent and will be "corrected"
52 C | later with a surface pressure gradient.
53 C | - Normally for tracers the result is the new field at time
54 C | level [n+1} *BUT* in the case of implicit diffusion the result
55 C | is also *only* a prediction.
56 C | - We denote "predictors" with an asterisk (*).
57 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
58 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
59 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
60 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
61 C | With implicit diffusion:
62 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
64 C | (1 + dt * K * d_zz) theta[n] = theta*
65 C | (1 + dt * K * d_zz) salt[n] = salt*
66 C |
67 C *==========================================================*
68 C \ev
69
70 C !USES:
71 IMPLICIT NONE
72 C == Global variables ===
73 #include "SIZE.h"
74 #include "EEPARAMS.h"
75 #include "PARAMS.h"
76 #include "DYNVARS.h"
77 #include "GRID.h"
78 #include "GAD.h"
79 #ifdef ALLOW_PASSIVE_TRACER
80 #include "TR1.h"
81 #endif
82 #ifdef ALLOW_AUTODIFF_TAMC
83 # include "tamc.h"
84 # include "tamc_keys.h"
85 # include "FFIELDS.h"
86 # ifdef ALLOW_KPP
87 # include "KPP.h"
88 # endif
89 # ifdef ALLOW_GMREDI
90 # include "GMREDI.h"
91 # endif
92 #endif /* ALLOW_AUTODIFF_TAMC */
93 #ifdef ALLOW_TIMEAVE
94 #include "TIMEAVE_STATV.h"
95 #endif
96
97 C !INPUT/OUTPUT PARAMETERS:
98 C == Routine arguments ==
99 C myTime - Current time in simulation
100 C myIter - Current iteration number in simulation
101 C myThid - Thread number for this instance of the routine.
102 _RL myTime
103 INTEGER myIter
104 INTEGER myThid
105
106 C !LOCAL VARIABLES:
107 C == Local variables
108 C xA, yA - Per block temporaries holding face areas
109 C uTrans, vTrans, rTrans - Per block temporaries holding flow
110 C transport
111 C o uTrans: Zonal transport
112 C o vTrans: Meridional transport
113 C o rTrans: Vertical transport
114 C maskUp o maskUp: land/water mask for W points
115 C fVer[STUV] o fVer: Vertical flux term - note fVer
116 C is "pipelined" in the vertical
117 C so we need an fVer for each
118 C variable.
119 C rhoK, rhoKM1 - Density at current level, and level above
120 C phiHyd - Hydrostatic part of the potential phiHydi.
121 C In z coords phiHydiHyd is the hydrostatic
122 C Potential (=pressure/rho0) anomaly
123 C In p coords phiHydiHyd is the geopotential
124 C surface height anomaly.
125 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
126 C phiSurfY or geopotentiel (atmos) in X and Y direction
127 C KappaRT, - Total diffusion in vertical for T and S.
128 C KappaRS (background + spatially varying, isopycnal term).
129 C iMin, iMax - Ranges and sub-block indices on which calculations
130 C jMin, jMax are applied.
131 C bi, bj
132 C k, kup, - Index for layer above and below. kup and kDown
133 C kDown, km1 are switched with layer to be the appropriate
134 C index into fVerTerm.
135 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
136 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
137 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
141 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
142 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
143 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
144 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
145 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
150 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
151 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
152 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
153 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
154 C This is currently used by IVDC and Diagnostics
155 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
156 INTEGER iMin, iMax
157 INTEGER jMin, jMax
158 INTEGER bi, bj
159 INTEGER i, j
160 INTEGER k, km1, kup, kDown
161
162 CEOP
163
164 #ifdef ALLOW_AUTODIFF_TAMC
165 C-- dummy statement to end declaration part
166 ikey = 1
167 #endif /* ALLOW_AUTODIFF_TAMC */
168
169 C-- Set up work arrays with valid (i.e. not NaN) values
170 C These inital values do not alter the numerical results. They
171 C just ensure that all memory references are to valid floating
172 C point numbers. This prevents spurious hardware signals due to
173 C uninitialised but inert locations.
174 DO j=1-OLy,sNy+OLy
175 DO i=1-OLx,sNx+OLx
176 xA(i,j) = 0. _d 0
177 yA(i,j) = 0. _d 0
178 uTrans(i,j) = 0. _d 0
179 vTrans(i,j) = 0. _d 0
180 rhok (i,j) = 0. _d 0
181 phiSurfX(i,j) = 0. _d 0
182 phiSurfY(i,j) = 0. _d 0
183 ENDDO
184 ENDDO
185
186
187 #ifdef ALLOW_AUTODIFF_TAMC
188 C-- HPF directive to help TAMC
189 CHPF$ INDEPENDENT
190 #endif /* ALLOW_AUTODIFF_TAMC */
191
192 DO bj=myByLo(myThid),myByHi(myThid)
193
194 #ifdef ALLOW_AUTODIFF_TAMC
195 C-- HPF directive to help TAMC
196 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
197 CHPF$& ,phiHyd,utrans,vtrans,xA,yA
198 CHPF$& ,KappaRT,KappaRS
199 CHPF$& )
200 #endif /* ALLOW_AUTODIFF_TAMC */
201
202 DO bi=myBxLo(myThid),myBxHi(myThid)
203
204 #ifdef ALLOW_AUTODIFF_TAMC
205 act1 = bi - myBxLo(myThid)
206 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
207 act2 = bj - myByLo(myThid)
208 max2 = myByHi(myThid) - myByLo(myThid) + 1
209 act3 = myThid - 1
210 max3 = nTx*nTy
211 act4 = ikey_dynamics - 1
212 ikey = (act1 + 1) + act2*max1
213 & + act3*max1*max2
214 & + act4*max1*max2*max3
215 #endif /* ALLOW_AUTODIFF_TAMC */
216
217 C-- Set up work arrays that need valid initial values
218 DO j=1-OLy,sNy+OLy
219 DO i=1-OLx,sNx+OLx
220 rTrans (i,j) = 0. _d 0
221 fVerT (i,j,1) = 0. _d 0
222 fVerT (i,j,2) = 0. _d 0
223 fVerS (i,j,1) = 0. _d 0
224 fVerS (i,j,2) = 0. _d 0
225 fVerTr1(i,j,1) = 0. _d 0
226 fVerTr1(i,j,2) = 0. _d 0
227 rhoKM1 (i,j) = 0. _d 0
228 ENDDO
229 ENDDO
230
231 DO k=1,Nr
232 DO j=1-OLy,sNy+OLy
233 DO i=1-OLx,sNx+OLx
234 C This is currently also used by IVDC and Diagnostics
235 phiHyd(i,j,k) = 0. _d 0
236 sigmaX(i,j,k) = 0. _d 0
237 sigmaY(i,j,k) = 0. _d 0
238 sigmaR(i,j,k) = 0. _d 0
239 ConvectCount(i,j,k) = 0.
240 KappaRT(i,j,k) = 0. _d 0
241 KappaRS(i,j,k) = 0. _d 0
242 #ifdef ALLOW_AUTODIFF_TAMC
243 gT(i,j,k,bi,bj) = 0. _d 0
244 gS(i,j,k,bi,bj) = 0. _d 0
245 #ifdef ALLOW_PASSIVE_TRACER
246 gTr1(i,j,k,bi,bj) = 0. _d 0
247 #endif
248 #ifdef ALLOW_GMREDI
249 Kwx(i,j,k,bi,bj) = 0. _d 0
250 Kwy(i,j,k,bi,bj) = 0. _d 0
251 Kwz(i,j,k,bi,bj) = 0. _d 0
252 #ifdef GM_NON_UNITY_DIAGONAL
253 Kux(i,j,k,bi,bj) = 0. _d 0
254 Kvy(i,j,k,bi,bj) = 0. _d 0
255 #endif
256 #endif /* ALLOW_GMREDI */
257 #endif
258 ENDDO
259 ENDDO
260 ENDDO
261
262 iMin = 1-OLx
263 iMax = sNx+OLx
264 jMin = 1-OLy
265 jMax = sNy+OLy
266
267
268 #ifdef ALLOW_AUTODIFF_TAMC
269 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
270 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
271 #ifdef ALLOW_KPP
272 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
273 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
274 #endif
275 #endif /* ALLOW_AUTODIFF_TAMC */
276
277 C-- Start of diagnostic loop
278 DO k=Nr,1,-1
279
280 #ifdef ALLOW_AUTODIFF_TAMC
281 C? Patrick, is this formula correct now that we change the loop range?
282 C? Do we still need this?
283 cph kkey formula corrected.
284 cph Needed for rhok, rhokm1, in the case useGMREDI.
285 kkey = (ikey-1)*Nr + k
286 CADJ STORE rhokm1(:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
287 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
288 #endif /* ALLOW_AUTODIFF_TAMC */
289
290 C-- Integrate continuity vertically for vertical velocity
291 CALL INTEGRATE_FOR_W(
292 I bi, bj, k, uVel, vVel,
293 O wVel,
294 I myThid )
295
296 #ifdef ALLOW_OBCS
297 #ifdef ALLOW_NONHYDROSTATIC
298 C-- Apply OBC to W if in N-H mode
299 IF (useOBCS.AND.nonHydrostatic) THEN
300 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
301 ENDIF
302 #endif /* ALLOW_NONHYDROSTATIC */
303 #endif /* ALLOW_OBCS */
304
305 C-- Calculate gradients of potential density for isoneutral
306 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
307 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
308 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
309 #ifdef ALLOW_AUTODIFF_TAMC
310 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
311 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
312 #endif /* ALLOW_AUTODIFF_TAMC */
313 CALL FIND_RHO(
314 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
315 I theta, salt,
316 O rhoK,
317 I myThid )
318 IF (k.GT.1) THEN
319 #ifdef ALLOW_AUTODIFF_TAMC
320 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
321 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
322 #endif /* ALLOW_AUTODIFF_TAMC */
323 CALL FIND_RHO(
324 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
325 I theta, salt,
326 O rhoKm1,
327 I myThid )
328 ENDIF
329 CALL GRAD_SIGMA(
330 I bi, bj, iMin, iMax, jMin, jMax, k,
331 I rhoK, rhoKm1, rhoK,
332 O sigmaX, sigmaY, sigmaR,
333 I myThid )
334 ENDIF
335
336 C-- Implicit Vertical Diffusion for Convection
337 c ==> should use sigmaR !!!
338 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
339 CALL CALC_IVDC(
340 I bi, bj, iMin, iMax, jMin, jMax, k,
341 I rhoKm1, rhoK,
342 U ConvectCount, KappaRT, KappaRS,
343 I myTime, myIter, myThid)
344 ENDIF
345
346 C-- end of diagnostic k loop (Nr:1)
347 ENDDO
348
349 #ifdef ALLOW_AUTODIFF_TAMC
350 cph avoids recomputation of integrate_for_w
351 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
352 #endif /* ALLOW_AUTODIFF_TAMC */
353
354 #ifdef ALLOW_OBCS
355 C-- Calculate future values on open boundaries
356 IF (useOBCS) THEN
357 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
358 I uVel, vVel, wVel, theta, salt,
359 I myThid )
360 ENDIF
361 #endif /* ALLOW_OBCS */
362
363 C-- Determines forcing terms based on external fields
364 C relaxation terms, etc.
365 CALL EXTERNAL_FORCING_SURF(
366 I bi, bj, iMin, iMax, jMin, jMax,
367 I myThid )
368 #ifdef ALLOW_AUTODIFF_TAMC
369 cph needed for KPP
370 CADJ STORE surfacetendencyU(:,:,bi,bj)
371 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
372 CADJ STORE surfacetendencyV(:,:,bi,bj)
373 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
374 CADJ STORE surfacetendencyS(:,:,bi,bj)
375 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
376 CADJ STORE surfacetendencyT(:,:,bi,bj)
377 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
378 #endif /* ALLOW_AUTODIFF_TAMC */
379
380 #ifdef ALLOW_GMREDI
381
382 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
383 IF (useGMRedi) THEN
384 CALL GMREDI_CALC_TENSOR(
385 I bi, bj, iMin, iMax, jMin, jMax,
386 I sigmaX, sigmaY, sigmaR,
387 I myThid )
388 #ifdef ALLOW_AUTODIFF_TAMC
389 ELSE
390 CALL GMREDI_CALC_TENSOR_DUMMY(
391 I bi, bj, iMin, iMax, jMin, jMax,
392 I sigmaX, sigmaY, sigmaR,
393 I myThid )
394 #endif /* ALLOW_AUTODIFF_TAMC */
395 ENDIF
396
397 #ifdef ALLOW_AUTODIFF_TAMC
398 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
399 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
400 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
401 #endif /* ALLOW_AUTODIFF_TAMC */
402
403 #endif /* ALLOW_GMREDI */
404
405 #ifdef ALLOW_KPP
406 C-- Compute KPP mixing coefficients
407 IF (useKPP) THEN
408 CALL KPP_CALC(
409 I bi, bj, myTime, myThid )
410 #ifdef ALLOW_AUTODIFF_TAMC
411 ELSE
412 CALL KPP_CALC_DUMMY(
413 I bi, bj, myTime, myThid )
414 #endif /* ALLOW_AUTODIFF_TAMC */
415 ENDIF
416
417 #ifdef ALLOW_AUTODIFF_TAMC
418 CADJ STORE KPPghat (:,:,:,bi,bj)
419 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
420 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
421 CADJ & , KPPfrac (:,: ,bi,bj)
422 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
423 #endif /* ALLOW_AUTODIFF_TAMC */
424
425 #endif /* ALLOW_KPP */
426
427 #ifdef ALLOW_AUTODIFF_TAMC
428 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key=ikey, byte=isbyte
429 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key=ikey, byte=isbyte
430 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
431 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
432 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
433 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
434 #ifdef ALLOW_PASSIVE_TRACER
435 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
436 #endif
437 #endif /* ALLOW_AUTODIFF_TAMC */
438
439 #ifdef ALLOW_AIM
440 C AIM - atmospheric intermediate model, physics package code.
441 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
442 IF ( useAIM ) THEN
443 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
444 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, bi, bj, myTime, myThid )
445 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
446 ENDIF
447 #endif /* ALLOW_AIM */
448
449 #ifdef ALLOW_TIMEAVE
450 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
451 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
452 I deltaTclock, bi, bj, myThid)
453 ENDIF
454 #endif /* ALLOW_TIMEAVE */
455
456 #ifndef DISABLE_MULTIDIM_ADVECTION
457 C-- Some advection schemes are better calculated using a multi-dimensional
458 C method in the absence of any other terms and, if used, is done here.
459 C
460 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
461 C The default is to use multi-dimensinal advection for non-linear advection
462 C schemes. However, for the sake of efficiency of the adjoint it is necessary
463 C to be able to exclude this scheme to avoid excessive storage and
464 C recomputation. It *is* differentiable, if you need it.
465 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
466 C disable this section of code.
467 IF (multiDimAdvection) THEN
468 IF (tempStepping .AND.
469 & tempAdvScheme.NE.ENUM_CENTERED_2ND .AND.
470 & tempAdvScheme.NE.ENUM_UPWIND_3RD .AND.
471 & tempAdvScheme.NE.ENUM_CENTERED_4TH ) THEN
472 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
473 U theta,gT,
474 I myTime,myIter,myThid)
475 ENDIF
476 IF (saltStepping .AND.
477 & saltAdvScheme.NE.ENUM_CENTERED_2ND .AND.
478 & saltAdvScheme.NE.ENUM_UPWIND_3RD .AND.
479 & saltAdvScheme.NE.ENUM_CENTERED_4TH ) THEN
480 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
481 U salt,gS,
482 I myTime,myIter,myThid)
483 ENDIF
484 ENDIF
485 C Since passive tracers are configurable separately from T,S we
486 C call the multi-dimensional method for PTRACERS regardless
487 C of whether multiDimAdvection is set or not.
488 #ifdef ALLOW_PTRACERS
489 IF ( usePTRACERS ) THEN
490 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
491 ENDIF
492 #endif /* ALLOW_PTRACERS */
493 #endif /* DISABLE_MULTIDIM_ADVECTION */
494
495 C-- Start of thermodynamics loop
496 DO k=Nr,1,-1
497 #ifdef ALLOW_AUTODIFF_TAMC
498 C? Patrick Is this formula correct?
499 cph Yes, but I rewrote it.
500 cph Also, the KappaR? need the index and subscript k!
501 kkey = (ikey-1)*Nr + k
502 #endif /* ALLOW_AUTODIFF_TAMC */
503
504 C-- km1 Points to level above k (=k-1)
505 C-- kup Cycles through 1,2 to point to layer above
506 C-- kDown Cycles through 2,1 to point to current layer
507
508 km1 = MAX(1,k-1)
509 kup = 1+MOD(k+1,2)
510 kDown= 1+MOD(k,2)
511
512 iMin = 1-OLx
513 iMax = sNx+OLx
514 jMin = 1-OLy
515 jMax = sNy+OLy
516
517 C-- Get temporary terms used by tendency routines
518 CALL CALC_COMMON_FACTORS (
519 I bi,bj,iMin,iMax,jMin,jMax,k,
520 O xA,yA,uTrans,vTrans,rTrans,maskUp,
521 I myThid)
522
523 #ifdef ALLOW_GMREDI
524 C-- Residual transp = Bolus transp + Eulerian transp
525 IF (useGMRedi) THEN
526 CALL GMREDI_CALC_UVFLOW(
527 & uTrans, vTrans, bi, bj, k, myThid)
528 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
529 & rTrans, bi, bj, k, myThid)
530 ENDIF
531 #endif /* ALLOW_GMREDI */
532
533 #ifdef ALLOW_AUTODIFF_TAMC
534 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
535 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
536 #endif /* ALLOW_AUTODIFF_TAMC */
537
538 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
539 C-- Calculate the total vertical diffusivity
540 CALL CALC_DIFFUSIVITY(
541 I bi,bj,iMin,iMax,jMin,jMax,k,
542 I maskUp,
543 O KappaRT,KappaRS,
544 I myThid)
545 #endif
546
547 iMin = 1-OLx+2
548 iMax = sNx+OLx-1
549 jMin = 1-OLy+2
550 jMax = sNy+OLy-1
551
552 C-- Calculate active tracer tendencies (gT,gS,...)
553 C and step forward storing result in gTnm1, gSnm1, etc.
554 IF ( tempStepping ) THEN
555 CALL CALC_GT(
556 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
557 I xA,yA,uTrans,vTrans,rTrans,maskUp,
558 I KappaRT,
559 U fVerT,
560 I myTime,myIter,myThid)
561 CALL TIMESTEP_TRACER(
562 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
563 I theta, gT,
564 I myIter, myThid)
565 ENDIF
566 IF ( saltStepping ) THEN
567 CALL CALC_GS(
568 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
569 I xA,yA,uTrans,vTrans,rTrans,maskUp,
570 I KappaRS,
571 U fVerS,
572 I myTime,myIter,myThid)
573 CALL TIMESTEP_TRACER(
574 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
575 I salt, gS,
576 I myIter, myThid)
577 ENDIF
578 #ifdef ALLOW_PASSIVE_TRACER
579 IF ( tr1Stepping ) THEN
580 CALL CALC_GTR1(
581 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
582 I xA,yA,uTrans,vTrans,rTrans,maskUp,
583 I KappaRT,
584 U fVerTr1,
585 I myTime,myIter,myThid)
586 CALL TIMESTEP_TRACER(
587 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
588 I Tr1, gTr1,
589 I myIter,myThid)
590 ENDIF
591 #endif
592 #ifdef ALLOW_PTRACERS
593 IF ( usePTRACERS ) THEN
594 CALL PTRACERS_INTEGERATE(
595 I bi,bj,k,
596 I xA,yA,uTrans,vTrans,rTrans,maskUp,
597 X KappaRS,
598 I myIter,myTime,myThid)
599 ENDIF
600 #endif /* ALLOW_PTRACERS */
601
602 #ifdef ALLOW_OBCS
603 C-- Apply open boundary conditions
604 IF (useOBCS) THEN
605 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
606 END IF
607 #endif /* ALLOW_OBCS */
608
609 C-- Freeze water
610 IF (allowFreezing) THEN
611 #ifdef ALLOW_AUTODIFF_TAMC
612 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
613 CADJ & , key = kkey, byte = isbyte
614 #endif /* ALLOW_AUTODIFF_TAMC */
615 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
616 END IF
617
618 C-- end of thermodynamic k loop (Nr:1)
619 ENDDO
620
621
622 #ifdef ALLOW_AUTODIFF_TAMC
623 C? Patrick? What about this one?
624 cph Keys iikey and idkey dont seem to be needed
625 cph since storing occurs on different tape for each
626 cph impldiff call anyways.
627 cph Thus, common block comlev1_impl isnt needed either.
628 cph Storing below needed in the case useGMREDI.
629 iikey = (ikey-1)*maximpl
630 #endif /* ALLOW_AUTODIFF_TAMC */
631
632 C-- Implicit diffusion
633 IF (implicitDiffusion) THEN
634
635 IF (tempStepping) THEN
636 #ifdef ALLOW_AUTODIFF_TAMC
637 idkey = iikey + 1
638 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
639 #endif /* ALLOW_AUTODIFF_TAMC */
640 CALL IMPLDIFF(
641 I bi, bj, iMin, iMax, jMin, jMax,
642 I deltaTtracer, KappaRT, recip_HFacC,
643 U gT,
644 I myThid )
645 ENDIF
646
647 IF (saltStepping) THEN
648 #ifdef ALLOW_AUTODIFF_TAMC
649 idkey = iikey + 2
650 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
651 #endif /* ALLOW_AUTODIFF_TAMC */
652 CALL IMPLDIFF(
653 I bi, bj, iMin, iMax, jMin, jMax,
654 I deltaTtracer, KappaRS, recip_HFacC,
655 U gS,
656 I myThid )
657 ENDIF
658
659 #ifdef ALLOW_PASSIVE_TRACER
660 IF (tr1Stepping) THEN
661 #ifdef ALLOW_AUTODIFF_TAMC
662 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
663 #endif /* ALLOW_AUTODIFF_TAMC */
664 CALL IMPLDIFF(
665 I bi, bj, iMin, iMax, jMin, jMax,
666 I deltaTtracer, KappaRT, recip_HFacC,
667 U gTr1,
668 I myThid )
669 ENDIF
670 #endif
671
672 #ifdef ALLOW_PTRACERS
673 C Vertical diffusion (implicit) for passive tracers
674 IF ( usePTRACERS ) THEN
675 CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
676 ENDIF
677 #endif /* ALLOW_PTRACERS */
678
679 #ifdef ALLOW_OBCS
680 C-- Apply open boundary conditions
681 IF (useOBCS) THEN
682 DO K=1,Nr
683 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
684 ENDDO
685 END IF
686 #endif /* ALLOW_OBCS */
687
688 C-- End If implicitDiffusion
689 ENDIF
690
691 Ccs-
692 ENDDO
693 ENDDO
694
695 #ifdef ALLOW_AIM
696 IF ( useAIM ) THEN
697 CALL AIM_AIM2DYN_EXCHANGES( myTime, myThid )
698 ENDIF
699 _EXCH_XYZ_R8(gT,myThid)
700 _EXCH_XYZ_R8(gS,myThid)
701 #else
702 IF (staggerTimeStep.AND.useCubedSphereExchange) THEN
703 _EXCH_XYZ_R8(gT,myThid)
704 _EXCH_XYZ_R8(gS,myThid)
705 ENDIF
706 #endif /* ALLOW_AIM */
707
708 #ifndef DISABLE_DEBUGMODE
709 If (debugMode) THEN
710 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
711 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
712 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
713 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
714 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
715 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
716 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
717 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
718 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
719 #ifdef ALLOW_PTRACERS
720 IF ( usePTRACERS ) THEN
721 CALL PTRACERS_DEBUG(myThid)
722 ENDIF
723 #endif /* ALLOW_PTRACERS */
724 ENDIF
725 #endif
726
727 RETURN
728 END

  ViewVC Help
Powered by ViewVC 1.1.22