/[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.14 by cnh, Mon Jun 8 21:43:01 1998 UTC revision 1.37 by cnh, Tue Nov 3 15:28:04 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(myTime, myIter, myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6  C     /==========================================================\  C     /==========================================================\
# Line 38  C     myThid - Thread number for this in Line 38  C     myThid - Thread number for this in
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     wVel                     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 wVel:   Vertical velocity at upper and lower  C                              o rVel:   Vertical velocity at upper and lower
46  C                                        cell faces.  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
# Line 59  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        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
88        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
98        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
99        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
101        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL K13     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111        _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _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...  C---    The algorithm...
125  C  C
# Line 115  C Line 134  C
134  C       "Calculation of Gs"  C       "Calculation of Gs"
135  C       ===================  C       ===================
136  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
137  C       physics, parameterizations etc...) are calculated  C       phiHydysics, parameterizations etc...) are calculated
138  C         w = sum_z ( div. u[n] )  C         rVel = sum_r ( div. u[n] )
139  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
140    C         b   = b(rho, theta)
141  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
142  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
143  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
144  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
145  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
146  C  C
147  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
148  C       ================================  C       ================================
# Line 164  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           KappaZT(i,j,k) = 0. _d 0           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          rhotmp(i,j)  = 0. _d 0          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    
206        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
207         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
208    
209  C--     Set up work arrays that need valid initial values  C--     Set up work arrays that need valid initial values
210          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
211           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
212            wTrans(i,j)  = 0. _d 0            rTrans(i,j)   = 0. _d 0
213            wVel  (i,j,1) = 0. _d 0            rVel  (i,j,1) = 0. _d 0
214            wVel  (i,j,2) = 0. _d 0            rVel  (i,j,2) = 0. _d 0
215            fVerT(i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
216            fVerT(i,j,2) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
217            fVerS(i,j,1) = 0. _d 0            fVerS (i,j,1) = 0. _d 0
218            fVerS(i,j,2) = 0. _d 0            fVerS (i,j,2) = 0. _d 0
219            fVerU(i,j,1) = 0. _d 0            fVerU (i,j,1) = 0. _d 0
220            fVerU(i,j,2) = 0. _d 0            fVerU (i,j,2) = 0. _d 0
221            fVerV(i,j,1) = 0. _d 0            fVerV (i,j,1) = 0. _d 0
222            fVerV(i,j,2) = 0. _d 0            fVerV (i,j,2) = 0. _d 0
223            pH(i,j,1) = 0. _d 0            phiHyd(i,j,1) = 0. _d 0
224            K13(i,j,1) = 0. _d 0            K13   (i,j,1) = 0. _d 0
225            K23(i,j,1) = 0. _d 0            K23   (i,j,1) = 0. _d 0
226            K33(i,j,1) = 0. _d 0            K33   (i,j,1) = 0. _d 0
227            KapGM(i,j) = 0. _d 0            KapGM (i,j)   = GMkbackground
228           ENDDO           ENDDO
229          ENDDO          ENDDO
230    
# Line 207  C--     Set up work arrays that need val Line 233  C--     Set up work arrays that need val
233          jMin = 1-OLy+1          jMin = 1-OLy+1
234          jMax = sNy+OLy          jMax = sNy+OLy
235    
236    
237            K = 1
238            BOTTOM_LAYER = K .EQ. Nr
239    
240  C--     Calculate gradient of surface pressure  C--     Calculate gradient of surface pressure
241          CALL GRAD_PSURF(          CALL CALC_GRAD_ETA_SURF(
242       I       bi,bj,iMin,iMax,jMin,jMax,       I       bi,bj,iMin,iMax,jMin,jMax,
243       O       pSurfX,pSurfY,       O       etaSurfX,etaSurfY,
244       I       myThid)       I       myThid)
   
245  C--     Update fields in top level according to tendency terms  C--     Update fields in top level according to tendency terms
246          CALL CORRECTION_STEP(          CALL CORRECTION_STEP(
247       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,K,
248         I       etaSurfX,etaSurfY,myTime,myThid)
249            IF ( .NOT. BOTTOM_LAYER ) THEN
250    C--      Update fields in layer below according to tendency terms
251             CALL CORRECTION_STEP(
252         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
253         I        etaSurfX,etaSurfY,myTime,myThid)
254            ENDIF
255  C--     Density of 1st level (below W(1)) reference to level 1  C--     Density of 1st level (below W(1)) reference to level 1
256          CALL FIND_RHO(          CALL FIND_RHO(
257       I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,       I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
258       O     rhoKm1,       O     rhoKm1,
259       I     myThid )       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,1,rhoKm1,rhoKm1,  
      U      pH,  
      I      myThid )  
         DO J=1-Oly,sNy+Oly  
          DO I=1-Olx,sNx+Olx  
           rhoKp1(I,J)=rhoKm1(I,J)  
          ENDDO  
         ENDDO  
260    
261          DO K=2,Nz          IF ( .NOT. BOTTOM_LAYER ) THEN
262  C--     Update fields in Kth level according to tendency terms  C--      Check static stability with layer below
263          CALL CORRECTION_STEP(  C--      and mix as needed.
264       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)           CALL FIND_RHO(
265  C--     Density of K-1 level (above W(K)) reference to K-1 level       I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
266  copt    CALL FIND_RHO(       O      rhoKp1,
267  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,       I      myThid )
268  copt O     rhoKm1,           CALL CONVECT(
269  copt I     myThid )       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
270  C       rhoKm1=rhoKp1       I       myTime,myIter,myThid)
271          DO J=1-Oly,sNy+Oly  C--      Recompute density after mixing
272           DO I=1-Olx,sNx+Olx           CALL FIND_RHO(
273            rhoKm1(I,J)=rhoKp1(I,J)       I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
274           ENDDO       O      rhoKm1,
275          ENDDO       I      myThid )
276  C--     Density of K level (below W(K)) reference to K level          ENDIF
277          CALL FIND_RHO(  C--     Calculate buoyancy
278       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,          CALL CALC_BUOYANCY(
279       O     rhoKp1,       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
280       I     myThid )       O      buoyKm1,
281  C--     Density of K-1 level (above W(K)) reference to K level       I      myThid )
282          CALL FIND_RHO(  C--     Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
283       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,          CALL CALC_PHI_HYD(
284       O     rhotmp,       I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
285       I     myThid )       U      phiHyd,
 C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
         CALL CALC_ISOSLOPES(  
      I            bi, bj, iMin, iMax, jMin, jMax, K,  
      I            rhoKm1, rhoKp1, rhotmp,  
      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,rhotmp,rhoKp1,  
      I      myTime,myIter,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, eosType,  
      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, eosType,  
      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,  
286       I      myThid )       I      myThid )
287    
288            DO K=2,Nr
289             BOTTOM_LAYER = K .EQ. Nr
290             IF ( .NOT. BOTTOM_LAYER ) THEN
291    C--       Update fields in layer below according to tendency terms
292              CALL CORRECTION_STEP(
293         I         bi,bj,iMin,iMax,jMin,jMax,K+1,
294         I         etaSurfX,etaSurfY,myTime,myThid)
295             ENDIF
296    C--      Density of K level (below W(K)) reference to K level
297             CALL FIND_RHO(
298         I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
299         O      rhoK,
300         I      myThid )
301             IF ( .NOT. BOTTOM_LAYER ) THEN
302    C--       Check static stability with layer below and mix as needed.
303    C--       Density of K+1 level (below W(K+1)) reference to K level.
304              CALL FIND_RHO(
305         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
306         O       rhoKp1,
307         I       myThid )
308              CALL CONVECT(
309         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
310         I        myTime,myIter,myThid)
311    C--       Recompute density after mixing
312              CALL FIND_RHO(
313         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
314         O       rhoK,
315         I       myThid )
316             ENDIF
317    C--      Calculate buoyancy
318             CALL CALC_BUOYANCY(
319         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
320         O       buoyK,
321         I       myThid )
322    C--      Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
323             CALL CALC_PHI_HYD(
324         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
325         U        phiHyd,
326         I        myThid )
327    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
328             CALL FIND_RHO(
329         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
330         O        rhoTmp,
331         I        myThid )
332             CALL CALC_ISOSLOPES(
333         I        bi, bj, iMin, iMax, jMin, jMax, K,
334         I        rhoKm1, rhoK, rhotmp,
335         O        K13, K23, K33, KapGM,
336         I        myThid )
337             DO J=jMin,jMax
338              DO I=iMin,iMax
339               rhoKm1 (I,J) = rhoK(I,J)
340               buoyKm1(I,J) = buoyK(I,J)
341              ENDDO
342             ENDDO
343          ENDDO ! K          ENDDO ! K
344    
345  C--     Initial boundary condition on barotropic divergence integral          DO K = Nr, 1, -1
         DO j=1-OLy,sNy+OLy  
          DO i=1-OLx,sNx+OLx  
           cg2d_b(i,j,bi,bj) = 0. _d 0  
          ENDDO  
         ENDDO  
346    
         DO K = Nz, 1, -1  
347           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
348           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
349           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 305  C--     Initial boundary condition on ba Line 355  C--     Initial boundary condition on ba
355  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
356           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
357       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
358       O        xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
359       I        myThid)       I        myThid)
   
360  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
361           CALL CALC_DIFFUSIVITY(           CALL CALC_DIFFUSIVITY(
362       I        bi,bj,iMin,iMax,jMin,jMax,K,       I        bi,bj,iMin,iMax,jMin,jMax,K,
363       I        maskC,maskUp,KapGM,K33,       I        maskC,maskUp,KapGM,K33,
364       O        KappaZT,       O        KappaRT,KappaRS,
365       I        myThid)       I        myThid)
   
366  C--      Calculate accelerations in the momentum equations  C--      Calculate accelerations in the momentum equations
367           IF ( momStepping ) THEN           IF ( momStepping ) THEN
368            CALL CALC_MOM_RHS(            CALL CALC_MOM_RHS(
369       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
370       I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,       I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
371       I         pH,       I         phiHyd,
372       U         aTerm,xTerm,cTerm,mTerm,pTerm,       U         aTerm,xTerm,cTerm,mTerm,pTerm,
373       U         fZon, fMer, fVerU, fVerV,       U         fZon, fMer, fVerU, fVerV,
374       I         myThid)       I         myThid)
375           ENDIF           ENDIF
   
376  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
377           IF ( tempStepping ) THEN           IF ( tempStepping ) THEN
378            CALL CALC_GT(            CALL CALC_GT(
379       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
380       I         xA,yA,uTrans,vTrans,wTrans,maskUp,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
381       I         K13,K23,KappaZT,KapGM,       I         K13,K23,KappaRT,KapGM,
382       U         aTerm,xTerm,fZon,fMer,fVerT,       U         aTerm,xTerm,fZon,fMer,fVerT,
383       I         myThid)       I         myTime, myThid)
384             ENDIF
385             IF ( saltStepping ) THEN
386              CALL CALC_GS(
387         I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
388         I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
389         I         K13,K23,KappaRS,KapGM,
390         U         aTerm,xTerm,fZon,fMer,fVerS,
391         I         myTime, myThid)
392           ENDIF           ENDIF
 Cdbg     CALL CALC_GS(  
 Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  
 Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  
 Cdbg I        K13,K23,K33,KapGM,  
 Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  
 Cdbg I        myThid)  
   
393  C--      Prediction step (step forward all model variables)  C--      Prediction step (step forward all model variables)
394           CALL TIMESTEP(           CALL TIMESTEP(
395       I       bi,bj,iMin,iMax,jMin,jMax,K,       I       bi,bj,iMin,iMax,jMin,jMax,K,
396       I       myThid)       I       myThid)
   
397  C--      Diagnose barotropic divergence of predicted fields  C--      Diagnose barotropic divergence of predicted fields
398           CALL DIV_G(           CALL CALC_DIV_GHAT(
399       I       bi,bj,iMin,iMax,jMin,jMax,K,       I       bi,bj,iMin,iMax,jMin,jMax,K,
400       I       xA,yA,       I       xA,yA,
401       I       myThid)       I       myThid)
402    
403    C--      Cumulative diagnostic calculations (ie. time-averaging)
404    #ifdef ALLOW_DIAGNOSTICS
405             IF (taveFreq.GT.0.) THEN
406              CALL DO_TIME_AVERAGES(
407         I                           myTime, myIter, bi, bj, K, kUp, kDown,
408         I                           K13, K23, rVel, KapGM,
409         I                           myThid )
410             ENDIF
411    #endif
412    
413          ENDDO ! K          ENDDO ! K
414    
415  C--     Implicit diffusion  C--     Implicit diffusion
416          IF (implicitDiffusion) THEN          IF (implicitDiffusion) THEN
417           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
418       I                  KappaZT,       I                  KappaRT,KappaRS,
419       I                  myThid )       I                  myThid )
420          ENDIF          ENDIF
421    
422         ENDDO         ENDDO
423        ENDDO        ENDDO
424    
425        write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
426        write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
427       &                         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.),
428        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.)
429       &                         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.),
430        write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
431       &                         maxval(K13(1:sNx,1:sNy,:))  C     write(0,*) 'dynamics: rVel(1) ',
432        write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
433       &                         maxval(K23(1:sNx,1:sNy,:))  C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
434        write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  C     write(0,*) 'dynamics: rVel(2) ',
435       &                         maxval(K33(1:sNx,1:sNy,:))  C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
436        write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),  C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
437       &                         maxval(gT(1:sNx,1:sNy,:,:,:))  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
438        write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),  cblk &                           maxval(K13(1:sNx,1:sNy,:))
439       &                         maxval(Theta(1:sNx,1:sNy,:,:,:))  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
440        write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),  cblk &                           maxval(K23(1:sNx,1:sNy,:))
441       &                          maxval(pH/(Gravity*Rhonil))  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
442    cblk &                           maxval(K33(1:sNx,1:sNy,:))
443    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
444    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
445    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
446    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
447    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
448    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
449    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
450    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
451    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
452    C    &                           maxval(phiHyd/(Gravity*Rhonil))
453    C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
454    C    &Nr, 1, myThid )
455    C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
456    C    &Nr, 1, myThid )
457    C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
458    C    &Nr, 1, myThid )
459    C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
460    C    &Nr, 1, myThid )
461    
462    
463        RETURN        RETURN
464        END        END

Legend:
Removed from v.1.14  
changed lines
  Added in v.1.37

  ViewVC Help
Powered by ViewVC 1.1.22