/[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.2 by cnh, Fri Apr 24 02:05:40 1998 UTC revision 1.46 by adcroft, Mon Aug 30 18:25:33 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"
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 51  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 phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111          _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    #ifdef INCLUDE_CONVECT_CALL
129          _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
130    #endif
131    
132        INTEGER iMin, iMax        INTEGER iMin, iMax
133        INTEGER jMin, jMax        INTEGER jMin, jMax
134        INTEGER bi, bj        INTEGER bi, bj
135        INTEGER i, j        INTEGER i, j
136        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
137          LOGICAL BOTTOM_LAYER
138    
139    C---    The algorithm...
140    C
141    C       "Correction Step"
142    C       =================
143    C       Here we update the horizontal velocities with the surface
144    C       pressure such that the resulting flow is either consistent
145    C       with the free-surface evolution or the rigid-lid:
146    C         U[n] = U* + dt x d/dx P
147    C         V[n] = V* + dt x d/dy P
148    C
149    C       "Calculation of Gs"
150    C       ===================
151    C       This is where all the accelerations and tendencies (ie.
152    C       phiHydysics, parameterizations etc...) are calculated
153    C         rVel = sum_r ( div. u[n] )
154    C         rho = rho ( theta[n], salt[n] )
155    C         b   = b(rho, theta)
156    C         K31 = K31 ( rho )
157    C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
158    C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
159    C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
160    C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
161    C
162    C       "Time-stepping" or "Prediction"
163    C       ================================
164    C       The models variables are stepped forward with the appropriate
165    C       time-stepping scheme (currently we use Adams-Bashforth II)
166    C       - For momentum, the result is always *only* a "prediction"
167    C       in that the flow may be divergent and will be "corrected"
168    C       later with a surface pressure gradient.
169    C       - Normally for tracers the result is the new field at time
170    C       level [n+1} *BUT* in the case of implicit diffusion the result
171    C       is also *only* a prediction.
172    C       - We denote "predictors" with an asterisk (*).
173    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
174    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
175    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
176    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
177    C       With implicit diffusion:
178    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
179    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
180    C         (1 + dt * K * d_zz) theta[n] = theta*
181    C         (1 + dt * K * d_zz) salt[n] = salt*
182    C---
183    
184  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
185  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
# Line 89  C     point numbers. This prevents spuri Line 188  C     point numbers. This prevents spuri
188  C     uninitialised but inert locations.  C     uninitialised but inert locations.
189        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
190         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
191          xA(i,j)      = 0.*1. _d 37          xA(i,j)      = 0. _d 0
192          yA(i,j)      = 0.*1. _d 37          yA(i,j)      = 0. _d 0
193          uTrans(i,j)  = 0.*1. _d 37          uTrans(i,j)  = 0. _d 0
194          vTrans(i,j)  = 0.*1. _d 37          vTrans(i,j)  = 0. _d 0
195          aTerm(i,j)   = 0.*1. _d 37          aTerm(i,j)   = 0. _d 0
196          xTerm(i,j)   = 0.*1. _d 37          xTerm(i,j)   = 0. _d 0
197          cTerm(i,j)   = 0.*1. _d 37          cTerm(i,j)   = 0. _d 0
198          mTerm(i,j)   = 0.*1. _d 37          mTerm(i,j)   = 0. _d 0
199          pTerm(i,j)   = 0.*1. _d 37          pTerm(i,j)   = 0. _d 0
200          fZon(i,j)    = 0.*1. _d 37          fZon(i,j)    = 0. _d 0
201          fMer(i,j)    = 0.*1. _d 37          fMer(i,j)    = 0. _d 0
202          DO K=1,nZ          DO K=1,Nr
203           pH (i,j,k)  = 0.*1. _d 37           phiHyd (i,j,k)  = 0. _d 0
204             K13(i,j,k)  = 0. _d 0
205             K23(i,j,k)  = 0. _d 0
206             K33(i,j,k)  = 0. _d 0
207             KappaRU(i,j,k) = 0. _d 0
208             KappaRV(i,j,k) = 0. _d 0
209          ENDDO          ENDDO
210            rhoKM1 (i,j) = 0. _d 0
211            rhok   (i,j) = 0. _d 0
212            rhoKP1 (i,j) = 0. _d 0
213            rhoTMP (i,j) = 0. _d 0
214            buoyKM1(i,j) = 0. _d 0
215            buoyK  (i,j) = 0. _d 0
216            maskC  (i,j) = 0. _d 0
217         ENDDO         ENDDO
218        ENDDO        ENDDO
219  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  
220    
221        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
222         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
223    
224    C--     Set up work arrays that need valid initial values
225            DO j=1-OLy,sNy+OLy
226             DO i=1-OLx,sNx+OLx
227              rTrans(i,j)   = 0. _d 0
228              rVel  (i,j,1) = 0. _d 0
229              rVel  (i,j,2) = 0. _d 0
230              fVerT (i,j,1) = 0. _d 0
231              fVerT (i,j,2) = 0. _d 0
232              fVerS (i,j,1) = 0. _d 0
233              fVerS (i,j,2) = 0. _d 0
234              fVerU (i,j,1) = 0. _d 0
235              fVerU (i,j,2) = 0. _d 0
236              fVerV (i,j,1) = 0. _d 0
237              fVerV (i,j,2) = 0. _d 0
238              phiHyd(i,j,1) = 0. _d 0
239              K13   (i,j,1) = 0. _d 0
240              K23   (i,j,1) = 0. _d 0
241              K33   (i,j,1) = 0. _d 0
242              KapGM (i,j)   = GMkbackground
243             ENDDO
244            ENDDO
245    
246            DO k=1,Nr
247             DO j=1-OLy,sNy+OLy
248              DO i=1-OLx,sNx+OLx
249    #ifdef INCLUDE_CONVECT_CALL
250               ConvectCount(i,j,k) = 0.
251    #endif
252               KappaRT(i,j,k) = 0. _d 0
253               KappaRS(i,j,k) = 0. _d 0
254              ENDDO
255             ENDDO
256            ENDDO
257    
258          iMin = 1-OLx+1          iMin = 1-OLx+1
259          iMax = sNx+OLx          iMax = sNx+OLx
260          jMin = 1-OLy+1          jMin = 1-OLy+1
261          jMax = sNy+OLy          jMax = sNy+OLy
262    
 C--     Update fields according to tendency terms  
         CALL TIMESTEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
263    
264  C--     Calculate rho with the appropriate equation of state          K = 1
265          CALL FIND_RHO(          BOTTOM_LAYER = K .EQ. Nr
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
   
 C--     Calculate static stability and mix where convectively unstable  
         CALL CONVECT(  
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
266    
267  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #ifdef DO_PIPELINED_CORRECTION_STEP
268          CALL CALC_PH(  C--     Calculate gradient of surface pressure
269            CALL CALC_GRAD_ETA_SURF(
270       I       bi,bj,iMin,iMax,jMin,jMax,       I       bi,bj,iMin,iMax,jMin,jMax,
271       O       pH,       O       etaSurfX,etaSurfY,
272         I       myThid)
273    C--     Update fields in top level according to tendency terms
274            CALL CORRECTION_STEP(
275         I       bi,bj,iMin,iMax,jMin,jMax,K,
276         I       etaSurfX,etaSurfY,myTime,myThid)
277    #ifdef ALLOW_OBCS
278            IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K, myThid )
279    #endif
280            IF ( .NOT. BOTTOM_LAYER ) THEN
281    C--      Update fields in layer below according to tendency terms
282             CALL CORRECTION_STEP(
283         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
284         I        etaSurfX,etaSurfY,myTime,myThid)
285    #ifdef ALLOW_OBCS
286             IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
287    #endif
288            ENDIF
289    #endif
290    C--     Density of 1st level (below W(1)) reference to level 1
291    #ifdef  INCLUDE_FIND_RHO_CALL
292            CALL FIND_RHO(
293         I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
294         O     rhoKm1,
295         I     myThid )
296    #endif
297    
298            IF (       (.NOT. BOTTOM_LAYER)
299    #ifdef ALLOW_KPP
300         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
301    #endif
302         &     ) THEN
303    C--      Check static stability with layer below
304    C--      and mix as needed.
305    #ifdef  INCLUDE_FIND_RHO_CALL
306             CALL FIND_RHO(
307         I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
308         O      rhoKp1,
309         I      myThid )
310    #endif
311    #ifdef  INCLUDE_CONVECT_CALL
312             CALL CONVECT(
313         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
314         U       ConvectCount,
315         I       myTime,myIter,myThid)
316    #endif
317    C--      Implicit Vertical Diffusion for Convection
318             IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
319         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
320         U       ConvectCount, KappaRT, KappaRS,
321         I       myTime,myIter,myThid)
322    C--      Recompute density after mixing
323    #ifdef  INCLUDE_FIND_RHO_CALL
324             CALL FIND_RHO(
325         I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
326         O      rhoKm1,
327         I      myThid )
328    #endif
329            ENDIF
330    C--     Calculate buoyancy
331            CALL CALC_BUOYANCY(
332         I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
333         O      buoyKm1,
334         I      myThid )
335    C--     Integrate hydrostatic balance for phiHyd with BC of
336    C--     phiHyd(z=0)=0
337            CALL CALC_PHI_HYD(
338         I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
339         U      phiHyd,
340         I      myThid )
341    
342            DO K=2,Nr
343             BOTTOM_LAYER = K .EQ. Nr
344    #ifdef DO_PIPELINED_CORRECTION_STEP
345             IF ( .NOT. BOTTOM_LAYER ) THEN
346    C--       Update fields in layer below according to tendency terms
347              CALL CORRECTION_STEP(
348         I         bi,bj,iMin,iMax,jMin,jMax,K+1,
349         I         etaSurfX,etaSurfY,myTime,myThid)
350    #ifdef ALLOW_OBCS
351              IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
352    #endif
353             ENDIF
354    #endif
355    C--      Density of K level (below W(K)) reference to K level
356    #ifdef  INCLUDE_FIND_RHO_CALL
357             CALL FIND_RHO(
358         I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
359         O      rhoK,
360         I      myThid )
361    #endif
362             IF (       (.NOT. BOTTOM_LAYER)
363    #ifdef ALLOW_KPP
364         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
365    #endif
366         &      ) THEN
367    C--       Check static stability with layer below and mix as needed.
368    C--       Density of K+1 level (below W(K+1)) reference to K level.
369    #ifdef  INCLUDE_FIND_RHO_CALL
370              CALL FIND_RHO(
371         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
372         O       rhoKp1,
373       I       myThid )       I       myThid )
374    #endif
375    #ifdef  INCLUDE_CONVECT_CALL
376              CALL CONVECT(
377         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
378         U        ConvectCount,
379         I        myTime,myIter,myThid)
380    #endif
381    C--      Implicit Vertical Diffusion for Convection
382             IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
383         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
384         U       ConvectCount, KappaRT, KappaRS,
385         I       myTime,myIter,myThid)
386    C--       Recompute density after mixing
387    #ifdef  INCLUDE_FIND_RHO_CALL
388              CALL FIND_RHO(
389         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
390         O       rhoK,
391         I       myThid )
392    #endif
393             ENDIF
394    C--      Calculate buoyancy
395             CALL CALC_BUOYANCY(
396         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
397         O       buoyK,
398         I       myThid )
399    C--      Integrate hydrostatic balance for phiHyd with BC of
400    C--      phiHyd(z=0)=0
401             CALL CALC_PHI_HYD(
402         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
403         U        phiHyd,
404         I        myThid )
405    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
406    #ifdef  INCLUDE_FIND_RHO_CALL
407             CALL FIND_RHO(
408         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
409         O        rhoTmp,
410         I        myThid )
411    #endif
412    #ifdef  INCLUDE_CALC_ISOSLOPES_CALL
413             CALL CALC_ISOSLOPES(
414         I        bi, bj, iMin, iMax, jMin, jMax, K,
415         I        rhoKm1, rhoK, rhotmp,
416         O        K13, K23, K33, KapGM,
417         I        myThid )
418    #endif
419             DO J=jMin,jMax
420              DO I=iMin,iMax
421    #ifdef  INCLUDE_FIND_RHO_CALL
422               rhoKm1 (I,J) = rhoK(I,J)
423    #endif
424               buoyKm1(I,J) = buoyK(I,J)
425              ENDDO
426             ENDDO
427            ENDDO ! K
428    
429    #ifdef ALLOW_KPP
430    C--     Compute KPP mixing coefficients
431            IF (usingKPPmixing) THEN
432             CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
433         I          , myThid)
434             CALL KVMIX(
435         I               bi, bj, myTime, myThid )
436             CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
437         I        , myThid)
438            ENDIF
439    #endif
440    
441            DO K = Nr, 1, -1
442    
         DO K = Nz, 1, -1  
443           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
444           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
445           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 158  C--     Integrate hydrostatic balance fo Line 451  C--     Integrate hydrostatic balance fo
451  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
452           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
453       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
454       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
455       I        myThid)       I        myThid)
456    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
457  C--      Calculate accelerations in the momentum equations  C--      Calculate the total vertical diffusivity
458           CALL CALC_MOM_RHS(           CALL CALC_DIFFUSIVITY(
459       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,K,
460       I        xA,yA,uTrans,vTrans,wTrans,maskC,       I        maskC,maskUp,KapGM,K33,
461       I        pH,       O        KappaRT,KappaRS,KappaRU,KappaRV,
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
462       I        myThid)       I        myThid)
463    #endif
464    C--      Calculate accelerations in the momentum equations
465             IF ( momStepping ) THEN
466              CALL CALC_MOM_RHS(
467         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
468         I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
469         I         phiHyd,KappaRU,KappaRV,
470         U         aTerm,xTerm,cTerm,mTerm,pTerm,
471         U         fZon, fMer, fVerU, fVerV,
472         I         myTime, myThid)
473             ENDIF
474  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
475           CALL CALC_GT(           IF ( tempStepping ) THEN
476       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
477       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
478       U        aTerm,xTerm,fZon,fMer,fVerT,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
479       I        myThid)       I         K13,K23,KappaRT,KapGM,
480  Cdbg     CALL CALC_GS(       U         aTerm,xTerm,fZon,fMer,fVerT,
481  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         myTime, myThid)
482  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,           ENDIF
483  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,           IF ( saltStepping ) THEN
484  Cdbg I        myThid)            CALL CALC_GS(
485         I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
486          ENDDO       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
487         I         K13,K23,KappaRS,KapGM,
488         U         aTerm,xTerm,fZon,fMer,fVerS,
489         I         myTime, myThid)
490             ENDIF
491    C--      Prediction step (step forward all model variables)
492             CALL TIMESTEP(
493         I       bi,bj,iMin,iMax,jMin,jMax,K,
494         I       myIter, myThid)
495    #ifdef ALLOW_OBCS
496    C--      Apply open boundary conditions
497             IF (openBoundaries) CALL APPLY_OBCS2( bi, bj, K, myThid )
498    #endif
499    C--      Freeze water
500             IF (allowFreezing)
501         &   CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
502    C--      Diagnose barotropic divergence of predicted fields
503             CALL CALC_DIV_GHAT(
504         I       bi,bj,iMin,iMax,jMin,jMax,K,
505         I       xA,yA,
506         I       myThid)
507    
508    C--      Cumulative diagnostic calculations (ie. time-averaging)
509    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
510             IF (taveFreq.GT.0.) THEN
511              CALL DO_TIME_AVERAGES(
512         I                           myTime, myIter, bi, bj, K, kUp, kDown,
513         I                           K13, K23, rVel, KapGM, ConvectCount,
514         I                           myThid )
515             ENDIF
516    #endif
517    
518    
519            ENDDO ! K
520    
521    C--     Implicit diffusion
522            IF (implicitDiffusion) THEN
523             IF (tempStepping) CALL IMPLDIFF(
524         I         bi, bj, iMin, iMax, jMin, jMax,
525         I         deltaTtracer, KappaRT,recip_HFacC,
526         U         gTNm1,
527         I         myThid )
528             IF (saltStepping) CALL IMPLDIFF(
529         I         bi, bj, iMin, iMax, jMin, jMax,
530         I         deltaTtracer, KappaRS,recip_HFacC,
531         U         gSNm1,
532         I         myThid )
533            ENDIF ! implicitDiffusion
534    C--     Implicit viscosity
535            IF (implicitViscosity) THEN
536             IF (momStepping) THEN
537              CALL IMPLDIFF(
538         I         bi, bj, iMin, iMax, jMin, jMax,
539         I         deltaTmom, KappaRU,recip_HFacW,
540         U         gUNm1,
541         I         myThid )
542              CALL IMPLDIFF(
543         I         bi, bj, iMin, iMax, jMin, jMax,
544         I         deltaTmom, KappaRV,recip_HFacS,
545         U         gVNm1,
546         I         myThid )
547    #ifdef INCLUDE_CD_CODE
548              CALL IMPLDIFF(
549         I         bi, bj, iMin, iMax, jMin, jMax,
550         I         deltaTmom, KappaRU,recip_HFacW,
551         U         vVelD,
552         I         myThid )
553              CALL IMPLDIFF(
554         I         bi, bj, iMin, iMax, jMin, jMax,
555         I         deltaTmom, KappaRV,recip_HFacS,
556         U         uVelD,
557         I         myThid )
558    #endif
559             ENDIF ! momStepping
560            ENDIF ! implicitViscosity
561    
562         ENDDO         ENDDO
563        ENDDO        ENDDO
564    
565    C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
566    C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
567    C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
568    C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
569    C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
570    C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
571    C     write(0,*) 'dynamics: rVel(1) ',
572    C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
573    C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
574    C     write(0,*) 'dynamics: rVel(2) ',
575    C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
576    C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
577    cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
578    cblk &                           maxval(K13(1:sNx,1:sNy,:))
579    cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
580    cblk &                           maxval(K23(1:sNx,1:sNy,:))
581    cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
582    cblk &                           maxval(K33(1:sNx,1:sNy,:))
583    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
584    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
585    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
586    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
587    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
588    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
589    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
590    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
591    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
592    C    &                           maxval(phiHyd/(Gravity*Rhonil))
593    C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
594    C    &Nr, 1, myThid )
595    C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
596    C    &Nr, 1, myThid )
597    C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
598    C    &Nr, 1, myThid )
599    C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
600    C    &Nr, 1, myThid )
601    C     CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
602    C    &Nr, 1, myThid )
603    
604    
605        RETURN        RETURN
606        END        END

Legend:
Removed from v.1.2  
changed lines
  Added in v.1.46

  ViewVC Help
Powered by ViewVC 1.1.22