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

Legend:
Removed from v.1.47  
changed lines
  Added in v.1.167

  ViewVC Help
Powered by ViewVC 1.1.22