/[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.7 by cnh, Mon May 25 16:17:36 1998 UTC revision 1.31 by cnh, Sat Aug 22 17:51:08 1998 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 29  C     == Global variables === Line 29  C     == Global variables ===
29  #include "DYNVARS.h"  #include "DYNVARS.h"
30    
31  C     == Routine arguments ==  C     == Routine arguments ==
32    C     myTime - Current time in simulation
33    C     myIter - Current iteration number in simulation
34  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
35        INTEGER myThid        INTEGER myThid
36          _RL myTime
37          INTEGER myIter
38    
39  C     == Local variables  C     == Local variables
40  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
41  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C     uTrans, vTrans, rTrans - Per block temporaries holding flow transport
42  C                              o uTrans: Zonal transport  C     rVel                     o uTrans: Zonal transport
43  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
44  C                              o wTrans: Vertical transport  C                              o rTrans: Vertical transport
45    C                              o rVel:   Vertical velocity at upper and lower
46    C                                        cell faces.
47  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
48  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
49  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in
# Line 53  C                              o fVer: V Line 59  C                              o fVer: V
59  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
60  C                                      so we need an fVer for each  C                                      so we need an fVer for each
61  C                                      variable.  C                                      variable.
62  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     rhoK, rhoKM1   - Density at current level, level above and level below.
63  C     jMin, jMax   are applied.  C     rhoKP1                                                                  
64    C     buoyK, buoyKM1 - Buoyancy at current level and level above.
65    C     phiHyd         - Hydrostatic part of the potential phiHydi.
66    C                      In z coords phiHydiHyd is the hydrostatic pressure anomaly
67    C                      In p coords phiHydiHyd is the geopotential surface height
68    C                      anomaly.
69    C     etaSurfX,      - Holds surface elevation gradient in X and Y.
70    C     etaSurfY
71    C     K13, K23, K33  - Non-zero elements of small-angle approximation
72    C                      diffusion tensor.
73    C     KapGM          - Spatially varying Visbeck et. al mixing coeff.
74    C     KappaRT,       - Total diffusion in vertical for T and S.
75    C     KappaRS          ( background + spatially varying, isopycnal term).
76    C     iMin, iMax     - Ranges and sub-block indices on which calculations
77    C     jMin, jMax       are applied.
78  C     bi, bj  C     bi, bj
79  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
80  C                          are switched with layer to be the appropriate index  C     kDown, kM1       are switched with layer to be the appropriate index
81  C                          into fVerTerm  C                      into fVerTerm
82        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
88        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
98        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
99        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
101        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109          _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110          _RL K13     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111          _RL K23     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
112          _RL K33     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
113          _RL KapGM   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114          _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
115          _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
116    
117        INTEGER iMin, iMax        INTEGER iMin, iMax
118        INTEGER jMin, jMax        INTEGER jMin, jMax
119        INTEGER bi, bj        INTEGER bi, bj
120        INTEGER i, j        INTEGER i, j
121        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
122          LOGICAL BOTTOM_LAYER
123    
124    C---    The algorithm...
125    C
126    C       "Correction Step"
127    C       =================
128    C       Here we update the horizontal velocities with the surface
129    C       pressure such that the resulting flow is either consistent
130    C       with the free-surface evolution or the rigid-lid:
131    C         U[n] = U* + dt x d/dx P
132    C         V[n] = V* + dt x d/dy P
133    C
134    C       "Calculation of Gs"
135    C       ===================
136    C       This is where all the accelerations and tendencies (ie.
137    C       phiHydysics, parameterizations etc...) are calculated
138    C         rVel = sum_r ( div. u[n] )
139    C         rho = rho ( theta[n], salt[n] )
140    C         b   = b(rho, theta)
141    C         K31 = K31 ( rho )
142    C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
143    C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
144    C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
145    C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
146    C
147    C       "Time-stepping" or "Prediction"
148    C       ================================
149    C       The models variables are stepped forward with the appropriate
150    C       time-stepping scheme (currently we use Adams-Bashforth II)
151    C       - For momentum, the result is always *only* a "prediction"
152    C       in that the flow may be divergent and will be "corrected"
153    C       later with a surface pressure gradient.
154    C       - Normally for tracers the result is the new field at time
155    C       level [n+1} *BUT* in the case of implicit diffusion the result
156    C       is also *only* a prediction.
157    C       - We denote "predictors" with an asterisk (*).
158    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
159    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
160    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
161    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
162    C       With implicit diffusion:
163    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
164    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
165    C         (1 + dt * K * d_zz) theta[n] = theta*
166    C         (1 + dt * K * d_zz) salt[n] = salt*
167    C---
168    
169  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
170  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
# Line 110  C     uninitialised but inert locations. Line 184  C     uninitialised but inert locations.
184          pTerm(i,j)   = 0. _d 0          pTerm(i,j)   = 0. _d 0
185          fZon(i,j)    = 0. _d 0          fZon(i,j)    = 0. _d 0
186          fMer(i,j)    = 0. _d 0          fMer(i,j)    = 0. _d 0
187          DO K=1,nZ          DO K=1,Nr
188           pH (i,j,k)  = 0. _d 0           phiHyd (i,j,k)  = 0. _d 0
189           K13(i,j,k) = 0. _d 0           K13(i,j,k)  = 0. _d 0
190           K23(i,j,k) = 0. _d 0           K23(i,j,k)  = 0. _d 0
191           K33(i,j,k) = 0. _d 0           K33(i,j,k)  = 0. _d 0
192             KappaRT(i,j,k) = 0. _d 0
193             KappaRS(i,j,k) = 0. _d 0
194          ENDDO          ENDDO
195          rhokm1(i,j)  = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
196          rhokp1(i,j)  = 0. _d 0          rhok   (i,j) = 0. _d 0
197            rhoKP1 (i,j) = 0. _d 0
198            rhoTMP (i,j) = 0. _d 0
199            buoyKM1(i,j) = 0. _d 0
200            buoyK  (i,j) = 0. _d 0
201            maskC  (i,j) = 0. _d 0
202         ENDDO         ENDDO
203        ENDDO        ENDDO
204    
205        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
206         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
207    
 C--     Boundary condition on hydrostatic pressure is pH(z=0)=0  
         DO j=1-OLy,sNy+OLy  
          DO i=1-OLx,sNx+OLx  
           pH(i,j,1) = 0. _d 0  
           K13(i,j,1) = 0. _d 0  
           K23(i,j,1) = 0. _d 0  
           K33(i,j,1) = 0. _d 0  
           KapGM(i,j) = 0. _d 0  
          ENDDO  
         ENDDO  
   
208  C--     Set up work arrays that need valid initial values  C--     Set up work arrays that need valid initial values
209          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
210           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
211            wTrans(i,j)  = 0. _d 0            rTrans(i,j)   = 0. _d 0
212            fVerT(i,j,1) = 0. _d 0            rVel  (i,j,1) = 0. _d 0
213            fVerT(i,j,2) = 0. _d 0            rVel  (i,j,2) = 0. _d 0
214            fVerS(i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
215            fVerS(i,j,2) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
216            fVerU(i,j,1) = 0. _d 0            fVerS (i,j,1) = 0. _d 0
217            fVerU(i,j,2) = 0. _d 0            fVerS (i,j,2) = 0. _d 0
218            fVerV(i,j,1) = 0. _d 0            fVerU (i,j,1) = 0. _d 0
219            fVerV(i,j,2) = 0. _d 0            fVerU (i,j,2) = 0. _d 0
220              fVerV (i,j,1) = 0. _d 0
221              fVerV (i,j,2) = 0. _d 0
222              phiHyd(i,j,1) = 0. _d 0
223              K13   (i,j,1) = 0. _d 0
224              K23   (i,j,1) = 0. _d 0
225              K33   (i,j,1) = 0. _d 0
226              KapGM (i,j)   = GMkbackground
227           ENDDO           ENDDO
228          ENDDO          ENDDO
229    
# Line 155  C--     Set up work arrays that need val Line 232  C--     Set up work arrays that need val
232          jMin = 1-OLy+1          jMin = 1-OLy+1
233          jMax = sNy+OLy          jMax = sNy+OLy
234    
235            K = 1
236            BOTTOM_LAYER = K .EQ. Nr
237    
238  C--     Calculate gradient of surface pressure  C--     Calculate gradient of surface pressure
239          CALL GRAD_PSURF(          CALL CALC_GRAD_ETA_SURF(
240       I       bi,bj,iMin,iMax,jMin,jMax,       I       bi,bj,iMin,iMax,jMin,jMax,
241       O       pSurfX,pSurfY,       O       etaSurfX,etaSurfY,
242       I       myThid)       I       myThid)
   
243  C--     Update fields in top level according to tendency terms  C--     Update fields in top level according to tendency terms
244          CALL TIMESTEP(          CALL CORRECTION_STEP(
245       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,K,
246         I       etaSurfX,etaSurfY,myTime,myThid)
247            IF ( .NOT. BOTTOM_LAYER ) THEN
248    C--      Update fields in layer below according to tendency terms
249             CALL CORRECTION_STEP(
250         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
251         I        etaSurfX,etaSurfY,myTime,myThid)
252            ENDIF
253  C--     Density of 1st level (below W(1)) reference to level 1  C--     Density of 1st level (below W(1)) reference to level 1
254          CALL FIND_RHO(          CALL FIND_RHO(
255       I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, 'LINEAR',       I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
256       O     rhoKm1,       O     rhoKm1,
257       I     myThid )       I     myThid )
258  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
259          CALL CALC_PH(          IF ( .NOT. BOTTOM_LAYER ) THEN
260       I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,  C--      Check static stability with layer below
261       U      pH,  C--      and mix as needed.
262             CALL FIND_RHO(
263         I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
264         O      rhoKp1,
265         I      myThid )
266             CALL CONVECT(
267         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
268         I       myTime,myIter,myThid)
269    C--      Recompute density after mixing
270             CALL FIND_RHO(
271         I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
272         O      rhoKm1,
273         I      myThid )
274            ENDIF
275    C--     Calculate buoyancy
276            CALL CALC_BUOY(
277         I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
278         O      buoyKm1,
279         I      myThid )
280    C--     Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
281            CALL CALC_PHI_HYD(
282         I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
283         U      phiHyd,
284       I      myThid )       I      myThid )
285    
286          DO K=2,Nz          DO K=2,Nr
287  C--     Update fields in Kth level according to tendency terms           BOTTOM_LAYER = K .EQ. Nr
288          CALL TIMESTEP(           IF ( .NOT. BOTTOM_LAYER ) THEN
289       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)  C--       Update fields in layer below according to tendency terms
290  C--     Density of K-1 level (above W(K)) reference to K level            CALL CORRECTION_STEP(
291          CALL FIND_RHO(       I         bi,bj,iMin,iMax,jMin,jMax,K+1,
292       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',       I         etaSurfX,etaSurfY,myTime,myThid)
293       O     rhoKm1,           ENDIF
294       I     myThid )  C--      Density of K level (below W(K)) reference to K level
295  C--     Density of K level (below W(K)) reference to K level           CALL FIND_RHO(
296          CALL FIND_RHO(       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
297       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',       O      rhoK,
      O     rhoKp1,  
      I     myThid )  
 C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
         CALL CALC_ISOSLOPES(  
      I            bi, bj, iMin, iMax, jMin, jMax, K,  
      I            rhoKm1, rhoKp1,  
      O            K13, K23, K33, KapGM,  
      I            myThid )  
 C--     Calculate static stability and mix where convectively unstable  
         CALL CONVECT(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,myThid)  
 C--     Density of K-1 level (above W(K)) reference to K-1 level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, 'LINEAR',  
      O     rhoKm1,  
      I     myThid )  
 C--     Density of K level (below W(K)) referenced to K level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',  
      O     rhoKp1,  
      I     myThid )  
 C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
         CALL CALC_PH(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,  
      U      pH,  
298       I      myThid )       I      myThid )
299             IF ( .NOT. BOTTOM_LAYER ) THEN
300    C--       Check static stability with layer below and mix as needed.
301    C--       Density of K+1 level (below W(K+1)) reference to K level.
302              CALL FIND_RHO(
303         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
304         O       rhoKp1,
305         I       myThid )
306              CALL CONVECT(
307         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
308         I        myTime,myIter,myThid)
309    C--       Recompute density after mixing
310              CALL FIND_RHO(
311         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
312         O       rhoK,
313         I       myThid )
314             ENDIF
315    C--      Calculate buoyancy
316             CALL CALC_BUOY(
317         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
318         O       buoyK,
319         I       myThid )
320    C--      Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
321             CALL CALC_PHI_HYD(
322         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
323         U        phiHyd,
324         I        myThid )
325    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
326             CALL FIND_RHO(
327         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
328         O        rhoTmp,
329         I        myThid )
330             CALL CALC_ISOSLOPES(
331         I        bi, bj, iMin, iMax, jMin, jMax, K,
332         I        rhoKm1, rhoK, rhotmp,
333         O        K13, K23, K33, KapGM,
334         I        myThid )
335             DO J=jMin,jMax
336              DO I=iMin,iMax
337               rhoKm1 (I,J) = rhoK(I,J)
338               buoyKm1(I,J) = buoyK(I,J)
339              ENDDO
340             ENDDO
341            ENDDO ! K
342    
343          ENDDO          DO K = Nr, 1, -1
344    
         DO K = Nz, 1, -1  
345           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
346           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
347           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 229  C--     Integrate hydrostatic balance fo Line 353  C--     Integrate hydrostatic balance fo
353  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
354           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
355       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
356       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
357       I        myThid)       I        myThid)
358    C--      Calculate the total vertical diffusivity
359  C--      Calculate accelerations in the momentum equations           CALL CALC_DIFFUSIVITY(
360           CALL CALC_MOM_RHS(       I        bi,bj,iMin,iMax,jMin,jMax,K,
361       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        maskC,maskUp,KapGM,K33,
362       I        xA,yA,uTrans,vTrans,wTrans,maskC,       O        KappaRT,KappaRS,
      I        pH,  
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
363       I        myThid)       I        myThid)
364    C--      Calculate accelerations in the momentum equations
365             IF ( momStepping ) THEN
366              CALL CALC_MOM_RHS(
367         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
368         I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
369         I         phiHyd,
370         U         aTerm,xTerm,cTerm,mTerm,pTerm,
371         U         fZon, fMer, fVerU, fVerV,
372         I         myThid)
373             ENDIF
374  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
375           CALL CALC_GT(           IF ( tempStepping ) THEN
376       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
377       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
378       I        K13,K23,K33,KapGM,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
379       U        aTerm,xTerm,fZon,fMer,fVerT,       I         K13,K23,KappaRT,KapGM,
380       I        myThid)       U         aTerm,xTerm,fZon,fMer,fVerT,
381  Cdbg     CALL CALC_GS(       I         myThid)
382  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,           ENDIF
383  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,           IF ( saltStepping ) THEN
384  Cdbg I        K13,K23,K33,KapGM,            CALL CALC_GS(
385  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
386  Cdbg I        myThid)       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
387         I         K13,K23,KappaRS,KapGM,
388         U         aTerm,xTerm,fZon,fMer,fVerS,
389         I         myThid)
390             ENDIF
391    C--      Prediction step (step forward all model variables)
392             CALL TIMESTEP(
393         I       bi,bj,iMin,iMax,jMin,jMax,K,
394         I       myThid)
395    C--      Diagnose barotropic divergence of predicted fields
396             CALL CALC_DIV_GHAT(
397         I       bi,bj,iMin,iMax,jMin,jMax,K,
398         I       xA,yA,
399         I       myThid)
400    
401          ENDDO  C--      Cumulative diagnostic calculations (ie. time-averaging)
402    #ifdef ALLOW_DIAGNOSTICS
403             IF (taveFreq.GT.0.) THEN
404              CALL DO_TIME_AVERAGES(
405         I                           myTime, myIter, bi, bj, K, kUp, kDown,
406         I                           K13, K23, rVel, KapGM,
407         I                           myThid )
408             ENDIF
409    #endif
410    
411            ENDDO ! K
412    
413    C--     Implicit diffusion
414            IF (implicitDiffusion) THEN
415             CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
416         I                  KappaRT,KappaRS,
417         I                  myThid )
418            ENDIF
419    
420         ENDDO         ENDDO
421        ENDDO        ENDDO
422    
423  !dbg  write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
424  !dbg  write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
425  !dbg &                         maxval(uVel(1:sNx,1:sNy,:,:,:))  C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
426  !dbg  write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
427  !dbg &                         maxval(vVel(1:sNx,1:sNy,:,:,:))  C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
428  !dbg  write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
429  !dbg &                         maxval(gT(1:sNx,1:sNy,:,:,:))  C     write(0,*) 'dynamics: rVel(1) ',
430  !dbg  write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),  C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
431  !dbg &                         maxval(Theta(1:sNx,1:sNy,:,:,:))  C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
432  !dbg  write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),  C     write(0,*) 'dynamics: rVel(2) ',
433  !dbg &                          maxval(pH/(Gravity*Rhonil))  C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
434    C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
435    cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
436    cblk &                           maxval(K13(1:sNx,1:sNy,:))
437    cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
438    cblk &                           maxval(K23(1:sNx,1:sNy,:))
439    cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
440    cblk &                           maxval(K33(1:sNx,1:sNy,:))
441    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
442    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
443    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
444    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
445    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
446    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
447    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
448    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
449    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
450    C    &                           maxval(phiHyd/(Gravity*Rhonil))
451    
452        RETURN        RETURN
453        END        END

Legend:
Removed from v.1.7  
changed lines
  Added in v.1.31

  ViewVC Help
Powered by ViewVC 1.1.22