/[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.110 by jmc, Wed Nov 10 03:02:00 2004 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
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    #ifdef ALLOW_CD_CODE
77    #include "CD_CODE_VARS.h"
78    #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 TIMEAVE_CUMUL_1T
112    C      |-- CALL DEBUG_STATS_RL
113    
114    C     !INPUT/OUTPUT PARAMETERS:
115  C     == Routine arguments ==  C     == Routine arguments ==
116  C     myTime - Current time in simulation  C     myTime - Current time in simulation
117  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
118  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
119        _RL myTime        _RL myTime
120        INTEGER myIter        INTEGER myIter
121          INTEGER myThid
122    
123    C     !LOCAL VARIABLES:
124  C     == Local variables  C     == Local variables
125  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[STUV]               o fVer: Vertical flux term - note fVer
 C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  
 C     wVel                     o uTrans: Zonal transport  
 C                              o vTrans: Meridional transport  
 C                              o wTrans: Vertical transport  
 C                              o wVel:   Vertical velocity at upper and lower  
 C                                        cell faces.  
 C     maskC,maskUp             o maskC: land/water mask for tracer cells  
 C                              o maskUp: land/water mask for W points  
 C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  
 C     mTerm, pTerm,            tendency equations.  
 C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  
 C                              o xTerm: Mixing term  
 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  
126  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
127  C                                      so we need an fVer for each  C                                      so we need an fVer for each
128  C                                      variable.  C                                      variable.
129  C     rhoK, rhoKM1   - Density at current level, level above and level below.  C     phiHydC    :: hydrostatic potential anomaly at cell center
130  C     rhoKP1                                                                    C                   In z coords phiHyd is the hydrostatic potential
131  C     buoyK, buoyKM1 - Buoyancy at current level and level above.  C                      (=pressure/rho0) anomaly
132  C     phiHyd         - Hydrostatic part of the potential phi.  C                   In p coords phiHyd is the geopotential height anomaly.
133  C                      In z coords phiHyd is the hydrostatic pressure anomaly  C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
134  C                      In p coords phiHyd is the geopotential surface height anomaly.  C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
135  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
136  C     jMin, jMax   are applied.  C     phiSurfY             or geopotential (atmos) in X and Y direction
137    C     guDissip   :: dissipation tendency (all explicit terms), u component
138    C     gvDissip   :: dissipation tendency (all explicit terms), v component
139    C     iMin, iMax     - Ranges and sub-block indices on which calculations
140    C     jMin, jMax       are applied.
141  C     bi, bj  C     bi, bj
142  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
143  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
144  C                          into fVerTerm  C                      index into fVerTerm.
145        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
147        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
151        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156        _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)  
157    
158        INTEGER iMin, iMax        INTEGER iMin, iMax
159        INTEGER jMin, jMax        INTEGER jMin, jMax
160        INTEGER bi, bj        INTEGER bi, bj
161        INTEGER i, j        INTEGER i, j
162        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kp1, kup, kDown
       LOGICAL BOTTOM_LAYER  
163    
164          LOGICAL  DIFFERENT_MULTIPLE
165          EXTERNAL DIFFERENT_MULTIPLE
166    
167  C---    The algorithm...  C---    The algorithm...
168  C  C
169  C       "Correction Step"  C       "Correction Step"
# Line 127  C       "Calculation of Gs" Line 178  C       "Calculation of Gs"
178  C       ===================  C       ===================
179  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
180  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         w = sum_z ( div. u[n] )  
181  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
182    C         b   = b(rho, theta)
183  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
184  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
185  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
186  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
187  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
188  C  C
189  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
190  C       ================================  C       ================================
# Line 156  C         salt* = salt[n] + dt x ( 3/2 G Line 207  C         salt* = salt[n] + dt x ( 3/2 G
207  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
208  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
209  C---  C---
210    CEOP
211    
212    C-- Call to routine for calculation of
213    C   Eliassen-Palm-flux-forced U-tendency,
214    C   if desired:
215    #ifdef INCLUDE_EP_FORCING_CODE
216          CALL CALC_EP_FORCING(myThid)
217    #endif
218    
219    #ifdef ALLOW_AUTODIFF_TAMC
220    C--   HPF directive to help TAMC
221    CHPF$ INDEPENDENT
222    #endif /* ALLOW_AUTODIFF_TAMC */
223    
224          DO bj=myByLo(myThid),myByHi(myThid)
225    
226    #ifdef ALLOW_AUTODIFF_TAMC
227    C--    HPF directive to help TAMC
228    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
229    CHPF$&                  ,phiHydF
230    CHPF$&                  ,KappaRU,KappaRV
231    CHPF$&                  )
232    #endif /* ALLOW_AUTODIFF_TAMC */
233    
234           DO bi=myBxLo(myThid),myBxHi(myThid)
235    
236    #ifdef ALLOW_AUTODIFF_TAMC
237              act1 = bi - myBxLo(myThid)
238              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
239              act2 = bj - myByLo(myThid)
240              max2 = myByHi(myThid) - myByLo(myThid) + 1
241              act3 = myThid - 1
242              max3 = nTx*nTy
243              act4 = ikey_dynamics - 1
244              idynkey = (act1 + 1) + act2*max1
245         &                      + act3*max1*max2
246         &                      + act4*max1*max2*max3
247    #endif /* ALLOW_AUTODIFF_TAMC */
248    
249  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
250  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
251  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
252  C     point numbers. This prevents spurious hardware signals due to  C     point numbers. This prevents spurious hardware signals due to
253  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  
254    
255        DO bj=myByLo(myThid),myByHi(myThid)          DO k=1,Nr
256         DO bi=myBxLo(myThid),myBxHi(myThid)           DO j=1-OLy,sNy+OLy
257              DO i=1-OLx,sNx+OLx
258  C--     Set up work arrays that need valid initial values             KappaRU(i,j,k) = 0. _d 0
259               KappaRV(i,j,k) = 0. _d 0
260    #ifdef ALLOW_AUTODIFF_TAMC
261    cph(
262    c--   need some re-initialisation here to break dependencies
263    c--   totphihyd is assumed zero from ini_pressure, i.e.
264    c--   avoiding iterate pressure p = integral of (g*rho(p)*dz)
265    cph)
266               totPhiHyd(i,j,k,bi,bj) = 0. _d 0
267               gu(i,j,k,bi,bj) = 0. _d 0
268               gv(i,j,k,bi,bj) = 0. _d 0
269    #endif
270              ENDDO
271             ENDDO
272            ENDDO
273          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
274           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
275            wTrans(i,j)  = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
276            wVel  (i,j,1) = 0. _d 0            fVerU  (i,j,2) = 0. _d 0
277            wVel  (i,j,2) = 0. _d 0            fVerV  (i,j,1) = 0. _d 0
278            fVerT(i,j,1) = 0. _d 0            fVerV  (i,j,2) = 0. _d 0
279            fVerT(i,j,2) = 0. _d 0            phiHydF (i,j)  = 0. _d 0
280            fVerS(i,j,1) = 0. _d 0            phiHydC (i,j)  = 0. _d 0
281            fVerS(i,j,2) = 0. _d 0            dPhiHydX(i,j)  = 0. _d 0
282            fVerU(i,j,1) = 0. _d 0            dPhiHydY(i,j)  = 0. _d 0
283            fVerU(i,j,2) = 0. _d 0            phiSurfX(i,j)  = 0. _d 0
284            fVerV(i,j,1) = 0. _d 0            phiSurfY(i,j)  = 0. _d 0
285            fVerV(i,j,2) = 0. _d 0            guDissip(i,j)  = 0. _d 0
286            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  
287           ENDDO           ENDDO
288          ENDDO          ENDDO
289    
290          iMin = 1-OLx+1  C--     Start computation of dynamics
291          iMax = sNx+OLx          iMin = 0
292          jMin = 1-OLy+1          iMax = sNx+1
293          jMax = sNy+OLy          jMin = 0
294            jMax = sNy+1
295          K = 1  
296          BOTTOM_LAYER = K .EQ. Nz  #ifdef ALLOW_AUTODIFF_TAMC
297    CADJ STORE wvel (:,:,:,bi,bj) =
298  C--     Calculate gradient of surface pressure  CADJ &     comlev1_bibj, key = idynkey, byte = isbyte
299          CALL GRAD_PSURF(  #endif /* ALLOW_AUTODIFF_TAMC */
300       I       bi,bj,iMin,iMax,jMin,jMax,  
301       O       pSurfX,pSurfY,  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
302       I       myThid)  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
303            IF (implicSurfPress.NE.1.) THEN
304  C--     Update fields in top level according to tendency terms            CALL CALC_GRAD_PHI_SURF(
305          CALL CORRECTION_STEP(       I         bi,bj,iMin,iMax,jMin,jMax,
306       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid)       I         etaN,
307         O         phiSurfX,phiSurfY,
308          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)  
309          ENDIF          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 )  
         ENDIF  
   
 C--     Calculate buoyancy  
         CALL CALC_BUOY(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,  
      O      buoyKm1,  
      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,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  
310    
311          DO K = Nz, 1, -1  #ifdef ALLOW_AUTODIFF_TAMC
312           kM1  =max(1,k-1)   ! Points to level above k (=k-1)  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
313           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
314           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer  #ifdef ALLOW_KPP
315           iMin = 1-OLx+2  CADJ STORE KPPviscAz (:,:,:,bi,bj)
316           iMax = sNx+OLx-1  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
317           jMin = 1-OLy+2  #endif /* ALLOW_KPP */
318           jMax = sNy+OLy-1  #endif /* ALLOW_AUTODIFF_TAMC */
   
 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)  
319    
320    #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        KappaZT,KappaZS,       O        KappaRU,KappaRV,
326       I        myThid)       I        myThid)
327           ENDDO
328    #endif
329    
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--      Calculate accelerations in the momentum equations  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,wTrans,wVel,maskC,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
375       I         phiHyd,       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         myThid)  #endif
379           ENDIF  #ifdef ALLOW_MOM_VECINV
380               IF (vectorInvariantMomentum) CALL MOM_VECINV(
381         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
382         I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
383         U         fVerU, fVerV,
384         O         guDissip, gvDissip,
385         I         myTime, myIter, myThid)
386    #endif
387               CALL TIMESTEP(
388         I         bi,bj,iMin,iMax,jMin,jMax,k,
389         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
390         I         guDissip, gvDissip,
391         I         myTime, myIter, myThid)
392    
393    #ifdef   ALLOW_OBCS
394    C--      Apply open boundary conditions
395               IF (useOBCS) THEN
396                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
397               ENDIF
398    #endif   /* ALLOW_OBCS */
399    
 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,maskC,  
      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,maskC,  
      I         K13,K23,KappaZS,KapGM,  
      U         aTerm,xTerm,fZon,fMer,fVerS,  
      I         myThid)  
400           ENDIF           ENDIF
401    
 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)  
   
 C--      Cumulative diagnostic calculations (ie. time-averaging)  
 #ifdef ALLOW_DIAGNOSTICS  
          IF (taveFreq.GT.0.) THEN  
           CALL DO_TIME_AVERAGES(  
      I                           myTime, myIter, bi, bj, K, kUp, kDown,  
      I                           K13, K23, wVel, KapGM,  
      I                           myThid )  
          ENDIF  
 #endif  
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           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,          IF ( momImplVertAdv ) THEN
409       I                  KappaZT,KappaZS,            CALL MOM_U_IMPLICIT_R( kappaRU,
410       I                  myThid )       I                           bi, bj, myTime, myIter, myThid )
411              CALL MOM_V_IMPLICIT_R( kappaRV,
412         I                           bi, bj, myTime, myIter, myThid )
413            ELSEIF ( implicitViscosity ) THEN
414    #else /* INCLUDE_IMPLVERTADV_CODE */
415            IF     ( implicitViscosity ) THEN
416    #endif /* INCLUDE_IMPLVERTADV_CODE */
417    #ifdef    ALLOW_AUTODIFF_TAMC
418    CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
419    CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
420    #endif    /* ALLOW_AUTODIFF_TAMC */
421              CALL IMPLDIFF(
422         I         bi, bj, iMin, iMax, jMin, jMax,
423         I         deltaTmom, KappaRU,recip_HFacW,
424         U         gU,
425         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(
431         I         bi, bj, iMin, iMax, jMin, jMax,
432         I         deltaTmom, KappaRV,recip_HFacS,
433         U         gV,
434         I         myThid )
435          ENDIF          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(
452         I         bi, bj, iMin, iMax, jMin, jMax,
453         I         deltaTmom, KappaRU,recip_HFacW,
454         U         vVelD,
455         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(
460         I         bi, bj, iMin, iMax, jMin, jMax,
461         I         deltaTmom, KappaRV,recip_HFacS,
462         U         uVelD,
463         I         myThid )
464            ENDIF
465    #endif    /* ALLOW_CD_CODE */
466    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: wVel(1) ',  Cml(
478  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
479  C    &            maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.)  C     run with etaH of a z/p-coordinate run the drift of phiHydLow
480  C     write(0,*) 'dynamics: wVel(2) ',  C     has to be removed by something like the following subroutine:
481  C    &            minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.),  C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
482  C    &            maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.)  C     &                'phiHydLow', myThid )
483  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  Cml)
484  cblk &                           maxval(K13(1:sNx,1:sNy,:))  
485  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  #ifdef ALLOW_DEBUG
486  cblk &                           maxval(K23(1:sNx,1:sNy,:))        If ( debugLevel .GE. debLevB ) THEN
487  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),         CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
488  cblk &                           maxval(K33(1:sNx,1:sNy,:))         CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
489  C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
490  C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
491  C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
492  C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
493  C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
494  C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
495  C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
496  C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
497  C     write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.),         CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
498  C    &                           maxval(pH/(Gravity*Rhonil))         CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
499           CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
500           CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
501          ENDIF
502    #endif
503    
504        RETURN        RETURN
505        END        END

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

  ViewVC Help
Powered by ViewVC 1.1.22