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

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

  ViewVC Help
Powered by ViewVC 1.1.22