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

Diff of /MITgcm/model/src/dynamics.F

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

revision 1.43 by adcroft, Mon May 24 15:42:23 1999 UTC revision 1.170 by jmc, Fri Apr 4 20:54:11 2014 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4    #include "PACKAGES_CONFIG.h"
5  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
6    #ifdef ALLOW_AUTODIFF
7    # include "AUTODIFF_OPTIONS.h"
8    #endif
9    #ifdef ALLOW_MOM_COMMON
10    # include "MOM_COMMON_OPTIONS.h"
11    #endif
12    #ifdef ALLOW_OBCS
13    # include "OBCS_OPTIONS.h"
14    #endif
15    
16    #undef DYNAMICS_GUGV_EXCH_CHECK
17    
18    CBOP
19    C     !ROUTINE: DYNAMICS
20    C     !INTERFACE:
21        SUBROUTINE DYNAMICS(myTime, myIter, myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
22  C     /==========================================================\  C     !DESCRIPTION: \bv
23  C     | SUBROUTINE DYNAMICS                                      |  C     *==========================================================*
24  C     | o Controlling routine for the explicit part of the model |  C     | SUBROUTINE DYNAMICS
25  C     |   dynamics.                                              |  C     | o Controlling routine for the explicit part of the model
26  C     |==========================================================|  C     |   dynamics.
27  C     | This routine evaluates the "dynamics" terms for each     |  C     *==========================================================*
28  C     | block of ocean in turn. Because the blocks of ocean have |  C     | This routine evaluates the "dynamics" terms for each
29  C     | overlap regions they are independent of one another.     |  C     | block of ocean in turn. Because the blocks of ocean have
30  C     | If terms involving lateral integrals are needed in this  |  C     | overlap regions they are independent of one another.
31  C     | routine care will be needed. Similarly finite-difference |  C     | If terms involving lateral integrals are needed in this
32  C     | operations with stencils wider than the overlap region   |  C     | routine care will be needed. Similarly finite-difference
33  C     | require special consideration.                           |  C     | operations with stencils wider than the overlap region
34  C     | Notes                                                    |  C     | require special consideration.
35  C     | =====                                                    |  C     | The algorithm...
36  C     | C*P* comments indicating place holders for which code is |  C     |
37  C     |      presently being developed.                          |  C     | "Correction Step"
38  C     \==========================================================/  C     | =================
39    C     | Here we update the horizontal velocities with the surface
40    C     | pressure such that the resulting flow is either consistent
41    C     | with the free-surface evolution or the rigid-lid:
42    C     |   U[n] = U* + dt x d/dx P
43    C     |   V[n] = V* + dt x d/dy P
44    C     |   W[n] = W* + dt x d/dz P  (NH mode)
45    C     |
46    C     | "Calculation of Gs"
47    C     | ===================
48    C     | This is where all the accelerations and tendencies (ie.
49    C     | physics, parameterizations etc...) are calculated
50    C     |   rho = rho ( theta[n], salt[n] )
51    C     |   b   = b(rho, theta)
52    C     |   K31 = K31 ( rho )
53    C     |   Gu[n] = Gu( u[n], v[n], wVel, b, ... )
54    C     |   Gv[n] = Gv( u[n], v[n], wVel, b, ... )
55    C     |   Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
56    C     |   Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
57    C     |
58    C     | "Time-stepping" or "Prediction"
59    C     | ================================
60    C     | The models variables are stepped forward with the appropriate
61    C     | time-stepping scheme (currently we use Adams-Bashforth II)
62    C     | - For momentum, the result is always *only* a "prediction"
63    C     | in that the flow may be divergent and will be "corrected"
64    C     | later with a surface pressure gradient.
65    C     | - Normally for tracers the result is the new field at time
66    C     | level [n+1} *BUT* in the case of implicit diffusion the result
67    C     | is also *only* a prediction.
68    C     | - We denote "predictors" with an asterisk (*).
69    C     |   U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
70    C     |   V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
71    C     |   theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
72    C     |   salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
73    C     | With implicit diffusion:
74    C     |   theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
75    C     |   salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
76    C     |   (1 + dt * K * d_zz) theta[n] = theta*
77    C     |   (1 + dt * K * d_zz) salt[n] = salt*
78    C     |
79    C     *==========================================================*
80    C     \ev
81    C     !USES:
82        IMPLICIT NONE        IMPLICIT NONE
   
83  C     == Global variables ===  C     == Global variables ===
84  #include "SIZE.h"  #include "SIZE.h"
85  #include "EEPARAMS.h"  #include "EEPARAMS.h"
 #include "CG2D.h"  
86  #include "PARAMS.h"  #include "PARAMS.h"
 #include "DYNVARS.h"  
87  #include "GRID.h"  #include "GRID.h"
88  #ifdef ALLOW_KPP  #include "DYNVARS.h"
89  #include "KPPMIX.h"  #ifdef ALLOW_MOM_COMMON
90    # include "MOM_VISC.h"
91    #endif
92    #ifdef ALLOW_CD_CODE
93    # include "CD_CODE_VARS.h"
94  #endif  #endif
95    #ifdef ALLOW_AUTODIFF
96    # include "tamc.h"
97    # include "tamc_keys.h"
98    # include "FFIELDS.h"
99    # include "EOS.h"
100    # ifdef ALLOW_KPP
101    #  include "KPP.h"
102    # endif
103    # ifdef ALLOW_PTRACERS
104    #  include "PTRACERS_SIZE.h"
105    #  include "PTRACERS_FIELDS.h"
106    # endif
107    # ifdef ALLOW_OBCS
108    #  include "OBCS_PARAMS.h"
109    #  include "OBCS_FIELDS.h"
110    #  ifdef ALLOW_PTRACERS
111    #   include "OBCS_PTRACERS.h"
112    #  endif
113    # endif
114    # ifdef ALLOW_MOM_FLUXFORM
115    #  include "MOM_FLUXFORM.h"
116    # endif
117    #endif /* ALLOW_AUTODIFF */
118    
119    C     !CALLING SEQUENCE:
120    C     DYNAMICS()
121    C      |
122    C      |-- CALC_EP_FORCING
123    C      |
124    C      |-- CALC_GRAD_PHI_SURF
125    C      |
126    C      |-- CALC_VISCOSITY
127    C      |
128    C      |-- MOM_CALC_3D_STRAIN
129    C      |
130    C      |-- CALC_EDDY_STRESS
131    C      |
132    C      |-- CALC_PHI_HYD
133    C      |
134    C      |-- MOM_FLUXFORM
135    C      |
136    C      |-- MOM_VECINV
137    C      |
138    C      |-- MOM_CALC_SMAG_3D
139    C      |-- MOM_UV_SMAG_3D
140    C      |
141    C      |-- TIMESTEP
142    C      |
143    C      |-- MOM_U_IMPLICIT_R
144    C      |-- MOM_V_IMPLICIT_R
145    C      |
146    C      |-- IMPLDIFF
147    C      |
148    C      |-- OBCS_APPLY_UV
149    C      |
150    C      |-- CALC_GW
151    C      |
152    C      |-- DIAGNOSTICS_FILL
153    C      |-- DEBUG_STATS_RL
154    
155    C     !INPUT/OUTPUT PARAMETERS:
156  C     == Routine arguments ==  C     == Routine arguments ==
157  C     myTime - Current time in simulation  C     myTime :: Current time in simulation
158  C     myIter - Current iteration number in simulation  C     myIter :: Current iteration number in simulation
159  C     myThid - Thread number for this instance of the routine.  C     myThid :: Thread number for this instance of the routine.
       INTEGER myThid  
160        _RL myTime        _RL myTime
161        INTEGER myIter        INTEGER myIter
162          INTEGER myThid
163    
164    C     !FUNCTIONS:
165    #ifdef ALLOW_DIAGNOSTICS
166          LOGICAL  DIAGNOSTICS_IS_ON
167          EXTERNAL DIAGNOSTICS_IS_ON
168    #endif
169    
170    C     !LOCAL VARIABLES:
171  C     == Local variables  C     == Local variables
172  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[UV]               o fVer: Vertical flux term - note fVer
173  C     uTrans, vTrans, rTrans - Per block temporaries holding flow  C                                    is "pipelined" in the vertical
174  C                              transport  C                                    so we need an fVer for each
175  C     rVel                     o uTrans: Zonal transport  C                                    variable.
176  C                              o vTrans: Meridional transport  C     phiHydC    :: hydrostatic potential anomaly at cell center
177  C                              o rTrans: Vertical transport  C                   In z coords phiHyd is the hydrostatic potential
178  C                              o rVel:   Vertical velocity at upper and  C                      (=pressure/rho0) anomaly
179  C                                        lower cell faces.  C                   In p coords phiHyd is the geopotential height anomaly.
180  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
181  C                              o maskUp: land/water mask for W points  C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
182  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
183  C     mTerm, pTerm,            tendency equations.  C     phiSurfY             or geopotential (atmos) in X and Y direction
184  C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  C     guDissip   :: dissipation tendency (all explicit terms), u component
185  C                              o xTerm: Mixing term  C     gvDissip   :: dissipation tendency (all explicit terms), v component
186  C                              o cTerm: Coriolis term  C     KappaRU    :: vertical viscosity for velocity U-component
187  C                              o mTerm: Metric term  C     KappaRV    :: vertical viscosity for velocity V-component
188  C                              o pTerm: Pressure term  C     iMin, iMax :: Ranges and sub-block indices on which calculations
189  C                              o fZon: Zonal flux term  C     jMin, jMax    are applied.
190  C                              o fMer: Meridional flux term  C     bi, bj     :: tile indices
191  C                              o fVer: Vertical flux term - note fVer  C     k          :: current level index
192  C                                      is "pipelined" in the vertical  C     km1, kp1   :: index of level above (k-1) and below (k+1)
193  C                                      so we need an fVer for each  C     kUp, kDown :: Index for interface above and below. kUp and kDown are
194  C                                      variable.  C                   are switched with k to be the appropriate index into fVerU,V
 C     rhoK, rhoKM1   - Density at current level, level above and level  
 C                      below.  
 C     rhoKP1                                                                    
 C     buoyK, buoyKM1 - Buoyancy at current level and level above.  
 C     phiHyd         - Hydrostatic part of the potential phiHydi.  
 C                      In z coords phiHydiHyd is the hydrostatic  
 C                      pressure anomaly  
 C                      In p coords phiHydiHyd is the geopotential  
 C                      surface height  
 C                      anomaly.  
 C     etaSurfX,      - Holds surface elevation gradient in X and Y.  
 C     etaSurfY  
 C     K13, K23, K33  - Non-zero elements of small-angle approximation  
 C                      diffusion tensor.  
 C     KapGM          - Spatially varying Visbeck et. al mixing coeff.  
 C     KappaRT,       - Total diffusion in vertical for T and S.  
 C     KappaRS          (background + spatially varying, isopycnal term).  
 C     iMin, iMax     - Ranges and sub-block indices on which calculations  
 C     jMin, jMax       are applied.  
 C     bi, bj  
 C     k, kUp,        - Index for layer above and below. kUp and kDown  
 C     kDown, kM1       are switched with layer to be the appropriate  
 C                      index into fVerTerm.  
       _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
195        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
196        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
197        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
198        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
199        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
200        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
201        _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
202        _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
203        _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
204        _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
205        _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
206        _RL K13     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaRV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
207        _RL K23     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  #ifdef ALLOW_SMAG_3D
208        _RL K33     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  C     str11       :: strain component Vxx @ grid-cell center
209        _RL KapGM   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     str22       :: strain component Vyy @ grid-cell center
210        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)  C     str33       :: strain component Vzz @ grid-cell center
211        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)  C     str12       :: strain component Vxy @ grid-cell corner
212        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)  C     str13       :: strain component Vxz @ above uVel
213        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)  C     str23       :: strain component Vyz @ above vVel
214    C     viscAh3d_00 :: Smagorinsky viscosity @ grid-cell center
215    C     viscAh3d_12 :: Smagorinsky viscosity @ grid-cell corner
216    C     viscAh3d_13 :: Smagorinsky viscosity @ above uVel
217    C     viscAh3d_23 :: Smagorinsky viscosity @ above vVel
218    C     addDissU    :: zonal momentum tendency from 3-D Smag. viscosity
219    C     addDissV    :: merid momentum tendency from 3-D Smag. viscosity
220          _RL str11(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
221          _RL str22(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
222          _RL str33(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
223          _RL str12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
224          _RL str13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
225          _RL str23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
226          _RL viscAh3d_00(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
227          _RL viscAh3d_12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
228          _RL viscAh3d_13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
229          _RL viscAh3d_23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
230          _RL addDissU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
231          _RL addDissV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
232    #elif ( defined ALLOW_NONHYDROSTATIC )
233          _RL str13(1), str23(1), str33(1)
234          _RL viscAh3d_00(1), viscAh3d_13(1), viscAh3d_23(1)
235    #endif
236    
       INTEGER iMin, iMax  
       INTEGER jMin, jMax  
237        INTEGER bi, bj        INTEGER bi, bj
238        INTEGER i, j        INTEGER i, j
239        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kp1, kUp, kDown
240        LOGICAL BOTTOM_LAYER        INTEGER iMin, iMax
241          INTEGER jMin, jMax
242          PARAMETER( iMin = 0 , iMax = sNx+1 )
243          PARAMETER( jMin = 0 , jMax = sNy+1 )
244    
245    #ifdef ALLOW_DIAGNOSTICS
246          LOGICAL dPhiHydDiagIsOn
247          _RL tmpFac
248    #endif /* ALLOW_DIAGNOSTICS */
249    
250  C---    The algorithm...  C---    The algorithm...
251  C  C
# Line 145  C Line 260  C
260  C       "Calculation of Gs"  C       "Calculation of Gs"
261  C       ===================  C       ===================
262  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
263  C       phiHydysics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         rVel = sum_r ( div. u[n] )  
264  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
265  C         b   = b(rho, theta)  C         b   = b(rho, theta)
266  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
267  C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
268  C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
269  C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
270  C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
271  C  C
272  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
273  C       ================================  C       ================================
# Line 176  C         salt* = salt[n] + dt x ( 3/2 G Line 290  C         salt* = salt[n] + dt x ( 3/2 G
290  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
291  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
292  C---  C---
293    CEOP
294    
295  C--   Set up work arrays with valid (i.e. not NaN) values  #ifdef ALLOW_DEBUG
296  C     These inital values do not alter the numerical results. They        IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
297  C     just ensure that all memory references are to valid floating  #endif
298  C     point numbers. This prevents spurious hardware signals due to  
299  C     uninitialised but inert locations.  #ifdef ALLOW_DIAGNOSTICS
300        DO j=1-OLy,sNy+OLy        dPhiHydDiagIsOn = .FALSE.
301         DO i=1-OLx,sNx+OLx        IF ( useDiagnostics )
302          xA(i,j)      = 0. _d 0       &  dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
303          yA(i,j)      = 0. _d 0       &               .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
304          uTrans(i,j)  = 0. _d 0  #endif
         vTrans(i,j)  = 0. _d 0  
         aTerm(i,j)   = 0. _d 0  
         xTerm(i,j)   = 0. _d 0  
         cTerm(i,j)   = 0. _d 0  
         mTerm(i,j)   = 0. _d 0  
         pTerm(i,j)   = 0. _d 0  
         fZon(i,j)    = 0. _d 0  
         fMer(i,j)    = 0. _d 0  
         DO K=1,Nr  
          phiHyd (i,j,k)  = 0. _d 0  
          K13(i,j,k)  = 0. _d 0  
          K23(i,j,k)  = 0. _d 0  
          K33(i,j,k)  = 0. _d 0  
          KappaRT(i,j,k) = 0. _d 0  
          KappaRS(i,j,k) = 0. _d 0  
         ENDDO  
         rhoKM1 (i,j) = 0. _d 0  
         rhok   (i,j) = 0. _d 0  
         rhoKP1 (i,j) = 0. _d 0  
         rhoTMP (i,j) = 0. _d 0  
         buoyKM1(i,j) = 0. _d 0  
         buoyK  (i,j) = 0. _d 0  
         maskC  (i,j) = 0. _d 0  
        ENDDO  
       ENDDO  
305    
306    C-- Call to routine for calculation of Eliassen-Palm-flux-forced
307    C    U-tendency, if desired:
308    #ifdef INCLUDE_EP_FORCING_CODE
309          CALL CALC_EP_FORCING(myThid)
310    #endif
311    
312    #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
313          CALL DUMMY_IN_DYNAMICS( myTime, myIter, myThid )
314    #endif
315    
316    #ifdef ALLOW_AUTODIFF_TAMC
317    C--   HPF directive to help TAMC
318    CHPF$ INDEPENDENT
319    #endif /* ALLOW_AUTODIFF_TAMC */
320    
321        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
322    
323    #ifdef ALLOW_AUTODIFF_TAMC
324    C--    HPF directive to help TAMC
325    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
326    CHPF$&                  ,phiHydF
327    CHPF$&                  ,KappaRU,KappaRV
328    CHPF$&                  )
329    #endif /* ALLOW_AUTODIFF_TAMC */
330    
331         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
332    
333  C--     Set up work arrays that need valid initial values  #ifdef ALLOW_AUTODIFF_TAMC
334              act1 = bi - myBxLo(myThid)
335              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
336              act2 = bj - myByLo(myThid)
337              max2 = myByHi(myThid) - myByLo(myThid) + 1
338              act3 = myThid - 1
339              max3 = nTx*nTy
340              act4 = ikey_dynamics - 1
341              idynkey = (act1 + 1) + act2*max1
342         &                      + act3*max1*max2
343         &                      + act4*max1*max2*max3
344    #endif /* ALLOW_AUTODIFF_TAMC */
345    
346    C--   Set up work arrays with valid (i.e. not NaN) values
347    C     These initial values do not alter the numerical results. They
348    C     just ensure that all memory references are to valid floating
349    C     point numbers. This prevents spurious hardware signals due to
350    C     uninitialised but inert locations.
351    
352    #ifdef ALLOW_AUTODIFF
353            DO k=1,Nr
354             DO j=1-OLy,sNy+OLy
355              DO i=1-OLx,sNx+OLx
356    c--   need some re-initialisation here to break dependencies
357               gU(i,j,k,bi,bj) = 0. _d 0
358               gV(i,j,k,bi,bj) = 0. _d 0
359              ENDDO
360             ENDDO
361            ENDDO
362    #endif /* ALLOW_AUTODIFF */
363          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
364           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
365            rTrans(i,j)   = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
366            rVel  (i,j,1) = 0. _d 0            fVerU  (i,j,2) = 0. _d 0
367            rVel  (i,j,2) = 0. _d 0            fVerV  (i,j,1) = 0. _d 0
368            fVerT (i,j,1) = 0. _d 0            fVerV  (i,j,2) = 0. _d 0
369            fVerT (i,j,2) = 0. _d 0            phiHydF (i,j)  = 0. _d 0
370            fVerS (i,j,1) = 0. _d 0            phiHydC (i,j)  = 0. _d 0
371            fVerS (i,j,2) = 0. _d 0  #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
372            fVerU (i,j,1) = 0. _d 0            dPhiHydX(i,j)  = 0. _d 0
373            fVerU (i,j,2) = 0. _d 0            dPhiHydY(i,j)  = 0. _d 0
374            fVerV (i,j,1) = 0. _d 0  #endif
375            fVerV (i,j,2) = 0. _d 0            phiSurfX(i,j)  = 0. _d 0
376            phiHyd(i,j,1) = 0. _d 0            phiSurfY(i,j)  = 0. _d 0
377            K13   (i,j,1) = 0. _d 0            guDissip(i,j)  = 0. _d 0
378            K23   (i,j,1) = 0. _d 0            gvDissip(i,j)  = 0. _d 0
379            K33   (i,j,1) = 0. _d 0  #ifdef ALLOW_AUTODIFF
380            KapGM (i,j)   = GMkbackground            phiHydLow(i,j,bi,bj) = 0. _d 0
381    # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
382    #  ifndef DISABLE_RSTAR_CODE
383              dWtransC(i,j,bi,bj) = 0. _d 0
384              dWtransU(i,j,bi,bj) = 0. _d 0
385              dWtransV(i,j,bi,bj) = 0. _d 0
386    #  endif
387    # endif
388    #endif /* ALLOW_AUTODIFF */
389           ENDDO           ENDDO
390          ENDDO          ENDDO
391    
392          iMin = 1-OLx+1  C--     Start computation of dynamics
393          iMax = sNx+OLx  
394          jMin = 1-OLy+1  #ifdef ALLOW_AUTODIFF_TAMC
395          jMax = sNy+OLy  CADJ STORE wVel (:,:,:,bi,bj) =
396    CADJ &     comlev1_bibj, key=idynkey, byte=isbyte
397    #endif /* ALLOW_AUTODIFF_TAMC */
398          K = 1  
399          BOTTOM_LAYER = K .EQ. Nr  C--     Explicit part of the Surface Potential Gradient (add in TIMESTEP)
400    C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
401  #ifdef DO_PIPELINED_CORRECTION_STEP          IF (implicSurfPress.NE.1.) THEN
402  C--     Calculate gradient of surface pressure            CALL CALC_GRAD_PHI_SURF(
403          CALL CALC_GRAD_ETA_SURF(       I         bi,bj,iMin,iMax,jMin,jMax,
404       I       bi,bj,iMin,iMax,jMin,jMax,       I         etaN,
405       O       etaSurfX,etaSurfY,       O         phiSurfX,phiSurfY,
406       I       myThid)       I         myThid )
 C--     Update fields in top level according to tendency terms  
         CALL CORRECTION_STEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,  
      I       etaSurfX,etaSurfY,myTime,myThid)  
 #ifdef ALLOW_OBCS  
         IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K, myThid )  
 #endif  
         IF ( .NOT. BOTTOM_LAYER ) THEN  
 C--      Update fields in layer below according to tendency terms  
          CALL CORRECTION_STEP(  
      I        bi,bj,iMin,iMax,jMin,jMax,K+1,  
      I        etaSurfX,etaSurfY,myTime,myThid)  
 #ifdef ALLOW_OBCS  
          IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )  
 #endif  
407          ENDIF          ENDIF
 #endif  
 C--     Density of 1st level (below W(1)) reference to level 1  
 #ifdef  INCLUDE_FIND_RHO_CALL  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  
      O     rhoKm1,  
      I     myThid )  
 #endif  
408    
409          IF (       (.NOT. BOTTOM_LAYER)  #ifdef ALLOW_AUTODIFF_TAMC
410    CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
411    CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
412  #ifdef ALLOW_KPP  #ifdef ALLOW_KPP
413       &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing  CADJ STORE KPPviscAz (:,:,:,bi,bj)
414  #endif  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
415       &     ) THEN  #endif /* ALLOW_KPP */
416  C--      Check static stability with layer below  #endif /* ALLOW_AUTODIFF_TAMC */
417  C--      and mix as needed.  
418  #ifdef  INCLUDE_FIND_RHO_CALL  #ifndef ALLOW_AUTODIFF
419           CALL FIND_RHO(          IF ( .NOT.momViscosity ) THEN
420       I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,  #endif
421       O      rhoKp1,            DO k=1,Nr
422       I      myThid )             DO j=1-OLy,sNy+OLy
423  #endif              DO i=1-OLx,sNx+OLx
424  #ifdef  INCLUDE_CONVECT_CALL               KappaRU(i,j,k) = 0. _d 0
425           CALL CONVECT(               KappaRV(i,j,k) = 0. _d 0
426       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,              ENDDO
427       I       myTime,myIter,myThid)             ENDDO
428  #endif            ENDDO
429  C--      Recompute density after mixing  #ifndef ALLOW_AUTODIFF
 #ifdef  INCLUDE_FIND_RHO_CALL  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  
      O      rhoKm1,  
      I      myThid )  
 #endif  
430          ENDIF          ENDIF
 C--     Calculate buoyancy  
         CALL CALC_BUOYANCY(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,  
      O      buoyKm1,  
      I      myThid )  
 C--     Integrate hydrostatic balance for phiHyd with BC of  
 C--     phiHyd(z=0)=0  
         CALL CALC_PHI_HYD(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,  
      U      phiHyd,  
      I      myThid )  
   
         DO K=2,Nr  
          BOTTOM_LAYER = K .EQ. Nr  
 #ifdef DO_PIPELINED_CORRECTION_STEP  
          IF ( .NOT. BOTTOM_LAYER ) THEN  
 C--       Update fields in layer below according to tendency terms  
           CALL CORRECTION_STEP(  
      I         bi,bj,iMin,iMax,jMin,jMax,K+1,  
      I         etaSurfX,etaSurfY,myTime,myThid)  
 #ifdef ALLOW_OBCS  
           IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )  
 #endif  
          ENDIF  
 #endif  
 C--      Density of K level (below W(K)) reference to K level  
 #ifdef  INCLUDE_FIND_RHO_CALL  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
      O      rhoK,  
      I      myThid )  
 #endif  
          IF (       (.NOT. BOTTOM_LAYER)  
 #ifdef ALLOW_KPP  
      &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing  
 #endif  
      &      ) THEN  
 C--       Check static stability with layer below and mix as needed.  
 C--       Density of K+1 level (below W(K+1)) reference to K level.  
 #ifdef  INCLUDE_FIND_RHO_CALL  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,  
      O       rhoKp1,  
      I       myThid )  
 #endif  
 #ifdef  INCLUDE_CONVECT_CALL  
           CALL CONVECT(  
      I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,  
      I        myTime,myIter,myThid)  
 #endif  
 C--       Recompute density after mixing  
 #ifdef  INCLUDE_FIND_RHO_CALL  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  
      O       rhoK,  
      I       myThid )  
 #endif  
          ENDIF  
 C--      Calculate buoyancy  
          CALL CALC_BUOYANCY(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,  
      O       buoyK,  
      I       myThid )  
 C--      Integrate hydrostatic balance for phiHyd with BC of  
 C--      phiHyd(z=0)=0  
          CALL CALC_PHI_HYD(  
      I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,  
      U        phiHyd,  
      I        myThid )  
 C--      Calculate iso-neutral slopes for the GM/Redi parameterisation  
 #ifdef  INCLUDE_FIND_RHO_CALL  
          CALL FIND_RHO(  
      I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,  
      O        rhoTmp,  
      I        myThid )  
 #endif  
 #ifdef  INCLUDE_CALC_ISOSLOPES_CALL  
          CALL CALC_ISOSLOPES(  
      I        bi, bj, iMin, iMax, jMin, jMax, K,  
      I        rhoKm1, rhoK, rhotmp,  
      O        K13, K23, K33, KapGM,  
      I        myThid )  
 #endif  
          DO J=jMin,jMax  
           DO I=iMin,iMax  
 #ifdef  INCLUDE_FIND_RHO_CALL  
            rhoKm1 (I,J) = rhoK(I,J)  
431  #endif  #endif
432             buoyKm1(I,J) = buoyK(I,J)  #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
433            ENDDO  C--     Calculate the total vertical viscosity
434           ENDDO          IF ( momViscosity ) THEN
435          ENDDO ! K            CALL CALC_VISCOSITY(
436         I            bi,bj, iMin,iMax,jMin,jMax,
437         O            KappaRU, KappaRV,
438         I            myThid )
439            ENDIF
440    #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
441    
442  #ifdef ALLOW_KPP  #ifdef ALLOW_SMAG_3D
443  C--     Compute KPP mixing coefficients          IF ( useSmag3D ) THEN
444          IF (usingKPPmixing) THEN            CALL MOM_CALC_3D_STRAIN(
445           CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'       O         str11, str22, str33, str12, str13, str23,
446       I          , myThid)       I         bi, bj, myThid )
          CALL KVMIX(  
      I               bi, bj, myTime, myThid )  
          CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'  
      I        , myThid)  
447          ENDIF          ENDIF
448  #endif  #endif /* ALLOW_SMAG_3D */
449    
450          DO K = Nr, 1, -1  #ifdef ALLOW_AUTODIFF_TAMC
451    CADJ STORE KappaRU(:,:,:)
452    CADJ &     = comlev1_bibj, key=idynkey, byte=isbyte
453    CADJ STORE KappaRV(:,:,:)
454    CADJ &     = comlev1_bibj, key=idynkey, byte=isbyte
455    #endif /* ALLOW_AUTODIFF_TAMC */
456    
457           kM1  =max(1,k-1)   ! Points to level above k (=k-1)  #ifdef ALLOW_OBCS
458           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above  C--   For Stevens boundary conditions velocities need to be extrapolated
459           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer  C     (copied) to a narrow strip outside the domain
460           iMin = 1-OLx+2          IF ( useOBCS ) THEN
461           iMax = sNx+OLx-1            CALL OBCS_COPY_UV_N(
462           jMin = 1-OLy+2       U         uVel(1-OLx,1-OLy,1,bi,bj),
463           jMax = sNy+OLy-1       U         vVel(1-OLx,1-OLy,1,bi,bj),
464         I         Nr, bi, bj, myThid )
465  C--      Get temporary terms used by tendency routines          ENDIF
466           CALL CALC_COMMON_FACTORS (  #endif /* ALLOW_OBCS */
467       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
468       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,  #ifdef ALLOW_EDDYPSI
469       I        myThid)          CALL CALC_EDDY_STRESS(bi,bj,myThid)
 #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL  
 C--      Calculate the total vertical diffusivity  
          CALL CALC_DIFFUSIVITY(  
      I        bi,bj,iMin,iMax,jMin,jMax,K,  
      I        maskC,maskUp,KapGM,K33,  
      O        KappaRT,KappaRS,KappaRU,KappaRV,  
      I        myThid)  
470  #endif  #endif
471  C--      Calculate accelerations in the momentum equations  
472           IF ( momStepping ) THEN  C--     Start of dynamics loop
473            CALL CALC_MOM_RHS(          DO k=1,Nr
474       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
475       I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,  C--       km1    Points to level above k (=k-1)
476       I         phiHyd,KappaRU,KappaRV,  C--       kup    Cycles through 1,2 to point to layer above
477       U         aTerm,xTerm,cTerm,mTerm,pTerm,  C--       kDown  Cycles through 2,1 to point to current layer
478       U         fZon, fMer, fVerU, fVerV,  
479       I         myTime, myThid)            km1  = MAX(1,k-1)
480           ENDIF            kp1  = MIN(k+1,Nr)
481  C--      Calculate active tracer tendencies            kup  = 1+MOD(k+1,2)
482           IF ( tempStepping ) THEN            kDown= 1+MOD(k,2)
483            CALL CALC_GT(  
484       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  #ifdef ALLOW_AUTODIFF_TAMC
485       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,           kkey = (idynkey-1)*Nr + k
486       I         K13,K23,KappaRT,KapGM,  CADJ STORE totPhiHyd (:,:,k,bi,bj)
487       U         aTerm,xTerm,fZon,fMer,fVerT,  CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
488       I         myTime, myThid)  CADJ STORE phiHydLow (:,:,bi,bj)
489    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
490    CADJ STORE theta (:,:,k,bi,bj)
491    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
492    CADJ STORE salt  (:,:,k,bi,bj)
493    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
494    CADJ STORE gT(:,:,k,bi,bj)
495    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
496    CADJ STORE gS(:,:,k,bi,bj)
497    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
498    # ifdef NONLIN_FRSURF
499    cph-test
500    CADJ STORE  phiHydC (:,:)
501    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
502    CADJ STORE  phiHydF (:,:)
503    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
504    CADJ STORE gU(:,:,k,bi,bj)
505    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
506    CADJ STORE gV(:,:,k,bi,bj)
507    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
508    #  ifndef ALLOW_ADAMSBASHFORTH_3
509    CADJ STORE guNm1(:,:,k,bi,bj)
510    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
511    CADJ STORE gvNm1(:,:,k,bi,bj)
512    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
513    #  else
514    CADJ STORE guNm(:,:,k,bi,bj,1)
515    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
516    CADJ STORE guNm(:,:,k,bi,bj,2)
517    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
518    CADJ STORE gvNm(:,:,k,bi,bj,1)
519    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
520    CADJ STORE gvNm(:,:,k,bi,bj,2)
521    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
522    #  endif
523    #  ifdef ALLOW_CD_CODE
524    CADJ STORE uNM1(:,:,k,bi,bj)
525    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
526    CADJ STORE vNM1(:,:,k,bi,bj)
527    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
528    CADJ STORE uVelD(:,:,k,bi,bj)
529    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
530    CADJ STORE vVelD(:,:,k,bi,bj)
531    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
532    #  endif
533    # endif /* NONLIN_FRSURF */
534    #endif /* ALLOW_AUTODIFF_TAMC */
535    
536    C--      Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
537             IF ( implicitIntGravWave ) THEN
538               CALL CALC_PHI_HYD(
539         I        bi,bj,iMin,iMax,jMin,jMax,k,
540         I        gT, gS,
541         U        phiHydF,
542         O        phiHydC, dPhiHydX, dPhiHydY,
543         I        myTime, myIter, myThid )
544             ELSE
545               CALL CALC_PHI_HYD(
546         I        bi,bj,iMin,iMax,jMin,jMax,k,
547         I        theta, salt,
548         U        phiHydF,
549         O        phiHydC, dPhiHydX, dPhiHydY,
550         I        myTime, myIter, myThid )
551           ENDIF           ENDIF
552           IF ( saltStepping ) THEN  #ifdef ALLOW_DIAGNOSTICS
553            CALL CALC_GS(           IF ( dPhiHydDiagIsOn ) THEN
554       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,             tmpFac = -1. _d 0
555       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,             CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
556       I         K13,K23,KappaRS,KapGM,       &                           'Um_dPHdx', k, 1, 2, bi, bj, myThid )
557       U         aTerm,xTerm,fZon,fMer,fVerS,             CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
558       I         myTime, myThid)       &                           'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
559           ENDIF           ENDIF
560  C--      Prediction step (step forward all model variables)  #endif /* ALLOW_DIAGNOSTICS */
561           CALL TIMESTEP(  
562       I       bi,bj,iMin,iMax,jMin,jMax,K,  C--      Calculate accelerations in the momentum equations (gU, gV, ...)
563       I       myThid)  C        and step forward storing the result in gU, gV, etc...
564  #ifdef ALLOW_OBCS           IF ( momStepping ) THEN
565  C--      Apply open boundary conditions  #ifdef ALLOW_AUTODIFF
566           IF (openBoundaries) CALL APPLY_OBCS2( bi, bj, K, myThid )             DO j=1-OLy,sNy+OLy
567  #endif              DO i=1-OLx,sNx+OLx
568  C--      Freeze water                guDissip(i,j)  = 0. _d 0
569           IF (allowFreezing)                gvDissip(i,j)  = 0. _d 0
570       &   CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )              ENDDO
571  C--      Diagnose barotropic divergence of predicted fields             ENDDO
572           CALL CALC_DIV_GHAT(  #endif /* ALLOW_AUTODIFF */
573       I       bi,bj,iMin,iMax,jMin,jMax,K,  #ifdef ALLOW_AUTODIFF_TAMC
574       I       xA,yA,  # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
575       I       myThid)  #  ifndef DISABLE_RSTAR_CODE
576    CADJ STORE dWtransC(:,:,bi,bj)
577  C--      Cumulative diagnostic calculations (ie. time-averaging)  CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
578  #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE  CADJ STORE dWtransU(:,:,bi,bj)
579           IF (taveFreq.GT.0.) THEN  CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
580            CALL DO_TIME_AVERAGES(  CADJ STORE dWtransV(:,:,bi,bj)
581       I                           myTime, myIter, bi, bj, K, kUp, kDown,  CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
582       I                           K13, K23, rVel, KapGM,  #  endif
583       I                           myThid )  # endif /* NONLIN_FRSURF and ALLOW_MOM_FLUXFORM */
584    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
585    CADJ STORE fVerU(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
586    CADJ STORE fVerV(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
587    # endif
588    #endif /* ALLOW_AUTODIFF_TAMC */
589               IF (.NOT. vectorInvariantMomentum) THEN
590    #ifdef ALLOW_MOM_FLUXFORM
591                  CALL MOM_FLUXFORM(
592         I         bi,bj,k,iMin,iMax,jMin,jMax,
593         I         KappaRU, KappaRV,
594         U         fVerU(1-OLx,1-OLy,kUp),   fVerV(1-OLx,1-OLy,kUp),
595         O         fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
596         O         guDissip, gvDissip,
597         I         myTime, myIter, myThid)
598    #endif
599               ELSE
600    #ifdef ALLOW_MOM_VECINV
601                 CALL MOM_VECINV(
602         I         bi,bj,k,iMin,iMax,jMin,jMax,
603         I         KappaRU, KappaRV,
604         I         fVerU(1-OLx,1-OLy,kUp),   fVerV(1-OLx,1-OLy,kUp),
605         O         fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
606         O         guDissip, gvDissip,
607         I         myTime, myIter, myThid)
608    #endif
609               ENDIF
610    
611    #ifdef ALLOW_SMAG_3D
612               IF ( useSmag3D ) THEN
613                 CALL MOM_CALC_SMAG_3D(
614         I         str11, str22, str33, str12, str13, str23,
615         O         viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
616         I         smag3D_hLsC, smag3D_hLsW, smag3D_hLsS, smag3D_hLsZ,
617         I         k, bi, bj, myThid )
618                 CALL MOM_UV_SMAG_3D(
619         I         str11, str22, str12, str13, str23,
620         I         viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
621         O         addDissU, addDissV,
622         I         iMin,iMax,jMin,jMax, k, bi, bj, myThid )
623                 DO j= jMin,jMax
624                  DO i= iMin,iMax
625                   guDissip(i,j) = guDissip(i,j) + addDissU(i,j)
626                   gvDissip(i,j) = gvDissip(i,j) + addDissV(i,j)
627                  ENDDO
628                 ENDDO
629               ENDIF
630    #endif /* ALLOW_SMAG_3D */
631    
632               CALL TIMESTEP(
633         I         bi,bj,iMin,iMax,jMin,jMax,k,
634         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
635         I         guDissip, gvDissip,
636         I         myTime, myIter, myThid)
637    
638           ENDIF           ENDIF
 #endif  
639    
640          ENDDO ! K  C--     end of dynamics k loop (1:Nr)
641            ENDDO
642    
643  C--     Implicit diffusion  C--     Implicit Vertical advection & viscosity
644          IF (implicitDiffusion) THEN  #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
645           IF (tempStepping) CALL IMPLDIFF(       defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
646       I         bi, bj, iMin, iMax, jMin, jMax,          IF ( momImplVertAdv ) THEN
647       I         deltaTtracer, KappaRT,recip_HFacC,            CALL MOM_U_IMPLICIT_R( kappaRU,
648       U         gTNm1,       I                           bi, bj, myTime, myIter, myThid )
649       I         myThid )            CALL MOM_V_IMPLICIT_R( kappaRV,
650           IF (saltStepping) CALL IMPLDIFF(       I                           bi, bj, myTime, myIter, myThid )
651       I         bi, bj, iMin, iMax, jMin, jMax,          ELSEIF ( implicitViscosity ) THEN
652       I         deltaTtracer, KappaRS,recip_HFacC,  #else /* INCLUDE_IMPLVERTADV_CODE */
653       U         gSNm1,          IF     ( implicitViscosity ) THEN
654       I         myThid )  #endif /* INCLUDE_IMPLVERTADV_CODE */
655           IF (momStepping) THEN  #ifdef    ALLOW_AUTODIFF_TAMC
656    CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
657    #endif    /* ALLOW_AUTODIFF_TAMC */
658            CALL IMPLDIFF(            CALL IMPLDIFF(
659       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
660       I         deltaTmom, KappaRU,recip_HFacW,       I         -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
661       U         gUNm1,       U         gU,
662       I         myThid )       I         myThid )
663    #ifdef    ALLOW_AUTODIFF_TAMC
664    CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
665    #endif    /* ALLOW_AUTODIFF_TAMC */
666            CALL IMPLDIFF(            CALL IMPLDIFF(
667       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
668       I         deltaTmom, KappaRV,recip_HFacS,       I         -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
669       U         gVNm1,       U         gV,
670       I         myThid )       I         myThid )
671  #ifdef INCLUDE_CD_CODE          ENDIF
672    
673    #ifdef ALLOW_OBCS
674    C--      Apply open boundary conditions
675            IF ( useOBCS ) THEN
676    C--      but first save intermediate velocities to be used in the
677    C        next time step for the Stevens boundary conditions
678              CALL OBCS_SAVE_UV_N(
679         I        bi, bj, iMin, iMax, jMin, jMax, 0,
680         I        gU, gV, myThid )
681              CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
682            ENDIF
683    #endif /* ALLOW_OBCS */
684    
685    #ifdef    ALLOW_CD_CODE
686            IF (implicitViscosity.AND.useCDscheme) THEN
687    #ifdef    ALLOW_AUTODIFF_TAMC
688    CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
689    #endif    /* ALLOW_AUTODIFF_TAMC */
690            CALL IMPLDIFF(            CALL IMPLDIFF(
691       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
692       I         deltaTmom, KappaRU,recip_HFacW,       I         0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
693       U         vVelD,       U         vVelD,
694       I         myThid )       I         myThid )
695    #ifdef    ALLOW_AUTODIFF_TAMC
696    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
697    #endif    /* ALLOW_AUTODIFF_TAMC */
698            CALL IMPLDIFF(            CALL IMPLDIFF(
699       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
700       I         deltaTmom, KappaRV,recip_HFacS,       I         0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
701       U         uVelD,       U         uVelD,
702       I         myThid )       I         myThid )
703            ENDIF
704    #endif    /* ALLOW_CD_CODE */
705    C--     End implicit Vertical advection & viscosity
706    
707    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
708    
709    #ifdef ALLOW_NONHYDROSTATIC
710    C--   Step forward W field in N-H algorithm
711            IF ( nonHydrostatic ) THEN
712    #ifdef ALLOW_DEBUG
713             IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
714    #endif
715             CALL TIMER_START('CALC_GW          [DYNAMICS]',myThid)
716             CALL CALC_GW(
717         I                 bi,bj, KappaRU, KappaRV,
718         I                 str13, str23, str33,
719         I                 viscAh3d_00, viscAh3d_13, viscAh3d_23,
720         I                 myTime, myIter, myThid )
721            ENDIF
722            IF ( nonHydrostatic.OR.implicitIntGravWave )
723         &   CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
724            IF ( nonHydrostatic )
725         &   CALL TIMER_STOP ('CALC_GW          [DYNAMICS]',myThid)
726  #endif  #endif
727           ENDIF ! momStepping  
728          ENDIF ! implicitDiffusion  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
729    
730    C-    end of bi,bj loops
731         ENDDO         ENDDO
732        ENDDO        ENDDO
733    
734  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),  #ifdef ALLOW_OBCS
735  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))        IF (useOBCS) THEN
736  C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),          CALL OBCS_EXCHANGES( myThid )
737  C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)        ENDIF
738  C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),  #endif
739  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)  
740  C     write(0,*) 'dynamics: rVel(1) ',  Cml(
741  C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),  C     In order to compare the variance of phiHydLow of a p/z-coordinate
742  C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)  C     run with etaH of a z/p-coordinate run the drift of phiHydLow
743  C     write(0,*) 'dynamics: rVel(2) ',  C     has to be removed by something like the following subroutine:
744  C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),  C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
745  C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)  C     &                     'phiHydLow', myTime, myThid )
746  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  Cml)
747  cblk &                           maxval(K13(1:sNx,1:sNy,:))  
748  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  #ifdef ALLOW_DIAGNOSTICS
749  cblk &                           maxval(K23(1:sNx,1:sNy,:))        IF ( useDiagnostics ) THEN
750  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  
751  cblk &                           maxval(K33(1:sNx,1:sNy,:))         CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD  ',0,Nr,0,1,1,myThid)
752  C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),         CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT  ',0, 1,0,1,1,myThid)
753  C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))  
754  C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),         tmpFac = 1. _d 0
755  C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))         CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
756  C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),       &                                 'PHIHYDSQ',0,Nr,0,1,1,myThid)
757  C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))  
758  C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),         CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
759  C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))       &                                 'PHIBOTSQ',0, 1,0,1,1,myThid)
760  C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),  
761  C    &                           maxval(phiHyd/(Gravity*Rhonil))        ENDIF
762  C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,  #endif /* ALLOW_DIAGNOSTICS */
763  C    &Nr, 1, myThid )  
764  C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,  #ifdef ALLOW_DEBUG
765  C    &Nr, 1, myThid )        IF ( debugLevel .GE. debLevD ) THEN
766  C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,         CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
767  C    &Nr, 1, myThid )         CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
768  C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,         CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
769  C    &Nr, 1, myThid )         CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
770  C     CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,         CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
771  C    &Nr, 1, myThid )         CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
772           CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
773           CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
774           CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
775           CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
776    #ifndef ALLOW_ADAMSBASHFORTH_3
777           CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
778           CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
779           CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
780           CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
781    #endif
782          ENDIF
783    #endif
784    
785    #ifdef DYNAMICS_GUGV_EXCH_CHECK
786    C- jmc: For safety checking only: This Exchange here should not change
787    C       the solution. If solution changes, it means something is wrong,
788    C       but it does not mean that it is less wrong with this exchange.
789          IF ( debugLevel .GE. debLevE ) THEN
790           CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
791          ENDIF
792    #endif
793    
794    #ifdef ALLOW_DEBUG
795          IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
796    #endif
797    
798        RETURN        RETURN
799        END        END

Legend:
Removed from v.1.43  
changed lines
  Added in v.1.170

  ViewVC Help
Powered by ViewVC 1.1.22