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

Legend:
Removed from v.1.6  
changed lines
  Added in v.1.26

  ViewVC Help
Powered by ViewVC 1.1.22