/[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.54.2.1 by adcroft, Thu Jan 4 21:41:33 2001 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 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"
# Line 27  C     == Global variables === Line 28  C     == Global variables ===
28  #include "CG2D.h"  #include "CG2D.h"
29  #include "PARAMS.h"  #include "PARAMS.h"
30  #include "DYNVARS.h"  #include "DYNVARS.h"
31    #include "GRID.h"
32    
33    #ifdef ALLOW_AUTODIFF_TAMC
34    # include "tamc.h"
35    # include "tamc_keys.h"
36    #endif /* ALLOW_AUTODIFF_TAMC */
37    
38    #ifdef ALLOW_KPP
39    # include "KPP.h"
40    #endif
41    
42  C     == Routine arguments ==  C     == Routine arguments ==
43  C     myTime - Current time in simulation  C     myTime - Current time in simulation
44  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
45  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
46        _RL myTime        _RL myTime
47        INTEGER myIter        INTEGER myIter
48          INTEGER myThid
49    
50  C     == Local variables  C     == Local variables
51  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
52  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C     uTrans, vTrans, rTrans - Per block temporaries holding flow
53  C     wVel                     o uTrans: Zonal transport  C                              transport
54    C     rVel                     o uTrans: Zonal transport
55  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
56  C                              o wTrans: Vertical transport  C                              o rTrans: Vertical transport
57  C                              o wVel:   Vertical velocity at upper and lower  C                              o rVel:   Vertical velocity at upper and
58  C                                        cell faces.  C                                        lower cell faces.
59  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
60  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
61  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in
# Line 59  C                              o fVer: V Line 71  C                              o fVer: V
71  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
72  C                                      so we need an fVer for each  C                                      so we need an fVer for each
73  C                                      variable.  C                                      variable.
74  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     rhoK, rhoKM1   - Density at current level, level above and level
75  C     jMin, jMax   are applied.  C                      below.
76    C     rhoKP1                                                                  
77    C     buoyK, buoyKM1 - Buoyancy at current level and level above.
78    C     phiHyd         - Hydrostatic part of the potential phiHydi.
79    C                      In z coords phiHydiHyd is the hydrostatic
80    C                      pressure anomaly
81    C                      In p coords phiHydiHyd is the geopotential
82    C                      surface height
83    C                      anomaly.
84    C     etaSurfX,      - Holds surface elevation gradient in X and Y.
85    C     etaSurfY
86    C     KappaRT,       - Total diffusion in vertical for T and S.
87    C     KappaRS          (background + spatially varying, isopycnal term).
88    C     iMin, iMax     - Ranges and sub-block indices on which calculations
89    C     jMin, jMax       are applied.
90  C     bi, bj  C     bi, bj
91  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
92  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
93  C                          into fVerTerm  C                      index into fVerTerm.
94        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99        _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
110        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
111        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
112        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
113        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
114        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116        _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
123        _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
124          _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
125          _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
126          _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
127          _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
128          _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
129    
130    C This is currently also used by IVDC and Diagnostics
131    C #ifdef INCLUDE_CONVECT_CALL
132          _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
133    C #endif
134    
135        INTEGER iMin, iMax        INTEGER iMin, iMax
136        INTEGER jMin, jMax        INTEGER jMin, jMax
137        INTEGER bi, bj        INTEGER bi, bj
138        INTEGER i, j        INTEGER i, j
139        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kup, kDown
140          LOGICAL BOTTOM_LAYER
141    
142    #ifdef ALLOW_AUTODIFF_TAMC
143          INTEGER    isbyte
144          PARAMETER( isbyte = 4 )
145    
146          INTEGER act1, act2, act3, act4
147          INTEGER max1, max2, max3
148          INTEGER iikey, kkey
149          INTEGER maximpl
150    #endif /* ALLOW_AUTODIFF_TAMC */
151    
152  C---    The algorithm...  C---    The algorithm...
153  C  C
# Line 116  C       "Calculation of Gs" Line 163  C       "Calculation of Gs"
163  C       ===================  C       ===================
164  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
165  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
166  C         w = sum_z ( div. u[n] )  C         rVel = sum_r ( div. u[n] )
167  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
168    C         b   = b(rho, theta)
169  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
170  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
171  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
172  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
173  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
174  C  C
175  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
176  C       ================================  C       ================================
# Line 146  C         (1 + dt * K * d_zz) theta[n] = Line 194  C         (1 + dt * K * d_zz) theta[n] =
194  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
195  C---  C---
196    
197    #ifdef ALLOW_AUTODIFF_TAMC
198    C--   dummy statement to end declaration part
199          ikey = 1
200    #endif /* ALLOW_AUTODIFF_TAMC */
201    
202    
203  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
204  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
205  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
# Line 164  C     uninitialised but inert locations. Line 218  C     uninitialised but inert locations.
218          pTerm(i,j)   = 0. _d 0          pTerm(i,j)   = 0. _d 0
219          fZon(i,j)    = 0. _d 0          fZon(i,j)    = 0. _d 0
220          fMer(i,j)    = 0. _d 0          fMer(i,j)    = 0. _d 0
221          DO K=1,nZ          DO k=1,Nr
222           pH (i,j,k)  = 0. _d 0           phiHyd (i,j,k)  = 0. _d 0
223           K13(i,j,k) = 0. _d 0           KappaRU(i,j,k) = 0. _d 0
224           K23(i,j,k) = 0. _d 0           KappaRV(i,j,k) = 0. _d 0
225           K33(i,j,k) = 0. _d 0           sigmaX(i,j,k) = 0. _d 0
226           KappaZT(i,j,k) = 0. _d 0           sigmaY(i,j,k) = 0. _d 0
227             sigmaR(i,j,k) = 0. _d 0
228          ENDDO          ENDDO
229          rhokm1(i,j)  = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
230          rhokp1(i,j)  = 0. _d 0          rhok   (i,j) = 0. _d 0
231          rhotmp(i,j)  = 0. _d 0          rhoKP1 (i,j) = 0. _d 0
232            rhoTMP (i,j) = 0. _d 0
233            buoyKM1(i,j) = 0. _d 0
234            buoyK  (i,j) = 0. _d 0
235            maskC  (i,j) = 0. _d 0
236         ENDDO         ENDDO
237        ENDDO        ENDDO
238    
239    
240    #ifdef ALLOW_AUTODIFF_TAMC
241    C--   HPF directive to help TAMC
242    CHPF$ INDEPENDENT
243    #endif /* ALLOW_AUTODIFF_TAMC */
244    
245        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
246    
247    #ifdef ALLOW_AUTODIFF_TAMC
248    C--    HPF directive to help TAMC
249    CHPF$  INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
250    CHPF$&                  ,phiHyd,utrans,vtrans,maskc,xA,yA
251    CHPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV
252    CHPF$&                  )
253    #endif /* ALLOW_AUTODIFF_TAMC */
254    
255         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
256    
257    #ifdef ALLOW_AUTODIFF_TAMC
258              act1 = bi - myBxLo(myThid)
259              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
260    
261              act2 = bj - myByLo(myThid)
262              max2 = myByHi(myThid) - myByLo(myThid) + 1
263    
264              act3 = myThid - 1
265              max3 = nTx*nTy
266    
267              act4 = ikey_dynamics - 1
268    
269              ikey = (act1 + 1) + act2*max1
270         &                      + act3*max1*max2
271         &                      + act4*max1*max2*max3
272    #endif /* ALLOW_AUTODIFF_TAMC */
273    
274  C--     Set up work arrays that need valid initial values  C--     Set up work arrays that need valid initial values
275          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
276           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
277            wTrans(i,j)  = 0. _d 0            rTrans(i,j)   = 0. _d 0
278            wVel  (i,j,1) = 0. _d 0            rVel  (i,j,1) = 0. _d 0
279            wVel  (i,j,2) = 0. _d 0            rVel  (i,j,2) = 0. _d 0
280            fVerT(i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
281            fVerT(i,j,2) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
282            fVerS(i,j,1) = 0. _d 0            fVerS (i,j,1) = 0. _d 0
283            fVerS(i,j,2) = 0. _d 0            fVerS (i,j,2) = 0. _d 0
284            fVerU(i,j,1) = 0. _d 0            fVerU (i,j,1) = 0. _d 0
285            fVerU(i,j,2) = 0. _d 0            fVerU (i,j,2) = 0. _d 0
286            fVerV(i,j,1) = 0. _d 0            fVerV (i,j,1) = 0. _d 0
287            fVerV(i,j,2) = 0. _d 0            fVerV (i,j,2) = 0. _d 0
288            pH(i,j,1) = 0. _d 0            phiHyd(i,j,1) = 0. _d 0
289            K13(i,j,1) = 0. _d 0           ENDDO
290            K23(i,j,1) = 0. _d 0          ENDDO
291            K33(i,j,1) = 0. _d 0  
292            KapGM(i,j) = 0. _d 0          DO k=1,Nr
293             DO j=1-OLy,sNy+OLy
294              DO i=1-OLx,sNx+OLx
295    #ifdef INCLUDE_CONVECT_CALL
296               ConvectCount(i,j,k) = 0.
297    #endif
298               KappaRT(i,j,k) = 0. _d 0
299               KappaRS(i,j,k) = 0. _d 0
300              ENDDO
301           ENDDO           ENDDO
302          ENDDO          ENDDO
303    
# Line 207  C--     Set up work arrays that need val Line 306  C--     Set up work arrays that need val
306          jMin = 1-OLy+1          jMin = 1-OLy+1
307          jMax = sNy+OLy          jMax = sNy+OLy
308    
309  C--     Calculate gradient of surface pressure          k = 1
310          CALL GRAD_PSURF(          BOTTOM_LAYER = k .EQ. Nr
      I       bi,bj,iMin,iMax,jMin,jMax,  
      O       pSurfX,pSurfY,  
      I       myThid)  
   
 C--     Update fields in top level according to tendency terms  
         CALL CORRECTION_STEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)  
311    
312  C--     Density of 1st level (below W(1)) reference to level 1  C--     Density of 1st level (below W(1)) reference to level 1
313    #ifdef  INCLUDE_FIND_RHO_CALL
314    #ifdef ALLOW_AUTODIFF_TAMC
315    CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
316    CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
317    #endif /* ALLOW_AUTODIFF_TAMC */
318          CALL FIND_RHO(          CALL FIND_RHO(
319       I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,       I     bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
320       O     rhoKm1,       O     rhoKm1,
321       I     myThid )       I     myThid )
322  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #endif
323          CALL CALC_PH(  
324       I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,          IF (.NOT. BOTTOM_LAYER) THEN
325       U      pH,  
326    C--      Check static stability with layer below
327    C--      and mix as needed.
328    #ifdef  INCLUDE_FIND_RHO_CALL
329    #ifdef ALLOW_AUTODIFF_TAMC
330    CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj
331    CADJ &   , key = ikey, byte = isbyte
332    CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj
333    CADJ &   , key = ikey, byte = isbyte
334    #endif /* ALLOW_AUTODIFF_TAMC */
335             CALL FIND_RHO(
336         I      bi, bj, iMin, iMax, jMin, jMax, k+1, k, eosType,
337         O      rhoKp1,
338       I      myThid )       I      myThid )
339          DO J=1-Oly,sNy+Oly  #endif
          DO I=1-Olx,sNx+Olx  
           rhoKp1(I,J)=rhoKm1(I,J)  
          ENDDO  
         ENDDO  
340    
341          DO K=2,Nz  #ifdef ALLOW_AUTODIFF_TAMC
342  C--     Update fields in Kth level according to tendency terms  CADJ STORE rhoKm1(:,:)  = comlev1_bibj, key = ikey, byte = isbyte
343          CALL CORRECTION_STEP(  CADJ STORE rhoKp1(:,:)  = comlev1_bibj, key = ikey, byte = isbyte
344       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)  #endif /* ALLOW_AUTODIFF_TAMC */
345  C--     Density of K-1 level (above W(K)) reference to K-1 level  
346  copt    CALL FIND_RHO(  C--      Implicit Vertical Diffusion for Convection
347  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,           IF (ivdc_kappa.NE.0.) THEN
348  copt O     rhoKm1,              CALL CALC_IVDC(
349  copt I     myThid )       I       bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,
350  C       rhoKm1=rhoKp1       U       ConvectCount, KappaRT, KappaRS,
351          DO J=1-Oly,sNy+Oly       I       myTime,myIter,myThid)
352           DO I=1-Olx,sNx+Olx           ENDIF
353            rhoKm1(I,J)=rhoKp1(I,J)  
354           ENDDO          ENDIF
355          ENDDO  
356  C--     Density of K level (below W(K)) reference to K level  C--     Calculate buoyancy
357          CALL FIND_RHO(          CALL CALC_BUOYANCY(
358       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,       I      bi,bj,iMin,iMax,jMin,jMax,k,rhoKm1,
359       O     rhoKp1,       O      buoyKm1,
360       I     myThid )       I      myThid )
361  C--     Density of K-1 level (above W(K)) reference to K level  
362          CALL FIND_RHO(  C--     Integrate hydrostatic balance for phiHyd with BC of
363       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,  C--     phiHyd(z=0)=0
364       O     rhotmp,          CALL CALC_PHI_HYD(
365       I     myThid )       I      bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyKm1,
366  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation       U      phiHyd,
367          CALL CALC_ISOSLOPES(       I      myThid )
368       I            bi, bj, iMin, iMax, jMin, jMax, K,  
369       I            rhoKm1, rhoKp1, rhotmp,  #ifdef ALLOW_GMREDI
370       O            K13, K23, K33, KapGM,          IF ( useGMRedi ) THEN
371            CALL GRAD_SIGMA(
372         I            bi, bj, iMin, iMax, jMin, jMax, k,
373         I            rhoKm1, rhoKm1, rhoKm1,
374         O            sigmaX, sigmaY, sigmaR,
375       I            myThid )       I            myThid )
376  C--     Calculate static stability and mix where convectively unstable          ENDIF
377          CALL CONVECT(  #endif
378       I      bi,bj,iMin,iMax,jMin,jMax,K,rhotmp,rhoKp1,  
379       I      myTime,myIter,myThid)  C--     Start of downward loop
380  C--     Density of K-1 level (above W(K)) reference to K-1 level          DO k=2,Nr
381          CALL FIND_RHO(  
382       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,  #ifdef ALLOW_AUTODIFF_TAMC
383       O     rhoKm1,           kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
384       I     myThid )  #endif /* ALLOW_AUTODIFF_TAMC */
385  C--     Density of K level (below W(K)) referenced to K level  
386          CALL FIND_RHO(           BOTTOM_LAYER = k .EQ. Nr
387       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
388       O     rhoKp1,  #ifdef DO_PIPELINED_CORRECTION_STEP
389       I     myThid )           IF ( .NOT. BOTTOM_LAYER ) THEN
390  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  C--       Update fields in layer below according to tendency terms
391          CALL CALC_PH(  #ifdef ALLOW_OBCS
392       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,            IF (openBoundaries) THEN
393       U      pH,  #ifdef ALLOW_AUTODIFF_TAMC
394    CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj_k
395    CADJ &   , key = kkey, byte = isbyte
396    CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj_k
397    CADJ &   , key = kkey, byte = isbyte
398    CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k
399    CADJ &   , key = kkey, byte = isbyte
400    CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k
401    CADJ &   , key = kkey, byte = isbyte
402    #endif /* ALLOW_AUTODIFF_TAMC */
403              END IF
404    #endif
405             ENDIF
406    #endif /* DO_PIPELINED_CORRECTION_STEP */
407    
408    C--      Density of k level (below W(k)) reference to k level
409    #ifdef  INCLUDE_FIND_RHO_CALL
410    #ifdef ALLOW_AUTODIFF_TAMC
411    CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k
412    CADJ &   , key = kkey, byte = isbyte
413    CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k
414    CADJ &   , key = kkey, byte = isbyte
415    #endif /* ALLOW_AUTODIFF_TAMC */
416             CALL FIND_RHO(
417         I      bi, bj, iMin, iMax, jMin, jMax,  k, k, eosType,
418         O      rhoK,
419       I      myThid )       I      myThid )
420    
421          ENDDO ! K  #ifdef ALLOW_AUTODIFF_TAMC
422    cph(   storing not necessary
423    cphCADJ STORE rhoK(:,:)  = comlev1_bibj_k, key = kkey, byte = isbyte
424    cph)
425    #endif /* ALLOW_AUTODIFF_TAMC */
426    #endif
427    
428             IF (.NOT. BOTTOM_LAYER) THEN
429    
430    C--       Check static stability with layer below and mix as needed.
431    C--       Density of k+1 level (below W(k+1)) reference to k level.
432    #ifdef  INCLUDE_FIND_RHO_CALL
433    #ifdef ALLOW_AUTODIFF_TAMC
434    CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj_k
435    CADJ &   , key = kkey, byte = isbyte
436    CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj_k
437    CADJ &   , key = kkey, byte = isbyte
438    #endif /* ALLOW_AUTODIFF_TAMC */
439              CALL FIND_RHO(
440         I       bi, bj, iMin, iMax, jMin, jMax,  k+1, k, eosType,
441         O       rhoKp1,
442         I       myThid )
443    #ifdef ALLOW_AUTODIFF_TAMC
444    CADJ STORE rhoKp1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
445    #endif /* ALLOW_AUTODIFF_TAMC */
446    #endif
447    
448    C--      Implicit Vertical Diffusion for Convection
449             IF (ivdc_kappa.NE.0.) THEN
450    #ifdef ALLOW_AUTODIFF_TAMC
451    CADJ STORE rhoKm1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
452    #endif /* ALLOW_AUTODIFF_TAMC */
453                CALL CALC_IVDC(
454         I       bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,
455         U       ConvectCount, KappaRT, KappaRS,
456         I       myTime,myIter,myThid)
457             END IF
458    
459  C--     Initial boundary condition on barotropic divergence integral  C--            IF (.NOT. BOTTOM_LAYER) ends here
460          DO j=1-OLy,sNy+OLy           ENDIF
461           DO i=1-OLx,sNx+OLx  
462            cg2d_b(i,j,bi,bj) = 0. _d 0  C--      Calculate buoyancy
463             CALL CALC_BUOYANCY(
464         I       bi,bj,iMin,iMax,jMin,jMax,k,rhoK,
465         O       buoyK,
466         I       myThid )
467    
468    C--      Integrate hydrostatic balance for phiHyd with BC of
469    C--      phiHyd(z=0)=0
470             CALL CALC_PHI_HYD(
471         I        bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyK,
472         U        phiHyd,
473         I        myThid )
474    
475    #ifdef  INCLUDE_FIND_RHO_CALL
476    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
477    
478    #ifdef ALLOW_AUTODIFF_TAMC
479    CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k
480    CADJ &   , key = kkey, byte = isbyte
481    CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k
482    CADJ &   , key = kkey, byte = isbyte
483    #endif /* ALLOW_AUTODIFF_TAMC */
484    
485             CALL FIND_RHO(
486         I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
487         O        rhoTmp,
488         I        myThid )
489    #endif
490    
491    
492    #ifdef ALLOW_GMREDI
493             IF ( useGMRedi ) THEN
494             CALL GRAD_SIGMA(
495         I             bi, bj, iMin, iMax, jMin, jMax, k,
496         I             rhoK, rhotmp, rhoK,
497         O             sigmaX, sigmaY, sigmaR,
498         I             myThid )
499             ENDIF
500    #endif
501    
502             DO J=jMin,jMax
503              DO I=iMin,iMax
504    #ifdef  INCLUDE_FIND_RHO_CALL
505               rhoKm1 (I,J) = rhoK(I,J)
506    #endif
507               buoyKm1(I,J) = buoyK(I,J)
508              ENDDO
509           ENDDO           ENDDO
510    
511    C--     end of k loop
512          ENDDO          ENDDO
513    
514          DO K = Nz, 1, -1  C     Determines forcing terms based on external fields
515           kM1  =max(1,k-1)   ! Points to level above k (=k-1)  C     relaxation terms, etc.
516           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above        CALL EXTERNAL_FORCING_SURF(
517           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer       I             bi, bj, iMin, iMax, jMin, jMax,
518         I             myThid )
519    
520    #ifdef ALLOW_GMREDI
521            IF (useGMRedi) THEN
522              DO k=1, Nr
523                CALL GMREDI_CALC_TENSOR(
524         I             bi, bj, iMin, iMax, jMin, jMax, k,
525         I             sigmaX, sigmaY, sigmaR,
526         I             myThid )
527              ENDDO
528            ENDIF
529    #endif
530    
531    #ifdef ALLOW_AUTODIFF_TAMC
532    CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
533    CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
534    
535    CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
536    CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
537    CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
538    CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
539    
540    C--     dummy initialization to break data flow because
541    C--     calc_div_ghat has a condition for initialization
542            DO J=jMin,jMax
543               DO I=iMin,iMax
544                  cg2d_b(i,j,bi,bj) = 0.0
545               ENDDO
546            ENDDO
547    #endif /* ALLOW_AUTODIFF_TAMC */
548    
549    #ifdef ALLOW_KPP
550    C--   Compute KPP mixing coefficients
551            IF (useKPP) THEN
552    
553              CALL TIMER_START('KPP_CALC               [DYNAMICS]', myThid)
554              CALL KPP_CALC(
555         I                  bi, bj, myTime, myThid )
556              CALL TIMER_STOP ('KPP_CALC               [DYNAMICS]', myThid)
557    
558    #ifdef ALLOW_AUTODIFF_TAMC
559            ELSE
560              DO j=1-OLy,sNy+OLy
561                DO i=1-OLx,sNx+OLx
562                  KPPhbl (i,j,bi,bj) = 1.0
563                  KPPfrac(i,j,bi,bj) = 0.0
564                  DO k = 1,Nr
565                     KPPghat   (i,j,k,bi,bj) = 0.0
566                     KPPviscAz (i,j,k,bi,bj) = viscAz
567                     KPPdiffKzT(i,j,k,bi,bj) = diffKzT
568                     KPPdiffKzS(i,j,k,bi,bj) = diffKzS
569                  ENDDO
570                ENDDO
571              ENDDO
572    #endif /* ALLOW_AUTODIFF_TAMC */
573            ENDIF
574    
575    #ifdef ALLOW_AUTODIFF_TAMC
576    CADJ STORE KPPghat   (:,:,:,bi,bj)
577    CADJ &   , KPPviscAz (:,:,:,bi,bj)
578    CADJ &   , KPPdiffKzT(:,:,:,bi,bj)
579    CADJ &   , KPPdiffKzS(:,:,:,bi,bj)
580    CADJ &   , KPPfrac   (:,:  ,bi,bj)
581    CADJ &                 = comlev1_bibj, key = ikey, byte = isbyte
582    #endif /* ALLOW_AUTODIFF_TAMC */
583    
584    #endif /* ALLOW_KPP */
585    
586    C--     Start of upward loop
587            DO k = Nr, 1, -1
588    
589    C--      km1    Points to level above k (=k-1)
590    C--      kup    Cycles through 1,2 to point to layer above
591    C--      kDown  Cycles through 2,1 to point to current layer
592    
593             km1  =max(1,k-1)
594             kup  =1+MOD(k+1,2)
595             kDown=1+MOD(k,2)
596    
597           iMin = 1-OLx+2           iMin = 1-OLx+2
598           iMax = sNx+OLx-1           iMax = sNx+OLx-1
599           jMin = 1-OLy+2           jMin = 1-OLy+2
600           jMax = sNy+OLy-1           jMax = sNy+OLy-1
601    
602    #ifdef ALLOW_AUTODIFF_TAMC
603             kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
604    
605    CADJ STORE rvel  (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
606    CADJ STORE rTrans(:,:)       = comlev1_bibj_k, key = kkey, byte = isbyte
607    CADJ STORE KappaRT(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
608    CADJ STORE KappaRS(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
609    #endif /* ALLOW_AUTODIFF_TAMC */
610    
611  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
612           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
613       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
614       O        xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
615       I        myThid)       I        myThid)
616    
617    #ifdef ALLOW_OBCS
618            IF (openBoundaries) THEN
619             CALL APPLY_OBCS3( bi, bj, k, kup, rTrans, rVel, myThid )
620            ENDIF
621    #endif
622    
623    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
624  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
625           CALL CALC_DIFFUSIVITY(           CALL CALC_DIFFUSIVITY(
626       I        bi,bj,iMin,iMax,jMin,jMax,K,       I        bi,bj,iMin,iMax,jMin,jMax,k,
627       I        maskC,maskUp,KapGM,K33,       I        maskC,maskUp,
628       O        KappaZT,       O        KappaRT,KappaRS,KappaRU,KappaRV,
629       I        myThid)       I        myThid)
630    #endif
631  C--      Calculate accelerations in the momentum equations  C--      Calculate accelerations in the momentum equations
632           IF ( momStepping ) THEN           IF ( momStepping ) THEN
633            CALL CALC_MOM_RHS(            CALL CALC_MOM_RHS(
634       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
635       I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,       I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
636       I         pH,       I         phiHyd,KappaRU,KappaRV,
637       U         aTerm,xTerm,cTerm,mTerm,pTerm,       U         aTerm,xTerm,cTerm,mTerm,pTerm,
638       U         fZon, fMer, fVerU, fVerV,       U         fZon, fMer, fVerU, fVerV,
639       I         myThid)       I         myTime, myThid)
640    #ifdef ALLOW_AUTODIFF_TAMC
641    #ifdef INCLUDE_CD_CODE
642             ELSE
643                DO j=1-OLy,sNy+OLy
644                   DO i=1-OLx,sNx+OLx
645                      guCD(i,j,k,bi,bj) = 0.0
646                      gvCD(i,j,k,bi,bj) = 0.0
647                   END DO
648                END DO
649    #endif
650    #endif /* ALLOW_AUTODIFF_TAMC */
651           ENDIF           ENDIF
   
652  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
653           IF ( tempStepping ) THEN           IF ( tempStepping ) THEN
654            CALL CALC_GT(            CALL CALC_GT(
655       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
656       I         xA,yA,uTrans,vTrans,wTrans,maskUp,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
657       I         K13,K23,KappaZT,KapGM,       I         KappaRT,
658       U         aTerm,xTerm,fZon,fMer,fVerT,       U         aTerm,xTerm,fZon,fMer,fVerT,
659       I         myThid)       I         myTime, myThid)
660           ENDIF           ENDIF
661  Cdbg     CALL CALC_GS(           IF ( saltStepping ) THEN
662  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GS(
663  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
664  Cdbg I        K13,K23,K33,KapGM,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
665  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,       I         KappaRS,
666  Cdbg I        myThid)       U         aTerm,xTerm,fZon,fMer,fVerS,
667         I         myTime, myThid)
668             ENDIF
669    #ifdef ALLOW_OBCS
670    C--      Calculate future values on open boundaries
671             IF (openBoundaries) THEN
672    Caja      CALL CYCLE_OBCS( k, bi, bj, myThid )
673              CALL SET_OBCS( k, bi, bj, myTime+deltaTclock, myThid )
674             ENDIF
675    #endif
676  C--      Prediction step (step forward all model variables)  C--      Prediction step (step forward all model variables)
677           CALL TIMESTEP(           CALL TIMESTEP(
678       I       bi,bj,iMin,iMax,jMin,jMax,K,       I       bi,bj,iMin,iMax,jMin,jMax,k,
679       I       myThid)       I       myIter, myThid)
680    #ifdef ALLOW_OBCS
681    C--      Apply open boundary conditions
682             IF (openBoundaries) THEN
683    #ifdef ALLOW_AUTODIFF_TAMC
684    CADJ STORE gunm1(:,:,k,bi,bj) = comlev1_bibj_k
685    CADJ &   , key = kkey, byte = isbyte
686    CADJ STORE gvnm1(:,:,k,bi,bj) = comlev1_bibj_k
687    CADJ &   , key = kkey, byte = isbyte
688    CADJ STORE gwnm1(:,:,k,bi,bj) = comlev1_bibj_k
689    CADJ &   , key = kkey, byte = isbyte
690    #endif /* ALLOW_AUTODIFF_TAMC */
691    
692                CALL APPLY_OBCS2( bi, bj, k, myThid )
693             END IF
694    #endif
695    C--      Freeze water
696             IF (allowFreezing) THEN
697    #ifdef ALLOW_AUTODIFF_TAMC
698    CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
699    CADJ &   , key = kkey, byte = isbyte
700    #endif /* ALLOW_AUTODIFF_TAMC */
701                CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
702             END IF
703    
704    #ifdef DIVG_IN_DYNAMICS
705  C--      Diagnose barotropic divergence of predicted fields  C--      Diagnose barotropic divergence of predicted fields
706           CALL DIV_G(           CALL CALC_DIV_GHAT(
707       I       bi,bj,iMin,iMax,jMin,jMax,K,       I       bi,bj,iMin,iMax,jMin,jMax,k,
708       I       xA,yA,       I       xA,yA,
709       I       myThid)       I       myThid)
710    #endif /* DIVG_IN_DYNAMICS */
711    
712    C--      Cumulative diagnostic calculations (ie. time-averaging)
713    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
714             IF (taveFreq.GT.0.) THEN
715              CALL DO_TIME_AVERAGES(
716         I                           myTime, myIter, bi, bj, k, kup, kDown,
717         I                           rVel, ConvectCount,
718         I                           myThid )
719             ENDIF
720    #endif
721    
722    
723    C--     k loop
724            ENDDO
725    
726          ENDDO ! K  #ifdef ALLOW_AUTODIFF_TAMC
727               maximpl = 6
728               iikey = (ikey-1)*maximpl
729    #endif /* ALLOW_AUTODIFF_TAMC */
730    
731  C--     Implicit diffusion  C--     Implicit diffusion
732          IF (implicitDiffusion) THEN          IF (implicitDiffusion) THEN
733           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,  
734       I                  KappaZT,           IF (tempStepping) THEN
735       I                  myThid )  #ifdef ALLOW_AUTODIFF_TAMC
736                idkey = iikey + 1
737    #endif /* ALLOW_AUTODIFF_TAMC */
738                CALL IMPLDIFF(
739         I         bi, bj, iMin, iMax, jMin, jMax,
740         I         deltaTtracer, KappaRT,recip_HFacC,
741         U         gTNm1,
742         I         myThid )
743             END IF
744    
745             IF (saltStepping) THEN
746    #ifdef ALLOW_AUTODIFF_TAMC
747             idkey = iikey + 2
748    #endif /* ALLOW_AUTODIFF_TAMC */
749                CALL IMPLDIFF(
750         I         bi, bj, iMin, iMax, jMin, jMax,
751         I         deltaTtracer, KappaRS,recip_HFacC,
752         U         gSNm1,
753         I         myThid )
754             END IF
755    
756    C--     implicitDiffusion
757            ENDIF
758    
759    C--     Implicit viscosity
760            IF (implicitViscosity) THEN
761    
762             IF (momStepping) THEN
763    #ifdef ALLOW_AUTODIFF_TAMC
764             idkey = iikey + 3
765    #endif /* ALLOW_AUTODIFF_TAMC */
766              CALL IMPLDIFF(
767         I         bi, bj, iMin, iMax, jMin, jMax,
768         I         deltaTmom, KappaRU,recip_HFacW,
769         U         gUNm1,
770         I         myThid )
771    #ifdef ALLOW_AUTODIFF_TAMC
772             idkey = iikey + 4
773    #endif /* ALLOW_AUTODIFF_TAMC */
774              CALL IMPLDIFF(
775         I         bi, bj, iMin, iMax, jMin, jMax,
776         I         deltaTmom, KappaRV,recip_HFacS,
777         U         gVNm1,
778         I         myThid )
779    
780    #ifdef INCLUDE_CD_CODE
781    
782    #ifdef ALLOW_AUTODIFF_TAMC
783             idkey = iikey + 5
784    #endif /* ALLOW_AUTODIFF_TAMC */
785              CALL IMPLDIFF(
786         I         bi, bj, iMin, iMax, jMin, jMax,
787         I         deltaTmom, KappaRU,recip_HFacW,
788         U         vVelD,
789         I         myThid )
790    #ifdef ALLOW_AUTODIFF_TAMC
791            idkey = iikey + 6
792    #endif /* ALLOW_AUTODIFF_TAMC */
793              CALL IMPLDIFF(
794         I         bi, bj, iMin, iMax, jMin, jMax,
795         I         deltaTmom, KappaRV,recip_HFacS,
796         U         uVelD,
797         I         myThid )
798    
799    #endif
800    
801    C--      momStepping
802             ENDIF
803    
804    C--     implicitViscosity
805          ENDIF          ENDIF
806    
807         ENDDO         ENDDO
808        ENDDO        ENDDO
809    
       write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)  
       write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),  
      &                         maxval(uVel(1:sNx,1:sNy,:,:,:))  
       write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),  
      &                         maxval(vVel(1:sNx,1:sNy,:,:,:))  
       write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  
      &                         maxval(K13(1:sNx,1:sNy,:))  
       write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  
      &                         maxval(K23(1:sNx,1:sNy,:))  
       write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  
      &                         maxval(K33(1:sNx,1:sNy,:))  
       write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),  
      &                         maxval(gT(1:sNx,1:sNy,:,:,:))  
       write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),  
      &                         maxval(Theta(1:sNx,1:sNy,:,:,:))  
       write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),  
      &                          maxval(pH/(Gravity*Rhonil))  
   
810        RETURN        RETURN
811        END        END

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

  ViewVC Help
Powered by ViewVC 1.1.22