/[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.26 by cnh, Wed Aug 19 16:20:49 1998 UTC revision 1.129 by heimbach, Sat Feb 25 16:20:01 2006 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    #undef DYNAMICS_GUGV_EXCH_CHECK
7    
8    CBOP
9    C     !ROUTINE: DYNAMICS
10    C     !INTERFACE:
11        SUBROUTINE DYNAMICS(myTime, myIter, myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
12  C     /==========================================================\  C     !DESCRIPTION: \bv
13  C     | SUBROUTINE DYNAMICS                                      |  C     *==========================================================*
14  C     | o Controlling routine for the explicit part of the model |  C     | SUBROUTINE DYNAMICS                                      
15  C     |   dynamics.                                              |  C     | o Controlling routine for the explicit part of the model  
16  C     |==========================================================|  C     |   dynamics.                                              
17  C     | This routine evaluates the "dynamics" terms for each     |  C     *==========================================================*
18  C     | block of ocean in turn. Because the blocks of ocean have |  C     | This routine evaluates the "dynamics" terms for each      
19  C     | overlap regions they are independent of one another.     |  C     | block of ocean in turn. Because the blocks of ocean have  
20  C     | If terms involving lateral integrals are needed in this  |  C     | overlap regions they are independent of one another.      
21  C     | routine care will be needed. Similarly finite-difference |  C     | If terms involving lateral integrals are needed in this  
22  C     | operations with stencils wider than the overlap region   |  C     | routine care will be needed. Similarly finite-difference  
23  C     | require special consideration.                           |  C     | operations with stencils wider than the overlap region    
24  C     | Notes                                                    |  C     | require special consideration.                            
25  C     | =====                                                    |  C     | The algorithm...
26  C     | C*P* comments indicating place holders for which code is |  C     |
27  C     |      presently being developed.                          |  C     | "Correction Step"
28  C     \==========================================================/  C     | =================
29    C     | Here we update the horizontal velocities with the surface
30    C     | pressure such that the resulting flow is either consistent
31    C     | with the free-surface evolution or the rigid-lid:
32    C     |   U[n] = U* + dt x d/dx P
33    C     |   V[n] = V* + dt x d/dy P
34    C     |   W[n] = W* + dt x d/dz P  (NH mode)
35    C     |
36    C     | "Calculation of Gs"
37    C     | ===================
38    C     | This is where all the accelerations and tendencies (ie.
39    C     | physics, parameterizations etc...) are calculated
40    C     |   rho = rho ( theta[n], salt[n] )
41    C     |   b   = b(rho, theta)
42    C     |   K31 = K31 ( rho )
43    C     |   Gu[n] = Gu( u[n], v[n], wVel, b, ... )
44    C     |   Gv[n] = Gv( u[n], v[n], wVel, b, ... )
45    C     |   Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
46    C     |   Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
47    C     |
48    C     | "Time-stepping" or "Prediction"
49    C     | ================================
50    C     | The models variables are stepped forward with the appropriate
51    C     | time-stepping scheme (currently we use Adams-Bashforth II)
52    C     | - For momentum, the result is always *only* a "prediction"
53    C     | in that the flow may be divergent and will be "corrected"
54    C     | later with a surface pressure gradient.
55    C     | - Normally for tracers the result is the new field at time
56    C     | level [n+1} *BUT* in the case of implicit diffusion the result
57    C     | is also *only* a prediction.
58    C     | - We denote "predictors" with an asterisk (*).
59    C     |   U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
60    C     |   V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
61    C     |   theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62    C     |   salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63    C     | With implicit diffusion:
64    C     |   theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65    C     |   salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66    C     |   (1 + dt * K * d_zz) theta[n] = theta*
67    C     |   (1 + dt * K * d_zz) salt[n] = salt*
68    C     |
69    C     *==========================================================*
70    C     \ev
71    C     !USES:
72          IMPLICIT NONE
73  C     == Global variables ===  C     == Global variables ===
74  #include "SIZE.h"  #include "SIZE.h"
75  #include "EEPARAMS.h"  #include "EEPARAMS.h"
 #include "CG2D.h"  
76  #include "PARAMS.h"  #include "PARAMS.h"
77  #include "DYNVARS.h"  #include "DYNVARS.h"
78    #ifdef ALLOW_CD_CODE
79    #include "CD_CODE_VARS.h"
80    #endif
81    #include "GRID.h"
82    #ifdef ALLOW_AUTODIFF_TAMC
83    # include "tamc.h"
84    # include "tamc_keys.h"
85    # include "FFIELDS.h"
86    # include "EOS.h"
87    # ifdef ALLOW_KPP
88    #  include "KPP.h"
89    # endif
90    #endif /* ALLOW_AUTODIFF_TAMC */
91    
92    C     !CALLING SEQUENCE:
93    C     DYNAMICS()
94    C      |
95    C      |-- CALC_EP_FORCING
96    C      |
97    C      |-- CALC_GRAD_PHI_SURF
98    C      |
99    C      |-- CALC_VISCOSITY
100    C      |
101    C      |-- CALC_PHI_HYD  
102    C      |
103    C      |-- MOM_FLUXFORM  
104    C      |
105    C      |-- MOM_VECINV    
106    C      |
107    C      |-- TIMESTEP      
108    C      |
109    C      |-- OBCS_APPLY_UV
110    C      |
111    C      |-- MOM_U_IMPLICIT_R      
112    C      |-- MOM_V_IMPLICIT_R      
113    C      |
114    C      |-- IMPLDIFF      
115    C      |
116    C      |-- OBCS_APPLY_UV
117    C      |
118    C      |-- CALC_GW
119    C      |
120    C      |-- DIAGNOSTICS_FILL
121    C      |-- DEBUG_STATS_RL
122    
123    C     !INPUT/OUTPUT PARAMETERS:
124  C     == Routine arguments ==  C     == Routine arguments ==
125  C     myTime - Current time in simulation  C     myTime - Current time in simulation
126  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
127  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
128        _RL myTime        _RL myTime
129        INTEGER myIter        INTEGER myIter
130          INTEGER myThid
131    
132    C     !LOCAL VARIABLES:
133  C     == Local variables  C     == Local variables
134  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[UV]               o fVer: Vertical flux term - note fVer
135  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C                                    is "pipelined" in the vertical
136  C     wVel                     o uTrans: Zonal transport  C                                    so we need an fVer for each
137  C                              o vTrans: Meridional transport  C                                    variable.
138  C                              o wTrans: Vertical transport  C     phiHydC    :: hydrostatic potential anomaly at cell center
139  C                              o wVel:   Vertical velocity at upper and lower  C                   In z coords phiHyd is the hydrostatic potential
140  C                                        cell faces.  C                      (=pressure/rho0) anomaly
141  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C                   In p coords phiHyd is the geopotential height anomaly.
142  C                              o maskUp: land/water mask for W points  C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
143  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
144  C     mTerm, pTerm,            tendency equations.  C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
145  C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  C     phiSurfY             or geopotential (atmos) in X and Y direction
146  C                              o xTerm: Mixing term  C     guDissip   :: dissipation tendency (all explicit terms), u component
147  C                              o cTerm: Coriolis term  C     gvDissip   :: dissipation tendency (all explicit terms), v component
148  C                              o mTerm: Metric term  C     iMin, iMax     - Ranges and sub-block indices on which calculations
149  C                              o pTerm: Pressure term  C     jMin, jMax       are applied.
 C                              o fZon: Zonal flux term  
 C                              o fMer: Meridional flux term  
 C                              o fVer: Vertical flux term - note fVer  
 C                                      is "pipelined" in the vertical  
 C                                      so we need an fVer for each  
 C                                      variable.  
 C     rhoK, rhoKM1   - Density at current level, level above and level below.  
 C     rhoKP1                                                                    
 C     buoyK, buoyKM1 - Buoyancy at current level and level above.  
 C     phiHyd         - Hydrostatic part of the potential phi.  
 C                      In z coords phiHyd is the hydrostatic pressure anomaly  
 C                      In p coords phiHyd is the geopotential surface height anomaly.  
 C     iMin, iMax - Ranges and sub-block indices on which calculations  
 C     jMin, jMax   are applied.  
150  C     bi, bj  C     bi, bj
151  C     k, kUp, kDown, kM1 - Index for layer above and below. kUp and kDown  C     k, kup,        - Index for layer above and below. kup and kDown
152  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
153  C                          into fVerTerm  C                      index into fVerTerm.
154        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
155        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
156        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
158        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
159        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
160        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
163        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
164        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
165        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
       _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)  
       _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rhok  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL buoyKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)  
       _RL KappaZS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)  
166    
167        INTEGER iMin, iMax        INTEGER iMin, iMax
168        INTEGER jMin, jMax        INTEGER jMin, jMax
169        INTEGER bi, bj        INTEGER bi, bj
170        INTEGER i, j        INTEGER i, j
171        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kp1, kup, kDown
       LOGICAL BOTTOM_LAYER  
172    
173    #ifdef ALLOW_DIAGNOSTICS
174          _RL tmpFac
175    #endif /* ALLOW_DIAGNOSTICS */
176    
177    
178  C---    The algorithm...  C---    The algorithm...
179  C  C
180  C       "Correction Step"  C       "Correction Step"
# Line 127  C       "Calculation of Gs" Line 189  C       "Calculation of Gs"
189  C       ===================  C       ===================
190  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
191  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         w = sum_z ( div. u[n] )  
192  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
193    C         b   = b(rho, theta)
194  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
195  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
196  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
197  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
198  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
199  C  C
200  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
201  C       ================================  C       ================================
# Line 156  C         salt* = salt[n] + dt x ( 3/2 G Line 218  C         salt* = salt[n] + dt x ( 3/2 G
218  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
219  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
220  C---  C---
221    CEOP
222    
223    #ifdef ALLOW_DEBUG
224          IF ( debugLevel .GE. debLevB )
225         &   CALL DEBUG_ENTER( 'DYNAMICS', myThid )
226    #endif
227    
228    C-- Call to routine for calculation of
229    C   Eliassen-Palm-flux-forced U-tendency,
230    C   if desired:
231    #ifdef INCLUDE_EP_FORCING_CODE
232          CALL CALC_EP_FORCING(myThid)
233    #endif
234    
235    #ifdef ALLOW_AUTODIFF_TAMC
236    C--   HPF directive to help TAMC
237    CHPF$ INDEPENDENT
238    #endif /* ALLOW_AUTODIFF_TAMC */
239    
240          DO bj=myByLo(myThid),myByHi(myThid)
241    
242    #ifdef ALLOW_AUTODIFF_TAMC
243    C--    HPF directive to help TAMC
244    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
245    CHPF$&                  ,phiHydF
246    CHPF$&                  ,KappaRU,KappaRV
247    CHPF$&                  )
248    #endif /* ALLOW_AUTODIFF_TAMC */
249    
250           DO bi=myBxLo(myThid),myBxHi(myThid)
251    
252    #ifdef ALLOW_AUTODIFF_TAMC
253              act1 = bi - myBxLo(myThid)
254              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
255              act2 = bj - myByLo(myThid)
256              max2 = myByHi(myThid) - myByLo(myThid) + 1
257              act3 = myThid - 1
258              max3 = nTx*nTy
259              act4 = ikey_dynamics - 1
260              idynkey = (act1 + 1) + act2*max1
261         &                      + act3*max1*max2
262         &                      + act4*max1*max2*max3
263    #endif /* ALLOW_AUTODIFF_TAMC */
264    
265  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
266  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
267  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
268  C     point numbers. This prevents spurious hardware signals due to  C     point numbers. This prevents spurious hardware signals due to
269  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  
         buoyKM1(i,j) = 0. _d 0  
         buoyK  (i,j) = 0. _d 0  
         maskC (i,j)  = 0. _d 0  
        ENDDO  
       ENDDO  
270    
271        DO bj=myByLo(myThid),myByHi(myThid)          DO k=1,Nr
272         DO bi=myBxLo(myThid),myBxHi(myThid)           DO j=1-OLy,sNy+OLy
273              DO i=1-OLx,sNx+OLx
274  C--     Set up work arrays that need valid initial values             KappaRU(i,j,k) = 0. _d 0
275               KappaRV(i,j,k) = 0. _d 0
276    #ifdef ALLOW_AUTODIFF_TAMC
277    cph(
278    c--   need some re-initialisation here to break dependencies
279    cph)
280               gU(i,j,k,bi,bj) = 0. _d 0
281               gV(i,j,k,bi,bj) = 0. _d 0
282    #endif
283              ENDDO
284             ENDDO
285            ENDDO
286          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
287           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
288            wTrans(i,j)  = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
289            wVel  (i,j,1) = 0. _d 0            fVerU  (i,j,2) = 0. _d 0
290            wVel  (i,j,2) = 0. _d 0            fVerV  (i,j,1) = 0. _d 0
291            fVerT(i,j,1) = 0. _d 0            fVerV  (i,j,2) = 0. _d 0
292            fVerT(i,j,2) = 0. _d 0            phiHydF (i,j)  = 0. _d 0
293            fVerS(i,j,1) = 0. _d 0            phiHydC (i,j)  = 0. _d 0
294            fVerS(i,j,2) = 0. _d 0            dPhiHydX(i,j)  = 0. _d 0
295            fVerU(i,j,1) = 0. _d 0            dPhiHydY(i,j)  = 0. _d 0
296            fVerU(i,j,2) = 0. _d 0            phiSurfX(i,j)  = 0. _d 0
297            fVerV(i,j,1) = 0. _d 0            phiSurfY(i,j)  = 0. _d 0
298            fVerV(i,j,2) = 0. _d 0            guDissip(i,j)  = 0. _d 0
299            pH(i,j,1) = 0. _d 0            gvDissip(i,j)  = 0. _d 0
           K13(i,j,1) = 0. _d 0  
           K23(i,j,1) = 0. _d 0  
           K33(i,j,1) = 0. _d 0  
           KapGM(i,j) = GMkbackground  
300           ENDDO           ENDDO
301          ENDDO          ENDDO
302    
303          iMin = 1-OLx+1  C--     Start computation of dynamics
304          iMax = sNx+OLx          iMin = 0
305          jMin = 1-OLy+1          iMax = sNx+1
306          jMax = sNy+OLy          jMin = 0
307            jMax = sNy+1
308          K = 1  
309          BOTTOM_LAYER = K .EQ. Nz  #ifdef ALLOW_AUTODIFF_TAMC
310    CADJ STORE wvel (:,:,:,bi,bj) =
311  C--     Calculate gradient of surface pressure  CADJ &     comlev1_bibj, key = idynkey, byte = isbyte
312          CALL GRAD_PSURF(  #endif /* ALLOW_AUTODIFF_TAMC */
313       I       bi,bj,iMin,iMax,jMin,jMax,  
314       O       pSurfX,pSurfY,  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
315       I       myThid)  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
316            IF (implicSurfPress.NE.1.) THEN
317  C--     Update fields in top level according to tendency terms            CALL CALC_GRAD_PHI_SURF(
318          CALL CORRECTION_STEP(       I         bi,bj,iMin,iMax,jMin,jMax,
319       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid)       I         etaN,
320         O         phiSurfX,phiSurfY,
321          IF ( .NOT. BOTTOM_LAYER ) THEN       I         myThid )                        
 C--      Update fields in layer below according to tendency terms  
          CALL CORRECTION_STEP(  
      I        bi,bj,iMin,iMax,jMin,jMax,K+1,pSurfX,pSurfY,myTime,myThid)  
         ENDIF  
 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 )  
322          ENDIF          ENDIF
323    
324  C--     Calculate buoyancy  #ifdef ALLOW_AUTODIFF_TAMC
325          CALL CALC_BUOY(  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
326       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
327       O      buoyKm1,  #ifdef ALLOW_KPP
328       I      myThid )  CADJ STORE KPPviscAz (:,:,:,bi,bj)
329    CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
330  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #endif /* ALLOW_KPP */
331          CALL CALC_PHI_HYD(  #endif /* ALLOW_AUTODIFF_TAMC */
      I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,  
      U      phiHyd,  
      I      myThid )  
   
         DO K=2,Nz  
   
          BOTTOM_LAYER = K .EQ. Nz  
          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,pSurfX,pSurfY,myTime,myThid)  
          ENDIF  
 C--      Density of K level (below W(K)) reference to K level  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
      O      rhoK,  
      I      myThid )  
          IF ( .NOT. BOTTOM_LAYER ) THEN  
 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--      Calculate buoyancy  
          CALL CALC_BUOY(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,  
      O       buoyK,  
      I       myThid )  
   
 C--      Integrate hydrostatic balance for pH with BC of pH(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  
          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)  
            buoyKm1(I,J)=buoyK(I,J)  
           ENDDO  
          ENDDO  
   
         ENDDO ! K  
   
         DO K = Nz, 1, -1  
          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,wTrans,wVel,maskC,maskUp,  
      I        myThid)  
332    
333    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
334  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
335           CALL CALC_DIFFUSIVITY(          DO k=1,Nr
336       I        bi,bj,iMin,iMax,jMin,jMax,K,           CALL CALC_VISCOSITY(
337       I        maskC,maskUp,KapGM,K33,       I        bi,bj,iMin,iMax,jMin,jMax,k,
338       O        KappaZT,KappaZS,       O        KappaRU,KappaRV,
339       I        myThid)       I        myThid)
340           ENDDO
341    #endif
342    
343  C--      Calculate accelerations in the momentum equations  #ifdef ALLOW_AUTODIFF_TAMC
344           IF ( momStepping ) THEN  CADJ STORE KappaRU(:,:,:)
345            CALL CALC_MOM_RHS(  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
346       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  CADJ STORE KappaRV(:,:,:)
347       I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
348       I         phiHyd,  #endif /* ALLOW_AUTODIFF_TAMC */
349       U         aTerm,xTerm,cTerm,mTerm,pTerm,  
350       U         fZon, fMer, fVerU, fVerV,  C--     Start of dynamics loop
351       I         myThid)          DO k=1,Nr
352    
353    C--       km1    Points to level above k (=k-1)
354    C--       kup    Cycles through 1,2 to point to layer above
355    C--       kDown  Cycles through 2,1 to point to current layer
356    
357              km1  = MAX(1,k-1)
358              kp1  = MIN(k+1,Nr)
359              kup  = 1+MOD(k+1,2)
360              kDown= 1+MOD(k,2)
361    
362    #ifdef ALLOW_AUTODIFF_TAMC
363             kkey = (idynkey-1)*Nr + k
364    c
365    CADJ STORE totphihyd (:,:,k,bi,bj)
366    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
367    CADJ STORE theta (:,:,k,bi,bj)
368    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
369    CADJ STORE salt  (:,:,k,bi,bj)
370    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
371    CADJ STORE gt(:,:,k,bi,bj)
372    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
373    CADJ STORE gs(:,:,k,bi,bj)
374    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
375    # ifdef NONLIN_FRSURF
376    cph-test
377    CADJ STORE  phiHydC (:,:)
378    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
379    CADJ STORE  phiHydF (:,:)
380    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
381    CADJ STORE  gudissip (:,:)
382    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
383    CADJ STORE  gvdissip (:,:)
384    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
385    CADJ STORE  fVerU (:,:,:)
386    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
387    CADJ STORE  fVerV (:,:,:)
388    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
389    CADJ STORE gu(:,:,k,bi,bj)
390    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
391    CADJ STORE gv(:,:,k,bi,bj)
392    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
393    CADJ STORE gunm1(:,:,k,bi,bj)
394    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
395    CADJ STORE gvnm1(:,:,k,bi,bj)
396    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
397    #  ifdef ALLOW_CD_CODE
398    CADJ STORE unm1(:,:,k,bi,bj)
399    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
400    CADJ STORE vnm1(:,:,k,bi,bj)
401    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
402    CADJ STORE uVelD(:,:,k,bi,bj)
403    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
404    CADJ STORE vVelD(:,:,k,bi,bj)
405    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
406    #  endif
407    # endif
408    #endif /* ALLOW_AUTODIFF_TAMC */
409    
410    C--      Integrate hydrostatic balance for phiHyd with BC of
411    C        phiHyd(z=0)=0
412             IF ( implicitIntGravWave ) THEN
413               CALL CALC_PHI_HYD(
414         I        bi,bj,iMin,iMax,jMin,jMax,k,
415         I        gT, gS,
416         U        phiHydF,
417         O        phiHydC, dPhiHydX, dPhiHydY,
418         I        myTime, myIter, myThid )
419             ELSE
420               CALL CALC_PHI_HYD(
421         I        bi,bj,iMin,iMax,jMin,jMax,k,
422         I        theta, salt,
423         U        phiHydF,
424         O        phiHydC, dPhiHydX, dPhiHydY,
425         I        myTime, myIter, myThid )
426           ENDIF           ENDIF
427    
428  C--      Calculate active tracer tendencies  C--      Calculate accelerations in the momentum equations (gU, gV, ...)
429           IF ( tempStepping ) THEN  C        and step forward storing the result in gU, gV, etc...
430            CALL CALC_GT(           IF ( momStepping ) THEN
431       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  #ifdef ALLOW_MOM_FLUXFORM
432       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,             IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
433       I         K13,K23,KappaZT,KapGM,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
434       U         aTerm,xTerm,fZon,fMer,fVerT,       I         KappaRU, KappaRV,
435       I         myThid)       U         fVerU, fVerV,
436           ENDIF       O         guDissip, gvDissip,
437           IF ( saltStepping ) THEN       I         myTime, myIter, myThid)
438            CALL CALC_GS(  #endif
439       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  #ifdef ALLOW_MOM_VECINV
440       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,             IF (vectorInvariantMomentum) THEN
441       I         K13,K23,KappaZS,KapGM,  C
442       U         aTerm,xTerm,fZon,fMer,fVerS,  # ifdef ALLOW_AUTODIFF_TAMC
443       I         myThid)  #  ifdef NONLIN_FRSURF
444    CADJ STORE fVerU(:,:,:)
445    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
446    CADJ STORE fVerV(:,:,:)
447    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
448    #  endif
449    # endif /* ALLOW_AUTODIFF_TAMC */
450    C
451                 CALL MOM_VECINV(
452         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
453         I         KappaRU, KappaRV,
454         U         fVerU, fVerV,
455         O         guDissip, gvDissip,
456         I         myTime, myIter, myThid)
457               ENDIF
458    #endif
459               CALL TIMESTEP(
460         I         bi,bj,iMin,iMax,jMin,jMax,k,
461         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
462         I         guDissip, gvDissip,
463         I         myTime, myIter, myThid)
464    
465    #ifdef   ALLOW_OBCS
466    C--      Apply open boundary conditions
467               IF (useOBCS) THEN
468                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
469               ENDIF
470    #endif   /* ALLOW_OBCS */
471    
472           ENDIF           ENDIF
473    
 C--      Prediction step (step forward all model variables)  
          CALL TIMESTEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,  
      I       myThid)  
   
 C--      Diagnose barotropic divergence of predicted fields  
          CALL DIV_G(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,  
      I       xA,yA,  
      I       myThid)  
474    
475  C--      Cumulative diagnostic calculations (ie. time-averaging)  C--     end of dynamics k loop (1:Nr)
476  #ifdef ALLOW_DIAGNOSTICS          ENDDO
477           IF (taveFreq.GT.0.) THEN  
478            CALL DO_TIME_AVERAGES(  C--     Implicit Vertical advection & viscosity
479       I                           myTime, myIter, bi, bj, K, kUp, kDown,  #ifdef INCLUDE_IMPLVERTADV_CODE
480       I                           K13, K23, wVel, KapGM,          IF ( momImplVertAdv ) THEN
481       I                           myThid )            CALL MOM_U_IMPLICIT_R( kappaRU,
482           ENDIF       I                           bi, bj, myTime, myIter, myThid )
483  #endif            CALL MOM_V_IMPLICIT_R( kappaRV,
484         I                           bi, bj, myTime, myIter, myThid )
485            ELSEIF ( implicitViscosity ) THEN
486    #else /* INCLUDE_IMPLVERTADV_CODE */
487            IF     ( implicitViscosity ) THEN
488    #endif /* INCLUDE_IMPLVERTADV_CODE */
489    #ifdef    ALLOW_AUTODIFF_TAMC
490    CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
491    CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
492    #endif    /* ALLOW_AUTODIFF_TAMC */
493              CALL IMPLDIFF(
494         I         bi, bj, iMin, iMax, jMin, jMax,
495         I         -1, KappaRU,recip_HFacW,
496         U         gU,
497         I         myThid )
498    #ifdef    ALLOW_AUTODIFF_TAMC
499    CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
500    CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
501    #endif    /* ALLOW_AUTODIFF_TAMC */
502              CALL IMPLDIFF(
503         I         bi, bj, iMin, iMax, jMin, jMax,
504         I         -2, KappaRV,recip_HFacS,
505         U         gV,
506         I         myThid )
507            ENDIF
508    
509          ENDDO ! K  #ifdef   ALLOW_OBCS
510    C--      Apply open boundary conditions
511            IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
512               DO K=1,Nr
513                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
514               ENDDO
515            ENDIF
516    #endif   /* ALLOW_OBCS */
517    
518  C--     Implicit diffusion  #ifdef    ALLOW_CD_CODE
519          IF (implicitDiffusion) THEN          IF (implicitViscosity.AND.useCDscheme) THEN
520           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,  #ifdef    ALLOW_AUTODIFF_TAMC
521       I                  KappaZT,KappaZS,  CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
522       I                  myThid )  #endif    /* ALLOW_AUTODIFF_TAMC */
523              CALL IMPLDIFF(
524         I         bi, bj, iMin, iMax, jMin, jMax,
525         I         0, KappaRU,recip_HFacW,
526         U         vVelD,
527         I         myThid )
528    #ifdef    ALLOW_AUTODIFF_TAMC
529    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
530    #endif    /* ALLOW_AUTODIFF_TAMC */
531              CALL IMPLDIFF(
532         I         bi, bj, iMin, iMax, jMin, jMax,
533         I         0, KappaRV,recip_HFacS,
534         U         uVelD,
535         I         myThid )
536          ENDIF          ENDIF
537    #endif    /* ALLOW_CD_CODE */
538    C--     End implicit Vertical advection & viscosity
539    
540         ENDDO         ENDDO
541        ENDDO        ENDDO
542    
543  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),  #ifdef ALLOW_OBCS
544  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))        IF (useOBCS) THEN
545  C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),         CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
546  C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)        ENDIF
547  C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),  #endif
548  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)  
549  C     write(0,*) 'dynamics: wVel(1) ',  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
550  C    &            minval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.),  
551  C    &            maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.)  #ifdef ALLOW_NONHYDROSTATIC
552  C     write(0,*) 'dynamics: wVel(2) ',  C--   Step forward W field in N-H algorithm
553  C    &            minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.),        IF ( nonHydrostatic ) THEN
554  C    &            maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.)  #ifdef ALLOW_DEBUG
555  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),           IF ( debugLevel .GE. debLevB )
556  cblk &                           maxval(K13(1:sNx,1:sNy,:))       &     CALL DEBUG_CALL('CALC_GW', myThid )
557  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  #endif
558  cblk &                           maxval(K23(1:sNx,1:sNy,:))           CALL TIMER_START('CALC_GW          [DYNAMICS]',myThid)
559  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),           CALL CALC_GW( myTime, myIter, myThid )
560  cblk &                           maxval(K33(1:sNx,1:sNy,:))        ENDIF
561  C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),        IF ( nonHydrostatic.OR.implicitIntGravWave )
562  C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))       &   CALL TIMESTEP_WVEL( myTime, myIter, myThid )
563  C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),        IF ( nonHydrostatic )
564  C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))       &   CALL TIMER_STOP ('CALC_GW          [DYNAMICS]',myThid)
565  C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),  #endif
566  C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))  
567  C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
568  C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))  
569  C     write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.),  Cml(
570  C    &                           maxval(pH/(Gravity*Rhonil))  C     In order to compare the variance of phiHydLow of a p/z-coordinate
571    C     run with etaH of a z/p-coordinate run the drift of phiHydLow
572    C     has to be removed by something like the following subroutine:
573    C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
574    C     &                'phiHydLow', myThid )
575    Cml)
576    
577    #ifdef ALLOW_DIAGNOSTICS
578          IF ( usediagnostics ) THEN
579    
580           CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD  ',0,Nr,0,1,1,myThid)
581           CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT  ',0, 1,0,1,1,myThid)
582    
583           tmpFac = 1. _d 0
584           CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
585         &                                 'PHIHYDSQ',0,Nr,0,1,1,myThid)
586    
587           CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
588         &                                 'PHIBOTSQ',0, 1,0,1,1,myThid)
589    
590          ENDIF
591    #endif /* ALLOW_DIAGNOSTICS */
592          
593    #ifdef ALLOW_DEBUG
594          If ( debugLevel .GE. debLevB ) THEN
595           CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
596           CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
597           CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
598           CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
599           CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
600           CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
601           CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
602           CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
603           CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
604           CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
605    #ifndef ALLOW_ADAMSBASHFORTH_3
606           CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
607           CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
608           CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
609           CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
610    #endif
611          ENDIF
612    #endif
613    
614    #ifdef DYNAMICS_GUGV_EXCH_CHECK
615    C- jmc: For safety checking only: This Exchange here should not change
616    C       the solution. If solution changes, it means something is wrong,
617    C       but it does not mean that it is less wrong with this exchange.
618          IF ( debugLevel .GT. debLevB ) THEN
619           CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
620          ENDIF
621    #endif
622    
623    #ifdef ALLOW_DEBUG
624          IF ( debugLevel .GE. debLevB )
625         &   CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
626    #endif
627    
628        RETURN        RETURN
629        END        END

Legend:
Removed from v.1.26  
changed lines
  Added in v.1.129

  ViewVC Help
Powered by ViewVC 1.1.22