Back to home page

MITgcm

 
 

    


File indexing completed on 2023-05-28 05:09:52 UTC

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