Back to home page

MITgcm

 
 

    


File indexing completed on 2023-02-03 06:09:35 UTC

view on githubraw file Latest commit edb66560 on 2023-02-02 23:32:31 UTC
6d54cf9ca1 Ed H*0001 #include "PACKAGES_CONFIG.h"
b9377e56be Alis*0002 #include "CPP_OPTIONS.h"
517dbdc414 Jean*0003 #ifdef ALLOW_AUTODIFF
                0004 # include "AUTODIFF_OPTIONS.h"
                0005 #endif
aecc8b0f47 Mart*0006 #ifdef ALLOW_CTRL
                0007 # include "CTRL_OPTIONS.h"
                0008 #endif
efc81681f0 Jean*0009 #ifdef ALLOW_MOM_COMMON
                0010 # include "MOM_COMMON_OPTIONS.h"
                0011 #endif
65b69c6b17 Patr*0012 #ifdef ALLOW_OBCS
                0013 # include "OBCS_OPTIONS.h"
                0014 #endif
                0015 
2150052b9c Jean*0016 #undef DYNAMICS_GUGV_EXCH_CHECK
924557e60a Chri*0017 
9366854e02 Chri*0018 CBOP
                0019 C     !ROUTINE: DYNAMICS
                0020 C     !INTERFACE:
cf8488c0fd Chri*0021       SUBROUTINE DYNAMICS(myTime, myIter, myThid)
9366854e02 Chri*0022 C     !DESCRIPTION: \bv
                0023 C     *==========================================================*
bd9f9e3f69 Jean*0024 C     | SUBROUTINE DYNAMICS
                0025 C     | o Controlling routine for the explicit part of the model
                0026 C     |   dynamics.
9366854e02 Chri*0027 C     *==========================================================*
                0028 C     \ev
                0029 C     !USES:
9bef66796c Alis*0030       IMPLICIT NONE
924557e60a Chri*0031 C     == Global variables ===
                0032 #include "SIZE.h"
                0033 #include "EEPARAMS.h"
4e41eecc46 Alis*0034 #include "PARAMS.h"
ef4f8a46df Jean*0035 #include "GRID.h"
ab4f260161 Alis*0036 #include "DYNVARS.h"
df999eca2c Jean*0037 #include "FFIELDS.h"
efc81681f0 Jean*0038 #ifdef ALLOW_MOM_COMMON
                0039 # include "MOM_VISC.h"
                0040 #endif
a5bd9cbb61 Ed H*0041 #ifdef ALLOW_CD_CODE
ef4f8a46df Jean*0042 # include "CD_CODE_VARS.h"
a5bd9cbb61 Ed H*0043 #endif
517dbdc414 Jean*0044 #ifdef ALLOW_AUTODIFF
7c50f07931 Mart*0045 # ifdef ALLOW_AUTODIFF_TAMC
                0046 #  include "tamc.h"
                0047 # endif
7103bd8015 Patr*0048 # include "EOS.h"
b8452ee69a Patr*0049 # ifdef ALLOW_KPP
                0050 #  include "KPP.h"
                0051 # endif
65b69c6b17 Patr*0052 # ifdef ALLOW_PTRACERS
                0053 #  include "PTRACERS_SIZE.h"
85f77391e5 Jean*0054 #  include "PTRACERS_FIELDS.h"
65b69c6b17 Patr*0055 # endif
                0056 # ifdef ALLOW_OBCS
ef4f8a46df Jean*0057 #  include "OBCS_PARAMS.h"
a9eb030de2 Jean*0058 #  include "OBCS_FIELDS.h"
65b69c6b17 Patr*0059 #  ifdef ALLOW_PTRACERS
                0060 #   include "OBCS_PTRACERS.h"
                0061 #  endif
                0062 # endif
1b11d5157f Patr*0063 # ifdef ALLOW_MOM_FLUXFORM
                0064 #  include "MOM_FLUXFORM.h"
                0065 # endif
517dbdc414 Jean*0066 #endif /* ALLOW_AUTODIFF */
b91a70728d Jean*0067 
9366854e02 Chri*0068 C     !CALLING SEQUENCE:
                0069 C     DYNAMICS()
                0070 C      |
6d7fef84db Jean*0071 C      |-- CALC_EP_FORCING
                0072 C      |
9366854e02 Chri*0073 C      |-- CALC_GRAD_PHI_SURF
                0074 C      |
                0075 C      |-- CALC_VISCOSITY
                0076 C      |
efc81681f0 Jean*0077 C      |-- MOM_CALC_3D_STRAIN
                0078 C      |
bcc58d6972 Mich*0079 C      |-- CALC_EDDY_STRESS
                0080 C      |
aab550ebd0 Jean*0081 C      |-- CALC_PHI_HYD
9366854e02 Chri*0082 C      |
aab550ebd0 Jean*0083 C      |-- MOM_FLUXFORM
9366854e02 Chri*0084 C      |
aab550ebd0 Jean*0085 C      |-- MOM_VECINV
9366854e02 Chri*0086 C      |
efc81681f0 Jean*0087 C      |-- MOM_CALC_SMAG_3D
                0088 C      |-- MOM_UV_SMAG_3D
                0089 C      |
aab550ebd0 Jean*0090 C      |-- TIMESTEP
9366854e02 Chri*0091 C      |
aab550ebd0 Jean*0092 C      |-- MOM_U_IMPLICIT_R
                0093 C      |-- MOM_V_IMPLICIT_R
6d7fef84db Jean*0094 C      |
aab550ebd0 Jean*0095 C      |-- IMPLDIFF
9366854e02 Chri*0096 C      |
                0097 C      |-- OBCS_APPLY_UV
                0098 C      |
6d7fef84db Jean*0099 C      |-- CALC_GW
                0100 C      |
                0101 C      |-- DIAGNOSTICS_FILL
                0102 C      |-- DEBUG_STATS_RL
9366854e02 Chri*0103 
                0104 C     !INPUT/OUTPUT PARAMETERS:
924557e60a Chri*0105 C     == Routine arguments ==
64f8465da1 Jean*0106 C     myTime :: Current time in simulation
                0107 C     myIter :: Current iteration number in simulation
                0108 C     myThid :: Thread number for this instance of the routine.
cf8488c0fd Chri*0109       _RL myTime
                0110       INTEGER myIter
e5c33d9529 Alis*0111       INTEGER myThid
924557e60a Chri*0112 
1975511d27 Jean*0113 C     !FUNCTIONS:
                0114 #ifdef ALLOW_DIAGNOSTICS
df999eca2c Jean*0115 c     LOGICAL  DIAGNOSTICS_IS_ON
                0116 c     EXTERNAL DIAGNOSTICS_IS_ON
1975511d27 Jean*0117 #endif
                0118 
9366854e02 Chri*0119 C     !LOCAL VARIABLES:
924557e60a Chri*0120 C     == Local variables
b5c2f2589c Jean*0121 C     fVer[UV]               o fVer: Vertical flux term - note fVer
                0122 C                                    is "pipelined" in the vertical
                0123 C                                    so we need an fVer for each
                0124 C                                    variable.
6ac14281aa Jean*0125 C     phiHydC    :: hydrostatic potential anomaly at cell center
                0126 C                   In z coords phiHyd is the hydrostatic potential
                0127 C                      (=pressure/rho0) anomaly
                0128 C                   In p coords phiHyd is the geopotential height anomaly.
                0129 C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
                0130 C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
                0131 C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
7e1bd93d4f Jean*0132 C     phiSurfY             or geopotential (atmos) in X and Y direction
4667e4066d Jean*0133 C     guDissip   :: dissipation tendency (all explicit terms), u component
                0134 C     gvDissip   :: dissipation tendency (all explicit terms), v component
efd2c352d2 Jean*0135 C     kappaRU    :: vertical viscosity for velocity U-component
                0136 C     kappaRV    :: vertical viscosity for velocity V-component
07e17fa184 Jean*0137 C     iMin, iMax :: Ranges and sub-block indices on which calculations
                0138 C     jMin, jMax    are applied.
                0139 C     bi, bj     :: tile indices
                0140 C     k          :: current level index
                0141 C     km1, kp1   :: index of level above (k-1) and below (k+1)
                0142 C     kUp, kDown :: Index for interface above and below. kUp and kDown are
                0143 C                   are switched with k to be the appropriate index into fVerU,V
b8f0ad9e58 Chri*0144       _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
                0145       _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
6ac14281aa Jean*0146       _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0147       _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1022be604b Jean*0148       _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0149       _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
678051767f Jean*0150       _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0151       _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
4667e4066d Jean*0152       _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0153       _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
efd2c352d2 Jean*0154       _RL kappaRU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0155       _RL kappaRV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
efc81681f0 Jean*0156 #ifdef ALLOW_SMAG_3D
                0157 C     str11       :: strain component Vxx @ grid-cell center
                0158 C     str22       :: strain component Vyy @ grid-cell center
                0159 C     str33       :: strain component Vzz @ grid-cell center
                0160 C     str12       :: strain component Vxy @ grid-cell corner
                0161 C     str13       :: strain component Vxz @ above uVel
                0162 C     str23       :: strain component Vyz @ above vVel
                0163 C     viscAh3d_00 :: Smagorinsky viscosity @ grid-cell center
                0164 C     viscAh3d_12 :: Smagorinsky viscosity @ grid-cell corner
                0165 C     viscAh3d_13 :: Smagorinsky viscosity @ above uVel
                0166 C     viscAh3d_23 :: Smagorinsky viscosity @ above vVel
                0167 C     addDissU    :: zonal momentum tendency from 3-D Smag. viscosity
                0168 C     addDissV    :: merid momentum tendency from 3-D Smag. viscosity
                0169       _RL str11(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0170       _RL str22(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0171       _RL str33(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0172       _RL str12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0173       _RL str13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0174       _RL str23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0175       _RL viscAh3d_00(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0176       _RL viscAh3d_12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
                0177       _RL viscAh3d_13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0178       _RL viscAh3d_23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0179       _RL addDissU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0180       _RL addDissV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0181 #elif ( defined ALLOW_NONHYDROSTATIC )
                0182       _RL str13(1), str23(1), str33(1)
                0183       _RL viscAh3d_00(1), viscAh3d_13(1), viscAh3d_23(1)
                0184 #endif
779cd6d73d Alis*0185 
924557e60a Chri*0186       INTEGER bi, bj
                0187       INTEGER i, j
07e17fa184 Jean*0188       INTEGER k, km1, kp1, kUp, kDown
efc81681f0 Jean*0189       INTEGER iMin, iMax
                0190       INTEGER jMin, jMax
                0191       PARAMETER( iMin = 0 , iMax = sNx+1 )
                0192       PARAMETER( jMin = 0 , jMax = sNy+1 )
7c50f07931 Mart*0193 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0194 C     tkey :: tape key (tile dependent)
                0195 C     kkey :: tape key (level and tile dependent)
                0196       INTEGER tkey, kkey
7c50f07931 Mart*0197 #endif
924557e60a Chri*0198 
b5c2f2589c Jean*0199 #ifdef ALLOW_DIAGNOSTICS
df999eca2c Jean*0200 c     LOGICAL dPhiHydDiagIsOn
9ad6623f4d Jean*0201       _RL tmpFac
b5c2f2589c Jean*0202 #endif /* ALLOW_DIAGNOSTICS */
                0203 
0bb99fb476 Alis*0204 C---    The algorithm...
                0205 C
                0206 C       "Correction Step"
                0207 C       =================
                0208 C       Here we update the horizontal velocities with the surface
                0209 C       pressure such that the resulting flow is either consistent
                0210 C       with the free-surface evolution or the rigid-lid:
                0211 C         U[n] = U* + dt x d/dx P
                0212 C         V[n] = V* + dt x d/dy P
                0213 C
                0214 C       "Calculation of Gs"
                0215 C       ===================
                0216 C       This is where all the accelerations and tendencies (ie.
4e30ef22d3 Patr*0217 C       physics, parameterizations etc...) are calculated
0bb99fb476 Alis*0218 C         rho = rho ( theta[n], salt[n] )
f1a86f1057 Chri*0219 C         b   = b(rho, theta)
0bb99fb476 Alis*0220 C         K31 = K31 ( rho )
f15b75110d Jean*0221 C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
                0222 C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
                0223 C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
                0224 C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
0bb99fb476 Alis*0225 C
779cd6d73d Alis*0226 C       "Time-stepping" or "Prediction"
0bb99fb476 Alis*0227 C       ================================
                0228 C       The models variables are stepped forward with the appropriate
                0229 C       time-stepping scheme (currently we use Adams-Bashforth II)
                0230 C       - For momentum, the result is always *only* a "prediction"
                0231 C       in that the flow may be divergent and will be "corrected"
                0232 C       later with a surface pressure gradient.
                0233 C       - Normally for tracers the result is the new field at time
                0234 C       level [n+1} *BUT* in the case of implicit diffusion the result
                0235 C       is also *only* a prediction.
                0236 C       - We denote "predictors" with an asterisk (*).
                0237 C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
                0238 C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
                0239 C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
                0240 C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
779cd6d73d Alis*0241 C       With implicit diffusion:
0bb99fb476 Alis*0242 C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
                0243 C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
779cd6d73d Alis*0244 C         (1 + dt * K * d_zz) theta[n] = theta*
                0245 C         (1 + dt * K * d_zz) salt[n] = salt*
0bb99fb476 Alis*0246 C---
9366854e02 Chri*0247 CEOP
0bb99fb476 Alis*0248 
2a77e09ef0 Jean*0249 #ifdef ALLOW_DEBUG
8440e8ae5d Jean*0250       IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
2a77e09ef0 Jean*0251 #endif
                0252 
1975511d27 Jean*0253 #ifdef ALLOW_DIAGNOSTICS
df999eca2c Jean*0254 c     dPhiHydDiagIsOn = .FALSE.
                0255 c     IF ( useDiagnostics )
                0256 c    &  dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
                0257 c    &               .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
1975511d27 Jean*0258 #endif
                0259 
ef4f8a46df Jean*0260 C-- Call to routine for calculation of Eliassen-Palm-flux-forced
                0261 C    U-tendency, if desired:
b275747e24 Patr*0262 #ifdef INCLUDE_EP_FORCING_CODE
                0263       CALL CALC_EP_FORCING(myThid)
                0264 #endif
                0265 
467a85a3d2 Patr*0266 #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
1022be604b Jean*0267       CALL DUMMY_IN_DYNAMICS( myTime, myIter, myThid )
467a85a3d2 Patr*0268 #endif
                0269 
5646dbad6f Patr*0270 #ifdef ALLOW_AUTODIFF_TAMC
                0271 C--   HPF directive to help TAMC
                0272 CHPF$ INDEPENDENT
                0273 #endif /* ALLOW_AUTODIFF_TAMC */
                0274 
924557e60a Chri*0275       DO bj=myByLo(myThid),myByHi(myThid)
5646dbad6f Patr*0276 
                0277 #ifdef ALLOW_AUTODIFF_TAMC
                0278 C--    HPF directive to help TAMC
                0279 CHPF$  INDEPENDENT, NEW (fVerU,fVerV
6ac14281aa Jean*0280 CHPF$&                  ,phiHydF
efd2c352d2 Jean*0281 CHPF$&                  ,kappaRU,kappaRV
5646dbad6f Patr*0282 CHPF$&                  )
                0283 #endif /* ALLOW_AUTODIFF_TAMC */
                0284 
924557e60a Chri*0285        DO bi=myBxLo(myThid),myBxHi(myThid)
5646dbad6f Patr*0286 
                0287 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0288         tkey = bi + (bj-1)*nSx + (ikey_dynamics-1)*nSx*nSy
5646dbad6f Patr*0289 #endif /* ALLOW_AUTODIFF_TAMC */
                0290 
df4ecb0ca5 Patr*0291 C--   Set up work arrays with valid (i.e. not NaN) values
1022be604b Jean*0292 C     These initial values do not alter the numerical results. They
df4ecb0ca5 Patr*0293 C     just ensure that all memory references are to valid floating
                0294 C     point numbers. This prevents spurious hardware signals due to
                0295 C     uninitialised but inert locations.
                0296 
b1d6b00360 Jean*0297 #ifdef ALLOW_AUTODIFF
6ac14281aa Jean*0298         DO k=1,Nr
                0299          DO j=1-OLy,sNy+OLy
                0300           DO i=1-OLx,sNx+OLx
df4ecb0ca5 Patr*0301 c--   need some re-initialisation here to break dependencies
6d7fef84db Jean*0302            gU(i,j,k,bi,bj) = 0. _d 0
                0303            gV(i,j,k,bi,bj) = 0. _d 0
c4549b0208 Patr*0304           ENDDO
6ac14281aa Jean*0305          ENDDO
                0306         ENDDO
b1d6b00360 Jean*0307 #endif /* ALLOW_AUTODIFF */
6ac14281aa Jean*0308         DO j=1-OLy,sNy+OLy
                0309          DO i=1-OLx,sNx+OLx
5646dbad6f Patr*0310           fVerU  (i,j,1) = 0. _d 0
                0311           fVerU  (i,j,2) = 0. _d 0
                0312           fVerV  (i,j,1) = 0. _d 0
                0313           fVerV  (i,j,2) = 0. _d 0
aab550ebd0 Jean*0314           phiHydF (i,j)  = 0. _d 0
                0315           phiHydC (i,j)  = 0. _d 0
d7fe92a587 Jean*0316 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
7e1bd93d4f Jean*0317           dPhiHydX(i,j)  = 0. _d 0
aab550ebd0 Jean*0318           dPhiHydY(i,j)  = 0. _d 0
d7fe92a587 Jean*0319 #endif
df4ecb0ca5 Patr*0320           phiSurfX(i,j)  = 0. _d 0
                0321           phiSurfY(i,j)  = 0. _d 0
4667e4066d Jean*0322           guDissip(i,j)  = 0. _d 0
                0323           gvDissip(i,j)  = 0. _d 0
b1d6b00360 Jean*0324 #ifdef ALLOW_AUTODIFF
62c877cddf Gael*0325           phiHydLow(i,j,bi,bj) = 0. _d 0
3631fe78d1 Jean*0326 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
81b43af3f0 Patr*0327 #  ifndef DISABLE_RSTAR_CODE
                0328           dWtransC(i,j,bi,bj) = 0. _d 0
                0329           dWtransU(i,j,bi,bj) = 0. _d 0
                0330           dWtransV(i,j,bi,bj) = 0. _d 0
                0331 #  endif
                0332 # endif
b1d6b00360 Jean*0333 #endif /* ALLOW_AUTODIFF */
5646dbad6f Patr*0334          ENDDO
                0335         ENDDO
df999eca2c Jean*0336         IF ( useDiagnostics ) THEN
                0337          DO j=1-OLy,sNy+OLy
                0338           DO i=1-OLx,sNx+OLx
                0339            botDragU(i,j,bi,bj)  = 0. _d 0
                0340            botDragV(i,j,bi,bj)  = 0. _d 0
                0341           ENDDO
                0342          ENDDO
                0343         ENDIF
8adf9f02ba Patr*0344 
678051767f Jean*0345 C--     Start computation of dynamics
                0346 
5646dbad6f Patr*0347 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0348 CADJ STORE wVel (:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
5646dbad6f Patr*0349 #endif /* ALLOW_AUTODIFF_TAMC */
                0350 
1022be604b Jean*0351 C--     Explicit part of the Surface Potential Gradient (add in TIMESTEP)
678051767f Jean*0352 C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
                0353         IF (implicSurfPress.NE.1.) THEN
b864f02e61 Jean*0354           CALL CALC_GRAD_PHI_SURF(
                0355      I         bi,bj,iMin,iMax,jMin,jMax,
                0356      I         etaN,
                0357      O         phiSurfX,phiSurfY,
aab550ebd0 Jean*0358      I         myThid )
678051767f Jean*0359         ENDIF
fb481a83c2 Alis*0360 
67a1e439d8 Patr*0361 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0362 CADJ STORE uVel     (:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
                0363 CADJ STORE vVel     (:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
67a1e439d8 Patr*0364 #ifdef ALLOW_KPP
edb6656069 Mart*0365 CADJ STORE KPPviscAz(:,:,:,bi,bj) = comlev1_bibj, key=tkey, byte=isbyte
67a1e439d8 Patr*0366 #endif /* ALLOW_KPP */
                0367 #endif /* ALLOW_AUTODIFF_TAMC */
                0368 
b1d6b00360 Jean*0369 #ifndef ALLOW_AUTODIFF
ef4f8a46df Jean*0370         IF ( .NOT.momViscosity ) THEN
b1d6b00360 Jean*0371 #endif
efd2c352d2 Jean*0372           DO k=1,Nr+1
ef4f8a46df Jean*0373            DO j=1-OLy,sNy+OLy
                0374             DO i=1-OLx,sNx+OLx
efd2c352d2 Jean*0375              kappaRU(i,j,k) = 0. _d 0
                0376              kappaRV(i,j,k) = 0. _d 0
ef4f8a46df Jean*0377             ENDDO
                0378            ENDDO
                0379           ENDDO
b1d6b00360 Jean*0380 #ifndef ALLOW_AUTODIFF
                0381         ENDIF
                0382 #endif
ef4f8a46df Jean*0383 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
64f8465da1 Jean*0384 C--     Calculate the total vertical viscosity
ef4f8a46df Jean*0385         IF ( momViscosity ) THEN
                0386           CALL CALC_VISCOSITY(
64f8465da1 Jean*0387      I            bi,bj, iMin,iMax,jMin,jMax,
efd2c352d2 Jean*0388      O            kappaRU, kappaRV,
64f8465da1 Jean*0389      I            myThid )
ef4f8a46df Jean*0390         ENDIF
                0391 #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
1dcbfc76b6 Patr*0392 
efc81681f0 Jean*0393 #ifdef ALLOW_SMAG_3D
                0394         IF ( useSmag3D ) THEN
                0395           CALL MOM_CALC_3D_STRAIN(
                0396      O         str11, str22, str33, str12, str13, str23,
                0397      I         bi, bj, myThid )
                0398         ENDIF
                0399 #endif /* ALLOW_SMAG_3D */
                0400 
056b8e541a Patr*0401 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0402 CADJ STORE kappaRU(:,:,:) = comlev1_bibj, key=tkey, byte=isbyte
                0403 CADJ STORE kappaRV(:,:,:) = comlev1_bibj, key=tkey, byte=isbyte
056b8e541a Patr*0404 #endif /* ALLOW_AUTODIFF_TAMC */
                0405 
04605efaca Mart*0406 #ifdef ALLOW_OBCS
                0407 C--   For Stevens boundary conditions velocities need to be extrapolated
                0408 C     (copied) to a narrow strip outside the domain
ef4f8a46df Jean*0409         IF ( useOBCS ) THEN
d64c4d306c Jean*0410           CALL OBCS_COPY_UV_N(
1022be604b Jean*0411      U         uVel(1-OLx,1-OLy,1,bi,bj),
                0412      U         vVel(1-OLx,1-OLy,1,bi,bj),
04605efaca Mart*0413      I         Nr, bi, bj, myThid )
ef4f8a46df Jean*0414         ENDIF
04605efaca Mart*0415 #endif /* ALLOW_OBCS */
                0416 
bcc58d6972 Mich*0417 #ifdef ALLOW_EDDYPSI
                0418         CALL CALC_EDDY_STRESS(bi,bj,myThid)
                0419 #endif
                0420 
fb481a83c2 Alis*0421 C--     Start of dynamics loop
                0422         DO k=1,Nr
                0423 
                0424 C--       km1    Points to level above k (=k-1)
                0425 C--       kup    Cycles through 1,2 to point to layer above
                0426 C--       kDown  Cycles through 2,1 to point to current layer
                0427 
                0428           km1  = MAX(1,k-1)
1dcbfc76b6 Patr*0429           kp1  = MIN(k+1,Nr)
fb481a83c2 Alis*0430           kup  = 1+MOD(k+1,2)
                0431           kDown= 1+MOD(k,2)
                0432 
bd9f9e3f69 Jean*0433 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0434          kkey = k + (tkey-1)*Nr
1022be604b Jean*0435 CADJ STORE totPhiHyd (:,:,k,bi,bj)
1f46e82020 Patr*0436 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
af26a7fa45 Gael*0437 CADJ STORE phiHydLow (:,:,bi,bj)
517dbdc414 Jean*0438 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0439 CADJ STORE theta (:,:,k,bi,bj)
edc3e287c7 Patr*0440 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0441 CADJ STORE salt  (:,:,k,bi,bj)
edc3e287c7 Patr*0442 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
cdc9f269ae Patr*0443 # ifdef NONLIN_FRSURF
                0444 cph-test
3631fe78d1 Jean*0445 CADJ STORE  phiHydC (:,:)
cdc9f269ae Patr*0446 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0447 CADJ STORE  phiHydF (:,:)
cdc9f269ae Patr*0448 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0449 CADJ STORE gU(:,:,k,bi,bj)
cdc9f269ae Patr*0450 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0451 CADJ STORE gV(:,:,k,bi,bj)
cdc9f269ae Patr*0452 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
b7b55ab99c Patr*0453 #  ifndef ALLOW_ADAMSBASHFORTH_3
1022be604b Jean*0454 CADJ STORE guNm1(:,:,k,bi,bj)
cdc9f269ae Patr*0455 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0456 CADJ STORE gvNm1(:,:,k,bi,bj)
cdc9f269ae Patr*0457 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
b7b55ab99c Patr*0458 #  else
1022be604b Jean*0459 CADJ STORE guNm(:,:,k,bi,bj,1)
b7b55ab99c Patr*0460 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0461 CADJ STORE guNm(:,:,k,bi,bj,2)
b7b55ab99c Patr*0462 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0463 CADJ STORE gvNm(:,:,k,bi,bj,1)
b7b55ab99c Patr*0464 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0465 CADJ STORE gvNm(:,:,k,bi,bj,2)
b7b55ab99c Patr*0466 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
                0467 #  endif
cdc9f269ae Patr*0468 #  ifdef ALLOW_CD_CODE
1022be604b Jean*0469 CADJ STORE uNM1(:,:,k,bi,bj)
cdc9f269ae Patr*0470 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
1022be604b Jean*0471 CADJ STORE vNM1(:,:,k,bi,bj)
cdc9f269ae Patr*0472 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0473 CADJ STORE uVelD(:,:,k,bi,bj)
cdc9f269ae Patr*0474 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0475 CADJ STORE vVelD(:,:,k,bi,bj)
cdc9f269ae Patr*0476 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
                0477 #  endif
b1d6b00360 Jean*0478 # endif /* NONLIN_FRSURF */
5646dbad6f Patr*0479 #endif /* ALLOW_AUTODIFF_TAMC */
                0480 
ef4f8a46df Jean*0481 C--      Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
bd27360393 Jean*0482          CALL CALC_PHI_HYD(
fb481a83c2 Alis*0483      I        bi,bj,iMin,iMax,jMin,jMax,k,
6ac14281aa Jean*0484      U        phiHydF,
                0485      O        phiHydC, dPhiHydX, dPhiHydY,
7e1bd93d4f Jean*0486      I        myTime, myIter, myThid )
fb481a83c2 Alis*0487 
                0488 C--      Calculate accelerations in the momentum equations (gU, gV, ...)
c0911a93b9 Jean*0489 C        and step forward storing the result in gU, gV, etc...
fb481a83c2 Alis*0490          IF ( momStepping ) THEN
b1d6b00360 Jean*0491 #ifdef ALLOW_AUTODIFF
                0492            DO j=1-OLy,sNy+OLy
                0493             DO i=1-OLx,sNx+OLx
                0494               guDissip(i,j)  = 0. _d 0
                0495               gvDissip(i,j)  = 0. _d 0
                0496             ENDDO
                0497            ENDDO
                0498 #endif /* ALLOW_AUTODIFF */
81b43af3f0 Patr*0499 #ifdef ALLOW_AUTODIFF_TAMC
b1d6b00360 Jean*0500 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
                0501 #  ifndef DISABLE_RSTAR_CODE
3631fe78d1 Jean*0502 CADJ STORE dWtransC(:,:,bi,bj)
81b43af3f0 Patr*0503 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0504 CADJ STORE dWtransU(:,:,bi,bj)
81b43af3f0 Patr*0505 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
3631fe78d1 Jean*0506 CADJ STORE dWtransV(:,:,bi,bj)
81b43af3f0 Patr*0507 CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
                0508 #  endif
b1d6b00360 Jean*0509 # endif /* NONLIN_FRSURF and ALLOW_MOM_FLUXFORM */
                0510 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
                0511 CADJ STORE fVerU(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
                0512 CADJ STORE fVerV(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
                0513 # endif
07e17fa184 Jean*0514 #endif /* ALLOW_AUTODIFF_TAMC */
e2c5f2fe7c Patr*0515            IF (.NOT. vectorInvariantMomentum) THEN
99bc607d7a Ed H*0516 #ifdef ALLOW_MOM_FLUXFORM
e2c5f2fe7c Patr*0517               CALL MOM_FLUXFORM(
07e17fa184 Jean*0518      I         bi,bj,k,iMin,iMax,jMin,jMax,
efd2c352d2 Jean*0519      I         kappaRU, kappaRV,
07e17fa184 Jean*0520      U         fVerU(1-OLx,1-OLy,kUp),   fVerV(1-OLx,1-OLy,kUp),
                0521      O         fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
722a76e285 Jean*0522      O         guDissip, gvDissip,
3a279374db Alis*0523      I         myTime, myIter, myThid)
66d5e1666c Alis*0524 #endif
e2c5f2fe7c Patr*0525            ELSE
99bc607d7a Ed H*0526 #ifdef ALLOW_MOM_VECINV
cdc9f269ae Patr*0527              CALL MOM_VECINV(
07e17fa184 Jean*0528      I         bi,bj,k,iMin,iMax,jMin,jMax,
efd2c352d2 Jean*0529      I         kappaRU, kappaRV,
07e17fa184 Jean*0530      I         fVerU(1-OLx,1-OLy,kUp),   fVerV(1-OLx,1-OLy,kUp),
                0531      O         fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
4667e4066d Jean*0532      O         guDissip, gvDissip,
3a279374db Alis*0533      I         myTime, myIter, myThid)
66d5e1666c Alis*0534 #endif
e2c5f2fe7c Patr*0535            ENDIF
ef4f8a46df Jean*0536 
efc81681f0 Jean*0537 #ifdef ALLOW_SMAG_3D
                0538            IF ( useSmag3D ) THEN
                0539              CALL MOM_CALC_SMAG_3D(
                0540      I         str11, str22, str33, str12, str13, str23,
                0541      O         viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
                0542      I         smag3D_hLsC, smag3D_hLsW, smag3D_hLsS, smag3D_hLsZ,
                0543      I         k, bi, bj, myThid )
                0544              CALL MOM_UV_SMAG_3D(
                0545      I         str11, str22, str12, str13, str23,
                0546      I         viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
                0547      O         addDissU, addDissV,
                0548      I         iMin,iMax,jMin,jMax, k, bi, bj, myThid )
                0549              DO j= jMin,jMax
                0550               DO i= iMin,iMax
                0551                guDissip(i,j) = guDissip(i,j) + addDissU(i,j)
                0552                gvDissip(i,j) = gvDissip(i,j) + addDissV(i,j)
                0553               ENDDO
                0554              ENDDO
                0555            ENDIF
                0556 #endif /* ALLOW_SMAG_3D */
                0557 
fb481a83c2 Alis*0558            CALL TIMESTEP(
678051767f Jean*0559      I         bi,bj,iMin,iMax,jMin,jMax,k,
6ac14281aa Jean*0560      I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
4667e4066d Jean*0561      I         guDissip, gvDissip,
c0911a93b9 Jean*0562      I         myTime, myIter, myThid)
fb481a83c2 Alis*0563 
                0564          ENDIF
                0565 
                0566 C--     end of dynamics k loop (1:Nr)
                0567         ENDDO
                0568 
b34050da0f Jean*0569 C--     Implicit Vertical advection & viscosity
d666da8b50 Gael*0570 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
5dc0884b3e Jean*0571      defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF))
855d57fc61 Jean*0572         IF ( momImplVertAdv .OR. implicitViscosity
                0573      &                      .OR. selectImplicitDrag.GE.1 ) THEN
5685a5c914 Jean*0574 C      to recover older (prior to 2016-10-05) results:
                0575 c       IF ( momImplVertAdv ) THEN
aab550ebd0 Jean*0576           CALL MOM_U_IMPLICIT_R( kappaRU,
b34050da0f Jean*0577      I                           bi, bj, myTime, myIter, myThid )
                0578           CALL MOM_V_IMPLICIT_R( kappaRV,
                0579      I                           bi, bj, myTime, myIter, myThid )
                0580         ELSEIF ( implicitViscosity ) THEN
                0581 #else /* INCLUDE_IMPLVERTADV_CODE */
                0582         IF     ( implicitViscosity ) THEN
                0583 #endif /* INCLUDE_IMPLVERTADV_CODE */
fb481a83c2 Alis*0584 #ifdef    ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0585 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=tkey, byte=isbyte
fb481a83c2 Alis*0586 #endif    /* ALLOW_AUTODIFF_TAMC */
fb3dc7d949 Alis*0587           CALL IMPLDIFF(
                0588      I         bi, bj, iMin, iMax, jMin, jMax,
efd2c352d2 Jean*0589      I         -1, kappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
75a9b2b4ba Jean*0590      U         gU(1-OLx,1-OLy,1,bi,bj),
fb3dc7d949 Alis*0591      I         myThid )
fb481a83c2 Alis*0592 #ifdef    ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0593 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=tkey, byte=isbyte
fb481a83c2 Alis*0594 #endif    /* ALLOW_AUTODIFF_TAMC */
fb3dc7d949 Alis*0595           CALL IMPLDIFF(
                0596      I         bi, bj, iMin, iMax, jMin, jMax,
efd2c352d2 Jean*0597      I         -2, kappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
75a9b2b4ba Jean*0598      U         gV(1-OLx,1-OLy,1,bi,bj),
fb3dc7d949 Alis*0599      I         myThid )
b34050da0f Jean*0600         ENDIF
8adf9f02ba Patr*0601 
04605efaca Mart*0602 #ifdef ALLOW_OBCS
fb481a83c2 Alis*0603 C--      Apply open boundary conditions
f549345f7f Jean*0604         IF ( useOBCS ) THEN
d64c4d306c Jean*0605 C--      but first save intermediate velocities to be used in the
04605efaca Mart*0606 C        next time step for the Stevens boundary conditions
d64c4d306c Jean*0607           CALL OBCS_SAVE_UV_N(
                0608      I        bi, bj, iMin, iMax, jMin, jMax, 0,
04605efaca Mart*0609      I        gU, gV, myThid )
f549345f7f Jean*0610           CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
b34050da0f Jean*0611         ENDIF
04605efaca Mart*0612 #endif /* ALLOW_OBCS */
8adf9f02ba Patr*0613 
138482fdf6 Ed H*0614 #ifdef    ALLOW_CD_CODE
b34050da0f Jean*0615         IF (implicitViscosity.AND.useCDscheme) THEN
fb481a83c2 Alis*0616 #ifdef    ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0617 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=tkey, byte=isbyte
fb481a83c2 Alis*0618 #endif    /* ALLOW_AUTODIFF_TAMC */
fb3dc7d949 Alis*0619           CALL IMPLDIFF(
                0620      I         bi, bj, iMin, iMax, jMin, jMax,
efd2c352d2 Jean*0621      I         0, kappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
75a9b2b4ba Jean*0622      U         vVelD(1-OLx,1-OLy,1,bi,bj),
fb3dc7d949 Alis*0623      I         myThid )
fb481a83c2 Alis*0624 #ifdef    ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0625 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=tkey, byte=isbyte
fb481a83c2 Alis*0626 #endif    /* ALLOW_AUTODIFF_TAMC */
fb3dc7d949 Alis*0627           CALL IMPLDIFF(
                0628      I         bi, bj, iMin, iMax, jMin, jMax,
efd2c352d2 Jean*0629      I         0, kappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
75a9b2b4ba Jean*0630      U         uVelD(1-OLx,1-OLy,1,bi,bj),
fb3dc7d949 Alis*0631      I         myThid )
4e30ef22d3 Patr*0632         ENDIF
b34050da0f Jean*0633 #endif    /* ALLOW_CD_CODE */
                0634 C--     End implicit Vertical advection & viscosity
64e4c8fd27 Patr*0635 
b5c2f2589c Jean*0636 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0637 
6d7fef84db Jean*0638 #ifdef ALLOW_NONHYDROSTATIC
                0639 C--   Step forward W field in N-H algorithm
aab550ebd0 Jean*0640         IF ( nonHydrostatic ) THEN
6d7fef84db Jean*0641 #ifdef ALLOW_DEBUG
8440e8ae5d Jean*0642          IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
6d7fef84db Jean*0643 #endif
                0644          CALL TIMER_START('CALC_GW          [DYNAMICS]',myThid)
ecb1e63329 Bayl*0645          CALL CALC_GW(
efd2c352d2 Jean*0646      I                 bi,bj, kappaRU, kappaRV,
efc81681f0 Jean*0647      I                 str13, str23, str33,
                0648      I                 viscAh3d_00, viscAh3d_13, viscAh3d_23,
aab550ebd0 Jean*0649      I                 myTime, myIter, myThid )
                0650         ENDIF
                0651         IF ( nonHydrostatic.OR.implicitIntGravWave )
                0652      &   CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
                0653         IF ( nonHydrostatic )
cb7fa97db9 Jean*0654      &   CALL TIMER_STOP ('CALC_GW          [DYNAMICS]',myThid)
6d7fef84db Jean*0655 #endif
                0656 
                0657 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0658 
aab550ebd0 Jean*0659 C-    end of bi,bj loops
                0660        ENDDO
                0661       ENDDO
                0662 
                0663 #ifdef ALLOW_OBCS
14e8dd3ad1 Jean*0664       IF (useOBCS) THEN
2150416ce8 Jean*0665         CALL OBCS_EXCHANGES( myThid )
14e8dd3ad1 Jean*0666       ENDIF
aab550ebd0 Jean*0667 #endif
                0668 
60c223928f Mart*0669 Cml(
                0670 C     In order to compare the variance of phiHydLow of a p/z-coordinate
                0671 C     run with etaH of a z/p-coordinate run the drift of phiHydLow
                0672 C     has to be removed by something like the following subroutine:
bd9f9e3f69 Jean*0673 C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
375581a429 Jean*0674 C     &                     'phiHydLow', myTime, myIter, myThid )
60c223928f Mart*0675 Cml)
                0676 
fb2c52cfd4 aver*0677 #ifdef INCLUDE_SOUNDSPEED_CALC_CODE
                0678       CALL DIAGS_SOUND_SPEED( myThid )
                0679 #endif
                0680 
b5c2f2589c Jean*0681 #ifdef ALLOW_DIAGNOSTICS
d2454c5ee6 Jean*0682       IF ( useDiagnostics ) THEN
b5c2f2589c Jean*0683 
                0684        CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD  ',0,Nr,0,1,1,myThid)
9ad6623f4d Jean*0685        CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT  ',0, 1,0,1,1,myThid)
                0686        tmpFac = 1. _d 0
                0687        CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
                0688      &                                 'PHIHYDSQ',0,Nr,0,1,1,myThid)
725dc592ee Andr*0689 
9ad6623f4d Jean*0690        CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
                0691      &                                 'PHIBOTSQ',0, 1,0,1,1,myThid)
b5c2f2589c Jean*0692 
df999eca2c Jean*0693        IF ( selectImplicitDrag.EQ.0 .AND.
                0694      &      (  no_slip_bottom
                0695      &    .OR. selectBotDragQuadr.GE.0
                0696      &    .OR. bottomDragLinear.NE.0. ) ) THEN
a54ff488da Jean*0697         CALL DIAGNOSTICS_FILL_RS( botDragU, 'botTauX ',
87f9394d43 Jean*0698      &                            0, 1, 0, 1, 1, myThid )
a54ff488da Jean*0699         CALL DIAGNOSTICS_FILL_RS( botDragV, 'botTauY ',
87f9394d43 Jean*0700      &                            0, 1, 0, 1, 1, myThid )
df999eca2c Jean*0701        ENDIF
                0702 
b5c2f2589c Jean*0703       ENDIF
                0704 #endif /* ALLOW_DIAGNOSTICS */
aab550ebd0 Jean*0705 
49e3578e36 Ed H*0706 #ifdef ALLOW_DEBUG
23d1f65433 Jean*0707       IF ( debugLevel .GE. debLevD ) THEN
f1ad9ab213 Alis*0708        CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
031876176a Alis*0709        CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
f1ad9ab213 Alis*0710        CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
                0711        CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
                0712        CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
                0713        CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
53d94a28a7 Jean*0714        CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
                0715        CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
                0716 #ifndef ALLOW_ADAMSBASHFORTH_3
                0717        CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
                0718        CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
                0719        CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
                0720        CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
                0721 #endif
edc5b80aef Alis*0722       ENDIF
f1ad9ab213 Alis*0723 #endif
                0724 
2150052b9c Jean*0725 #ifdef DYNAMICS_GUGV_EXCH_CHECK
bd9f9e3f69 Jean*0726 C- jmc: For safety checking only: This Exchange here should not change
                0727 C       the solution. If solution changes, it means something is wrong,
2150052b9c Jean*0728 C       but it does not mean that it is less wrong with this exchange.
23d1f65433 Jean*0729       IF ( debugLevel .GE. debLevE ) THEN
2150052b9c Jean*0730        CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
                0731       ENDIF
                0732 #endif
                0733 
2a77e09ef0 Jean*0734 #ifdef ALLOW_DEBUG
8440e8ae5d Jean*0735       IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
2a77e09ef0 Jean*0736 #endif
                0737 
924557e60a Chri*0738       RETURN
                0739       END