/[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.14 by cnh, Mon Jun 8 21:43:01 1998 UTC revision 1.120 by jmc, Mon Jul 11 19:30:42 2005 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_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
117  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
118        _RL myTime        _RL myTime
119        INTEGER myIter        INTEGER myIter
120          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, wTrans - Per block temporaries holding flow transport  C                                    is "pipelined" in the vertical
126  C     wVel                     o uTrans: Zonal transport  C                                    so we need an fVer for each
127  C                              o vTrans: Meridional transport  C                                    variable.
128  C                              o wTrans: Vertical transport  C     phiHydC    :: hydrostatic potential anomaly at cell center
129  C                              o wVel:   Vertical velocity at upper and lower  C                   In z coords phiHyd is the hydrostatic potential
130  C                                        cell faces.  C                      (=pressure/rho0) anomaly
131  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C                   In p coords phiHyd is the geopotential height anomaly.
132  C                              o maskUp: land/water mask for W points  C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
133  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
134  C     mTerm, pTerm,            tendency equations.  C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
135  C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  C     phiSurfY             or geopotential (atmos) in X and Y direction
136  C                              o xTerm: Mixing term  C     guDissip   :: dissipation tendency (all explicit terms), u component
137  C                              o cTerm: Coriolis term  C     gvDissip   :: dissipation tendency (all explicit terms), v component
138  C                              o mTerm: Metric term  C     iMin, iMax     - Ranges and sub-block indices on which calculations
139  C                              o pTerm: Pressure term  C     jMin, jMax       are applied.
 C                              o fZon: Zonal flux term  
 C                              o fMer: Meridional flux term  
 C                              o fVer: Vertical flux term - note fVer  
 C                                      is "pipelined" in the vertical  
 C                                      so we need an fVer for each  
 C                                      variable.  
 C     iMin, iMax - Ranges and sub-block indices on which calculations  
 C     jMin, jMax   are applied.  
140  C     bi, bj  C     bi, bj
141  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
142  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
143  C                          into fVerTerm  C                      index into fVerTerm.
144        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
145        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
149        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
155        _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 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 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)  
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
162    
163    #ifdef ALLOW_DIAGNOSTICS
164          _RL tmpFac
165    #endif /* ALLOW_DIAGNOSTICS */
166    
167    
168  C---    The algorithm...  C---    The algorithm...
169  C  C
170  C       "Correction Step"  C       "Correction Step"
# Line 116  C       "Calculation of Gs" Line 179  C       "Calculation of Gs"
179  C       ===================  C       ===================
180  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
181  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         w = sum_z ( div. u[n] )  
182  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
183    C         b   = b(rho, theta)
184  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
185  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
186  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
187  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
188  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
189  C  C
190  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
191  C       ================================  C       ================================
# Line 145  C         salt* = salt[n] + dt x ( 3/2 G Line 208  C         salt* = salt[n] + dt x ( 3/2 G
208  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
209  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
210  C---  C---
211    CEOP
212    
213  C--   Set up work arrays with valid (i.e. not NaN) values  C-- Call to routine for calculation of
214  C     These inital values do not alter the numerical results. They  C   Eliassen-Palm-flux-forced U-tendency,
215  C     just ensure that all memory references are to valid floating  C   if desired:
216  C     point numbers. This prevents spurious hardware signals due to  #ifdef INCLUDE_EP_FORCING_CODE
217  C     uninitialised but inert locations.        CALL CALC_EP_FORCING(myThid)
218        DO j=1-OLy,sNy+OLy  #endif
219         DO i=1-OLx,sNx+OLx  
220          xA(i,j)      = 0. _d 0  #ifdef ALLOW_AUTODIFF_TAMC
221          yA(i,j)      = 0. _d 0  C--   HPF directive to help TAMC
222          uTrans(i,j)  = 0. _d 0  CHPF$ INDEPENDENT
223          vTrans(i,j)  = 0. _d 0  #endif /* ALLOW_AUTODIFF_TAMC */
         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  
         rhokp1(i,j)  = 0. _d 0  
         rhotmp(i,j)  = 0. _d 0  
        ENDDO  
       ENDDO  
224    
225        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
226    
227    #ifdef ALLOW_AUTODIFF_TAMC
228    C--    HPF directive to help TAMC
229    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
230    CHPF$&                  ,phiHydF
231    CHPF$&                  ,KappaRU,KappaRV
232    CHPF$&                  )
233    #endif /* ALLOW_AUTODIFF_TAMC */
234    
235         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
236    
237  C--     Set up work arrays that need valid initial values  #ifdef ALLOW_AUTODIFF_TAMC
238          DO j=1-OLy,sNy+OLy            act1 = bi - myBxLo(myThid)
239           DO i=1-OLx,sNx+OLx            max1 = myBxHi(myThid) - myBxLo(myThid) + 1
240            wTrans(i,j)  = 0. _d 0            act2 = bj - myByLo(myThid)
241            wVel  (i,j,1) = 0. _d 0            max2 = myByHi(myThid) - myByLo(myThid) + 1
242            wVel  (i,j,2) = 0. _d 0            act3 = myThid - 1
243            fVerT(i,j,1) = 0. _d 0            max3 = nTx*nTy
244            fVerT(i,j,2) = 0. _d 0            act4 = ikey_dynamics - 1
245            fVerS(i,j,1) = 0. _d 0            idynkey = (act1 + 1) + act2*max1
246            fVerS(i,j,2) = 0. _d 0       &                      + act3*max1*max2
247            fVerU(i,j,1) = 0. _d 0       &                      + act4*max1*max2*max3
248            fVerU(i,j,2) = 0. _d 0  #endif /* ALLOW_AUTODIFF_TAMC */
           fVerV(i,j,1) = 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) = 0. _d 0  
          ENDDO  
         ENDDO  
249    
250          iMin = 1-OLx+1  C--   Set up work arrays with valid (i.e. not NaN) values
251          iMax = sNx+OLx  C     These inital values do not alter the numerical results. They
252          jMin = 1-OLy+1  C     just ensure that all memory references are to valid floating
253          jMax = sNy+OLy  C     point numbers. This prevents spurious hardware signals due to
254    C     uninitialised but inert locations.
 C--     Calculate gradient of surface pressure  
         CALL GRAD_PSURF(  
      I       bi,bj,iMin,iMax,jMin,jMax,  
      O       pSurfX,pSurfY,  
      I       myThid)  
   
 C--     Update fields in top level according to tendency terms  
         CALL CORRECTION_STEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)  
   
 C--     Density of 1st level (below W(1)) reference to level 1  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,  
      O     rhoKm1,  
      I     myThid )  
 C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
         CALL CALC_PH(  
      I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,  
      U      pH,  
      I      myThid )  
         DO J=1-Oly,sNy+Oly  
          DO I=1-Olx,sNx+Olx  
           rhoKp1(I,J)=rhoKm1(I,J)  
          ENDDO  
         ENDDO  
255    
256          DO K=2,Nz          DO k=1,Nr
257  C--     Update fields in Kth level according to tendency terms           DO j=1-OLy,sNy+OLy
258          CALL CORRECTION_STEP(            DO i=1-OLx,sNx+OLx
259       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)             KappaRU(i,j,k) = 0. _d 0
260  C--     Density of K-1 level (above W(K)) reference to K-1 level             KappaRV(i,j,k) = 0. _d 0
261  copt    CALL FIND_RHO(  #ifdef ALLOW_AUTODIFF_TAMC
262  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,  cph(
263  copt O     rhoKm1,  c--   need some re-initialisation here to break dependencies
264  copt I     myThid )  cph)
265  C       rhoKm1=rhoKp1             gu(i,j,k,bi,bj) = 0. _d 0
266          DO J=1-Oly,sNy+Oly             gv(i,j,k,bi,bj) = 0. _d 0
267           DO I=1-Olx,sNx+Olx  #endif
268            rhoKm1(I,J)=rhoKp1(I,J)            ENDDO
269           ENDDO           ENDDO
270          ENDDO          ENDDO
 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     rhoKp1,  
      I     myThid )  
 C--     Density of K-1 level (above W(K)) reference to K level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,  
      O     rhotmp,  
      I     myThid )  
 C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
         CALL CALC_ISOSLOPES(  
      I            bi, bj, iMin, iMax, jMin, jMax, K,  
      I            rhoKm1, rhoKp1, rhotmp,  
      O            K13, K23, K33, KapGM,  
      I            myThid )  
 C--     Calculate static stability and mix where convectively unstable  
         CALL CONVECT(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhotmp,rhoKp1,  
      I      myTime,myIter,myThid)  
 C--     Density of K-1 level (above W(K)) reference to K-1 level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,  
      O     rhoKm1,  
      I     myThid )  
 C--     Density of K level (below W(K)) referenced to K level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
      O     rhoKp1,  
      I     myThid )  
 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,rhoKp1,  
      U      pH,  
      I      myThid )  
   
         ENDDO ! K  
   
 C--     Initial boundary condition on barotropic divergence integral  
271          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
272           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
273            cg2d_b(i,j,bi,bj) = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
274              fVerU  (i,j,2) = 0. _d 0
275              fVerV  (i,j,1) = 0. _d 0
276              fVerV  (i,j,2) = 0. _d 0
277              phiHydF (i,j)  = 0. _d 0
278              phiHydC (i,j)  = 0. _d 0
279              dPhiHydX(i,j)  = 0. _d 0
280              dPhiHydY(i,j)  = 0. _d 0
281              phiSurfX(i,j)  = 0. _d 0
282              phiSurfY(i,j)  = 0. _d 0
283              guDissip(i,j)  = 0. _d 0
284              gvDissip(i,j)  = 0. _d 0
285           ENDDO           ENDDO
286          ENDDO          ENDDO
287    
288          DO K = Nz, 1, -1  C--     Start computation of dynamics
289           kM1  =max(1,k-1)   ! Points to level above k (=k-1)          iMin = 0
290           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above          iMax = sNx+1
291           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer          jMin = 0
292           iMin = 1-OLx+2          jMax = sNy+1
293           iMax = sNx+OLx-1  
294           jMin = 1-OLy+2  #ifdef ALLOW_AUTODIFF_TAMC
295           jMax = sNy+OLy-1  CADJ STORE wvel (:,:,:,bi,bj) =
296    CADJ &     comlev1_bibj, key = idynkey, byte = isbyte
297  C--      Get temporary terms used by tendency routines  #endif /* ALLOW_AUTODIFF_TAMC */
298           CALL CALC_COMMON_FACTORS (  
299       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
300       O        xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp,  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
301       I        myThid)          IF (implicSurfPress.NE.1.) THEN
302              CALL CALC_GRAD_PHI_SURF(
303         I         bi,bj,iMin,iMax,jMin,jMax,
304         I         etaN,
305         O         phiSurfX,phiSurfY,
306         I         myThid )                        
307            ENDIF
308    
309    #ifdef ALLOW_AUTODIFF_TAMC
310    CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
311    CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
312    #ifdef ALLOW_KPP
313    CADJ STORE KPPviscAz (:,:,:,bi,bj)
314    CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
315    #endif /* ALLOW_KPP */
316    #endif /* ALLOW_AUTODIFF_TAMC */
317    
318    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
319  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
320           CALL CALC_DIFFUSIVITY(          DO k=1,Nr
321       I        bi,bj,iMin,iMax,jMin,jMax,K,           CALL CALC_VISCOSITY(
322       I        maskC,maskUp,KapGM,K33,       I        bi,bj,iMin,iMax,jMin,jMax,k,
323       O        KappaZT,       O        KappaRU,KappaRV,
324       I        myThid)       I        myThid)
325           ENDDO
326    #endif
327    
328  C--      Calculate accelerations in the momentum equations  #ifdef ALLOW_AUTODIFF_TAMC
329    CADJ STORE KappaRU(:,:,:)
330    CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
331    CADJ STORE KappaRV(:,:,:)
332    CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
333    #endif /* ALLOW_AUTODIFF_TAMC */
334    
335    C--     Start of dynamics loop
336            DO k=1,Nr
337    
338    C--       km1    Points to level above k (=k-1)
339    C--       kup    Cycles through 1,2 to point to layer above
340    C--       kDown  Cycles through 2,1 to point to current layer
341    
342              km1  = MAX(1,k-1)
343              kp1  = MIN(k+1,Nr)
344              kup  = 1+MOD(k+1,2)
345              kDown= 1+MOD(k,2)
346    
347    #ifdef ALLOW_AUTODIFF_TAMC
348             kkey = (idynkey-1)*Nr + k
349    c
350    CADJ STORE totphihyd (:,:,k,bi,bj)
351    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
352    CADJ STORE theta (:,:,k,bi,bj)
353    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
354    CADJ STORE salt  (:,:,k,bi,bj)
355    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
356    #endif /* ALLOW_AUTODIFF_TAMC */
357    
358    C--      Integrate hydrostatic balance for phiHyd with BC of
359    C        phiHyd(z=0)=0
360             CALL CALC_PHI_HYD(
361         I        bi,bj,iMin,iMax,jMin,jMax,k,
362         I        theta, salt,
363         U        phiHydF,
364         O        phiHydC, dPhiHydX, dPhiHydY,
365         I        myTime, myIter, myThid )
366    
367    C--      Calculate accelerations in the momentum equations (gU, gV, ...)
368    C        and step forward storing the result in gU, gV, etc...
369           IF ( momStepping ) THEN           IF ( momStepping ) THEN
370            CALL CALC_MOM_RHS(  #ifdef ALLOW_MOM_FLUXFORM
371       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,             IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
372       I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
373       I         pH,       I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
374       U         aTerm,xTerm,cTerm,mTerm,pTerm,       U         fVerU, fVerV,
375       U         fZon, fMer, fVerU, fVerV,       I         myTime, myIter, myThid)
376       I         myThid)  #endif
377           ENDIF  #ifdef ALLOW_MOM_VECINV
378               IF (vectorInvariantMomentum) CALL MOM_VECINV(
379         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
380         I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
381         U         fVerU, fVerV,
382         O         guDissip, gvDissip,
383         I         myTime, myIter, myThid)
384    #endif
385               CALL TIMESTEP(
386         I         bi,bj,iMin,iMax,jMin,jMax,k,
387         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
388         I         guDissip, gvDissip,
389         I         myTime, myIter, myThid)
390    
391    #ifdef   ALLOW_OBCS
392    C--      Apply open boundary conditions
393               IF (useOBCS) THEN
394                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
395               ENDIF
396    #endif   /* ALLOW_OBCS */
397    
 C--      Calculate active tracer tendencies  
          IF ( tempStepping ) THEN  
           CALL CALC_GT(  
      I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  
      I         xA,yA,uTrans,vTrans,wTrans,maskUp,  
      I         K13,K23,KappaZT,KapGM,  
      U         aTerm,xTerm,fZon,fMer,fVerT,  
      I         myThid)  
398           ENDIF           ENDIF
399  Cdbg     CALL CALC_GS(  
400  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  
401  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  C--     end of dynamics k loop (1:Nr)
402  Cdbg I        K13,K23,K33,KapGM,          ENDDO
403  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  
404  Cdbg I        myThid)  C--     Implicit Vertical advection & viscosity
405    #ifdef INCLUDE_IMPLVERTADV_CODE
406  C--      Prediction step (step forward all model variables)          IF ( momImplVertAdv ) THEN
407           CALL TIMESTEP(            CALL MOM_U_IMPLICIT_R( kappaRU,
408       I       bi,bj,iMin,iMax,jMin,jMax,K,       I                           bi, bj, myTime, myIter, myThid )
409       I       myThid)            CALL MOM_V_IMPLICIT_R( kappaRV,
410         I                           bi, bj, myTime, myIter, myThid )
411  C--      Diagnose barotropic divergence of predicted fields          ELSEIF ( implicitViscosity ) THEN
412           CALL DIV_G(  #else /* INCLUDE_IMPLVERTADV_CODE */
413       I       bi,bj,iMin,iMax,jMin,jMax,K,          IF     ( implicitViscosity ) THEN
414       I       xA,yA,  #endif /* INCLUDE_IMPLVERTADV_CODE */
415       I       myThid)  #ifdef    ALLOW_AUTODIFF_TAMC
416    CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
417          ENDDO ! K  CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
418    #endif    /* ALLOW_AUTODIFF_TAMC */
419  C--     Implicit diffusion            CALL IMPLDIFF(
420          IF (implicitDiffusion) THEN       I         bi, bj, iMin, iMax, jMin, jMax,
421           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,       I         0, KappaRU,recip_HFacW,
422       I                  KappaZT,       U         gU,
423       I                  myThid )       I         myThid )
424    #ifdef    ALLOW_AUTODIFF_TAMC
425    CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
426    CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
427    #endif    /* ALLOW_AUTODIFF_TAMC */
428              CALL IMPLDIFF(
429         I         bi, bj, iMin, iMax, jMin, jMax,
430         I         0, KappaRV,recip_HFacS,
431         U         gV,
432         I         myThid )
433            ENDIF
434    
435    #ifdef   ALLOW_OBCS
436    C--      Apply open boundary conditions
437            IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
438               DO K=1,Nr
439                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
440               ENDDO
441            ENDIF
442    #endif   /* ALLOW_OBCS */
443    
444    #ifdef    ALLOW_CD_CODE
445            IF (implicitViscosity.AND.useCDscheme) THEN
446    #ifdef    ALLOW_AUTODIFF_TAMC
447    CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
448    #endif    /* ALLOW_AUTODIFF_TAMC */
449              CALL IMPLDIFF(
450         I         bi, bj, iMin, iMax, jMin, jMax,
451         I         0, KappaRU,recip_HFacW,
452         U         vVelD,
453         I         myThid )
454    #ifdef    ALLOW_AUTODIFF_TAMC
455    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
456    #endif    /* ALLOW_AUTODIFF_TAMC */
457              CALL IMPLDIFF(
458         I         bi, bj, iMin, iMax, jMin, jMax,
459         I         0, KappaRV,recip_HFacS,
460         U         uVelD,
461         I         myThid )
462          ENDIF          ENDIF
463    #endif    /* ALLOW_CD_CODE */
464    C--     End implicit Vertical advection & viscosity
465    
466         ENDDO         ENDDO
467        ENDDO        ENDDO
468    
469        write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)  #ifdef ALLOW_OBCS
470        write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),        IF (useOBCS) THEN
471       &                         maxval(uVel(1:sNx,1:sNy,:,:,:))         CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
472        write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),        ENDIF
473       &                         maxval(vVel(1:sNx,1:sNy,:,:,:))  #endif
474        write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  
475       &                         maxval(K13(1:sNx,1:sNy,:))  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
476        write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  
477       &                         maxval(K23(1:sNx,1:sNy,:))  Cml(
478        write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  C     In order to compare the variance of phiHydLow of a p/z-coordinate
479       &                         maxval(K33(1:sNx,1:sNy,:))  C     run with etaH of a z/p-coordinate run the drift of phiHydLow
480        write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),  C     has to be removed by something like the following subroutine:
481       &                         maxval(gT(1:sNx,1:sNy,:,:,:))  C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
482        write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),  C     &                'phiHydLow', myThid )
483       &                         maxval(Theta(1:sNx,1:sNy,:,:,:))  Cml)
484        write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),  
485       &                          maxval(pH/(Gravity*Rhonil))  #ifdef ALLOW_DIAGNOSTICS
486          IF ( usediagnostics ) THEN
487    
488           CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD  ',0,Nr,0,1,1,myThid)
489           CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT  ',0, 1,0,1,1,myThid)
490    
491           tmpFac = 1. _d 0
492           CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
493         &                                 'PHIHYDSQ',0,Nr,0,1,1,myThid)
494    
495           CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
496         &                                 'PHIBOTSQ',0, 1,0,1,1,myThid)
497    
498          ENDIF
499    #endif /* ALLOW_DIAGNOSTICS */
500          
501    #ifdef ALLOW_DEBUG
502          If ( debugLevel .GE. debLevB ) THEN
503           CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
504           CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
505           CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
506           CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
507           CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
508           CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
509           CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
510           CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
511           CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
512           CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
513    #ifndef ALLOW_ADAMSBASHFORTH_3
514           CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
515           CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
516           CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
517           CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
518    #endif
519          ENDIF
520    #endif
521    
522        RETURN        RETURN
523        END        END

Legend:
Removed from v.1.14  
changed lines
  Added in v.1.120

  ViewVC Help
Powered by ViewVC 1.1.22