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

Legend:
Removed from v.1.13  
changed lines
  Added in v.1.125

  ViewVC Help
Powered by ViewVC 1.1.22