/[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.139 - (show annotations) (download)
Thu Dec 1 14:25:41 2011 UTC (12 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63g, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k
Changes since 1.138: +78 -20 lines
compute local recip_hFacC (no bi,bj) to use in vertical implicit solver

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.138 2011/11/07 20:21:05 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6 #ifdef ALLOW_GENERIC_ADVDIFF
7 # include "GAD_OPTIONS.h"
8 #endif
9 #ifdef ALLOW_LONGSTEP
10 #include "LONGSTEP_OPTIONS.h"
11 #endif
12
13 #ifdef ALLOW_AUTODIFF_TAMC
14 # ifdef ALLOW_GMREDI
15 # include "GMREDI_OPTIONS.h"
16 # endif
17 # ifdef ALLOW_KPP
18 # include "KPP_OPTIONS.h"
19 # endif
20 #endif /* ALLOW_AUTODIFF_TAMC */
21
22 CBOP
23 C !ROUTINE: THERMODYNAMICS
24 C !INTERFACE:
25 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
26 C !DESCRIPTION: \bv
27 C *==========================================================*
28 C | SUBROUTINE THERMODYNAMICS
29 C | o Controlling routine for the prognostic part of the
30 C | thermo-dynamics.
31 C *===========================================================
32 C | The algorithm...
33 C |
34 C | "Correction Step"
35 C | =================
36 C | Here we update the horizontal velocities with the surface
37 C | pressure such that the resulting flow is either consistent
38 C | with the free-surface evolution or the rigid-lid:
39 C | U[n] = U* + dt x d/dx P
40 C | V[n] = V* + dt x d/dy P
41 C |
42 C | "Calculation of Gs"
43 C | ===================
44 C | This is where all the accelerations and tendencies (ie.
45 C | physics, parameterizations etc...) are calculated
46 C | rho = rho ( theta[n], salt[n] )
47 C | b = b(rho, theta)
48 C | K31 = K31 ( rho )
49 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
50 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
51 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
52 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
53 C |
54 C | "Time-stepping" or "Prediction"
55 C | ================================
56 C | The models variables are stepped forward with the appropriate
57 C | time-stepping scheme (currently we use Adams-Bashforth II)
58 C | - For momentum, the result is always *only* a "prediction"
59 C | in that the flow may be divergent and will be "corrected"
60 C | later with a surface pressure gradient.
61 C | - Normally for tracers the result is the new field at time
62 C | level [n+1} *BUT* in the case of implicit diffusion the result
63 C | is also *only* a prediction.
64 C | - We denote "predictors" with an asterisk (*).
65 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
66 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
67 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
68 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69 C | With implicit diffusion:
70 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
71 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
72 C | (1 + dt * K * d_zz) theta[n] = theta*
73 C | (1 + dt * K * d_zz) salt[n] = salt*
74 C |
75 C *==========================================================*
76 C \ev
77
78 C !USES:
79 IMPLICIT NONE
80 C == Global variables ===
81 #include "SIZE.h"
82 #include "EEPARAMS.h"
83 #include "PARAMS.h"
84 #include "RESTART.h"
85 #include "DYNVARS.h"
86 #include "GRID.h"
87 #include "SURFACE.h"
88 #ifdef ALLOW_GENERIC_ADVDIFF
89 # include "GAD.h"
90 # include "GAD_SOM_VARS.h"
91 #endif
92 #ifdef ALLOW_PTRACERS
93 #ifndef ALLOW_LONGSTEP
94 # include "PTRACERS_SIZE.h"
95 # include "PTRACERS_PARAMS.h"
96 # include "PTRACERS_FIELDS.h"
97 #endif
98 #endif
99 #ifdef ALLOW_TIMEAVE
100 # include "TIMEAVE_STATV.h"
101 #endif
102
103 #ifdef ALLOW_AUTODIFF_TAMC
104 # include "tamc.h"
105 # include "tamc_keys.h"
106 # include "FFIELDS.h"
107 # include "EOS.h"
108 # ifdef ALLOW_KPP
109 # include "KPP.h"
110 # endif
111 # ifdef ALLOW_GMREDI
112 # include "GMREDI.h"
113 # endif
114 # ifdef ALLOW_EBM
115 # include "EBM.h"
116 # endif
117 # ifdef ALLOW_SALT_PLUME
118 # include "SALT_PLUME.h"
119 # endif
120 #endif /* ALLOW_AUTODIFF_TAMC */
121
122 C !INPUT/OUTPUT PARAMETERS:
123 C == Routine arguments ==
124 C myTime - Current time in simulation
125 C myIter - Current iteration number in simulation
126 C myThid - Thread number for this instance of the routine.
127 _RL myTime
128 INTEGER myIter
129 INTEGER myThid
130
131 #ifdef ALLOW_GENERIC_ADVDIFF
132 C !LOCAL VARIABLES:
133 C == Local variables
134 C xA, yA - Per block temporaries holding face areas
135 C uFld, vFld, wFld - Local copy of velocity field (3 components)
136 C uTrans, vTrans, rTrans - Per block temporaries holding flow transport
137 C o uTrans: Zonal transport
138 C o vTrans: Meridional transport
139 C o rTrans: Vertical transport
140 C rTransKp1 o vertical volume transp. at interface k+1
141 C maskUp o maskUp: land/water mask for W points
142 C fVer[STUV] o fVer: Vertical flux term - note fVer
143 C is "pipelined" in the vertical
144 C so we need an fVer for each
145 C variable.
146 C kappaRT, - Total diffusion in vertical at level k, for T and S
147 C kappaRS (background + spatially varying, isopycnal term).
148 C kappaRTr - Total diffusion in vertical at level k,
149 C for each passive Tracer
150 C kappaRk - Total diffusion in vertical, all levels, 1 tracer
151 C useVariableK = T when vertical diffusion is not constant
152 C iMin, iMax - Ranges and sub-block indices on which calculations
153 C jMin, jMax are applied.
154 C bi, bj
155 C k, kup, - Index for layer above and below. kup and kDown
156 C kDown, km1 are switched with layer to be the appropriate
157 C index into fVerTerm.
158 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
159 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
160 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
163 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
164 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
165 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
166 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
167 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
168 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
169 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
170 _RL kappaRT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
171 _RL kappaRS (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
172 #ifdef ALLOW_PTRACERS
173 #ifndef ALLOW_LONGSTEP
174 _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
175 _RL kappaRTr(1-OLx:sNx+OLx,1-OLy:sNy+OLy,PTRACERS_num)
176 #endif
177 #endif
178 _RL kappaRk (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
179 _RS recip_hFacNew(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
180 INTEGER iMin, iMax
181 INTEGER jMin, jMax
182 INTEGER bi, bj
183 INTEGER i, j
184 INTEGER k, km1, kup, kDown
185 #ifdef ALLOW_ADAMSBASHFORTH_3
186 INTEGER iterNb, m1, m2
187 _RL tmpFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
188 #endif
189 #ifdef ALLOW_TIMEAVE
190 LOGICAL useVariableK
191 #endif
192 #ifdef ALLOW_PTRACERS
193 #ifndef ALLOW_LONGSTEP
194 INTEGER iTracer, ip
195 #endif
196 #endif
197
198 CEOP
199
200 #ifdef ALLOW_DEBUG
201 IF (debugMode) CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
202 #endif
203
204 #ifdef ALLOW_AUTODIFF_TAMC
205 C-- dummy statement to end declaration part
206 ikey = 1
207 itdkey = 1
208 #endif /* ALLOW_AUTODIFF_TAMC */
209
210 #ifdef ALLOW_AUTODIFF_TAMC
211 C-- HPF directive to help TAMC
212 CHPF$ INDEPENDENT
213 #endif /* ALLOW_AUTODIFF_TAMC */
214
215 C-- Compute correction at the surface for Lin Free Surf.
216 #ifdef ALLOW_AUTODIFF_TAMC
217 TsurfCor = 0. _d 0
218 SsurfCor = 0. _d 0
219 #endif
220 IF (linFSConserveTr) THEN
221 #ifdef ALLOW_AUTODIFF_TAMC
222 CADJ STORE theta,salt,wvel = comlev1, key = ikey_dynamics, byte=isbyte
223 #endif
224 CALL CALC_WSURF_TR(theta,salt,wVel,
225 & myTime,myIter,myThid)
226 ENDIF
227
228 DO bj=myByLo(myThid),myByHi(myThid)
229
230 #ifdef ALLOW_AUTODIFF_TAMC
231 C-- HPF directive to help TAMC
232 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
233 CHPF$& ,utrans,vtrans,xA,yA
234 CHPF$& ,kappaRT,kappaRS
235 CHPF$& )
236 # ifdef ALLOW_PTRACERS
237 # ifndef ALLOW_LONGSTEP
238 CHPF$ INDEPENDENT, NEW (fVerP,kappaRTr)
239 # endif
240 # endif
241 #endif /* ALLOW_AUTODIFF_TAMC */
242
243 DO bi=myBxLo(myThid),myBxHi(myThid)
244
245 #ifdef ALLOW_AUTODIFF_TAMC
246 act1 = bi - myBxLo(myThid)
247 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
248 act2 = bj - myByLo(myThid)
249 max2 = myByHi(myThid) - myByLo(myThid) + 1
250 act3 = myThid - 1
251 max3 = nTx*nTy
252 act4 = ikey_dynamics - 1
253 itdkey = (act1 + 1) + act2*max1
254 & + act3*max1*max2
255 & + act4*max1*max2*max3
256 #endif /* ALLOW_AUTODIFF_TAMC */
257
258 C-- Set up work arrays with valid (i.e. not NaN) values
259 C These inital values do not alter the numerical results. They
260 C just ensure that all memory references are to valid floating
261 C point numbers. This prevents spurious hardware signals due to
262 C uninitialised but inert locations.
263
264 DO j=1-OLy,sNy+OLy
265 DO i=1-OLx,sNx+OLx
266 xA(i,j) = 0. _d 0
267 yA(i,j) = 0. _d 0
268 uTrans(i,j) = 0. _d 0
269 vTrans(i,j) = 0. _d 0
270 rTrans (i,j) = 0. _d 0
271 rTransKp1(i,j) = 0. _d 0
272 fVerT (i,j,1) = 0. _d 0
273 fVerT (i,j,2) = 0. _d 0
274 fVerS (i,j,1) = 0. _d 0
275 fVerS (i,j,2) = 0. _d 0
276 kappaRT(i,j) = 0. _d 0
277 kappaRS(i,j) = 0. _d 0
278 ENDDO
279 ENDDO
280
281 DO k=1,Nr
282 DO j=1-OLy,sNy+OLy
283 DO i=1-OLx,sNx+OLx
284 C This is currently also used by IVDC and Diagnostics
285 kappaRk(i,j,k) = 0. _d 0
286 C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
287 gT(i,j,k,bi,bj) = 0. _d 0
288 gS(i,j,k,bi,bj) = 0. _d 0
289 ENDDO
290 ENDDO
291 ENDDO
292
293 #ifdef ALLOW_PTRACERS
294 #ifndef ALLOW_LONGSTEP
295 IF ( usePTRACERS ) THEN
296 DO ip=1,PTRACERS_num
297 DO j=1-OLy,sNy+OLy
298 DO i=1-OLx,sNx+OLx
299 fVerP (i,j,1,ip) = 0. _d 0
300 fVerP (i,j,2,ip) = 0. _d 0
301 kappaRTr(i,j,ip) = 0. _d 0
302 ENDDO
303 ENDDO
304 ENDDO
305 C- set tracer tendency to zero:
306 DO iTracer=1,PTRACERS_num
307 DO k=1,Nr
308 DO j=1-OLy,sNy+OLy
309 DO i=1-OLx,sNx+OLx
310 gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
311 ENDDO
312 ENDDO
313 ENDDO
314 ENDDO
315 ENDIF
316 #endif
317 #endif
318
319 #ifdef ALLOW_ADAMSBASHFORTH_3
320 C- Apply AB on T,S :
321 iterNb = myIter
322 IF (staggerTimeStep) iterNb = myIter - 1
323 m1 = 1 + MOD(iterNb+1,2)
324 m2 = 1 + MOD( iterNb ,2)
325 C compute T^n+1/2 (stored in gtNm) extrapolating T forward in time
326 IF ( AdamsBashforth_T ) CALL ADAMS_BASHFORTH3(
327 I bi, bj, 0,
328 U theta, gtNm, tmpFld,
329 I tempStartAB, iterNb, myThid )
330 C compute S^n+1/2 (stored in gsNm) extrapolating S forward in time
331 IF ( AdamsBashforth_S ) CALL ADAMS_BASHFORTH3(
332 I bi, bj, 0,
333 U salt, gsNm, tmpFld,
334 I saltStartAB, iterNb, myThid )
335 #endif /* ALLOW_ADAMSBASHFORTH_3 */
336
337 c iMin = 1-OLx
338 c iMax = sNx+OLx
339 c jMin = 1-OLy
340 c jMax = sNy+OLy
341
342 #ifdef ALLOW_AUTODIFF_TAMC
343 cph avoids recomputation of integrate_for_w
344 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
345 #endif /* ALLOW_AUTODIFF_TAMC */
346
347 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
348 C-- MOST of THERMODYNAMICS will be disabled
349 #ifndef SINGLE_LAYER_MODE
350
351 #ifdef ALLOW_AUTODIFF_TAMC
352 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
353 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
354 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
355 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
356 # if ((defined ALLOW_DEPTH_CONTROL) || (defined NONLIN_FRSURF))
357 # ifndef ALLOW_ADAMSBASHFORTH_3
358 CADJ STORE gtnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
359 CADJ STORE gsnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
360 # endif
361 # endif
362 #endif /* ALLOW_AUTODIFF_TAMC */
363
364 #ifndef DISABLE_MULTIDIM_ADVECTION
365 C-- Some advection schemes are better calculated using a multi-dimensional
366 C method in the absence of any other terms and, if used, is done here.
367 C
368 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
369 C The default is to use multi-dimensinal advection for non-linear advection
370 C schemes. However, for the sake of efficiency of the adjoint it is necessary
371 C to be able to exclude this scheme to avoid excessive storage and
372 C recomputation. It *is* differentiable, if you need it.
373 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
374 C disable this section of code.
375 #ifdef GAD_ALLOW_TS_SOM_ADV
376 IF ( tempSOM_Advection ) THEN
377 #ifdef ALLOW_DEBUG
378 IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
379 #endif
380 CALL GAD_SOM_ADVECT(
381 I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
382 I GAD_TEMPERATURE, dTtracerLev,
383 I uVel, vVel, wVel, theta,
384 U som_T,
385 O gT,
386 I bi,bj,myTime,myIter,myThid)
387 ELSEIF (tempMultiDimAdvec) THEN
388 #else /* GAD_ALLOW_TS_SOM_ADV */
389 IF (tempMultiDimAdvec) THEN
390 #endif /* GAD_ALLOW_TS_SOM_ADV */
391 #ifdef ALLOW_DEBUG
392 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
393 #endif
394 CALL GAD_ADVECTION(
395 I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
396 I GAD_TEMPERATURE, dTtracerLev,
397 I uVel, vVel, wVel, theta,
398 O gT,
399 I bi,bj,myTime,myIter,myThid)
400 ENDIF
401 #ifdef GAD_ALLOW_TS_SOM_ADV
402 IF ( saltSOM_Advection ) THEN
403 #ifdef ALLOW_DEBUG
404 IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
405 #endif
406 CALL GAD_SOM_ADVECT(
407 I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
408 I GAD_SALINITY, dTtracerLev,
409 I uVel, vVel, wVel, salt,
410 U som_S,
411 O gS,
412 I bi,bj,myTime,myIter,myThid)
413 ELSEIF (saltMultiDimAdvec) THEN
414 #else /* GAD_ALLOW_TS_SOM_ADV */
415 IF (saltMultiDimAdvec) THEN
416 #endif /* GAD_ALLOW_TS_SOM_ADV */
417 #ifdef ALLOW_DEBUG
418 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
419 #endif
420 CALL GAD_ADVECTION(
421 I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
422 I GAD_SALINITY, dTtracerLev,
423 I uVel, vVel, wVel, salt,
424 O gS,
425 I bi,bj,myTime,myIter,myThid)
426 ENDIF
427
428 C Since passive tracers are configurable separately from T,S we
429 C call the multi-dimensional method for PTRACERS regardless
430 C of whether multiDimAdvection is set or not.
431 #ifdef ALLOW_PTRACERS
432 #ifndef ALLOW_LONGSTEP
433 IF ( usePTRACERS ) THEN
434 #ifdef ALLOW_DEBUG
435 IF (debugMode) CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
436 #endif
437 CALL PTRACERS_ADVECTION( bi,bj,myTime,myIter,myThid )
438 ENDIF
439 #endif /* ALLOW_LONGSTEP */
440 #endif /* ALLOW_PTRACERS */
441 #endif /* DISABLE_MULTIDIM_ADVECTION */
442
443 #ifdef ALLOW_DEBUG
444 IF (debugMode)
445 & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
446 #endif
447
448 #ifdef ALLOW_AUTODIFF_TAMC
449 # ifdef ALLOW_SALT_PLUME
450 CADJ STORE saltPlumeFlux(:,:,bi,bj) =
451 CADJ & comlev1_bibj, key=itdkey,kind = isbyte
452 CADJ STORE saltPlumeDepth(:,:,bi,bj) =
453 CADJ & comlev1_bibj, key=itdkey,kind = isbyte
454 # endif
455 #endif /* ALLOW_AUTODIFF_TAMC */
456
457 C-- Start of thermodynamics loop
458 DO k=Nr,1,-1
459 #ifdef ALLOW_AUTODIFF_TAMC
460 C? Patrick Is this formula correct?
461 cph Yes, but I rewrote it.
462 cph Also, the kappaR? need the index and subscript k!
463 kkey = (itdkey-1)*Nr + k
464 #endif /* ALLOW_AUTODIFF_TAMC */
465
466 C-- km1 Points to level above k (=k-1)
467 C-- kup Cycles through 1,2 to point to layer above
468 C-- kDown Cycles through 2,1 to point to current layer
469
470 km1 = MAX(1,k-1)
471 kup = 1+MOD(k+1,2)
472 kDown= 1+MOD(k,2)
473
474 iMin = 1-OLx
475 iMax = sNx+OLx
476 jMin = 1-OLy
477 jMax = sNy+OLy
478
479 IF (k.EQ.Nr) THEN
480 DO j=1-OLy,sNy+OLy
481 DO i=1-OLx,sNx+OLx
482 rTransKp1(i,j) = 0. _d 0
483 ENDDO
484 ENDDO
485 ELSE
486 DO j=1-OLy,sNy+OLy
487 DO i=1-OLx,sNx+OLx
488 rTransKp1(i,j) = rTrans(i,j)
489 ENDDO
490 ENDDO
491 ENDIF
492 #ifdef ALLOW_AUTODIFF_TAMC
493 CADJ STORE rTransKp1(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
494 #endif
495
496 C-- Get temporary terms used by tendency routines :
497 C- Calculate horizontal "volume transport" through tracer cell face
498 C anelastic: uTrans,vTrans are scaled by rhoFacC (~ mass transport)
499 CALL CALC_COMMON_FACTORS (
500 I uVel, vVel,
501 O uFld, vFld, uTrans, vTrans, xA, yA,
502 I k,bi,bj, myThid )
503
504 C- Calculate vertical "volume transport" through tracer cell face
505 IF (k.EQ.1) THEN
506 C- Surface interface :
507 DO j=1-OLy,sNy+OLy
508 DO i=1-OLx,sNx+OLx
509 wFld(i,j) = 0. _d 0
510 maskUp(i,j) = 0. _d 0
511 rTrans(i,j) = 0. _d 0
512 ENDDO
513 ENDDO
514 ELSE
515 C- Interior interface :
516 C anelastic: rTrans is scaled by rhoFacF (~ mass transport)
517 DO j=1-OLy,sNy+OLy
518 DO i=1-OLx,sNx+OLx
519 wFld(i,j) = wVel(i,j,k,bi,bj)
520 maskUp(i,j) = maskC(i,j,k-1,bi,bj)*maskC(i,j,k,bi,bj)
521 rTrans(i,j) = wFld(i,j)*rA(i,j,bi,bj)*maskUp(i,j)
522 & *deepFac2F(k)*rhoFacF(k)
523 ENDDO
524 ENDDO
525 ENDIF
526
527 #ifdef ALLOW_GMREDI
528 C-- Residual transp = Bolus transp + Eulerian transp
529 IF (useGMRedi) THEN
530 CALL GMREDI_CALC_UVFLOW(
531 U uFld, vFld, uTrans, vTrans,
532 I k, bi, bj, myThid )
533 IF (K.GE.2) THEN
534 CALL GMREDI_CALC_WFLOW(
535 U wFld, rTrans,
536 I k, bi, bj, myThid )
537 ENDIF
538 ENDIF
539 # ifdef ALLOW_AUTODIFF_TAMC
540 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
541 CADJ STORE wfld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
542 # ifdef GM_BOLUS_ADVEC
543 CADJ STORE ufld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
544 CADJ STORE vfld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
545 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
546 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
547 # endif
548 # endif /* ALLOW_AUTODIFF_TAMC */
549 #endif /* ALLOW_GMREDI */
550
551 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
552 C-- Calculate the total vertical diffusivity
553 IF ( .NOT.implicitDiffusion ) THEN
554 CALL CALC_DIFFUSIVITY(
555 I bi,bj,iMin,iMax,jMin,jMax,k,
556 I maskUp,
557 O kappaRT,kappaRS,
558 I myThid)
559 ENDIF
560 # ifdef ALLOW_AUTODIFF_TAMC
561 CADJ STORE kappaRT(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
562 CADJ STORE kappaRS(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
563 # endif /* ALLOW_AUTODIFF_TAMC */
564 #endif
565
566 iMin = 1-OLx+2
567 iMax = sNx+OLx-1
568 jMin = 1-OLy+2
569 jMax = sNy+OLy-1
570
571 C-- Calculate active tracer tendencies (gT,gS,...)
572 C and step forward storing result in gT, gS, etc.
573 C--
574 # ifdef ALLOW_AUTODIFF_TAMC
575 # if ((defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)) && (defined ALLOW_GMREDI)
576 # ifdef GM_NON_UNITY_DIAGONAL
577 CADJ STORE kux(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
578 CADJ STORE kvy(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
579 # endif
580 # ifdef GM_EXTRA_DIAGONAL
581 CADJ STORE kuz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
582 CADJ STORE kvz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
583 # endif
584 # endif
585 # endif /* ALLOW_AUTODIFF_TAMC */
586 C
587 #ifdef ALLOW_AUTODIFF_TAMC
588 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
589 cph-test
590 CADJ STORE uFld(:,:), vFld(:,:), wFld(:,:)
591 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
592 CADJ STORE uTrans(:,:), vTrans(:,:)
593 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
594 CADJ STORE xA(:,:), yA(:,:)
595 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
596 # ifdef ALLOW_ADAMSBASHFORTH_3
597 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
598 CADJ STORE gS(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
599 CADJ STORE gSnm(:,:,k,bi,bj,1)= comlev1_bibj_k, key=kkey, byte=isbyte
600 CADJ STORE gSnm(:,:,k,bi,bj,2)= comlev1_bibj_k, key=kkey, byte=isbyte
601 CADJ STORE gTnm(:,:,k,bi,bj,1)= comlev1_bibj_k, key=kkey, byte=isbyte
602 CADJ STORE gTnm(:,:,k,bi,bj,2)= comlev1_bibj_k, key=kkey, byte=isbyte
603 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
604 CADJ STORE salt(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
605 CADJ STORE fvert(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
606 CADJ STORE fvers(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
607 # endif /* ALLOW_ADAMSBASHFORTH_3 */
608 # endif /* NONLIN_FRSURF */
609 #endif /* ALLOW_AUTODIFF_TAMC */
610 C
611 IF ( tempStepping ) THEN
612 #ifdef ALLOW_AUTODIFF_TAMC
613 # ifndef ALLOW_ADAMSBASHFORTH_3
614 CADJ STORE gTnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
615 # else
616 # ifndef NONLIN_FRSURF
617 CADJ STORE gTnm(:,:,k,bi,bj,1)= comlev1_bibj_k, key=kkey, byte=isbyte
618 CADJ STORE gTnm(:,:,k,bi,bj,2)= comlev1_bibj_k, key=kkey, byte=isbyte
619 # endif /* ndef NONLIN_FRSURF */
620 # endif /* ndef ALLOW_ADAMSBASHFORTH_3 */
621 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
622 CADJ STORE gt(:,:,:,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
623 CADJ STORE fvert(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
624 # endif
625 #endif /* ALLOW_AUTODIFF_TAMC */
626 CALL CALC_GT(
627 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
628 I xA, yA, maskUp, uFld, vFld, wFld,
629 I uTrans, vTrans, rTrans, rTransKp1,
630 I kappaRT,
631 U fVerT,
632 I myTime,myIter,myThid)
633 #ifdef ALLOW_ADAMSBASHFORTH_3
634 IF ( AdamsBashforth_T ) THEN
635 CALL TIMESTEP_TRACER(
636 I bi, bj, k, dTtracerLev(k),
637 I gtNm(1-OLx,1-OLy,1,1,1,m2),
638 U gT,
639 I myIter, myThid )
640 ELSE
641 #endif
642 CALL TIMESTEP_TRACER(
643 I bi, bj, k, dTtracerLev(k),
644 I theta,
645 U gT,
646 I myIter, myThid )
647 #ifdef ALLOW_ADAMSBASHFORTH_3
648 ENDIF
649 #endif
650 ENDIF
651
652 #ifdef ALLOW_AUTODIFF_TAMC
653 # if (defined NONLIN_FRSURF) && (defined ALLOW_ADAMSBASHFORTH_3)
654 CADJ STORE gTnm(:,:,k,bi,bj,1)= comlev1_bibj_k, key=kkey, byte=isbyte
655 CADJ STORE gTnm(:,:,k,bi,bj,2)= comlev1_bibj_k, key=kkey, byte=isbyte
656 CADJ STORE fvert(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
657 # endif
658 #endif
659
660 IF ( saltStepping ) THEN
661 #ifdef ALLOW_AUTODIFF_TAMC
662 # ifndef ALLOW_ADAMSBASHFORTH_3
663 CADJ STORE gSnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
664 # else
665 # ifndef NONLIN_FRSURF
666 CADJ STORE gSnm(:,:,k,bi,bj,1)= comlev1_bibj_k, key=kkey, byte=isbyte
667 CADJ STORE gSnm(:,:,k,bi,bj,2)= comlev1_bibj_k, key=kkey, byte=isbyte
668 # endif /* ndef NONLIN_FRSURF */
669 # endif /* ndef ALLOW_ADAMSBASHFORTH_3 */
670 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
671 CADJ STORE gs(:,:,:,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
672 CADJ STORE fvers(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
673 # endif
674 #endif /* ALLOW_AUTODIFF_TAMC */
675
676 CALL CALC_GS(
677 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
678 I xA, yA, maskUp, uFld, vFld, wFld,
679 I uTrans, vTrans, rTrans, rTransKp1,
680 I kappaRS,
681 U fVerS,
682 I myTime,myIter,myThid)
683 #ifdef ALLOW_ADAMSBASHFORTH_3
684 IF ( AdamsBashforth_S ) THEN
685 CALL TIMESTEP_TRACER(
686 I bi, bj, k, dTtracerLev(k),
687 I gsNm(1-OLx,1-OLy,1,1,1,m2),
688 U gS,
689 I myIter, myThid )
690 ELSE
691 #endif
692 CALL TIMESTEP_TRACER(
693 I bi, bj, k, dTtracerLev(k),
694 I salt,
695 U gS,
696 I myIter, myThid )
697 #ifdef ALLOW_ADAMSBASHFORTH_3
698 ENDIF
699 #endif
700 ENDIF
701
702 #ifdef ALLOW_PTRACERS
703 #ifndef ALLOW_LONGSTEP
704 IF ( usePTRACERS ) THEN
705 IF ( .NOT.implicitDiffusion ) THEN
706 CALL PTRACERS_CALC_DIFF(
707 I bi,bj,iMin,iMax,jMin,jMax,k,
708 I maskUp,
709 O kappaRTr,
710 I myThid)
711 ENDIF
712 # ifdef ALLOW_AUTODIFF_TAMC
713 CADJ STORE kappaRTr(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
714 # endif /* ALLOW_AUTODIFF_TAMC */
715 CALL PTRACERS_INTEGRATE(
716 I bi,bj,k,
717 I xA, yA, maskUp, uFld, vFld, wFld,
718 I uTrans, vTrans, rTrans, rTransKp1,
719 I kappaRTr,
720 U fVerP,
721 I myTime,myIter,myThid)
722 ENDIF
723 #endif /*ALLOW_LONGSTEP */
724 #endif /* ALLOW_PTRACERS */
725
726 C-- Freeze water
727 C this bit of code is left here for backward compatibility.
728 C freezing at surface level has been moved to DO_OCEANIC_PHYS
729 IF ( useOldFreezing .AND. .NOT. useSEAICE
730 & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
731 #ifdef ALLOW_AUTODIFF_TAMC
732 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
733 CADJ & , key = kkey, byte = isbyte
734 #endif /* ALLOW_AUTODIFF_TAMC */
735 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
736 ENDIF
737
738 C-- end of thermodynamic k loop (Nr:1)
739 ENDDO
740
741 C-- Compute new reciprocal hFac for implicit calculation
742 #ifdef NONLIN_FRSURF
743 IF ( nonlinFreeSurf.GT.0 ) THEN
744 IF ( select_rStar.GT.0 ) THEN
745 # ifndef DISABLE_RSTAR_CODE
746 DO k=1,Nr
747 DO j=1-OLy,sNy+OLy
748 DO i=1-OLx,sNx+OLx
749 recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
750 & / rStarExpC(i,j,bi,bj)
751 ENDDO
752 ENDDO
753 ENDDO
754 # endif /* DISABLE_RSTAR_CODE */
755 ELSEIF ( selectSigmaCoord.NE.0 ) THEN
756 # ifndef DISABLE_SIGMA_CODE
757 DO k=1,Nr
758 DO j=1-OLy,sNy+OLy
759 DO i=1-OLx,sNx+OLx
760 recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
761 & /( 1. _d 0 + dEtaHdt(i,j,bi,bj)*deltaTfreesurf
762 & *dBHybSigF(k)*recip_drF(k)
763 & *recip_hFacC(i,j,k,bi,bj)
764 & )
765 ENDDO
766 ENDDO
767 ENDDO
768 # endif /* DISABLE_RSTAR_CODE */
769 ELSE
770 DO k=1,Nr
771 DO j=1-OLy,sNy+OLy
772 DO i=1-OLx,sNx+OLx
773 IF ( k.EQ.kSurfC(i,j,bi,bj) ) THEN
774 recip_hFacNew(i,j,k) = 1. _d 0 / hFac_surfC(i,j,bi,bj)
775 ELSE
776 recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
777 ENDIF
778 ENDDO
779 ENDDO
780 ENDDO
781 ENDIF
782 ELSE
783 #endif /* NONLIN_FRSURF */
784 DO k=1,Nr
785 DO j=1-OLy,sNy+OLy
786 DO i=1-OLx,sNx+OLx
787 recip_hFacNew(i,j,k) = _recip_hFacC(i,j,k,bi,bj)
788 ENDDO
789 ENDDO
790 ENDDO
791 #ifdef NONLIN_FRSURF
792 ENDIF
793 #endif /* NONLIN_FRSURF */
794
795 #ifdef ALLOW_DOWN_SLOPE
796 IF ( tempStepping .AND. useDOWN_SLOPE ) THEN
797 IF ( usingPCoords ) THEN
798 CALL DWNSLP_APPLY(
799 I GAD_TEMPERATURE, bi, bj, kSurfC,
800 I recip_drF, recip_hFacC, recip_rA,
801 I dTtracerLev,
802 I theta,
803 U gT,
804 I myTime, myIter, myThid )
805 ELSE
806 CALL DWNSLP_APPLY(
807 I GAD_TEMPERATURE, bi, bj, kLowC,
808 I recip_drF, recip_hFacC, recip_rA,
809 I dTtracerLev,
810 I theta,
811 U gT,
812 I myTime, myIter, myThid )
813 ENDIF
814 ENDIF
815 IF ( saltStepping .AND. useDOWN_SLOPE ) THEN
816 IF ( usingPCoords ) THEN
817 CALL DWNSLP_APPLY(
818 I GAD_SALINITY, bi, bj, kSurfC,
819 I recip_drF, recip_hFacC, recip_rA,
820 I dTtracerLev,
821 I salt,
822 U gS,
823 I myTime, myIter, myThid )
824 ELSE
825 CALL DWNSLP_APPLY(
826 I GAD_SALINITY, bi, bj, kLowC,
827 I recip_drF, recip_hFacC, recip_rA,
828 I dTtracerLev,
829 I salt,
830 U gS,
831 I myTime, myIter, myThid )
832 ENDIF
833 ENDIF
834 #ifdef ALLOW_PTRACERS
835 #ifndef ALLOW_LONGSTEP
836 IF ( usePTRACERS .AND. useDOWN_SLOPE ) THEN
837 CALL PTRACERS_DWNSLP_APPLY(
838 I bi, bj, myTime, myIter, myThid )
839 ENDIF
840 #endif /*ALLOW_LONGSTEP */
841 #endif /* ALLOW_PTRACERS */
842 #endif /* ALLOW_DOWN_SLOPE */
843
844 C All explicit advection/diffusion/sources should now be
845 C done. The updated tracer field is in gPtr. Accumalate
846 C explicit tendency and also reset gPtr to initial tracer
847 C field for implicit matrix calculation
848
849 #ifdef ALLOW_MATRIX
850 IF (useMATRIX)
851 & CALL MATRIX_STORE_TENDENCY_EXP(bi,bj, myTime,myIter,myThid)
852 #endif
853
854 iMin = 1
855 iMax = sNx
856 jMin = 1
857 jMax = sNy
858
859 C-- Implicit vertical advection & diffusion
860 #ifndef ALLOW_DEPTH_CONTROL
861 IF ( tempStepping .AND. implicitDiffusion ) THEN
862 CALL CALC_3D_DIFFUSIVITY(
863 I bi,bj,iMin,iMax,jMin,jMax,
864 I GAD_TEMPERATURE, useGMredi, useKPP,
865 O kappaRk,
866 I myThid)
867 ENDIF
868 #ifdef INCLUDE_IMPLVERTADV_CODE
869 IF ( tempImplVertAdv ) THEN
870 #ifdef ALLOW_AUTODIFF_TAMC
871 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
872 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
873 CADJ STORE wvel(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
874 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
875 #endif /* ALLOW_AUTODIFF_TAMC */
876 CALL GAD_IMPLICIT_R(
877 I tempImplVertAdv, tempVertAdvScheme, GAD_TEMPERATURE,
878 I dTtracerLev,
879 I kappaRk, recip_hFacNew, wVel, theta,
880 U gT,
881 I bi, bj, myTime, myIter, myThid )
882 ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN
883 #else /* INCLUDE_IMPLVERTADV_CODE */
884 IF ( tempStepping .AND. implicitDiffusion ) THEN
885 #endif /* INCLUDE_IMPLVERTADV_CODE */
886 #ifdef ALLOW_AUTODIFF_TAMC
887 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
888 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
889 #endif /* ALLOW_AUTODIFF_TAMC */
890 CALL IMPLDIFF(
891 I bi, bj, iMin, iMax, jMin, jMax,
892 I GAD_TEMPERATURE, kappaRk, recip_hFacNew,
893 U gT,
894 I myThid )
895 ENDIF
896
897 #ifdef ALLOW_TIMEAVE
898 useVariableK = useKPP .OR. usePP81 .OR. useMY82 .OR. useGGL90
899 & .OR. useGMredi .OR. ivdc_kappa.NE.0.
900 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
901 IF (implicitDiffusion) THEN
902 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRk,
903 I Nr, 3, deltaTclock, bi, bj, myThid)
904 c ELSE
905 c CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
906 c I Nr, 3, deltaTclock, bi, bj, myThid)
907 ENDIF
908 ENDIF
909 #endif /* ALLOW_TIMEAVE */
910
911 IF ( saltStepping .AND. implicitDiffusion ) THEN
912 CALL CALC_3D_DIFFUSIVITY(
913 I bi,bj,iMin,iMax,jMin,jMax,
914 I GAD_SALINITY, useGMredi, useKPP,
915 O kappaRk,
916 I myThid)
917 ENDIF
918
919 #ifdef INCLUDE_IMPLVERTADV_CODE
920 IF ( saltImplVertAdv ) THEN
921 #ifdef ALLOW_AUTODIFF_TAMC
922 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
923 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
924 CADJ STORE wvel(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
925 CADJ STORE salt(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
926 #endif /* ALLOW_AUTODIFF_TAMC */
927 CALL GAD_IMPLICIT_R(
928 I saltImplVertAdv, saltVertAdvScheme, GAD_SALINITY,
929 I dTtracerLev,
930 I kappaRk, recip_hFacNew, wVel, salt,
931 U gS,
932 I bi, bj, myTime, myIter, myThid )
933 ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN
934 #else /* INCLUDE_IMPLVERTADV_CODE */
935 IF ( saltStepping .AND. implicitDiffusion ) THEN
936 #endif /* INCLUDE_IMPLVERTADV_CODE */
937 #ifdef ALLOW_AUTODIFF_TAMC
938 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
939 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
940 #endif /* ALLOW_AUTODIFF_TAMC */
941 CALL IMPLDIFF(
942 I bi, bj, iMin, iMax, jMin, jMax,
943 I GAD_SALINITY, kappaRk, recip_hFacNew,
944 U gS,
945 I myThid )
946 ENDIF
947 #endif /* ndef ALLOW_DEPTH_CONTROL */
948
949 #ifdef ALLOW_PTRACERS
950 #ifndef ALLOW_LONGSTEP
951 IF ( usePTRACERS ) THEN
952 C-- Vertical advection/diffusion (implicit) for passive tracers
953 C Also apply open boundary conditions for each passive tracer
954 CALL PTRACERS_IMPLICIT(
955 U kappaRk,
956 I recip_hFacNew,
957 I bi, bj, myTime, myIter, myThid )
958 ENDIF
959 #endif /* ALLOW_LONGSTEP */
960 #endif /* ALLOW_PTRACERS */
961
962 #ifdef ALLOW_OBCS
963 C-- Apply open boundary conditions
964 IF ( useOBCS ) THEN
965 CALL OBCS_APPLY_TS( bi, bj, 0, gT, gS, myThid )
966 ENDIF
967 #endif /* ALLOW_OBCS */
968
969 #endif /* SINGLE_LAYER_MODE */
970
971 C-- end bi,bj loops.
972 ENDDO
973 ENDDO
974
975 #ifdef ALLOW_DEBUG
976 IF ( debugLevel.GE.debLevD ) THEN
977 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
978 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
979 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
980 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
981 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
982 CALL DEBUG_STATS_RL(Nr,gT,'Gt (THERMODYNAMICS)',myThid)
983 CALL DEBUG_STATS_RL(Nr,gS,'Gs (THERMODYNAMICS)',myThid)
984 #ifndef ALLOW_ADAMSBASHFORTH_3
985 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
986 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
987 #endif
988 #ifdef ALLOW_PTRACERS
989 #ifndef ALLOW_LONGSTEP
990 IF ( usePTRACERS ) THEN
991 CALL PTRACERS_DEBUG(myThid)
992 ENDIF
993 #endif /* ALLOW_LONGSTEP */
994 #endif /* ALLOW_PTRACERS */
995 ENDIF
996 #endif /* ALLOW_DEBUG */
997
998 #ifdef ALLOW_DEBUG
999 IF (debugMode) CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
1000 #endif
1001
1002 #endif /* ALLOW_GENERIC_ADVDIFF */
1003
1004 RETURN
1005 END

  ViewVC Help
Powered by ViewVC 1.1.22