/[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.23 by adcroft, Wed Jul 1 19:49:36 1998 UTC revision 1.106 by jmc, Sat Jan 3 01:01:34 2004 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CPP_EEOPTIONS.h"  #include "PACKAGES_CONFIG.h"
5    #include "CPP_OPTIONS.h"
6    
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_PASSIVE_TRACER
81    #include "TR1.h"
82    #endif
83    #ifdef ALLOW_AUTODIFF_TAMC
84    # include "tamc.h"
85    # include "tamc_keys.h"
86    # include "FFIELDS.h"
87    # include "EOS.h"
88    # ifdef ALLOW_KPP
89    #  include "KPP.h"
90    # endif
91    #endif /* ALLOW_AUTODIFF_TAMC */
92    
93    C     !CALLING SEQUENCE:
94    C     DYNAMICS()
95    C      |
96    C      |-- CALC_GRAD_PHI_SURF
97    C      |
98    C      |-- CALC_VISCOSITY
99    C      |
100    C      |-- CALC_PHI_HYD  
101    C      |
102    C      |-- MOM_FLUXFORM  
103    C      |
104    C      |-- MOM_VECINV    
105    C      |
106    C      |-- TIMESTEP      
107    C      |
108    C      |-- OBCS_APPLY_UV
109    C      |
110    C      |-- IMPLDIFF      
111    C      |
112    C      |-- OBCS_APPLY_UV
113    C      |
114    C      |-- CALL TIMEAVE_CUMUL_1T
115    C      |-- CALL DEBUG_STATS_RL
116    
117    C     !INPUT/OUTPUT PARAMETERS:
118  C     == Routine arguments ==  C     == Routine arguments ==
119  C     myTime - Current time in simulation  C     myTime - Current time in simulation
120  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
121  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
122        _RL myTime        _RL myTime
123        INTEGER myIter        INTEGER myIter
124          INTEGER myThid
125    
126    C     !LOCAL VARIABLES:
127  C     == Local variables  C     == Local variables
128  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  
129  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
130  C                                      so we need an fVer for each  C                                      so we need an fVer for each
131  C                                      variable.  C                                      variable.
132  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     phiHydC    :: hydrostatic potential anomaly at cell center
133  C     jMin, jMax   are applied.  C                   In z coords phiHyd is the hydrostatic potential
134    C                      (=pressure/rho0) anomaly
135    C                   In p coords phiHyd is the geopotential height anomaly.
136    C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
137    C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
138    C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
139    C     phiSurfY             or geopotential (atmos) in X and Y direction
140    C     iMin, iMax     - Ranges and sub-block indices on which calculations
141    C     jMin, jMax       are applied.
142  C     bi, bj  C     bi, bj
143  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
144  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
145  C                          into fVerTerm  C                      index into fVerTerm.
146        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
147        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
148        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
151        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
152        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
155        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
       _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)  
       _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL pH    (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 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)  
156    
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          LOGICAL  DIFFERENT_MULTIPLE
164          EXTERNAL DIFFERENT_MULTIPLE
165    
166  C---    The algorithm...  C---    The algorithm...
167  C  C
168  C       "Correction Step"  C       "Correction Step"
# Line 119  C       "Calculation of Gs" Line 177  C       "Calculation of Gs"
177  C       ===================  C       ===================
178  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
179  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         w = sum_z ( div. u[n] )  
180  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
181    C         b   = b(rho, theta)
182  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
183  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
184  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
185  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
186  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
187  C  C
188  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
189  C       ================================  C       ================================
# Line 148  C         salt* = salt[n] + dt x ( 3/2 G Line 206  C         salt* = salt[n] + dt x ( 3/2 G
206  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
207  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
208  C---  C---
209    CEOP
210    
211    C-- Call to routine for calculation of
212    C   Eliassen-Palm-flux-forced U-tendency,
213    C   if desired:
214    #ifdef INCLUDE_EP_FORCING_CODE
215          CALL CALC_EP_FORCING(myThid)
216    #endif
217    
218    #ifdef ALLOW_AUTODIFF_TAMC
219    C--   HPF directive to help TAMC
220    CHPF$ INDEPENDENT
221    #endif /* ALLOW_AUTODIFF_TAMC */
222    
223          DO bj=myByLo(myThid),myByHi(myThid)
224    
225    #ifdef ALLOW_AUTODIFF_TAMC
226    C--    HPF directive to help TAMC
227    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
228    CHPF$&                  ,phiHydF
229    CHPF$&                  ,KappaRU,KappaRV
230    CHPF$&                  )
231    #endif /* ALLOW_AUTODIFF_TAMC */
232    
233           DO bi=myBxLo(myThid),myBxHi(myThid)
234    
235    #ifdef ALLOW_AUTODIFF_TAMC
236              act1 = bi - myBxLo(myThid)
237              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
238              act2 = bj - myByLo(myThid)
239              max2 = myByHi(myThid) - myByLo(myThid) + 1
240              act3 = myThid - 1
241              max3 = nTx*nTy
242              act4 = ikey_dynamics - 1
243              idynkey = (act1 + 1) + act2*max1
244         &                      + act3*max1*max2
245         &                      + act4*max1*max2*max3
246    #endif /* ALLOW_AUTODIFF_TAMC */
247    
248  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
249  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
250  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
251  C     point numbers. This prevents spurious hardware signals due to  C     point numbers. This prevents spurious hardware signals due to
252  C     uninitialised but inert locations.  C     uninitialised but inert locations.
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         xA(i,j)      = 0. _d 0  
         yA(i,j)      = 0. _d 0  
         uTrans(i,j)  = 0. _d 0  
         vTrans(i,j)  = 0. _d 0  
         aTerm(i,j)   = 0. _d 0  
         xTerm(i,j)   = 0. _d 0  
         cTerm(i,j)   = 0. _d 0  
         mTerm(i,j)   = 0. _d 0  
         pTerm(i,j)   = 0. _d 0  
         fZon(i,j)    = 0. _d 0  
         fMer(i,j)    = 0. _d 0  
         DO K=1,nZ  
          pH (i,j,k)  = 0. _d 0  
          K13(i,j,k) = 0. _d 0  
          K23(i,j,k) = 0. _d 0  
          K33(i,j,k) = 0. _d 0  
          KappaZT(i,j,k) = 0. _d 0  
         ENDDO  
         rhokm1(i,j)  = 0. _d 0  
         rhok  (i,j)  = 0. _d 0  
         rhokp1(i,j)  = 0. _d 0  
         rhotmp(i,j)  = 0. _d 0  
         maskC (i,j)  = 0. _d 0  
        ENDDO  
       ENDDO  
   
       DO bj=myByLo(myThid),myByHi(myThid)  
        DO bi=myBxLo(myThid),myBxHi(myThid)  
253    
254  C--     Set up work arrays that need valid initial values          DO k=1,Nr
255             DO j=1-OLy,sNy+OLy
256              DO i=1-OLx,sNx+OLx
257               KappaRU(i,j,k) = 0. _d 0
258               KappaRV(i,j,k) = 0. _d 0
259    #ifdef ALLOW_AUTODIFF_TAMC
260    cph(
261    c--   need some re-initialisation here to break dependencies
262    c--   totphihyd is assumed zero from ini_pressure, i.e.
263    c--   avoiding iterate pressure p = integral of (g*rho(p)*dz)
264    cph)
265               totPhiHyd(i,j,k,bi,bj) = 0. _d 0
266               gu(i,j,k,bi,bj) = 0. _d 0
267               gv(i,j,k,bi,bj) = 0. _d 0
268    #endif
269              ENDDO
270             ENDDO
271            ENDDO
272          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
273           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
274            wTrans(i,j)  = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
275            wVel  (i,j,1) = 0. _d 0            fVerU  (i,j,2) = 0. _d 0
276            wVel  (i,j,2) = 0. _d 0            fVerV  (i,j,1) = 0. _d 0
277            fVerT(i,j,1) = 0. _d 0            fVerV  (i,j,2) = 0. _d 0
278            fVerT(i,j,2) = 0. _d 0            phiHydF (i,j)  = 0. _d 0
279            fVerS(i,j,1) = 0. _d 0            phiHydC (i,j)  = 0. _d 0
280            fVerS(i,j,2) = 0. _d 0            dPhiHydX(i,j)  = 0. _d 0
281            fVerU(i,j,1) = 0. _d 0            dPhiHydY(i,j)  = 0. _d 0
282            fVerU(i,j,2) = 0. _d 0            phiSurfX(i,j)  = 0. _d 0
283            fVerV(i,j,1) = 0. _d 0            phiSurfY(i,j)  = 0. _d 0
           fVerV(i,j,2) = 0. _d 0  
           pH(i,j,1) = 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  
284           ENDDO           ENDDO
285          ENDDO          ENDDO
286    
287          iMin = 1-OLx+1  C--     Start computation of dynamics
288          iMax = sNx+OLx          iMin = 0
289          jMin = 1-OLy+1          iMax = sNx+1
290          jMax = sNy+OLy          jMin = 0
291            jMax = sNy+1
292          K = 1  
293          BOTTOM_LAYER = K .EQ. Nz  #ifdef ALLOW_AUTODIFF_TAMC
294    CADJ STORE wvel (:,:,:,bi,bj) =
295  C--     Calculate gradient of surface pressure  CADJ &     comlev1_bibj, key = idynkey, byte = isbyte
296          CALL GRAD_PSURF(  #endif /* ALLOW_AUTODIFF_TAMC */
297       I       bi,bj,iMin,iMax,jMin,jMax,  
298       O       pSurfX,pSurfY,  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
299       I       myThid)  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
300            IF (implicSurfPress.NE.1.) THEN
301  C--     Update fields in top level according to tendency terms            CALL CALC_GRAD_PHI_SURF(
302          CALL CORRECTION_STEP(       I         bi,bj,iMin,iMax,jMin,jMax,
303       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid)       I         etaN,
304         O         phiSurfX,phiSurfY,
305          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)  
306          ENDIF          ENDIF
307    
308  C--     Density of 1st level (below W(1)) reference to level 1  #ifdef ALLOW_AUTODIFF_TAMC
309          CALL FIND_RHO(  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
310       I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
311       O     rhoKm1,  #ifdef ALLOW_KPP
312       I     myThid )  CADJ STORE KPPviscAz (:,:,:,bi,bj)
313    CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
314          IF ( .NOT. BOTTOM_LAYER ) THEN  #endif /* ALLOW_KPP */
315  C--      Check static stability with layer below  #endif /* ALLOW_AUTODIFF_TAMC */
 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--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
         CALL CALC_PH(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKm1,  
      U      pH,  
      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--      Update fields in layer below according to tendency terms  
 C        CALL CORRECTION_STEP(  
 C    I        bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid)  
   
 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--      Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
          CALL CALC_PH(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoK,  
      U       pH,  
      I       myThid )  
 C--      Calculate iso-neutral slopes for the GM/Redi parameterisation  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,  
      O      rhoTmp,  
      I      myThid )  
          CALL CALC_ISOSLOPES(  
      I             bi, bj, iMin, iMax, jMin, jMax, K,  
      I             rhoKm1, rhoK, rhotmp,  
      O             K13, K23, K33, KapGM,  
      I             myThid )  
          DO J=jMin,jMax  
           DO I=iMin,iMax  
            rhoKm1(I,J)=rhoK(I,J)  
           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)  
316    
317    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
318  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
319           CALL CALC_DIFFUSIVITY(          DO k=1,Nr
320       I        bi,bj,iMin,iMax,jMin,jMax,K,           CALL CALC_VISCOSITY(
321       I        maskC,maskUp,KapGM,K33,       I        bi,bj,iMin,iMax,jMin,jMax,k,
322       O        KappaZT,KappaZS,       O        KappaRU,KappaRV,
323       I        myThid)       I        myThid)
324           ENDDO
325    #endif
326    
327  C--      Calculate accelerations in the momentum equations  #ifdef ALLOW_AUTODIFF_TAMC
328           IF ( momStepping ) THEN  CADJ STORE KappaRU(:,:,:)
329            CALL CALC_MOM_RHS(  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
330       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  CADJ STORE KappaRV(:,:,:)
331       I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
332       I         pH,  #endif /* ALLOW_AUTODIFF_TAMC */
333       U         aTerm,xTerm,cTerm,mTerm,pTerm,  
334       U         fZon, fMer, fVerU, fVerV,  C--     Start of dynamics loop
335       I         myThid)          DO k=1,Nr
336    
337    C--       km1    Points to level above k (=k-1)
338    C--       kup    Cycles through 1,2 to point to layer above
339    C--       kDown  Cycles through 2,1 to point to current layer
340    
341              km1  = MAX(1,k-1)
342              kp1  = MIN(k+1,Nr)
343              kup  = 1+MOD(k+1,2)
344              kDown= 1+MOD(k,2)
345    
346    #ifdef ALLOW_AUTODIFF_TAMC
347             kkey = (idynkey-1)*Nr + k
348    c
349    CADJ STORE totphihyd (:,:,k,bi,bj)
350    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
351    CADJ STORE gt (:,:,k,bi,bj)
352    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
353    CADJ STORE gs (:,:,k,bi,bj)
354    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
355    CADJ STORE theta (:,:,k,bi,bj)
356    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
357    CADJ STORE salt  (:,:,k,bi,bj)
358    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
359    #endif /* ALLOW_AUTODIFF_TAMC */
360    
361    C--      Integrate hydrostatic balance for phiHyd with BC of
362    C        phiHyd(z=0)=0
363    C        distinguishe between Stagger and Non Stagger time stepping
364             IF (staggerTimeStep) THEN
365               CALL CALC_PHI_HYD(
366         I        bi,bj,iMin,iMax,jMin,jMax,k,
367         I        gT, gS,
368         U        phiHydF,
369         O        phiHydC, dPhiHydX, dPhiHydY,
370         I        myTime, myIter, myThid )
371             ELSE
372               CALL CALC_PHI_HYD(
373         I        bi,bj,iMin,iMax,jMin,jMax,k,
374         I        theta, salt,
375         U        phiHydF,
376         O        phiHydC, dPhiHydX, dPhiHydY,
377         I        myTime, myIter, myThid )
378           ENDIF           ENDIF
379    
380  C--      Calculate active tracer tendencies  C--      Calculate accelerations in the momentum equations (gU, gV, ...)
381           IF ( tempStepping ) THEN  C        and step forward storing the result in gU, gV, etc...
382            CALL CALC_GT(           IF ( momStepping ) THEN
383       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  #ifdef ALLOW_MOM_FLUXFORM
384       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,             IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
385       I         K13,K23,KappaZT,KapGM,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
386       U         aTerm,xTerm,fZon,fMer,fVerT,       I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
387       I         myThid)       U         fVerU, fVerV,
388           ENDIF       I         myTime, myIter, myThid)
389           IF ( saltStepping ) THEN  #endif
390            CALL CALC_GS(  #ifdef ALLOW_MOM_VECINV
391       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,             IF (vectorInvariantMomentum) CALL MOM_VECINV(
392       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
393       I         K13,K23,KappaZS,KapGM,       I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
394       U         aTerm,xTerm,fZon,fMer,fVerS,       U         fVerU, fVerV,
395       I         myThid)       I         myTime, myIter, myThid)
396           ENDIF  #endif
397               CALL TIMESTEP(
398         I         bi,bj,iMin,iMax,jMin,jMax,k,
399         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
400         I         myTime, myIter, myThid)
401    
402    #ifdef   ALLOW_OBCS
403    C--      Apply open boundary conditions
404               IF (useOBCS) THEN
405                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
406               ENDIF
407    #endif   /* ALLOW_OBCS */
408    
 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,  
      I                           myThid )  
409           ENDIF           ENDIF
 #endif  
410    
         ENDDO ! K  
411    
412  C--     Implicit diffusion  C--     end of dynamics k loop (1:Nr)
413          IF (implicitDiffusion) THEN          ENDDO
414           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,  
415       I                  KappaZT,KappaZS,  C--     Implicit Vertical advection & viscosity
416       I                  myThid )  #ifdef INCLUDE_IMPLVERTADV_CODE
417            IF ( momImplVertAdv ) THEN
418              CALL MOM_U_IMPLICIT_R( kappaRU,
419         I                           bi, bj, myTime, myIter, myThid )
420              CALL MOM_V_IMPLICIT_R( kappaRV,
421         I                           bi, bj, myTime, myIter, myThid )
422            ELSEIF ( implicitViscosity ) THEN
423    #else /* INCLUDE_IMPLVERTADV_CODE */
424            IF     ( implicitViscosity ) THEN
425    #endif /* INCLUDE_IMPLVERTADV_CODE */
426    #ifdef    ALLOW_AUTODIFF_TAMC
427    CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
428    CADJ STORE gU(:,:,:,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, KappaRU,recip_HFacW,
433         U         gU,
434         I         myThid )
435    #ifdef    ALLOW_AUTODIFF_TAMC
436    CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
437    CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
438    #endif    /* ALLOW_AUTODIFF_TAMC */
439              CALL IMPLDIFF(
440         I         bi, bj, iMin, iMax, jMin, jMax,
441         I         deltaTmom, KappaRV,recip_HFacS,
442         U         gV,
443         I         myThid )
444            ENDIF
445    
446    #ifdef   ALLOW_OBCS
447    C--      Apply open boundary conditions
448            IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
449               DO K=1,Nr
450                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
451               ENDDO
452          ENDIF          ENDIF
453    #endif   /* ALLOW_OBCS */
454    
455    #ifdef    ALLOW_CD_CODE
456            IF (implicitViscosity.AND.useCDscheme) THEN
457    #ifdef    ALLOW_AUTODIFF_TAMC
458    CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
459    #endif    /* ALLOW_AUTODIFF_TAMC */
460              CALL IMPLDIFF(
461         I         bi, bj, iMin, iMax, jMin, jMax,
462         I         deltaTmom, KappaRU,recip_HFacW,
463         U         vVelD,
464         I         myThid )
465    #ifdef    ALLOW_AUTODIFF_TAMC
466    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
467    #endif    /* ALLOW_AUTODIFF_TAMC */
468              CALL IMPLDIFF(
469         I         bi, bj, iMin, iMax, jMin, jMax,
470         I         deltaTmom, KappaRV,recip_HFacS,
471         U         uVelD,
472         I         myThid )
473            ENDIF
474    #endif    /* ALLOW_CD_CODE */
475    C--     End implicit Vertical advection & viscosity
476    
477         ENDDO         ENDDO
478        ENDDO        ENDDO
479    
480  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),  Cml(
481  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))  C     In order to compare the variance of phiHydLow of a p/z-coordinate
482  C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),  C     run with etaH of a z/p-coordinate run the drift of phiHydLow
483  C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)  C     has to be removed by something like the following subroutine:
484  C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),  C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
485  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)  C     &                'phiHydLow', myThid )
486  C     write(0,*) 'dynamics: wVel(1) ',  Cml)
487  C    &            minval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.),  
488  C    &            maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.)  #ifdef ALLOW_DEBUG
489  C     write(0,*) 'dynamics: wVel(2) ',        If ( debugLevel .GE. debLevB ) THEN
490  C    &            minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.),         CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
491  C    &            maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.)         CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
492  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),         CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
493  cblk &                           maxval(K13(1:sNx,1:sNy,:))         CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
494  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),         CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
495  cblk &                           maxval(K23(1:sNx,1:sNy,:))         CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
496  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),         CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
497  cblk &                           maxval(K33(1:sNx,1:sNy,:))         CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
498  C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
499  C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
500  C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
501  C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
502  C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),         CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
503  C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))         CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
504  C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),        ENDIF
505  C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))  #endif
 C     write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.),  
 C    &                           maxval(pH/(Gravity*Rhonil))  
506    
507        RETURN        RETURN
508        END        END

Legend:
Removed from v.1.23  
changed lines
  Added in v.1.106

  ViewVC Help
Powered by ViewVC 1.1.22