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

Legend:
Removed from v.1.48  
changed lines
  Added in v.1.117

  ViewVC Help
Powered by ViewVC 1.1.22