/[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.1 by cnh, Wed Apr 22 19:15:30 1998 UTC revision 1.44 by adcroft, Wed Jul 28 16:32:12 1999 UTC
# Line 1  Line 1 
1  C $Id$  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        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 89  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
207            rhok   (i,j) = 0. _d 0
208            rhoKP1 (i,j) = 0. _d 0
209            rhoTMP (i,j) = 0. _d 0
210            buoyKM1(i,j) = 0. _d 0
211            buoyK  (i,j) = 0. _d 0
212            maskC  (i,j) = 0. _d 0
213         ENDDO         ENDDO
214        ENDDO        ENDDO
215  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  
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--     Set up work arrays that need valid initial values
221            DO j=1-OLy,sNy+OLy
222             DO i=1-OLx,sNx+OLx
223              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
240            ENDDO
241    
242          iMin = 1-OLx+1          iMin = 1-OLx+1
243          iMax = sNx+OLx          iMax = sNx+OLx
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)  
   
 C--     Calculate rho with the appropriate equation of state  
         CALL FIND_RHO(  
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
247    
248  C--     Calculate static stability and mix where convectively unstable          K = 1
249          CALL CONVECT(          BOTTOM_LAYER = K .EQ. Nr
      I       bi,bj,iMin,iMax,jMin,jMax,myThid)  
250    
251  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #ifdef DO_PIPELINED_CORRECTION_STEP
252          CALL CALC_PH(  C--     Calculate gradient of surface pressure
253            CALL CALC_GRAD_ETA_SURF(
254       I       bi,bj,iMin,iMax,jMin,jMax,       I       bi,bj,iMin,iMax,jMin,jMax,
255       O       pH,       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    #ifdef ALLOW_OBCS
262            IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K, myThid )
263    #endif
264            IF ( .NOT. BOTTOM_LAYER ) THEN
265    C--      Update fields in layer below according to tendency terms
266             CALL CORRECTION_STEP(
267         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
268         I        etaSurfX,etaSurfY,myTime,myThid)
269    #ifdef ALLOW_OBCS
270             IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
271    #endif
272            ENDIF
273    #endif
274    C--     Density of 1st level (below W(1)) reference to level 1
275    #ifdef  INCLUDE_FIND_RHO_CALL
276            CALL FIND_RHO(
277         I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
278         O     rhoKm1,
279         I     myThid )
280    #endif
281    
282            IF (       (.NOT. BOTTOM_LAYER)
283    #ifdef ALLOW_KPP
284         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
285    #endif
286         &     ) THEN
287    C--      Check static stability with layer below
288    C--      and mix as needed.
289    #ifdef  INCLUDE_FIND_RHO_CALL
290             CALL FIND_RHO(
291         I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
292         O      rhoKp1,
293         I      myThid )
294    #endif
295    #ifdef  INCLUDE_CONVECT_CALL
296             CALL CONVECT(
297         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
298         I       myTime,myIter,myThid)
299    #endif
300    C--      Recompute density after mixing
301    #ifdef  INCLUDE_FIND_RHO_CALL
302             CALL FIND_RHO(
303         I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
304         O      rhoKm1,
305         I      myThid )
306    #endif
307            ENDIF
308    C--     Calculate buoyancy
309            CALL CALC_BUOYANCY(
310         I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
311         O      buoyKm1,
312         I      myThid )
313    C--     Integrate hydrostatic balance for phiHyd with BC of
314    C--     phiHyd(z=0)=0
315            CALL CALC_PHI_HYD(
316         I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
317         U      phiHyd,
318         I      myThid )
319    
320            DO K=2,Nr
321             BOTTOM_LAYER = K .EQ. Nr
322    #ifdef DO_PIPELINED_CORRECTION_STEP
323             IF ( .NOT. BOTTOM_LAYER ) THEN
324    C--       Update fields in layer below according to tendency terms
325              CALL CORRECTION_STEP(
326         I         bi,bj,iMin,iMax,jMin,jMax,K+1,
327         I         etaSurfX,etaSurfY,myTime,myThid)
328    #ifdef ALLOW_OBCS
329              IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
330    #endif
331             ENDIF
332    #endif
333    C--      Density of K level (below W(K)) reference to K level
334    #ifdef  INCLUDE_FIND_RHO_CALL
335             CALL FIND_RHO(
336         I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
337         O      rhoK,
338         I      myThid )
339    #endif
340             IF (       (.NOT. BOTTOM_LAYER)
341    #ifdef ALLOW_KPP
342         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
343    #endif
344         &      ) THEN
345    C--       Check static stability with layer below and mix as needed.
346    C--       Density of K+1 level (below W(K+1)) reference to K level.
347    #ifdef  INCLUDE_FIND_RHO_CALL
348              CALL FIND_RHO(
349         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
350         O       rhoKp1,
351         I       myThid )
352    #endif
353    #ifdef  INCLUDE_CONVECT_CALL
354              CALL CONVECT(
355         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
356         I        myTime,myIter,myThid)
357    #endif
358    C--       Recompute density after mixing
359    #ifdef  INCLUDE_FIND_RHO_CALL
360              CALL FIND_RHO(
361         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
362         O       rhoK,
363       I       myThid )       I       myThid )
364    #endif
365             ENDIF
366    C--      Calculate buoyancy
367             CALL CALC_BUOYANCY(
368         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
369         O       buoyK,
370         I       myThid )
371    C--      Integrate hydrostatic balance for phiHyd with BC of
372    C--      phiHyd(z=0)=0
373             CALL CALC_PHI_HYD(
374         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
375         U        phiHyd,
376         I        myThid )
377    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
378    #ifdef  INCLUDE_FIND_RHO_CALL
379             CALL FIND_RHO(
380         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
381         O        rhoTmp,
382         I        myThid )
383    #endif
384    #ifdef  INCLUDE_CALC_ISOSLOPES_CALL
385             CALL CALC_ISOSLOPES(
386         I        bi, bj, iMin, iMax, jMin, jMax, K,
387         I        rhoKm1, rhoK, rhotmp,
388         O        K13, K23, K33, KapGM,
389         I        myThid )
390    #endif
391             DO J=jMin,jMax
392              DO I=iMin,iMax
393    #ifdef  INCLUDE_FIND_RHO_CALL
394               rhoKm1 (I,J) = rhoK(I,J)
395    #endif
396               buoyKm1(I,J) = buoyK(I,J)
397              ENDDO
398             ENDDO
399            ENDDO ! K
400    
401    #ifdef ALLOW_KPP
402    C--     Compute KPP mixing coefficients
403            IF (usingKPPmixing) THEN
404             CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
405         I          , myThid)
406             CALL KVMIX(
407         I               bi, bj, myTime, myThid )
408             CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
409         I        , myThid)
410            ENDIF
411    #endif
412    
413            DO K = Nr, 1, -1
414    
         DO K = Nz, 1, -1  
415           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
416           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
417           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 423  C--     Integrate hydrostatic balance fo
423  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
424           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
425       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
426       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
427       I        myThid)       I        myThid)
428    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
429  C--      Calculate accelerations in the momentum equations  C--      Calculate the total vertical diffusivity
430           CALL CALC_MOM_RHS(           CALL CALC_DIFFUSIVITY(
431       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,K,
432       I        xA,yA,uTrans,vTrans,wTrans,maskC,       I        maskC,maskUp,KapGM,K33,
433       I        pH,       O        KappaRT,KappaRS,KappaRU,KappaRV,
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
434       I        myThid)       I        myThid)
435    #endif
436    C--      Calculate accelerations in the momentum equations
437             IF ( momStepping ) THEN
438              CALL CALC_MOM_RHS(
439         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
440         I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
441         I         phiHyd,KappaRU,KappaRV,
442         U         aTerm,xTerm,cTerm,mTerm,pTerm,
443         U         fZon, fMer, fVerU, fVerV,
444         I         myTime, myThid)
445             ENDIF
446  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
447           CALL CALC_GT(           IF ( tempStepping ) THEN
448       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
449       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
450       U        aTerm,xTerm,fZon,fMer,fVerT,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
451       I        myThid)       I         K13,K23,KappaRT,KapGM,
452  Cdbg     CALL CALC_GS(       U         aTerm,xTerm,fZon,fMer,fVerT,
453  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         myTime, myThid)
454  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,           ENDIF
455  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,           IF ( saltStepping ) THEN
456  Cdbg I        myThid)            CALL CALC_GS(
457         I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
458          ENDDO       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
459         I         K13,K23,KappaRS,KapGM,
460         U         aTerm,xTerm,fZon,fMer,fVerS,
461         I         myTime, myThid)
462             ENDIF
463    C--      Prediction step (step forward all model variables)
464             CALL TIMESTEP(
465         I       bi,bj,iMin,iMax,jMin,jMax,K,
466         I       myThid)
467    #ifdef ALLOW_OBCS
468    C--      Apply open boundary conditions
469             IF (openBoundaries) CALL APPLY_OBCS2( bi, bj, K, myThid )
470    #endif
471    C--      Freeze water
472             IF (allowFreezing)
473         &   CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
474    C--      Diagnose barotropic divergence of predicted fields
475             CALL CALC_DIV_GHAT(
476         I       bi,bj,iMin,iMax,jMin,jMax,K,
477         I       xA,yA,
478         I       myThid)
479    
480    C--      Cumulative diagnostic calculations (ie. time-averaging)
481    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
482             IF (taveFreq.GT.0.) THEN
483              CALL DO_TIME_AVERAGES(
484         I                           myTime, myIter, bi, bj, K, kUp, kDown,
485         I                           K13, K23, rVel, KapGM,
486         I                           myThid )
487             ENDIF
488    #endif
489    
490            ENDDO ! K
491    
492    C--     Implicit diffusion
493            IF (implicitDiffusion) THEN
494             IF (tempStepping) CALL IMPLDIFF(
495         I         bi, bj, iMin, iMax, jMin, jMax,
496         I         deltaTtracer, KappaRT,recip_HFacC,
497         U         gTNm1,
498         I         myThid )
499             IF (saltStepping) CALL IMPLDIFF(
500         I         bi, bj, iMin, iMax, jMin, jMax,
501         I         deltaTtracer, KappaRS,recip_HFacC,
502         U         gSNm1,
503         I         myThid )
504            ENDIF ! implicitDiffusion
505    C--     Implicit viscosity
506            IF (implicitViscosity) THEN
507             IF (momStepping) THEN
508              CALL IMPLDIFF(
509         I         bi, bj, iMin, iMax, jMin, jMax,
510         I         deltaTmom, KappaRU,recip_HFacW,
511         U         gUNm1,
512         I         myThid )
513              CALL IMPLDIFF(
514         I         bi, bj, iMin, iMax, jMin, jMax,
515         I         deltaTmom, KappaRV,recip_HFacS,
516         U         gVNm1,
517         I         myThid )
518    #ifdef INCLUDE_CD_CODE
519              CALL IMPLDIFF(
520         I         bi, bj, iMin, iMax, jMin, jMax,
521         I         deltaTmom, KappaRU,recip_HFacW,
522         U         vVelD,
523         I         myThid )
524              CALL IMPLDIFF(
525         I         bi, bj, iMin, iMax, jMin, jMax,
526         I         deltaTmom, KappaRV,recip_HFacS,
527         U         uVelD,
528         I         myThid )
529    #endif
530             ENDIF ! momStepping
531            ENDIF ! implicitViscosity
532    
533         ENDDO         ENDDO
534        ENDDO        ENDDO
535    
536    C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
537    C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
538    C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
539    C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
540    C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
541    C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
542    C     write(0,*) 'dynamics: rVel(1) ',
543    C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
544    C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
545    C     write(0,*) 'dynamics: rVel(2) ',
546    C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
547    C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
548    cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
549    cblk &                           maxval(K13(1:sNx,1:sNy,:))
550    cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
551    cblk &                           maxval(K23(1:sNx,1:sNy,:))
552    cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
553    cblk &                           maxval(K33(1:sNx,1:sNy,:))
554    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
555    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
556    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
557    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
558    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
559    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
560    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
561    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
562    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
563    C    &                           maxval(phiHyd/(Gravity*Rhonil))
564    C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
565    C    &Nr, 1, myThid )
566    C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
567    C    &Nr, 1, myThid )
568    C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
569    C    &Nr, 1, myThid )
570    C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
571    C    &Nr, 1, myThid )
572    C     CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
573    C    &Nr, 1, myThid )
574    
575    
576        RETURN        RETURN
577        END        END

Legend:
Removed from v.1.1  
changed lines
  Added in v.1.44

  ViewVC Help
Powered by ViewVC 1.1.22