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

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

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


Revision 1.57 - (hide annotations) (download)
Thu Nov 6 22:01:43 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52, checkpoint52a_pre, ecco_c52_e35, checkpoint51u_post
Changes since 1.56: +8 -1 lines
o merging from ecco-branch
o minor CPP options update

1 heimbach 1.57 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.56 2003/11/04 18:40:58 edhill Exp $
2 adcroft 1.1 C $Name: $
3    
4 edhill 1.51 #include "PACKAGES_CONFIG.h"
5 adcroft 1.1 #include "CPP_OPTIONS.h"
6 edhill 1.51
7 jmc 1.21 #ifdef ALLOW_AUTODIFF_TAMC
8     # ifdef ALLOW_GMREDI
9     # include "GMREDI_OPTIONS.h"
10     # endif
11     # ifdef ALLOW_KPP
12     # include "KPP_OPTIONS.h"
13 heimbach 1.42 # endif
14 dimitri 1.48 #ifdef ALLOW_PTRACERS
15     # include "PTRACERS_OPTIONS.h"
16     #endif
17 jmc 1.21 #endif /* ALLOW_AUTODIFF_TAMC */
18 adcroft 1.1
19 cnh 1.9 CBOP
20     C !ROUTINE: THERMODYNAMICS
21     C !INTERFACE:
22 adcroft 1.1 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
23 cnh 1.9 C !DESCRIPTION: \bv
24     C *==========================================================*
25     C | SUBROUTINE THERMODYNAMICS
26     C | o Controlling routine for the prognostic part of the
27     C | thermo-dynamics.
28     C *===========================================================
29     C | The algorithm...
30     C |
31     C | "Correction Step"
32     C | =================
33     C | Here we update the horizontal velocities with the surface
34     C | pressure such that the resulting flow is either consistent
35     C | with the free-surface evolution or the rigid-lid:
36     C | U[n] = U* + dt x d/dx P
37     C | V[n] = V* + dt x d/dy P
38     C |
39     C | "Calculation of Gs"
40     C | ===================
41     C | This is where all the accelerations and tendencies (ie.
42     C | physics, parameterizations etc...) are calculated
43     C | rho = rho ( theta[n], salt[n] )
44     C | b = b(rho, theta)
45     C | K31 = K31 ( rho )
46     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
47     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
48     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
49     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
50     C |
51     C | "Time-stepping" or "Prediction"
52     C | ================================
53     C | The models variables are stepped forward with the appropriate
54     C | time-stepping scheme (currently we use Adams-Bashforth II)
55     C | - For momentum, the result is always *only* a "prediction"
56     C | in that the flow may be divergent and will be "corrected"
57     C | later with a surface pressure gradient.
58     C | - Normally for tracers the result is the new field at time
59     C | level [n+1} *BUT* in the case of implicit diffusion the result
60     C | is also *only* a prediction.
61     C | - We denote "predictors" with an asterisk (*).
62     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
63     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
64     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66     C | With implicit diffusion:
67     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
68     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69     C | (1 + dt * K * d_zz) theta[n] = theta*
70     C | (1 + dt * K * d_zz) salt[n] = salt*
71     C |
72     C *==========================================================*
73     C \ev
74    
75     C !USES:
76 adcroft 1.1 IMPLICIT NONE
77     C == Global variables ===
78     #include "SIZE.h"
79     #include "EEPARAMS.h"
80     #include "PARAMS.h"
81     #include "DYNVARS.h"
82     #include "GRID.h"
83 adcroft 1.4 #include "GAD.h"
84 adcroft 1.1 #ifdef ALLOW_PASSIVE_TRACER
85     #include "TR1.h"
86     #endif
87 jmc 1.45 #ifdef ALLOW_PTRACERS
88     #include "PTRACERS.h"
89     #endif
90 heimbach 1.42 #ifdef ALLOW_TIMEAVE
91     #include "TIMEAVE_STATV.h"
92     #endif
93    
94 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
95     # include "tamc.h"
96     # include "tamc_keys.h"
97     # include "FFIELDS.h"
98 heimbach 1.30 # include "EOS.h"
99 adcroft 1.1 # ifdef ALLOW_KPP
100     # include "KPP.h"
101     # endif
102     # ifdef ALLOW_GMREDI
103     # include "GMREDI.h"
104     # endif
105 heimbach 1.57 # ifdef ALLOW_SEAICE
106     # include "SEAICE_EXTERNAL.h"
107     # endif
108 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
109    
110 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
111 adcroft 1.1 C == Routine arguments ==
112     C myTime - Current time in simulation
113     C myIter - Current iteration number in simulation
114     C myThid - Thread number for this instance of the routine.
115     _RL myTime
116     INTEGER myIter
117     INTEGER myThid
118    
119 cnh 1.9 C !LOCAL VARIABLES:
120 adcroft 1.1 C == Local variables
121     C xA, yA - Per block temporaries holding face areas
122     C uTrans, vTrans, rTrans - Per block temporaries holding flow
123     C transport
124     C o uTrans: Zonal transport
125     C o vTrans: Meridional transport
126     C o rTrans: Vertical transport
127     C maskUp o maskUp: land/water mask for W points
128     C fVer[STUV] o fVer: Vertical flux term - note fVer
129     C is "pipelined" in the vertical
130     C so we need an fVer for each
131     C variable.
132     C rhoK, rhoKM1 - Density at current level, and level above
133     C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
134     C phiSurfY or geopotentiel (atmos) in X and Y direction
135     C KappaRT, - Total diffusion in vertical for T and S.
136     C KappaRS (background + spatially varying, isopycnal term).
137 jmc 1.39 C useVariableK = T when vertical diffusion is not constant
138 adcroft 1.1 C iMin, iMax - Ranges and sub-block indices on which calculations
139     C jMin, jMax are applied.
140     C bi, bj
141     C k, kup, - Index for layer above and below. kup and kDown
142     C kDown, km1 are switched with layer to be the appropriate
143     C index into fVerTerm.
144     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
151     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
152 heimbach 1.55 #ifdef ALLOW_PASSIVE_TRACER
153 adcroft 1.1 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
154 heimbach 1.55 #endif
155     #ifdef ALLOW_PTRACERS
156     _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
157     #endif
158 adcroft 1.1 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
159     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
160     _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162     _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
163     _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
164     _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
165     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
166     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
167 cnh 1.9 C This is currently used by IVDC and Diagnostics
168 adcroft 1.1 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
169 jmc 1.39 LOGICAL useVariableK
170 adcroft 1.1 INTEGER iMin, iMax
171     INTEGER jMin, jMax
172     INTEGER bi, bj
173     INTEGER i, j
174     INTEGER k, km1, kup, kDown
175 heimbach 1.55 INTEGER iTracer, ip
176 adcroft 1.1
177 cnh 1.9 CEOP
178 adcroft 1.40
179 edhill 1.56 #ifdef ALLOW_DEBUG
180 heimbach 1.43 IF ( debugLevel .GE. debLevB )
181     & CALL DEBUG_ENTER('FORWARD_STEP',myThid)
182 adcroft 1.40 #endif
183 adcroft 1.1
184     #ifdef ALLOW_AUTODIFF_TAMC
185     C-- dummy statement to end declaration part
186     ikey = 1
187 heimbach 1.30 itdkey = 1
188 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
189    
190     #ifdef ALLOW_AUTODIFF_TAMC
191     C-- HPF directive to help TAMC
192     CHPF$ INDEPENDENT
193     #endif /* ALLOW_AUTODIFF_TAMC */
194    
195     DO bj=myByLo(myThid),myByHi(myThid)
196    
197     #ifdef ALLOW_AUTODIFF_TAMC
198     C-- HPF directive to help TAMC
199 heimbach 1.2 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
200 jmc 1.37 CHPF$& ,utrans,vtrans,xA,yA
201 heimbach 1.2 CHPF$& ,KappaRT,KappaRS
202 adcroft 1.1 CHPF$& )
203     #endif /* ALLOW_AUTODIFF_TAMC */
204    
205     DO bi=myBxLo(myThid),myBxHi(myThid)
206    
207     #ifdef ALLOW_AUTODIFF_TAMC
208     act1 = bi - myBxLo(myThid)
209     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
210     act2 = bj - myByLo(myThid)
211     max2 = myByHi(myThid) - myByLo(myThid) + 1
212     act3 = myThid - 1
213     max3 = nTx*nTy
214     act4 = ikey_dynamics - 1
215 heimbach 1.30 itdkey = (act1 + 1) + act2*max1
216 adcroft 1.1 & + act3*max1*max2
217     & + act4*max1*max2*max3
218     #endif /* ALLOW_AUTODIFF_TAMC */
219    
220 heimbach 1.41 C-- Set up work arrays with valid (i.e. not NaN) values
221     C These inital values do not alter the numerical results. They
222     C just ensure that all memory references are to valid floating
223     C point numbers. This prevents spurious hardware signals due to
224     C uninitialised but inert locations.
225    
226 adcroft 1.1 DO j=1-OLy,sNy+OLy
227     DO i=1-OLx,sNx+OLx
228 heimbach 1.41 xA(i,j) = 0. _d 0
229     yA(i,j) = 0. _d 0
230     uTrans(i,j) = 0. _d 0
231     vTrans(i,j) = 0. _d 0
232     rhok (i,j) = 0. _d 0
233     rhoKM1 (i,j) = 0. _d 0
234     phiSurfX(i,j) = 0. _d 0
235     phiSurfY(i,j) = 0. _d 0
236 adcroft 1.1 rTrans (i,j) = 0. _d 0
237     fVerT (i,j,1) = 0. _d 0
238     fVerT (i,j,2) = 0. _d 0
239     fVerS (i,j,1) = 0. _d 0
240     fVerS (i,j,2) = 0. _d 0
241 heimbach 1.55 #ifdef ALLOW_PASSIVE_TRACER
242 adcroft 1.1 fVerTr1(i,j,1) = 0. _d 0
243     fVerTr1(i,j,2) = 0. _d 0
244 heimbach 1.55 #endif
245     #ifdef ALLOW_PTRACERS
246     DO ip=1,PTRACERS_num
247     fVerP (i,j,1,ip) = 0. _d 0
248     fVerP (i,j,2,ip) = 0. _d 0
249     ENDDO
250     #endif
251 adcroft 1.1 ENDDO
252     ENDDO
253    
254     DO k=1,Nr
255     DO j=1-OLy,sNy+OLy
256     DO i=1-OLx,sNx+OLx
257     C This is currently also used by IVDC and Diagnostics
258 heimbach 1.20 sigmaX(i,j,k) = 0. _d 0
259     sigmaY(i,j,k) = 0. _d 0
260     sigmaR(i,j,k) = 0. _d 0
261 adcroft 1.1 ConvectCount(i,j,k) = 0.
262 heimbach 1.30 KappaRT(i,j,k) = 0. _d 0
263     KappaRS(i,j,k) = 0. _d 0
264 jmc 1.45 C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
265 heimbach 1.30 gT(i,j,k,bi,bj) = 0. _d 0
266     gS(i,j,k,bi,bj) = 0. _d 0
267     # ifdef ALLOW_PASSIVE_TRACER
268 edhill 1.51 ceh3 needs an IF ( use PASSIVE_TRACER) THEN
269 heimbach 1.5 gTr1(i,j,k,bi,bj) = 0. _d 0
270 heimbach 1.30 # endif
271 heimbach 1.42 # ifdef ALLOW_PTRACERS
272 edhill 1.51 ceh3 this should have an IF ( usePTRACERS ) THEN
273 heimbach 1.42 DO iTracer=1,PTRACERS_numInUse
274     gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
275     ENDDO
276     # endif
277 jmc 1.45 #ifdef ALLOW_AUTODIFF_TAMC
278     cph all the following init. are necessary for TAF
279     cph although some of these are re-initialised later.
280 heimbach 1.30 # ifdef ALLOW_GMREDI
281     Kwx(i,j,k,bi,bj) = 0. _d 0
282     Kwy(i,j,k,bi,bj) = 0. _d 0
283     Kwz(i,j,k,bi,bj) = 0. _d 0
284     # ifdef GM_NON_UNITY_DIAGONAL
285     Kux(i,j,k,bi,bj) = 0. _d 0
286     Kvy(i,j,k,bi,bj) = 0. _d 0
287     # endif
288     # ifdef GM_EXTRA_DIAGONAL
289     Kuz(i,j,k,bi,bj) = 0. _d 0
290     Kvz(i,j,k,bi,bj) = 0. _d 0
291     # endif
292     # ifdef GM_BOLUS_ADVEC
293     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
294     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
295 heimbach 1.36 # endif
296     # ifdef GM_VISBECK_VARIABLE_K
297     VisbeckK(i,j,bi,bj) = 0. _d 0
298 heimbach 1.30 # endif
299     # endif /* ALLOW_GMREDI */
300     #endif /* ALLOW_AUTODIFF_TAMC */
301 adcroft 1.1 ENDDO
302     ENDDO
303     ENDDO
304    
305 jmc 1.21 iMin = 1-OLx
306 adcroft 1.1 iMax = sNx+OLx
307 jmc 1.21 jMin = 1-OLy
308 adcroft 1.1 jMax = sNy+OLy
309    
310     #ifdef ALLOW_AUTODIFF_TAMC
311 heimbach 1.30 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
312     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
313 heimbach 1.38 CADJ STORE totphihyd
314     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
315 heimbach 1.11 #ifdef ALLOW_KPP
316 heimbach 1.30 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
317     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
318 heimbach 1.11 #endif
319 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
320    
321 edhill 1.56 #ifdef ALLOW_DEBUG
322 heimbach 1.43 IF ( debugLevel .GE. debLevB )
323     & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
324 adcroft 1.40 #endif
325    
326 adcroft 1.1 C-- Start of diagnostic loop
327     DO k=Nr,1,-1
328    
329     #ifdef ALLOW_AUTODIFF_TAMC
330     C? Patrick, is this formula correct now that we change the loop range?
331     C? Do we still need this?
332     cph kkey formula corrected.
333     cph Needed for rhok, rhokm1, in the case useGMREDI.
334 heimbach 1.30 kkey = (itdkey-1)*Nr + k
335 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
336    
337     C-- Integrate continuity vertically for vertical velocity
338 jmc 1.27 c CALL INTEGRATE_FOR_W(
339     c I bi, bj, k, uVel, vVel,
340     c O wVel,
341     c I myThid )
342 adcroft 1.1
343     #ifdef ALLOW_OBCS
344     #ifdef ALLOW_NONHYDROSTATIC
345     C-- Apply OBC to W if in N-H mode
346 jmc 1.27 c IF (useOBCS.AND.nonHydrostatic) THEN
347     c CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
348     c ENDIF
349 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
350     #endif /* ALLOW_OBCS */
351    
352 heimbach 1.22 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
353     C-- MOST of THERMODYNAMICS will be disabled
354     #ifndef SINGLE_LAYER_MODE
355    
356 adcroft 1.1 C-- Calculate gradients of potential density for isoneutral
357     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
358     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
359     IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
360 edhill 1.56 #ifdef ALLOW_DEBUG
361 heimbach 1.43 IF ( debugLevel .GE. debLevB )
362     & CALL DEBUG_CALL('FIND_RHO',myThid)
363 adcroft 1.40 #endif
364 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
365     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
366     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
367     #endif /* ALLOW_AUTODIFF_TAMC */
368     CALL FIND_RHO(
369 mlosch 1.26 I bi, bj, iMin, iMax, jMin, jMax, k, k,
370 adcroft 1.1 I theta, salt,
371     O rhoK,
372     I myThid )
373 heimbach 1.30
374 adcroft 1.1 IF (k.GT.1) THEN
375     #ifdef ALLOW_AUTODIFF_TAMC
376     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
377     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
378     #endif /* ALLOW_AUTODIFF_TAMC */
379     CALL FIND_RHO(
380 mlosch 1.26 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
381 adcroft 1.1 I theta, salt,
382     O rhoKm1,
383     I myThid )
384     ENDIF
385 edhill 1.56 #ifdef ALLOW_DEBUG
386 heimbach 1.43 IF ( debugLevel .GE. debLevB )
387     & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
388 adcroft 1.40 #endif
389 adcroft 1.1 CALL GRAD_SIGMA(
390     I bi, bj, iMin, iMax, jMin, jMax, k,
391     I rhoK, rhoKm1, rhoK,
392     O sigmaX, sigmaY, sigmaR,
393     I myThid )
394     ENDIF
395    
396 heimbach 1.30 #ifdef ALLOW_AUTODIFF_TAMC
397     CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
398     CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
399     #endif /* ALLOW_AUTODIFF_TAMC */
400 adcroft 1.1 C-- Implicit Vertical Diffusion for Convection
401     c ==> should use sigmaR !!!
402     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
403 edhill 1.56 #ifdef ALLOW_DEBUG
404 heimbach 1.43 IF ( debugLevel .GE. debLevB )
405     & CALL DEBUG_CALL('CALC_IVDC',myThid)
406 adcroft 1.40 #endif
407 adcroft 1.1 CALL CALC_IVDC(
408     I bi, bj, iMin, iMax, jMin, jMax, k,
409     I rhoKm1, rhoK,
410     U ConvectCount, KappaRT, KappaRS,
411     I myTime, myIter, myThid)
412     ENDIF
413    
414 heimbach 1.22 #endif /* SINGLE_LAYER_MODE */
415    
416 adcroft 1.1 C-- end of diagnostic k loop (Nr:1)
417     ENDDO
418    
419     #ifdef ALLOW_AUTODIFF_TAMC
420     cph avoids recomputation of integrate_for_w
421 heimbach 1.30 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
422 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
423    
424     #ifdef ALLOW_OBCS
425     C-- Calculate future values on open boundaries
426     IF (useOBCS) THEN
427 edhill 1.56 #ifdef ALLOW_DEBUG
428 heimbach 1.43 IF ( debugLevel .GE. debLevB )
429     & CALL DEBUG_CALL('OBCS_CALC',myThid)
430 adcroft 1.40 #endif
431 jmc 1.16 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
432 adcroft 1.1 I uVel, vVel, wVel, theta, salt,
433     I myThid )
434     ENDIF
435     #endif /* ALLOW_OBCS */
436    
437 edhill 1.54
438 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
439 jmc 1.44 IF (useThermSeaIce) THEN
440 edhill 1.56 #ifdef ALLOW_DEBUG
441 heimbach 1.43 IF ( debugLevel .GE. debLevB )
442     & CALL DEBUG_CALL('ICE_FORCING',myThid)
443 adcroft 1.40 #endif
444 cheisey 1.31 C-- Determines forcing terms based on external fields
445 jmc 1.44 C including effects from ice
446 cheisey 1.31 CALL ICE_FORCING(
447     I bi, bj, iMin, iMax, jMin, jMax,
448     I myThid )
449 edhill 1.54 ELSE
450     #endif /* ALLOW_THERM_SEAICE */
451    
452     C-- Determines forcing terms based on external fields
453     C relaxation terms, etc.
454 edhill 1.56 #ifdef ALLOW_DEBUG
455 edhill 1.54 IF ( debugLevel .GE. debLevB )
456     & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
457     #endif
458     CALL EXTERNAL_FORCING_SURF(
459     I bi, bj, iMin, iMax, jMin, jMax,
460     I myTime, myIter, myThid )
461    
462     #ifdef ALLOW_THERM_SEAICE
463     C-- end of if/else block useThermSeaIce --
464 jmc 1.50 ENDIF
465     #endif /* ALLOW_THERM_SEAICE */
466 cheisey 1.31
467 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
468     cph needed for KPP
469     CADJ STORE surfacetendencyU(:,:,bi,bj)
470 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
471 adcroft 1.1 CADJ STORE surfacetendencyV(:,:,bi,bj)
472 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
473 adcroft 1.1 CADJ STORE surfacetendencyS(:,:,bi,bj)
474 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
475 adcroft 1.1 CADJ STORE surfacetendencyT(:,:,bi,bj)
476 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
477 heimbach 1.57 # ifdef ALLOW_SEAICE
478     CADJ STORE surfacetendencyTice(:,:,bi,bj)
479     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
480     # endif
481 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
482    
483 heimbach 1.22 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
484     C-- MOST of THERMODYNAMICS will be disabled
485     #ifndef SINGLE_LAYER_MODE
486    
487 adcroft 1.1 #ifdef ALLOW_GMREDI
488    
489 heimbach 1.22 #ifdef ALLOW_AUTODIFF_TAMC
490 heimbach 1.34 cph storing here is needed only for one GMREDI_OPTIONS:
491     cph define GM_BOLUS_ADVEC
492     cph but I've avoided the #ifdef for now, in case more things change
493     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
494     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
495     CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
496 heimbach 1.22 #endif /* ALLOW_AUTODIFF_TAMC */
497 heimbach 1.35
498 adcroft 1.1 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
499     IF (useGMRedi) THEN
500 edhill 1.56 #ifdef ALLOW_DEBUG
501 heimbach 1.43 IF ( debugLevel .GE. debLevB )
502     & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
503 adcroft 1.40 #endif
504 jmc 1.15 CALL GMREDI_CALC_TENSOR(
505     I bi, bj, iMin, iMax, jMin, jMax,
506 adcroft 1.1 I sigmaX, sigmaY, sigmaR,
507     I myThid )
508     #ifdef ALLOW_AUTODIFF_TAMC
509     ELSE
510 jmc 1.15 CALL GMREDI_CALC_TENSOR_DUMMY(
511     I bi, bj, iMin, iMax, jMin, jMax,
512 adcroft 1.1 I sigmaX, sigmaY, sigmaR,
513     I myThid )
514     #endif /* ALLOW_AUTODIFF_TAMC */
515     ENDIF
516    
517     #ifdef ALLOW_AUTODIFF_TAMC
518 heimbach 1.30 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
519     CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
520     CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
521 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
522    
523     #endif /* ALLOW_GMREDI */
524    
525     #ifdef ALLOW_KPP
526     C-- Compute KPP mixing coefficients
527     IF (useKPP) THEN
528 edhill 1.56 #ifdef ALLOW_DEBUG
529 heimbach 1.43 IF ( debugLevel .GE. debLevB )
530     & CALL DEBUG_CALL('KPP_CALC',myThid)
531 adcroft 1.40 #endif
532 adcroft 1.1 CALL KPP_CALC(
533     I bi, bj, myTime, myThid )
534     #ifdef ALLOW_AUTODIFF_TAMC
535     ELSE
536     CALL KPP_CALC_DUMMY(
537     I bi, bj, myTime, myThid )
538     #endif /* ALLOW_AUTODIFF_TAMC */
539     ENDIF
540    
541     #ifdef ALLOW_AUTODIFF_TAMC
542     CADJ STORE KPPghat (:,:,:,bi,bj)
543     CADJ & , KPPdiffKzT(:,:,:,bi,bj)
544     CADJ & , KPPdiffKzS(:,:,:,bi,bj)
545     CADJ & , KPPfrac (:,: ,bi,bj)
546 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
547 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
548    
549     #endif /* ALLOW_KPP */
550    
551     #ifdef ALLOW_AUTODIFF_TAMC
552 heimbach 1.30 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
553     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
554     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
555     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
556 adcroft 1.1 #ifdef ALLOW_PASSIVE_TRACER
557 heimbach 1.30 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
558 adcroft 1.1 #endif
559 heimbach 1.42 #ifdef ALLOW_PTRACERS
560     cph-- moved to forward_step to avoid key computation
561     cphCADJ STORE ptracer(:,:,:,bi,bj,itracer) = comlev1_bibj,
562     cphCADJ & key=itdkey, byte=isbyte
563     #endif
564 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
565    
566     #ifdef ALLOW_AIM
567     C AIM - atmospheric intermediate model, physics package code.
568     IF ( useAIM ) THEN
569 edhill 1.56 #ifdef ALLOW_DEBUG
570 heimbach 1.43 IF ( debugLevel .GE. debLevB )
571     & CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
572 adcroft 1.40 #endif
573 jmc 1.33 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
574     CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
575     CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
576 adcroft 1.1 ENDIF
577     #endif /* ALLOW_AIM */
578 jmc 1.14
579 adcroft 1.12 #ifndef DISABLE_MULTIDIM_ADVECTION
580 adcroft 1.4 C-- Some advection schemes are better calculated using a multi-dimensional
581     C method in the absence of any other terms and, if used, is done here.
582 adcroft 1.13 C
583     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
584     C The default is to use multi-dimensinal advection for non-linear advection
585     C schemes. However, for the sake of efficiency of the adjoint it is necessary
586     C to be able to exclude this scheme to avoid excessive storage and
587     C recomputation. It *is* differentiable, if you need it.
588     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
589     C disable this section of code.
590 jmc 1.24 IF (tempMultiDimAdvec) THEN
591 edhill 1.56 #ifdef ALLOW_DEBUG
592 heimbach 1.43 IF ( debugLevel .GE. debLevB )
593     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
594 adcroft 1.40 #endif
595 heimbach 1.11 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
596     U theta,gT,
597 adcroft 1.6 I myTime,myIter,myThid)
598 jmc 1.23 ENDIF
599 jmc 1.24 IF (saltMultiDimAdvec) THEN
600 edhill 1.56 #ifdef ALLOW_DEBUG
601 heimbach 1.43 IF ( debugLevel .GE. debLevB )
602     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
603 adcroft 1.40 #endif
604 heimbach 1.11 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
605     U salt,gS,
606 adcroft 1.6 I myTime,myIter,myThid)
607     ENDIF
608 adcroft 1.17 C Since passive tracers are configurable separately from T,S we
609     C call the multi-dimensional method for PTRACERS regardless
610     C of whether multiDimAdvection is set or not.
611     #ifdef ALLOW_PTRACERS
612     IF ( usePTRACERS ) THEN
613 edhill 1.56 #ifdef ALLOW_DEBUG
614 heimbach 1.43 IF ( debugLevel .GE. debLevB )
615     & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
616 adcroft 1.40 #endif
617 adcroft 1.17 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
618     ENDIF
619     #endif /* ALLOW_PTRACERS */
620 adcroft 1.12 #endif /* DISABLE_MULTIDIM_ADVECTION */
621 adcroft 1.1
622 edhill 1.56 #ifdef ALLOW_DEBUG
623 heimbach 1.43 IF ( debugLevel .GE. debLevB )
624     & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
625 adcroft 1.40 #endif
626    
627 adcroft 1.1 C-- Start of thermodynamics loop
628     DO k=Nr,1,-1
629     #ifdef ALLOW_AUTODIFF_TAMC
630     C? Patrick Is this formula correct?
631     cph Yes, but I rewrote it.
632     cph Also, the KappaR? need the index and subscript k!
633 heimbach 1.30 kkey = (itdkey-1)*Nr + k
634 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
635    
636     C-- km1 Points to level above k (=k-1)
637     C-- kup Cycles through 1,2 to point to layer above
638     C-- kDown Cycles through 2,1 to point to current layer
639    
640     km1 = MAX(1,k-1)
641     kup = 1+MOD(k+1,2)
642     kDown= 1+MOD(k,2)
643    
644     iMin = 1-OLx
645     iMax = sNx+OLx
646     jMin = 1-OLy
647     jMax = sNy+OLy
648    
649     C-- Get temporary terms used by tendency routines
650     CALL CALC_COMMON_FACTORS (
651     I bi,bj,iMin,iMax,jMin,jMax,k,
652     O xA,yA,uTrans,vTrans,rTrans,maskUp,
653     I myThid)
654 jmc 1.19
655     #ifdef ALLOW_GMREDI
656 heimbach 1.35
657 jmc 1.19 C-- Residual transp = Bolus transp + Eulerian transp
658     IF (useGMRedi) THEN
659     CALL GMREDI_CALC_UVFLOW(
660     & uTrans, vTrans, bi, bj, k, myThid)
661     IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
662     & rTrans, bi, bj, k, myThid)
663     ENDIF
664 heimbach 1.35
665     #ifdef ALLOW_AUTODIFF_TAMC
666     #ifdef GM_BOLUS_ADVEC
667     CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
668     CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
669     CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
670     #endif
671     #endif /* ALLOW_AUTODIFF_TAMC */
672    
673 jmc 1.19 #endif /* ALLOW_GMREDI */
674 adcroft 1.1
675     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
676     C-- Calculate the total vertical diffusivity
677     CALL CALC_DIFFUSIVITY(
678     I bi,bj,iMin,iMax,jMin,jMax,k,
679     I maskUp,
680 heimbach 1.2 O KappaRT,KappaRS,
681 adcroft 1.1 I myThid)
682 heimbach 1.52 # ifdef ALLOW_AUTODIFF_TAMC
683     CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
684     CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
685     # endif /* ALLOW_AUTODIFF_TAMC */
686 adcroft 1.1 #endif
687    
688     iMin = 1-OLx+2
689     iMax = sNx+OLx-1
690     jMin = 1-OLy+2
691     jMax = sNy+OLy-1
692    
693     C-- Calculate active tracer tendencies (gT,gS,...)
694     C and step forward storing result in gTnm1, gSnm1, etc.
695     IF ( tempStepping ) THEN
696     CALL CALC_GT(
697     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
698     I xA,yA,uTrans,vTrans,rTrans,maskUp,
699     I KappaRT,
700     U fVerT,
701 adcroft 1.7 I myTime,myIter,myThid)
702 adcroft 1.1 CALL TIMESTEP_TRACER(
703 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
704 adcroft 1.1 I theta, gT,
705     I myIter, myThid)
706     ENDIF
707 jmc 1.44
708 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
709 jmc 1.44 IF (useThermSeaIce .AND. k.EQ.1) THEN
710     CALL ICE_FREEZE( bi,bj, iMin,iMax,jMin,jMax, myThid )
711     ENDIF
712 cheisey 1.32 #endif
713 jmc 1.44
714 adcroft 1.1 IF ( saltStepping ) THEN
715     CALL CALC_GS(
716     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
717     I xA,yA,uTrans,vTrans,rTrans,maskUp,
718     I KappaRS,
719     U fVerS,
720 adcroft 1.7 I myTime,myIter,myThid)
721 adcroft 1.1 CALL TIMESTEP_TRACER(
722 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
723 adcroft 1.1 I salt, gS,
724     I myIter, myThid)
725     ENDIF
726     #ifdef ALLOW_PASSIVE_TRACER
727 edhill 1.51 ceh3 needs an IF ( usePASSIVE_TRACER ) THEN
728 adcroft 1.1 IF ( tr1Stepping ) THEN
729     CALL CALC_GTR1(
730     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
731     I xA,yA,uTrans,vTrans,rTrans,maskUp,
732     I KappaRT,
733     U fVerTr1,
734 heimbach 1.8 I myTime,myIter,myThid)
735 adcroft 1.1 CALL TIMESTEP_TRACER(
736 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
737 adcroft 1.1 I Tr1, gTr1,
738 heimbach 1.8 I myIter,myThid)
739 adcroft 1.1 ENDIF
740     #endif
741 adcroft 1.17 #ifdef ALLOW_PTRACERS
742     IF ( usePTRACERS ) THEN
743 heimbach 1.42 CALL PTRACERS_INTEGRATE(
744 adcroft 1.17 I bi,bj,k,
745     I xA,yA,uTrans,vTrans,rTrans,maskUp,
746 heimbach 1.55 X fVerP, KappaRS,
747 adcroft 1.17 I myIter,myTime,myThid)
748     ENDIF
749     #endif /* ALLOW_PTRACERS */
750 adcroft 1.1
751     #ifdef ALLOW_OBCS
752     C-- Apply open boundary conditions
753     IF (useOBCS) THEN
754 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
755 adcroft 1.1 END IF
756     #endif /* ALLOW_OBCS */
757 edhill 1.54
758     C-- Freeze water
759     IF ( allowFreezing .AND. .NOT. useSEAICE
760     & .AND. .NOT.(useThermSeaIce.AND.k.EQ.1) ) THEN
761     #ifdef ALLOW_AUTODIFF_TAMC
762     CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
763     CADJ & , key = kkey, byte = isbyte
764     #endif /* ALLOW_AUTODIFF_TAMC */
765     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
766     ENDIF
767 adcroft 1.1
768     C-- end of thermodynamic k loop (Nr:1)
769     ENDDO
770 cheisey 1.31
771     cswdice -- add ---
772 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
773 cheisey 1.31 c timeaveraging for ice model values
774 edhill 1.51 ceh3 This should be wrapped in an IF ( useThermSeaIce ) THEN
775 cheisey 1.31 CALL ICE_AVE(bi,bj,iMin,iMax,jMin,jMax,myThid )
776     #endif
777     cswdice --- end add ---
778    
779    
780    
781 adcroft 1.1
782     C-- Implicit diffusion
783     IF (implicitDiffusion) THEN
784    
785     IF (tempStepping) THEN
786     #ifdef ALLOW_AUTODIFF_TAMC
787 heimbach 1.52 CADJ STORE KappaRT(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
788 heimbach 1.30 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
789 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
790     CALL IMPLDIFF(
791     I bi, bj, iMin, iMax, jMin, jMax,
792     I deltaTtracer, KappaRT, recip_HFacC,
793 adcroft 1.7 U gT,
794 adcroft 1.1 I myThid )
795     ENDIF
796    
797     IF (saltStepping) THEN
798     #ifdef ALLOW_AUTODIFF_TAMC
799 heimbach 1.52 CADJ STORE KappaRS(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
800 heimbach 1.30 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
801 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
802     CALL IMPLDIFF(
803     I bi, bj, iMin, iMax, jMin, jMax,
804     I deltaTtracer, KappaRS, recip_HFacC,
805 adcroft 1.7 U gS,
806 adcroft 1.1 I myThid )
807     ENDIF
808    
809     #ifdef ALLOW_PASSIVE_TRACER
810     IF (tr1Stepping) THEN
811     #ifdef ALLOW_AUTODIFF_TAMC
812 heimbach 1.30 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
813 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
814     CALL IMPLDIFF(
815     I bi, bj, iMin, iMax, jMin, jMax,
816     I deltaTtracer, KappaRT, recip_HFacC,
817 adcroft 1.7 U gTr1,
818 adcroft 1.1 I myThid )
819     ENDIF
820     #endif
821    
822 adcroft 1.17 #ifdef ALLOW_PTRACERS
823     C Vertical diffusion (implicit) for passive tracers
824     IF ( usePTRACERS ) THEN
825     CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
826     ENDIF
827     #endif /* ALLOW_PTRACERS */
828    
829 adcroft 1.1 #ifdef ALLOW_OBCS
830     C-- Apply open boundary conditions
831     IF (useOBCS) THEN
832     DO K=1,Nr
833 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
834 adcroft 1.1 ENDDO
835     END IF
836     #endif /* ALLOW_OBCS */
837    
838     C-- End If implicitDiffusion
839     ENDIF
840 heimbach 1.22
841 jmc 1.39 #ifdef ALLOW_TIMEAVE
842 edhill 1.51 ceh3 needs an IF ( useTIMEAVE ) THEN
843 jmc 1.39 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
844     CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
845     I Nr, deltaTclock, bi, bj, myThid)
846     ENDIF
847     useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
848     IF (taveFreq.GT.0. .AND. useVariableK ) THEN
849     IF (implicitDiffusion) THEN
850     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
851     I Nr, 3, deltaTclock, bi, bj, myThid)
852     ELSE
853     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
854     I Nr, 3, deltaTclock, bi, bj, myThid)
855     ENDIF
856     ENDIF
857     #endif /* ALLOW_TIMEAVE */
858    
859 heimbach 1.22 #endif /* SINGLE_LAYER_MODE */
860 adcroft 1.1
861 jmc 1.39 C-- end bi,bj loops.
862 adcroft 1.1 ENDDO
863     ENDDO
864 adcroft 1.17
865 edhill 1.56 #ifdef ALLOW_DEBUG
866 adcroft 1.17 If (debugMode) THEN
867     CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
868     CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
869     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
870     CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
871     CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
872     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
873     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
874     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
875     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
876 adcroft 1.18 #ifdef ALLOW_PTRACERS
877     IF ( usePTRACERS ) THEN
878     CALL PTRACERS_DEBUG(myThid)
879     ENDIF
880     #endif /* ALLOW_PTRACERS */
881 adcroft 1.17 ENDIF
882 adcroft 1.40 #endif
883    
884 edhill 1.56 #ifdef ALLOW_DEBUG
885 heimbach 1.43 IF ( debugLevel .GE. debLevB )
886     & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
887 adcroft 1.17 #endif
888 adcroft 1.1
889     RETURN
890     END

  ViewVC Help
Powered by ViewVC 1.1.22