/[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.3 by adcroft, Wed Apr 29 21:31:09 1998 UTC revision 1.42 by adcroft, Tue May 18 18:01:12 1999 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    
3  #include "CPP_EEOPTIONS.h"  #include "CPP_OPTIONS.h"
4    
5        SUBROUTINE DYNAMICS(myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6  C     /==========================================================\  C     /==========================================================\
7  C     | SUBROUTINE DYNAMICS                                      |  C     | SUBROUTINE DYNAMICS                                      |
8  C     | o Controlling routine for the explicit part of the model |  C     | o Controlling routine for the explicit part of the model |
# Line 20  C     | ===== Line 20  C     | =====
20  C     | C*P* comments indicating place holders for which code is |  C     | C*P* comments indicating place holders for which code is |
21  C     |      presently being developed.                          |  C     |      presently being developed.                          |
22  C     \==========================================================/  C     \==========================================================/
23          IMPLICIT NONE
24    
25  C     == Global variables ===  C     == Global variables ===
26  #include "SIZE.h"  #include "SIZE.h"
27  #include "EEPARAMS.h"  #include "EEPARAMS.h"
28  #include "CG2D.h"  #include "CG2D.h"
29    #include "PARAMS.h"
30  #include "DYNVARS.h"  #include "DYNVARS.h"
31    #include "GRID.h"
32    #ifdef ALLOW_KPP
33    #include "KPPMIX.h"
34    #endif
35    
36  C     == Routine arguments ==  C     == Routine arguments ==
37    C     myTime - Current time in simulation
38    C     myIter - Current iteration number in simulation
39  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
40        INTEGER myThid        INTEGER myThid
41          _RL myTime
42          INTEGER myIter
43    
44  C     == Local variables  C     == Local variables
45  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
46  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C     uTrans, vTrans, rTrans - Per block temporaries holding flow
47  C                              o uTrans: Zonal transport  C                              transport
48    C     rVel                     o uTrans: Zonal transport
49  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
50  C                              o wTrans: Vertical transport  C                              o rTrans: Vertical transport
51    C                              o rVel:   Vertical velocity at upper and
52    C                                        lower cell faces.
53  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
54  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
55  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in
# Line 52  C                              o fVer: V Line 65  C                              o fVer: V
65  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
66  C                                      so we need an fVer for each  C                                      so we need an fVer for each
67  C                                      variable.  C                                      variable.
68  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     rhoK, rhoKM1   - Density at current level, level above and level
69  C     jMin, jMax   are applied.  C                      below.
70    C     rhoKP1                                                                  
71    C     buoyK, buoyKM1 - Buoyancy at current level and level above.
72    C     phiHyd         - Hydrostatic part of the potential phiHydi.
73    C                      In z coords phiHydiHyd is the hydrostatic
74    C                      pressure anomaly
75    C                      In p coords phiHydiHyd is the geopotential
76    C                      surface height
77    C                      anomaly.
78    C     etaSurfX,      - Holds surface elevation gradient in X and Y.
79    C     etaSurfY
80    C     K13, K23, K33  - Non-zero elements of small-angle approximation
81    C                      diffusion tensor.
82    C     KapGM          - Spatially varying Visbeck et. al mixing coeff.
83    C     KappaRT,       - Total diffusion in vertical for T and S.
84    C     KappaRS          (background + spatially varying, isopycnal term).
85    C     iMin, iMax     - Ranges and sub-block indices on which calculations
86    C     jMin, jMax       are applied.
87  C     bi, bj  C     bi, bj
88  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
89  C                          are switched with layer to be the appropriate index  C     kDown, kM1       are switched with layer to be the appropriate
90  C                          into fVerTerm  C                      index into fVerTerm.
91        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
97        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
107        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
108        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
109        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
110        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112          _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113          _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114          _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115          _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116          _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117          _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118          _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119          _RL K13     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
120          _RL K23     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
121          _RL K33     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
122          _RL KapGM   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123          _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
124          _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
125          _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
126          _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
127    
128        INTEGER iMin, iMax        INTEGER iMin, iMax
129        INTEGER jMin, jMax        INTEGER jMin, jMax
130        INTEGER bi, bj        INTEGER bi, bj
131        INTEGER i, j        INTEGER i, j
132        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
133          LOGICAL BOTTOM_LAYER
134    
135    C---    The algorithm...
136    C
137    C       "Correction Step"
138    C       =================
139    C       Here we update the horizontal velocities with the surface
140    C       pressure such that the resulting flow is either consistent
141    C       with the free-surface evolution or the rigid-lid:
142    C         U[n] = U* + dt x d/dx P
143    C         V[n] = V* + dt x d/dy P
144    C
145    C       "Calculation of Gs"
146    C       ===================
147    C       This is where all the accelerations and tendencies (ie.
148    C       phiHydysics, parameterizations etc...) are calculated
149    C         rVel = sum_r ( div. u[n] )
150    C         rho = rho ( theta[n], salt[n] )
151    C         b   = b(rho, theta)
152    C         K31 = K31 ( rho )
153    C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
154    C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
155    C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
156    C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
157    C
158    C       "Time-stepping" or "Prediction"
159    C       ================================
160    C       The models variables are stepped forward with the appropriate
161    C       time-stepping scheme (currently we use Adams-Bashforth II)
162    C       - For momentum, the result is always *only* a "prediction"
163    C       in that the flow may be divergent and will be "corrected"
164    C       later with a surface pressure gradient.
165    C       - Normally for tracers the result is the new field at time
166    C       level [n+1} *BUT* in the case of implicit diffusion the result
167    C       is also *only* a prediction.
168    C       - We denote "predictors" with an asterisk (*).
169    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
170    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
171    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
172    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
173    C       With implicit diffusion:
174    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
175    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
176    C         (1 + dt * K * d_zz) theta[n] = theta*
177    C         (1 + dt * K * d_zz) salt[n] = salt*
178    C---
179    
180  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
181  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
# Line 92  C     point numbers. This prevents spuri Line 184  C     point numbers. This prevents spuri
184  C     uninitialised but inert locations.  C     uninitialised but inert locations.
185        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
186         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
187          xA(i,j)      = 0.*1. _d 37          xA(i,j)      = 0. _d 0
188          yA(i,j)      = 0.*1. _d 37          yA(i,j)      = 0. _d 0
189          uTrans(i,j)  = 0.*1. _d 37          uTrans(i,j)  = 0. _d 0
190          vTrans(i,j)  = 0.*1. _d 37          vTrans(i,j)  = 0. _d 0
191          aTerm(i,j)   = 0.*1. _d 37          aTerm(i,j)   = 0. _d 0
192          xTerm(i,j)   = 0.*1. _d 37          xTerm(i,j)   = 0. _d 0
193          cTerm(i,j)   = 0.*1. _d 37          cTerm(i,j)   = 0. _d 0
194          mTerm(i,j)   = 0.*1. _d 37          mTerm(i,j)   = 0. _d 0
195          pTerm(i,j)   = 0.*1. _d 37          pTerm(i,j)   = 0. _d 0
196          fZon(i,j)    = 0.*1. _d 37          fZon(i,j)    = 0. _d 0
197          fMer(i,j)    = 0.*1. _d 37          fMer(i,j)    = 0. _d 0
198          DO K=1,nZ          DO K=1,Nr
199           pH (i,j,k)  = 0.*1. _d 37           phiHyd (i,j,k)  = 0. _d 0
200             K13(i,j,k)  = 0. _d 0
201             K23(i,j,k)  = 0. _d 0
202             K33(i,j,k)  = 0. _d 0
203             KappaRT(i,j,k) = 0. _d 0
204             KappaRS(i,j,k) = 0. _d 0
205          ENDDO          ENDDO
206          rhokm1(i,j)    = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
207          rhokp1(i,j)    = 0. _d 0          rhok   (i,j) = 0. _d 0
208         ENDDO          rhoKP1 (i,j) = 0. _d 0
209        ENDDO          rhoTMP (i,j) = 0. _d 0
210  C--   Set up work arrays that need valid initial values          buoyKM1(i,j) = 0. _d 0
211        DO j=1-OLy,sNy+OLy          buoyK  (i,j) = 0. _d 0
212         DO i=1-OLx,sNx+OLx          maskC  (i,j) = 0. _d 0
         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  
213         ENDDO         ENDDO
214        ENDDO        ENDDO
215    
216    
217        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
218         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
219    
220  C--   Boundary condition on hydrostatic pressure is pH(z=0)=0  C--     Set up work arrays that need valid initial values
221          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
222           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
223            pH(i,j,1) = 0. _d 0            rTrans(i,j)   = 0. _d 0
224              rVel  (i,j,1) = 0. _d 0
225              rVel  (i,j,2) = 0. _d 0
226              fVerT (i,j,1) = 0. _d 0
227              fVerT (i,j,2) = 0. _d 0
228              fVerS (i,j,1) = 0. _d 0
229              fVerS (i,j,2) = 0. _d 0
230              fVerU (i,j,1) = 0. _d 0
231              fVerU (i,j,2) = 0. _d 0
232              fVerV (i,j,1) = 0. _d 0
233              fVerV (i,j,2) = 0. _d 0
234              phiHyd(i,j,1) = 0. _d 0
235              K13   (i,j,1) = 0. _d 0
236              K23   (i,j,1) = 0. _d 0
237              K33   (i,j,1) = 0. _d 0
238              KapGM (i,j)   = GMkbackground
239           ENDDO           ENDDO
240          ENDDO          ENDDO
241    
# Line 140  C--   Boundary condition on hydrostatic Line 244  C--   Boundary condition on hydrostatic
244          jMin = 1-OLy+1          jMin = 1-OLy+1
245          jMax = sNy+OLy          jMax = sNy+OLy
246    
 C--     Update fields according to tendency terms  
         CALL TIMESTEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
247    
248          DO K=2,Nz          K = 1
249  C Density of K-1 level (above W(K)) reference to K level          BOTTOM_LAYER = K .EQ. Nr
250           CALL FIND_RHO(  
251       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',  #ifdef DO_PIPELINED_CORRECTION_STEP
252       O      rhoKm1,  C--     Calculate gradient of surface pressure
253       I      myThid )          CALL CALC_GRAD_ETA_SURF(
254  C Density of K level (below W(K)) reference to K level       I       bi,bj,iMin,iMax,jMin,jMax,
255         O       etaSurfX,etaSurfY,
256         I       myThid)
257    C--     Update fields in top level according to tendency terms
258            CALL CORRECTION_STEP(
259         I       bi,bj,iMin,iMax,jMin,jMax,K,
260         I       etaSurfX,etaSurfY,myTime,myThid)
261            IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K, myThid )
262            IF ( .NOT. BOTTOM_LAYER ) THEN
263    C--      Update fields in layer below according to tendency terms
264             CALL CORRECTION_STEP(
265         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
266         I        etaSurfX,etaSurfY,myTime,myThid)
267             IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
268            ENDIF
269    #endif
270    C--     Density of 1st level (below W(1)) reference to level 1
271    #ifdef  INCLUDE_FIND_RHO_CALL
272            CALL FIND_RHO(
273         I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
274         O     rhoKm1,
275         I     myThid )
276    #endif
277    
278            IF (       (.NOT. BOTTOM_LAYER)
279    #ifdef ALLOW_KPP
280         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
281    #endif
282         &     ) THEN
283    C--      Check static stability with layer below
284    C--      and mix as needed.
285    #ifdef  INCLUDE_FIND_RHO_CALL
286           CALL FIND_RHO(           CALL FIND_RHO(
287       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',       I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
288       O      rhoKp1,       O      rhoKp1,
289       I      myThid )       I      myThid )
290  C--     Calculate static stability and mix where convectively unstable  #endif
291    #ifdef  INCLUDE_CONVECT_CALL
292           CALL CONVECT(           CALL CONVECT(
293       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
294  C Density of K-1 level (above W(K)) reference to K-1 level       I       myTime,myIter,myThid)
295    #endif
296    C--      Recompute density after mixing
297    #ifdef  INCLUDE_FIND_RHO_CALL
298           CALL FIND_RHO(           CALL FIND_RHO(
299       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, 'LINEAR',       I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
300       O      rhoKm1,       O      rhoKm1,
301       I      myThid )       I      myThid )
302  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #endif
303           CALL CALC_PH(          ENDIF
304       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,  C--     Calculate buoyancy
305       U       pH,          CALL CALC_BUOYANCY(
306       I       myThid )       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
307          ENDDO ! K       O      buoyKm1,
308         I      myThid )
309    C--     Integrate hydrostatic balance for phiHyd with BC of
310    C--     phiHyd(z=0)=0
311            CALL CALC_PHI_HYD(
312         I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
313         U      phiHyd,
314         I      myThid )
315    
316  C Density of Nz level (bottom level) reference to Nz level          DO K=2,Nr
317             BOTTOM_LAYER = K .EQ. Nr
318    #ifdef DO_PIPELINED_CORRECTION_STEP
319             IF ( .NOT. BOTTOM_LAYER ) THEN
320    C--       Update fields in layer below according to tendency terms
321              CALL CORRECTION_STEP(
322         I         bi,bj,iMin,iMax,jMin,jMax,K+1,
323         I         etaSurfX,etaSurfY,myTime,myThid)
324              IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
325             ENDIF
326    #endif
327    C--      Density of K level (below W(K)) reference to K level
328    #ifdef  INCLUDE_FIND_RHO_CALL
329           CALL FIND_RHO(           CALL FIND_RHO(
330       I      bi, bj, iMin, iMax, jMin, jMax,  Nz, Nz, 'LINEAR',       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
331       O      rhoKm1,       O      rhoK,
332       I      myThid )       I      myThid )
333  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #endif
334           CALL CALC_PH(           IF (       (.NOT. BOTTOM_LAYER)
335       I       bi,bj,iMin,iMax,jMin,jMax,Nz+1,rhoKm1,  #ifdef ALLOW_KPP
336       U       pH,       &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
337    #endif
338         &      ) THEN
339    C--       Check static stability with layer below and mix as needed.
340    C--       Density of K+1 level (below W(K+1)) reference to K level.
341    #ifdef  INCLUDE_FIND_RHO_CALL
342              CALL FIND_RHO(
343         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
344         O       rhoKp1,
345         I       myThid )
346    #endif
347    #ifdef  INCLUDE_CONVECT_CALL
348              CALL CONVECT(
349         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
350         I        myTime,myIter,myThid)
351    #endif
352    C--       Recompute density after mixing
353    #ifdef  INCLUDE_FIND_RHO_CALL
354              CALL FIND_RHO(
355         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
356         O       rhoK,
357       I       myThid )       I       myThid )
358    #endif
359             ENDIF
360    C--      Calculate buoyancy
361             CALL CALC_BUOYANCY(
362         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
363         O       buoyK,
364         I       myThid )
365    C--      Integrate hydrostatic balance for phiHyd with BC of
366    C--      phiHyd(z=0)=0
367             CALL CALC_PHI_HYD(
368         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
369         U        phiHyd,
370         I        myThid )
371    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
372    #ifdef  INCLUDE_FIND_RHO_CALL
373             CALL FIND_RHO(
374         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
375         O        rhoTmp,
376         I        myThid )
377    #endif
378    #ifdef  INCLUDE_CALC_ISOSLOPES_CALL
379             CALL CALC_ISOSLOPES(
380         I        bi, bj, iMin, iMax, jMin, jMax, K,
381         I        rhoKm1, rhoK, rhotmp,
382         O        K13, K23, K33, KapGM,
383         I        myThid )
384    #endif
385             DO J=jMin,jMax
386              DO I=iMin,iMax
387    #ifdef  INCLUDE_FIND_RHO_CALL
388               rhoKm1 (I,J) = rhoK(I,J)
389    #endif
390               buoyKm1(I,J) = buoyK(I,J)
391              ENDDO
392             ENDDO
393            ENDDO ! K
394    
395    #ifdef ALLOW_KPP
396    C--     Compute KPP mixing coefficients
397            IF (usingKPPmixing) THEN
398             CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
399         I          , myThid)
400             CALL KVMIX(
401         I               bi, bj, myTime, myThid )
402             CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
403         I        , myThid)
404            ENDIF
405    #endif
406    
407            DO K = Nr, 1, -1
408    
         DO K = Nz, 1, -1  
409           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
410           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
411           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer
# Line 193  C--     Integrate hydrostatic balance fo Line 417  C--     Integrate hydrostatic balance fo
417  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
418           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
419       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
420       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
421       I        myThid)       I        myThid)
422    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
423  C--      Calculate accelerations in the momentum equations  C--      Calculate the total vertical diffusivity
424           CALL CALC_MOM_RHS(           CALL CALC_DIFFUSIVITY(
425       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,K,
426       I        xA,yA,uTrans,vTrans,wTrans,maskC,       I        maskC,maskUp,KapGM,K33,
427       I        pH,       O        KappaRT,KappaRS,KappaRU,KappaRV,
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
428       I        myThid)       I        myThid)
429    #endif
430    C--      Calculate accelerations in the momentum equations
431             IF ( momStepping ) THEN
432              CALL CALC_MOM_RHS(
433         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
434         I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
435         I         phiHyd,KappaRU,KappaRV,
436         U         aTerm,xTerm,cTerm,mTerm,pTerm,
437         U         fZon, fMer, fVerU, fVerV,
438         I         myTime, myThid)
439             ENDIF
440  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
441           CALL CALC_GT(           IF ( tempStepping ) THEN
442       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
443       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
444       U        aTerm,xTerm,fZon,fMer,fVerT,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
445       I        myThid)       I         K13,K23,KappaRT,KapGM,
446  Cdbg     CALL CALC_GS(       U         aTerm,xTerm,fZon,fMer,fVerT,
447  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         myTime, myThid)
448  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,           ENDIF
449  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,           IF ( saltStepping ) THEN
450  Cdbg I        myThid)            CALL CALC_GS(
451         I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
452         I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
453         I         K13,K23,KappaRS,KapGM,
454         U         aTerm,xTerm,fZon,fMer,fVerS,
455         I         myTime, myThid)
456             ENDIF
457    C--      Prediction step (step forward all model variables)
458             CALL TIMESTEP(
459         I       bi,bj,iMin,iMax,jMin,jMax,K,
460         I       myThid)
461    C--      Apply open boundary conditions
462             IF (openBoundaries) CALL APPLY_OBCS2( bi, bj, K, myThid )
463    C--      Freeze water
464             IF (allowFreezing)
465         &   CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
466    C--      Diagnose barotropic divergence of predicted fields
467             CALL CALC_DIV_GHAT(
468         I       bi,bj,iMin,iMax,jMin,jMax,K,
469         I       xA,yA,
470         I       myThid)
471    
472    C--      Cumulative diagnostic calculations (ie. time-averaging)
473    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
474             IF (taveFreq.GT.0.) THEN
475              CALL DO_TIME_AVERAGES(
476         I                           myTime, myIter, bi, bj, K, kUp, kDown,
477         I                           K13, K23, rVel, KapGM,
478         I                           myThid )
479             ENDIF
480    #endif
481    
482          ENDDO          ENDDO ! K
483    
484    C--     Implicit diffusion
485            IF (implicitDiffusion) THEN
486             IF (tempStepping) CALL IMPLDIFF(
487         I         bi, bj, iMin, iMax, jMin, jMax,
488         I         deltaTtracer, KappaRT,recip_HFacC,
489         U         gTNm1,
490         I         myThid )
491             IF (saltStepping) CALL IMPLDIFF(
492         I         bi, bj, iMin, iMax, jMin, jMax,
493         I         deltaTtracer, KappaRS,recip_HFacC,
494         U         gSNm1,
495         I         myThid )
496             IF (momStepping) THEN
497              CALL IMPLDIFF(
498         I         bi, bj, iMin, iMax, jMin, jMax,
499         I         deltaTmom, KappaRU,recip_HFacW,
500         U         gUNm1,
501         I         myThid )
502              CALL IMPLDIFF(
503         I         bi, bj, iMin, iMax, jMin, jMax,
504         I         deltaTmom, KappaRV,recip_HFacS,
505         U         gVNm1,
506         I         myThid )
507    #ifdef INCLUDE_CD_CODE
508              CALL IMPLDIFF(
509         I         bi, bj, iMin, iMax, jMin, jMax,
510         I         deltaTmom, KappaRU,recip_HFacW,
511         U         vVelD,
512         I         myThid )
513              CALL IMPLDIFF(
514         I         bi, bj, iMin, iMax, jMin, jMax,
515         I         deltaTmom, KappaRV,recip_HFacS,
516         U         uVelD,
517         I         myThid )
518    #endif
519             ENDIF ! momStepping
520            ENDIF ! implicitDiffusion
521    
522         ENDDO         ENDDO
523        ENDDO        ENDDO
524    
525    C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
526    C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
527    C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
528    C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
529    C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
530    C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
531    C     write(0,*) 'dynamics: rVel(1) ',
532    C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
533    C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
534    C     write(0,*) 'dynamics: rVel(2) ',
535    C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
536    C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
537    cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
538    cblk &                           maxval(K13(1:sNx,1:sNy,:))
539    cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
540    cblk &                           maxval(K23(1:sNx,1:sNy,:))
541    cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
542    cblk &                           maxval(K33(1:sNx,1:sNy,:))
543    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
544    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
545    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
546    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
547    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
548    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
549    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
550    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
551    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
552    C    &                           maxval(phiHyd/(Gravity*Rhonil))
553    C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
554    C    &Nr, 1, myThid )
555    C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
556    C    &Nr, 1, myThid )
557    C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
558    C    &Nr, 1, myThid )
559    C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
560    C    &Nr, 1, myThid )
561    C     CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
562    C    &Nr, 1, myThid )
563    
564    
565        RETURN        RETURN
566        END        END

Legend:
Removed from v.1.3  
changed lines
  Added in v.1.42

  ViewVC Help
Powered by ViewVC 1.1.22