/[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.5 by adcroft, Mon May 4 16:32:10 1998 UTC revision 1.98.2.1 by edhill, Thu Oct 2 18:10:45 2003 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"
       SUBROUTINE DYNAMICS(myThid)  
 C     /==========================================================\  
 C     | SUBROUTINE DYNAMICS                                      |  
 C     | o Controlling routine for the explicit part of the model |  
 C     |   dynamics.                                              |  
 C     |==========================================================|  
 C     | This routine evaluates the "dynamics" terms for each     |  
 C     | block of ocean in turn. Because the blocks of ocean have |  
 C     | overlap regions they are independent of one another.     |  
 C     | If terms involving lateral integrals are needed in this  |  
 C     | routine care will be needed. Similarly finite-difference |  
 C     | operations with stencils wider than the overlap region   |  
 C     | require special consideration.                           |  
 C     | Notes                                                    |  
 C     | =====                                                    |  
 C     | C*P* comments indicating place holders for which code is |  
 C     |      presently being developed.                          |  
 C     \==========================================================/  
6    
7    CBOP
8    C     !ROUTINE: DYNAMICS
9    C     !INTERFACE:
10          SUBROUTINE DYNAMICS(myTime, myIter, myThid)
11    C     !DESCRIPTION: \bv
12    C     *==========================================================*
13    C     | SUBROUTINE DYNAMICS                                      
14    C     | o Controlling routine for the explicit part of the model  
15    C     |   dynamics.                                              
16    C     *==========================================================*
17    C     | This routine evaluates the "dynamics" terms for each      
18    C     | block of ocean in turn. Because the blocks of ocean have  
19    C     | overlap regions they are independent of one another.      
20    C     | If terms involving lateral integrals are needed in this  
21    C     | routine care will be needed. Similarly finite-difference  
22    C     | operations with stencils wider than the overlap region    
23    C     | require special consideration.                            
24    C     | The algorithm...
25    C     |
26    C     | "Correction Step"
27    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"
74  #include "CG2D.h"  #include "PARAMS.h"
75  #include "DYNVARS.h"  #include "DYNVARS.h"
76    #include "GRID.h"
77    #ifdef ALLOW_PASSIVE_TRACER
78    #include "TR1.h"
79    #endif
80    #ifdef ALLOW_AUTODIFF_TAMC
81    # include "tamc.h"
82    # include "tamc_keys.h"
83    # include "FFIELDS.h"
84    # include "EOS.h"
85    # ifdef ALLOW_KPP
86    #  include "KPP.h"
87    # endif
88    #endif /* ALLOW_AUTODIFF_TAMC */
89    
90    C     !CALLING SEQUENCE:
91    C     DYNAMICS()
92    C      |
93    C      |-- CALC_GRAD_PHI_SURF
94    C      |
95    C      |-- CALC_VISCOSITY
96    C      |
97    C      |-- CALC_PHI_HYD  
98    C      |
99    C      |-- MOM_FLUXFORM  
100    C      |
101    C      |-- MOM_VECINV    
102    C      |
103    C      |-- TIMESTEP      
104    C      |
105    C      |-- OBCS_APPLY_UV
106    C      |
107    C      |-- IMPLDIFF      
108    C      |
109    C      |-- OBCS_APPLY_UV
110    C      |
111    C      |-- CALL TIMEAVE_CUMUL_1T
112    C      |-- CALL DEBUG_STATS_RL
113    
114    C     !INPUT/OUTPUT PARAMETERS:
115  C     == Routine arguments ==  C     == Routine arguments ==
116    C     myTime - Current time in simulation
117    C     myIter - Current iteration number in simulation
118  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
119          _RL myTime
120          INTEGER myIter
121        INTEGER myThid        INTEGER myThid
122    
123    C     !LOCAL VARIABLES:
124  C     == Local variables  C     == Local variables
125  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[STUV]               o fVer: Vertical flux term - note fVer
 C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  
 C                              o uTrans: Zonal transport  
 C                              o vTrans: Meridional transport  
 C                              o wTrans: Vertical transport  
 C     maskC,maskUp             o maskC: land/water mask for tracer cells  
 C                              o maskUp: land/water mask for W points  
 C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  
 C     mTerm, pTerm,            tendency equations.  
 C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  
 C                              o xTerm: Mixing term  
 C                              o cTerm: Coriolis term  
 C                              o mTerm: Metric term  
 C                              o pTerm: Pressure term  
 C                              o fZon: Zonal flux term  
 C                              o fMer: Meridional flux term  
 C                              o fVer: Vertical flux term - note fVer  
126  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
127  C                                      so we need an fVer for each  C                                      so we need an fVer for each
128  C                                      variable.  C                                      variable.
129  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     phiHydC    :: hydrostatic potential anomaly at cell center
130  C     jMin, jMax   are applied.  C                   In z coords phiHyd is the hydrostatic potential
131    C                      (=pressure/rho0) anomaly
132    C                   In p coords phiHyd is the geopotential height anomaly.
133    C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
134    C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
135    C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
136    C     phiSurfY             or geopotential (atmos) in X and Y direction
137    C     iMin, iMax     - Ranges and sub-block indices on which calculations
138    C     jMin, jMax       are applied.
139  C     bi, bj  C     bi, bj
140  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
141  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
142  C                          into fVerTerm  C                      index into fVerTerm.
143        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
144        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
145        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
148        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
149        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
152        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
153        _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 pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
154        INTEGER iMin, iMax        INTEGER iMin, iMax
155        INTEGER jMin, jMax        INTEGER jMin, jMax
156        INTEGER bi, bj        INTEGER bi, bj
157        INTEGER i, j        INTEGER i, j
158        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kp1, kup, kDown
159    
160          LOGICAL  DIFFERENT_MULTIPLE
161          EXTERNAL DIFFERENT_MULTIPLE
162    
163    C---    The algorithm...
164    C
165    C       "Correction Step"
166    C       =================
167    C       Here we update the horizontal velocities with the surface
168    C       pressure such that the resulting flow is either consistent
169    C       with the free-surface evolution or the rigid-lid:
170    C         U[n] = U* + dt x d/dx P
171    C         V[n] = V* + dt x d/dy P
172    C
173    C       "Calculation of Gs"
174    C       ===================
175    C       This is where all the accelerations and tendencies (ie.
176    C       physics, parameterizations etc...) are calculated
177    C         rho = rho ( theta[n], salt[n] )
178    C         b   = b(rho, theta)
179    C         K31 = K31 ( rho )
180    C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
181    C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
182    C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
183    C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
184    C
185    C       "Time-stepping" or "Prediction"
186    C       ================================
187    C       The models variables are stepped forward with the appropriate
188    C       time-stepping scheme (currently we use Adams-Bashforth II)
189    C       - For momentum, the result is always *only* a "prediction"
190    C       in that the flow may be divergent and will be "corrected"
191    C       later with a surface pressure gradient.
192    C       - Normally for tracers the result is the new field at time
193    C       level [n+1} *BUT* in the case of implicit diffusion the result
194    C       is also *only* a prediction.
195    C       - We denote "predictors" with an asterisk (*).
196    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
197    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
198    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
199    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
200    C       With implicit diffusion:
201    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
202    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
203    C         (1 + dt * K * d_zz) theta[n] = theta*
204    C         (1 + dt * K * d_zz) salt[n] = salt*
205    C---
206    CEOP
207    
208    C-- Call to routine for calculation of
209    C   Eliassen-Palm-flux-forced U-tendency,
210    C   if desired:
211    #ifdef INCLUDE_EP_FORCING_CODE
212          CALL CALC_EP_FORCING(myThid)
213    #endif
214    
215    #ifdef ALLOW_AUTODIFF_TAMC
216    C--   HPF directive to help TAMC
217    CHPF$ INDEPENDENT
218    #endif /* ALLOW_AUTODIFF_TAMC */
219    
220          DO bj=myByLo(myThid),myByHi(myThid)
221    
222    #ifdef ALLOW_AUTODIFF_TAMC
223    C--    HPF directive to help TAMC
224    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
225    CHPF$&                  ,phiHydF
226    CHPF$&                  ,KappaRU,KappaRV
227    CHPF$&                  )
228    #endif /* ALLOW_AUTODIFF_TAMC */
229    
230           DO bi=myBxLo(myThid),myBxHi(myThid)
231    
232    #ifdef ALLOW_AUTODIFF_TAMC
233              act1 = bi - myBxLo(myThid)
234              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
235              act2 = bj - myByLo(myThid)
236              max2 = myByHi(myThid) - myByLo(myThid) + 1
237              act3 = myThid - 1
238              max3 = nTx*nTy
239              act4 = ikey_dynamics - 1
240              idynkey = (act1 + 1) + act2*max1
241         &                      + act3*max1*max2
242         &                      + act4*max1*max2*max3
243    #endif /* ALLOW_AUTODIFF_TAMC */
244    
245  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
246  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
247  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
248  C     point numbers. This prevents spurious hardware signals due to  C     point numbers. This prevents spurious hardware signals due to
249  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  
         ENDDO  
         rhokm1(i,j)  = 0. _d 0  
         rhokp1(i,j)  = 0. _d 0  
        ENDDO  
       ENDDO  
 C--   Set up work arrays that need valid initial values  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         wTrans(i,j)  = 0. _d 0  
         fVerT(i,j,1) = 0. _d 0  
         fVerT(i,j,2) = 0. _d 0  
         fVerS(i,j,1) = 0. _d 0  
         fVerS(i,j,2) = 0. _d 0  
         fVerU(i,j,1) = 0. _d 0  
         fVerU(i,j,2) = 0. _d 0  
         fVerV(i,j,1) = 0. _d 0  
         fVerV(i,j,2) = 0. _d 0  
        ENDDO  
       ENDDO  
   
       DO bj=myByLo(myThid),myByHi(myThid)  
        DO bi=myBxLo(myThid),myBxHi(myThid)  
250    
251  C--   Boundary condition on hydrostatic pressure is pH(z=0)=0          DO k=1,Nr
252             DO j=1-OLy,sNy+OLy
253              DO i=1-OLx,sNx+OLx
254               KappaRU(i,j,k) = 0. _d 0
255               KappaRV(i,j,k) = 0. _d 0
256    #ifdef ALLOW_AUTODIFF_TAMC
257    cph(
258    c--   need some re-initialisation here to break dependencies
259    c--   totphihyd is assumed zero from ini_pressure, i.e.
260    c--   avoiding iterate pressure p = integral of (g*rho(p)*dz)
261    cph)
262               totPhiHyd(i,j,k,bi,bj) = 0. _d 0
263               gu(i,j,k,bi,bj) = 0. _d 0
264               gv(i,j,k,bi,bj) = 0. _d 0
265    #endif
266              ENDDO
267             ENDDO
268            ENDDO
269          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
270           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
271            pH(i,j,1) = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
272              fVerU  (i,j,2) = 0. _d 0
273              fVerV  (i,j,1) = 0. _d 0
274              fVerV  (i,j,2) = 0. _d 0
275              phiHydF (i,j)  = 0. _d 0
276              phiHydC (i,j)  = 0. _d 0
277              dPhiHydX(i,j)  = 0. _d 0
278              dPhiHydY(i,j)  = 0. _d 0
279              phiSurfX(i,j)  = 0. _d 0
280              phiSurfY(i,j)  = 0. _d 0
281           ENDDO           ENDDO
282          ENDDO          ENDDO
283    
284          iMin = 1-OLx+1  C--     Start computation of dynamics
285          iMax = sNx+OLx          iMin = 0
286          jMin = 1-OLy+1          iMax = sNx+1
287          jMax = sNy+OLy          jMin = 0
288            jMax = sNy+1
289  C--     Calculate gradient of surface pressure  
290          CALL GRAD_PSURF(  #ifdef ALLOW_AUTODIFF_TAMC
291       I       bi,bj,iMin,iMax,jMin,jMax,  CADJ STORE wvel (:,:,:,bi,bj) =
292       O       pSurfX,pSurfY,  CADJ &     comlev1_bibj, key = idynkey, byte = isbyte
293       I       myThid)  #endif /* ALLOW_AUTODIFF_TAMC */
294    
295  C--     Update fields in top level according to tendency terms  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
296          CALL TIMESTEP(  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
297       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)          IF (implicSurfPress.NE.1.) THEN
298              CALL CALC_GRAD_PHI_SURF(
299  C Density of 1st level (below W(1)) reference to level 1       I         bi,bj,iMin,iMax,jMin,jMax,
300           CALL FIND_RHO(       I         etaN,
301       I      bi, bj, iMin, iMax, jMin, jMax, 1, 1, 'LINEAR',       O         phiSurfX,phiSurfY,
302       O      rhoKm1,       I         myThid )                        
303       I      myThid )          ENDIF
304  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
305           CALL CALC_PH(  #ifdef ALLOW_AUTODIFF_TAMC
306       I       bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
307       U       pH,  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
308       I       myThid )  #ifdef ALLOW_KPP
309    CADJ STORE KPPviscAz (:,:,:,bi,bj)
310          DO K=2,Nz  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
311  C--     Update fields in Kth level according to tendency terms  #endif /* ALLOW_KPP */
312          CALL TIMESTEP(  #endif /* ALLOW_AUTODIFF_TAMC */
313       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)  
314  C Density of K-1 level (above W(K)) reference to K level  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
315           CALL FIND_RHO(  C--      Calculate the total vertical diffusivity
316       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',          DO k=1,Nr
317       O      rhoKm1,           CALL CALC_VISCOSITY(
318       I      myThid )       I        bi,bj,iMin,iMax,jMin,jMax,k,
319  C Density of K level (below W(K)) reference to K level       O        KappaRU,KappaRV,
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',  
      O      rhoKp1,  
      I      myThid )  
 C--     Calculate static stability and mix where convectively unstable  
          CALL CONVECT(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,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, 'LINEAR',  
      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, 'LINEAR',  
      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  
   
         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,maskC,maskUp,  
320       I        myThid)       I        myThid)
321           ENDDO
322    #endif
323    
324  C--      Calculate accelerations in the momentum equations  C--     Start of dynamics loop
325           CALL CALC_MOM_RHS(          DO k=1,Nr
326       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
327       I        xA,yA,uTrans,vTrans,wTrans,maskC,  C--       km1    Points to level above k (=k-1)
328       I        pH,  C--       kup    Cycles through 1,2 to point to layer above
329       U        aTerm,xTerm,cTerm,mTerm,pTerm,  C--       kDown  Cycles through 2,1 to point to current layer
330       U        fZon, fMer, fVerU, fVerV,  
331       I        myThid)            km1  = MAX(1,k-1)
332              kp1  = MIN(k+1,Nr)
333              kup  = 1+MOD(k+1,2)
334              kDown= 1+MOD(k,2)
335    
336    #ifdef ALLOW_AUTODIFF_TAMC
337             kkey = (idynkey-1)*Nr + k
338    CADJ STORE totphihyd (:,:,k,bi,bj)
339    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
340    #endif /* ALLOW_AUTODIFF_TAMC */
341    
342    C--      Integrate hydrostatic balance for phiHyd with BC of
343    C        phiHyd(z=0)=0
344    C        distinguishe between Stagger and Non Stagger time stepping
345             IF (staggerTimeStep) THEN
346               CALL CALC_PHI_HYD(
347         I        bi,bj,iMin,iMax,jMin,jMax,k,
348         I        gT, gS,
349         U        phiHydF,
350         O        phiHydC, dPhiHydX, dPhiHydY,
351         I        myTime, myIter, myThid )
352             ELSE
353               CALL CALC_PHI_HYD(
354         I        bi,bj,iMin,iMax,jMin,jMax,k,
355         I        theta, salt,
356         U        phiHydF,
357         O        phiHydC, dPhiHydX, dPhiHydY,
358         I        myTime, myIter, myThid )
359             ENDIF
360    
361    C--      Calculate accelerations in the momentum equations (gU, gV, ...)
362    C        and step forward storing the result in gU, gV, etc...
363             IF ( momStepping ) THEN
364    #ifndef DISABLE_MOM_FLUXFORM
365               IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
366         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
367         I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
368         U         fVerU, fVerV,
369         I         myTime, myIter, myThid)
370    #endif
371    #ifndef DISABLE_MOM_VECINV
372               IF (vectorInvariantMomentum) CALL MOM_VECINV(
373         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
374         I         dPhiHydX,dPhiHydY,KappaRU,KappaRV,
375         U         fVerU, fVerV,
376         I         myTime, myIter, myThid)
377    #endif
378               CALL TIMESTEP(
379         I         bi,bj,iMin,iMax,jMin,jMax,k,
380         I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
381         I         myTime, myIter, myThid)
382    
383    #ifdef   ALLOW_OBCS
384    C--      Apply open boundary conditions
385               IF (useOBCS) THEN
386                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
387               ENDIF
388    #endif   /* ALLOW_OBCS */
389    
390             ENDIF
391    
 C--      Calculate active tracer tendencies  
          CALL CALC_GT(  
      I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  
      I        xA,yA,uTrans,vTrans,wTrans,maskUp,  
      U        aTerm,xTerm,fZon,fMer,fVerT,  
      I        myThid)  
 Cdbg     CALL CALC_GS(  
 Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  
 Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  
 Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  
 Cdbg I        myThid)  
392    
393          ENDDO ! K  C--     end of dynamics k loop (1:Nr)
394            ENDDO
395    
396    C--     Implicit viscosity
397            IF (implicitViscosity.AND.momStepping) THEN
398    #ifdef    ALLOW_AUTODIFF_TAMC
399    CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
400    #endif    /* ALLOW_AUTODIFF_TAMC */
401              CALL IMPLDIFF(
402         I         bi, bj, iMin, iMax, jMin, jMax,
403         I         deltaTmom, KappaRU,recip_HFacW,
404         U         gU,
405         I         myThid )
406    #ifdef    ALLOW_AUTODIFF_TAMC
407    CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
408    #endif    /* ALLOW_AUTODIFF_TAMC */
409              CALL IMPLDIFF(
410         I         bi, bj, iMin, iMax, jMin, jMax,
411         I         deltaTmom, KappaRV,recip_HFacS,
412         U         gV,
413         I         myThid )
414    
415    #ifdef   ALLOW_OBCS
416    C--      Apply open boundary conditions
417             IF (useOBCS) THEN
418               DO K=1,Nr
419                 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
420               ENDDO
421             END IF
422    #endif   /* ALLOW_OBCS */
423    
424    #ifdef    INCLUDE_CD_CODE
425    #ifdef    ALLOW_AUTODIFF_TAMC
426    CADJ STORE vVelD(:,:,:,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         deltaTmom, KappaRU,recip_HFacW,
431         U         vVelD,
432         I         myThid )
433    #ifdef    ALLOW_AUTODIFF_TAMC
434    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
435    #endif    /* ALLOW_AUTODIFF_TAMC */
436              CALL IMPLDIFF(
437         I         bi, bj, iMin, iMax, jMin, jMax,
438         I         deltaTmom, KappaRV,recip_HFacS,
439         U         uVelD,
440         I         myThid )
441    #endif    /* INCLUDE_CD_CODE */
442    C--     End If implicitViscosity.AND.momStepping
443            ENDIF
444    
445         ENDDO         ENDDO
446        ENDDO        ENDDO
447    
448    Cml(
449    C     In order to compare the variance of phiHydLow of a p/z-coordinate
450    C     run with etaH of a z/p-coordinate run the drift of phiHydLow
451    C     has to be removed by something like the following subroutine:
452    C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
453    C     &                'phiHydLow', myThid )
454    Cml)
455    
456    #ifndef DISABLE_DEBUGMODE
457          If ( debugLevel .GE. debLevB ) THEN
458           CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
459           CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
460           CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
461           CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
462           CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
463           CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
464           CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
465           CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
466           CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
467           CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
468           CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
469           CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
470           CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
471           CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
472          ENDIF
473    #endif
474    
475        RETURN        RETURN
476        END        END

Legend:
Removed from v.1.5  
changed lines
  Added in v.1.98.2.1

  ViewVC Help
Powered by ViewVC 1.1.22