Back to home page

MITgcm

 
 

    


File indexing completed on 2026-05-27 05:08:18 UTC

view on githubraw file Latest commit 24d87fa9 on 2026-05-26 09:10:47 UTC
6d54cf9ca1 Ed H*0001 #include "PACKAGES_CONFIG.h"
aea29c8517 Alis*0002 #include "CPP_OPTIONS.h"
517dbdc414 Jean*0003 #ifdef ALLOW_AUTODIFF
                0004 # include "AUTODIFF_OPTIONS.h"
                0005 #endif
7e62bca48a Jean*0006 #if (defined ALLOW_PTRACERS) && (!defined ALLOW_LONGSTEP)
                0007 # define DO_PTRACERS_HERE
4e66ab0b67 Oliv*0008 #endif
6d54cf9ca1 Ed H*0009 
517dbdc414 Jean*0010 #ifdef ALLOW_AUTODIFF
3fbf7b8b4a Jean*0011 # ifdef ALLOW_GMREDI
                0012 #  include "GMREDI_OPTIONS.h"
                0013 # endif
                0014 # ifdef ALLOW_KPP
                0015 #  include "KPP_OPTIONS.h"
                0016 # endif
80fd556159 Ou W*0017 # ifdef ALLOW_SALT_PLUME
                0018 #  include "SALT_PLUME_OPTIONS.h"
                0019 # endif
517dbdc414 Jean*0020 #endif /* ALLOW_AUTODIFF */
aecc8b0f47 Mart*0021 #ifdef ALLOW_CTRL
                0022 # include "CTRL_OPTIONS.h"
                0023 #endif
aea29c8517 Alis*0024 
9366854e02 Chri*0025 CBOP
                0026 C     !ROUTINE: THERMODYNAMICS
                0027 C     !INTERFACE:
7ffd8c9970 Alis*0028       SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
9366854e02 Chri*0029 C     !DESCRIPTION: \bv
                0030 C     *==========================================================*
a27159adf7 Jean*0031 C     | SUBROUTINE THERMODYNAMICS
                0032 C     | o Controlling routine for the prognostic part of the
                0033 C     |   thermo-dynamics.
9366854e02 Chri*0034 C     *===========================================================
                0035 C     \ev
aea29c8517 Alis*0036 
9366854e02 Chri*0037 C     !USES:
                0038       IMPLICIT NONE
aea29c8517 Alis*0039 C     == Global variables ===
                0040 #include "SIZE.h"
                0041 #include "EEPARAMS.h"
                0042 #include "PARAMS.h"
117ee419f5 Jean*0043 #include "RESTART.h"
aea29c8517 Alis*0044 #include "DYNVARS.h"
                0045 #include "GRID.h"
d0e9cdda89 Jean*0046 #include "SURFACE.h"
f4d71cd088 Jean*0047 #include "FFIELDS.h"
3ff07dd7e9 Jean*0048 #ifdef ALLOW_GENERIC_ADVDIFF
c4e839d3c4 Jean*0049 # include "GAD.h"
3ff07dd7e9 Jean*0050 #endif
7e62bca48a Jean*0051 #ifdef DO_PTRACERS_HERE
c4e839d3c4 Jean*0052 # include "PTRACERS_SIZE.h"
85f77391e5 Jean*0053 # include "PTRACERS_PARAMS.h"
                0054 # include "PTRACERS_FIELDS.h"
cf04370137 Jean*0055 #endif
9b89fcf692 antn*0056 #ifdef ALLOW_LAYERS
                0057 # include "LAYERS_P2SHARE.h"
                0058 #endif
cf2549e769 Patr*0059 
517dbdc414 Jean*0060 #ifdef ALLOW_AUTODIFF
7c50f07931 Mart*0061 # ifdef ALLOW_AUTODIFF_TAMC
                0062 #  include "tamc.h"
                0063 # endif
7103bd8015 Patr*0064 # include "EOS.h"
1389d71047 Chri*0065 # ifdef ALLOW_KPP
                0066 #  include "KPP.h"
                0067 # endif
                0068 # ifdef ALLOW_GMREDI
                0069 #  include "GMREDI.h"
                0070 # endif
b08554040b Patr*0071 # ifdef ALLOW_EBM
                0072 #  include "EBM.h"
d8206d87ee Patr*0073 # endif
6b89e1b973 Gael*0074 # ifdef ALLOW_SALT_PLUME
                0075 #  include "SALT_PLUME.h"
                0076 # endif
517dbdc414 Jean*0077 #endif /* ALLOW_AUTODIFF */
aea29c8517 Alis*0078 
b4daa24319 Shre*0079 #ifdef ALLOW_TAPENADE
                0080 # ifdef ALLOW_GENERIC_ADVDIFF
                0081 #  include "GAD_SOM_VARS.h"
                0082 # endif
                0083 #endif /* ALLOW_TAPENADE */
                0084 
9366854e02 Chri*0085 C     !INPUT/OUTPUT PARAMETERS:
aea29c8517 Alis*0086 C     == Routine arguments ==
4794324d0b Jean*0087 C     myTime :: Current time in simulation
                0088 C     myIter :: Current iteration number in simulation
                0089 C     myThid :: Thread number for this instance of the routine.
aea29c8517 Alis*0090       _RL myTime
                0091       INTEGER myIter
                0092       INTEGER myThid
                0093 
3ff07dd7e9 Jean*0094 #ifdef ALLOW_GENERIC_ADVDIFF
dbfd781425 Jean*0095 # ifdef ALLOW_MONITOR
                0096 C     !FUNCTIONS:
                0097       LOGICAL  DIFFERENT_MULTIPLE
                0098       EXTERNAL DIFFERENT_MULTIPLE
                0099 # endif /* ALLOW_MONITOR */
                0100 
9366854e02 Chri*0101 C     !LOCAL VARIABLES:
aea29c8517 Alis*0102 C     == Local variables
4794324d0b Jean*0103 C     uFld,vFld,wFld :: Local copy of velocity field (3 components)
                0104 C     kappaRk        :: Total diffusion in vertical, all levels, 1 tracer
0b191c5f5a Jean*0105 C     recip_hFacNew  :: reciprocal of futur time-step hFacC
4794324d0b Jean*0106 C     bi, bj         :: Tile indices
                0107 C     i, j, k        :: loop indices
                0108       _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
                0109       _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
                0110       _RL wFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
d0e9cdda89 Jean*0111       _RL kappaRk (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
                0112       _RS recip_hFacNew(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
aea29c8517 Alis*0113       INTEGER bi, bj
4794324d0b Jean*0114       INTEGER i, j, k
dbfd781425 Jean*0115 # ifdef ALLOW_MONITOR
                0116       LOGICAL monOutputCFL
                0117       _RL wrTime
                0118       _RL trAdvCFL(3,nSx,nSy)
                0119 # endif /* ALLOW_MONITOR */
7c50f07931 Mart*0120 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0121 C     tkey :: tape key (tile dependent)
                0122       INTEGER tkey
7c50f07931 Mart*0123 #endif
9366854e02 Chri*0124 CEOP
73ead277e0 Alis*0125 
49e3578e36 Ed H*0126 #ifdef ALLOW_DEBUG
1378a26402 Jean*0127       IF (debugMode) CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
73ead277e0 Alis*0128 #endif
a27159adf7 Jean*0129 
dbfd781425 Jean*0130 # ifdef ALLOW_MONITOR
                0131       monOutputCFL = .FALSE.
                0132       IF ( monitorSelect.GE.2 ) THEN
                0133         wrTime = myTime
                0134         IF ( .NOT.staggerTimeStep ) wrTime = myTime + deltaTClock
                0135         monOutputCFL =
                0136      &       DIFFERENT_MULTIPLE( monitorFreq, wrTime, deltaTClock )
                0137       ENDIF
                0138 # endif /* ALLOW_MONITOR */
ff1bd1eafd Jean*0139 
aea29c8517 Alis*0140 #ifdef ALLOW_AUTODIFF_TAMC
                0141 C--   dummy statement to end declaration part
edb6656069 Mart*0142       tkey = 1
aea29c8517 Alis*0143 
b7b61e618a Mart*0144 C--   HPF directive to help TAF
aea29c8517 Alis*0145 CHPF$ INDEPENDENT
                0146 #endif /* ALLOW_AUTODIFF_TAMC */
                0147 
f2d1ba7d38 Davi*0148 C-- Compute correction at the surface for Lin Free Surf.
517dbdc414 Jean*0149 #ifdef ALLOW_AUTODIFF
81b43af3f0 Patr*0150       TsurfCor = 0. _d 0
                0151       SsurfCor = 0. _d 0
                0152 #endif
4794324d0b Jean*0153       IF ( linFSConserveTr ) THEN
bef5cbe0ec Patr*0154 #ifdef ALLOW_AUTODIFF_TAMC
8a7ca847ff Patr*0155 CADJ STORE theta,salt,wvel = comlev1, key = ikey_dynamics, byte=isbyte
bef5cbe0ec Patr*0156 #endif
4794324d0b Jean*0157        CALL CALC_WSURF_TR( theta, salt, wVel,
                0158      &                     myTime, myIter, myThid )
50d8304171 Ryan*0159 #ifdef ALLOW_LAYERS
24d87fa921 Mart*0160        IF ( layers_useThermo ) THEN
28b880eddd Jean*0161         CALL LAYERS_WSURF_TR( theta, salt, wVel,
24d87fa921 Mart*0162      &                        myTime, myIter, myThid )
                0163        ENDIF
50d8304171 Ryan*0164 #endif /* ALLOW_LAYERS */
24d87fa921 Mart*0165       ENDIF
50d8304171 Ryan*0166 
4776222c09 Patr*0167 #ifdef DO_PTRACERS_HERE
28b880eddd Jean*0168 #ifdef ALLOW_AUTODIFF
cb262aefe7 Mart*0169       DO k=1,PTRACERS_numInUse
28b880eddd Jean*0170         meanSurfCorPTr(k) = 0.0 _d 0
cb262aefe7 Mart*0171       ENDDO
                0172 #endif
4da58558f5 Oliv*0173       IF ( PTRACERS_calcSurfCor ) THEN
cb262aefe7 Mart*0174 #ifdef ALLOW_AUTODIFF_TAMC
                0175 CADJ STORE ptracer = comlev1, key = ikey_dynamics, byte=isbyte
                0176 #endif
4da58558f5 Oliv*0177        CALL PTRACERS_CALC_WSURF_TR(wVel,myTime,myIter,myThid)
                0178       ENDIF
28b880eddd Jean*0179 #endif /* DO_PTRACERS_HERE */
4da58558f5 Oliv*0180 
aea29c8517 Alis*0181       DO bj=myByLo(myThid),myByHi(myThid)
                0182        DO bi=myBxLo(myThid),myBxHi(myThid)
                0183 
                0184 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0185         tkey = bi + (bj-1)*nSx + (ikey_dynamics-1)*nSx*nSy
aea29c8517 Alis*0186 #endif /* ALLOW_AUTODIFF_TAMC */
                0187 
df4ecb0ca5 Patr*0188 C--   Set up work arrays with valid (i.e. not NaN) values
                0189 C     These inital values do not alter the numerical results. They
                0190 C     just ensure that all memory references are to valid floating
                0191 C     point numbers. This prevents spurious hardware signals due to
                0192 C     uninitialised but inert locations.
                0193 
aea29c8517 Alis*0194         DO k=1,Nr
                0195          DO j=1-OLy,sNy+OLy
                0196           DO i=1-OLx,sNx+OLx
4794324d0b Jean*0197            recip_hFacNew(i,j,k) = 0. _d 0
aea29c8517 Alis*0198 C This is currently also used by IVDC and Diagnostics
a4e4b5f62b Jean*0199            kappaRk(i,j,k)    = 0. _d 0
aea29c8517 Alis*0200           ENDDO
                0201          ENDDO
                0202         ENDDO
                0203 
0b191c5f5a Jean*0204 C--     Compute new reciprocal hFac for implicit calculation
                0205 #ifdef NONLIN_FRSURF
                0206         IF ( nonlinFreeSurf.GT.0 ) THEN
                0207          IF ( select_rStar.GT.0 ) THEN
                0208 # ifndef DISABLE_RSTAR_CODE
                0209           DO k=1,Nr
                0210            DO j=1-OLy,sNy+OLy
                0211             DO i=1-OLx,sNx+OLx
                0212              recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
                0213      &                            / rStarExpC(i,j,bi,bj)
                0214             ENDDO
                0215            ENDDO
                0216           ENDDO
                0217 # endif /* DISABLE_RSTAR_CODE */
                0218          ELSEIF ( selectSigmaCoord.NE.0 ) THEN
                0219 # ifndef DISABLE_SIGMA_CODE
                0220           DO k=1,Nr
                0221            DO j=1-OLy,sNy+OLy
                0222             DO i=1-OLx,sNx+OLx
                0223              recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
                0224      &        /( 1. _d 0 + dEtaHdt(i,j,bi,bj)*deltaTFreeSurf
                0225      &                    *dBHybSigF(k)*recip_drF(k)
                0226      &                    *recip_hFacC(i,j,k,bi,bj)
                0227      &         )
                0228             ENDDO
                0229            ENDDO
                0230           ENDDO
                0231 # endif /* DISABLE_RSTAR_CODE */
                0232          ELSE
                0233           DO k=1,Nr
                0234            DO j=1-OLy,sNy+OLy
                0235             DO i=1-OLx,sNx+OLx
                0236              IF ( k.EQ.kSurfC(i,j,bi,bj) ) THEN
                0237               recip_hFacNew(i,j,k) = 1. _d 0 / hFac_surfC(i,j,bi,bj)
                0238              ELSE
                0239               recip_hFacNew(i,j,k) = recip_hFacC(i,j,k,bi,bj)
                0240              ENDIF
                0241             ENDDO
                0242            ENDDO
                0243           ENDDO
                0244          ENDIF
                0245         ELSE
                0246 #endif /* NONLIN_FRSURF */
                0247           DO k=1,Nr
                0248            DO j=1-OLy,sNy+OLy
                0249             DO i=1-OLx,sNx+OLx
                0250              recip_hFacNew(i,j,k) = _recip_hFacC(i,j,k,bi,bj)
                0251             ENDDO
                0252            ENDDO
                0253           ENDDO
                0254 #ifdef NONLIN_FRSURF
                0255         ENDIF
                0256 #endif /* NONLIN_FRSURF */
                0257 
4794324d0b Jean*0258 C--   Set up 3-D velocity field that we use to advect tracers:
                0259 C-    just do a local copy:
                0260         DO k=1,Nr
                0261          DO j=1-OLy,sNy+OLy
                0262           DO i=1-OLx,sNx+OLx
                0263            uFld(i,j,k) = uVel(i,j,k,bi,bj)
                0264            vFld(i,j,k) = vVel(i,j,k,bi,bj)
                0265            wFld(i,j,k) = wVel(i,j,k,bi,bj)
                0266           ENDDO
                0267          ENDDO
                0268         ENDDO
                0269 #ifdef ALLOW_GMREDI
                0270 C-    add Bolus velocity to Eulerian-mean velocity:
                0271         IF (useGMRedi) THEN
dbfd781425 Jean*0272           CALL GMREDI_RESIDUAL_FLOW(
4794324d0b Jean*0273      U                  uFld, vFld, wFld,
                0274      I                  bi, bj, myIter, myThid )
                0275         ENDIF
                0276 #endif /* ALLOW_GMREDI */
dbfd781425 Jean*0277 #ifdef ALLOW_MONITOR
                0278         IF ( monOutputCFL  ) THEN
                0279           CALL MON_CALC_ADVCFL_TILE( Nr, bi, bj,
                0280      I                         uFld, vFld, wFld, dTtracerLev,
                0281      O                         trAdvCFL(1,bi,bj),
                0282      I                         myIter, myThid )
                0283 c        WRITE(standardMessageUnit,'(A,I8,2I3,A,1P3E14.6)')
                0284 c    &     ' trAdv_CFL: it,bi,bj=', myIter,bi,bj,
                0285 c    &     ' , CFL =', (trAdvCFL(i,bi,bj),i=1,3)
                0286         ENDIF
                0287 #endif /* ALLOW_MONITOR */
1776f0dce8 Patr*0288 
ff1bd1eafd Jean*0289 C-    Apply AB on T,S : moved inside TEMP/SALT_INTEGRATE
                0290 
aea29c8517 Alis*0291 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0292 CADJ STORE recip_hFacNew(:,:,:) = comlev1_bibj, key=tkey, byte=isbyte
                0293 CADJ STORE uFld (:,:,:)         = comlev1_bibj, key=tkey, byte=isbyte
                0294 CADJ STORE vFld (:,:,:)         = comlev1_bibj, key=tkey, byte=isbyte
                0295 CADJ STORE wFld (:,:,:)         = comlev1_bibj, key=tkey, byte=isbyte
                0296 CADJ STORE theta(:,:,:,bi,bj)   = comlev1_bibj, key=tkey, byte=isbyte
                0297 CADJ STORE salt (:,:,:,bi,bj)   = comlev1_bibj, key=tkey, byte=isbyte
6b89e1b973 Gael*0298 # ifdef ALLOW_SALT_PLUME
edb6656069 Mart*0299 CADJ STORE saltPlumeFlux(:,:,bi,bj)  = comlev1_bibj,key=tkey,kind=isbyte
                0300 CADJ STORE saltPlumeDepth(:,:,bi,bj) = comlev1_bibj,key=tkey,kind=isbyte
6b89e1b973 Gael*0301 # endif
aecc8b0f47 Mart*0302 # if (defined NONLIN_FRSURF || defined ALLOW_DEPTH_CONTROL) \
                0303         && defined ALLOW_GMREDI
edb6656069 Mart*0304 CADJ STORE kux(:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
                0305 CADJ STORE kvy(:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
4794324d0b Jean*0306 #  ifdef GM_EXTRA_DIAGONAL
edb6656069 Mart*0307 CADJ STORE kuz(:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
                0308 CADJ STORE kvz(:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
0499e7ef48 Patr*0309 #  endif
d6e8c9ce95 Patr*0310 # endif
                0311 #endif /* ALLOW_AUTODIFF_TAMC */
f9de86979b Jean*0312 
0b6eaf7dba Jean*0313 C--     Calculate active tracer tendencies and step forward in time.
                0314 C       Active tracer arrays are updated while adjustments (filters,
                0315 C       conv.adjustment) are applied later in TRACERS_CORRECTION_STEP.
fc2298c7c4 Patr*0316 
fb8688a588 Jean*0317         IF ( tempStepping ) THEN
                0318 #ifdef ALLOW_DEBUG
                0319           IF (debugMode) CALL DEBUG_CALL('TEMP_INTEGRATE',myThid)
                0320 #endif
4794324d0b Jean*0321           CALL TEMP_INTEGRATE(
0b191c5f5a Jean*0322      I         bi, bj, recip_hFacNew,
fb8688a588 Jean*0323      I         uFld, vFld, wFld,
                0324      U         kappaRk,
4794324d0b Jean*0325      I         myTime, myIter, myThid )
                0326         ENDIF
6b89e1b973 Gael*0327 
4794324d0b Jean*0328         IF ( saltStepping ) THEN
fb8688a588 Jean*0329 #ifdef ALLOW_DEBUG
                0330           IF (debugMode) CALL DEBUG_CALL('SALT_INTEGRATE',myThid)
                0331 #endif
4794324d0b Jean*0332           CALL SALT_INTEGRATE(
0b191c5f5a Jean*0333      I         bi, bj, recip_hFacNew,
fb8688a588 Jean*0334      I         uFld, vFld, wFld,
                0335      U         kappaRk,
4794324d0b Jean*0336      I         myTime, myIter, myThid )
                0337         ENDIF
a27159adf7 Jean*0338 
7e62bca48a Jean*0339 #ifdef DO_PTRACERS_HERE
0b6eaf7dba Jean*0340 C--     Calculate passive tracer tendencies and step forward in time.
                0341 C       Passive tracer arrays are updated while adjustments (filters,
                0342 C       conv.adjustment) are applied later in TRACERS_CORRECTION_STEP.
0b191c5f5a Jean*0343 C       Also apply open boundary conditions for each passive tracer
4794324d0b Jean*0344         IF ( usePTRACERS ) THEN
fb8688a588 Jean*0345 #ifdef ALLOW_DEBUG
                0346           IF (debugMode) CALL DEBUG_CALL('PTRACERS_INTEGRATE',myThid)
                0347 #endif
f4d71cd088 Jean*0348           CALL PTRACERS_INTEGRATE(
0b191c5f5a Jean*0349      I          bi, bj, recip_hFacNew,
                0350      I          uFld, vFld, wFld,
4794324d0b Jean*0351      U          kappaRk,
                0352      I          myTime, myIter, myThid )
                0353         ENDIF
7e62bca48a Jean*0354 #endif /* DO_PTRACERS_HERE */
aea29c8517 Alis*0355 
                0356 #ifdef   ALLOW_OBCS
fb8688a588 Jean*0357 C--   Apply open boundary conditions
a135ccae2c Jean*0358         IF ( useOBCS ) THEN
f4d71cd088 Jean*0359           CALL OBCS_APPLY_TS( bi, bj, 0, theta, salt, myThid )
aea29c8517 Alis*0360         ENDIF
b34050da0f Jean*0361 #endif   /* ALLOW_OBCS */
aea29c8517 Alis*0362 
f4d71cd088 Jean*0363 #ifdef ALLOW_FRICTION_HEATING
                0364 #ifdef ALLOW_DIAGNOSTICS
                0365         IF ( addFrictionHeating .AND. useDiagnostics ) THEN
                0366           CALL DIAGNOSTICS_FILL_RS( frictionHeating, 'HeatDiss',
                0367      &                              0, Nr, 1, bi, bj, myThid )
                0368         ENDIF
                0369 #endif /* ALLOW_DIAGNOSTICS */
                0370 C-    Reset frictionHeating to zero
                0371         IF ( addFrictionHeating .AND. .NOT.staggerTimeStep ) THEN
                0372           DO k=1,Nr
                0373            DO j=1-OLy,sNy+OLy
                0374             DO i=1-OLx,sNx+OLx
                0375               frictionHeating(i,j,k,bi,bj) = 0. _d 0
                0376             ENDDO
                0377            ENDDO
                0378           ENDDO
                0379         ENDIF
                0380 #endif /* ALLOW_FRICTION_HEATING */
                0381 
408a59859d Jean*0382 C--   end bi,bj loops.
aea29c8517 Alis*0383        ENDDO
                0384       ENDDO
                0385 
dbfd781425 Jean*0386 #ifdef ALLOW_MONITOR
                0387       IF ( monOutputCFL ) THEN
                0388         CALL MON_CALC_ADVCFL_GLOB(
                0389      I                       trAdvCFL, myIter, myThid )
                0390       ENDIF
                0391 #endif /* ALLOW_MONITOR */
                0392 
49e3578e36 Ed H*0393 #ifdef ALLOW_DEBUG
23d1f65433 Jean*0394       IF ( debugLevel.GE.debLevD ) THEN
5c43c390b6 Alis*0395        CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
                0396        CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
                0397        CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
                0398        CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
                0399        CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
53d94a28a7 Jean*0400 #ifndef ALLOW_ADAMSBASHFORTH_3
                0401        CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
                0402        CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
                0403 #endif
7e62bca48a Jean*0404 #ifdef DO_PTRACERS_HERE
f128ea9284 Alis*0405        IF ( usePTRACERS ) THEN
                0406          CALL PTRACERS_DEBUG(myThid)
                0407        ENDIF
7e62bca48a Jean*0408 #endif /* DO_PTRACERS_HERE */
5c43c390b6 Alis*0409       ENDIF
53d94a28a7 Jean*0410 #endif /* ALLOW_DEBUG */
5c43c390b6 Alis*0411 
49e3578e36 Ed H*0412 #ifdef ALLOW_DEBUG
1378a26402 Jean*0413       IF (debugMode) CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
73ead277e0 Alis*0414 #endif
                0415 
3ff07dd7e9 Jean*0416 #endif /* ALLOW_GENERIC_ADVDIFF */
                0417 
aea29c8517 Alis*0418       RETURN
                0419       END