Back to home page

MITgcm

 
 

    


File indexing completed on 2025-03-03 06:10:47 UTC

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