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

Legend:
Removed from v.1.20  
changed lines
  Added in v.1.176

  ViewVC Help
Powered by ViewVC 1.1.22