/[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.26 by cnh, Wed Aug 19 16:20:49 1998 UTC revision 1.58 by adcroft, Fri Feb 2 21:04:48 2001 UTC
# 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     fVer[STUV]               o fVer: Vertical flux term - note fVer
 C     mTerm, pTerm,            tendency equations.  
 C     fZon, fMer, fVer[STUV]   o aTerm: Advection term  
 C                              o xTerm: Mixing term  
 C                              o cTerm: Coriolis term  
 C                              o mTerm: Metric term  
 C                              o pTerm: Pressure term  
 C                              o fZon: Zonal flux term  
 C                              o fMer: Meridional flux term  
 C                              o fVer: Vertical flux term - note fVer  
62  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
63  C                                      so we need an fVer for each  C                                      so we need an fVer for each
64  C                                      variable.  C                                      variable.
65  C     rhoK, rhoKM1   - Density at current level, level above and level below.  C     rhoK, rhoKM1   - Density at current level, and level above
66  C     rhoKP1                                                                    C     phiHyd         - Hydrostatic part of the potential phiHydi.
67  C     buoyK, buoyKM1 - Buoyancy at current level and level above.  C                      In z coords phiHydiHyd is the hydrostatic
68  C     phiHyd         - Hydrostatic part of the potential phi.  C                      pressure anomaly
69  C                      In z coords phiHyd is the hydrostatic pressure anomaly  C                      In p coords phiHydiHyd is the geopotential
70  C                      In p coords phiHyd is the geopotential surface height anomaly.  C                      surface height
71  C     iMin, iMax - Ranges and sub-block indices on which calculations  C                      anomaly.
72  C     jMin, jMax   are applied.  C     etaSurfX,      - Holds surface elevation gradient in X and Y.
73    C     etaSurfY
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
81  C                          into fVerTerm  C                      index 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 fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
91        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
92        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
93        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
94        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
98        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
99        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
100        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
101        _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
103        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
104        _RL rhok  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
105        _RL buoyKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C This is currently also used by IVDC and Diagnostics
106        _RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C #ifdef INCLUDE_CONVECT_CALL
107        _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
108        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C #endif
       _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)  
       _RL KappaZS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)  
109    
110        INTEGER iMin, iMax        INTEGER iMin, iMax
111        INTEGER jMin, jMax        INTEGER jMin, jMax
112        INTEGER bi, bj        INTEGER bi, bj
113        INTEGER i, j        INTEGER i, j
114        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kup, kDown
115        LOGICAL BOTTOM_LAYER  
116    #ifdef ALLOW_AUTODIFF_TAMC
117          INTEGER    isbyte
118          PARAMETER( isbyte = 4 )
119    
120          INTEGER act1, act2, act3, act4
121          INTEGER max1, max2, max3
122          INTEGER iikey, kkey
123          INTEGER maximpl
124    #endif /* ALLOW_AUTODIFF_TAMC */
125    
126  C---    The algorithm...  C---    The algorithm...
127  C  C
# Line 127  C       "Calculation of Gs" Line 137  C       "Calculation of Gs"
137  C       ===================  C       ===================
138  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
139  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
140  C         w = sum_z ( div. u[n] )  C         rVel = sum_r ( div. u[n] )
141  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
142    C         b   = b(rho, theta)
143  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
144  C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )  C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
145  C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )  C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
146  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
147  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
148  C  C
149  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
150  C       ================================  C       ================================
# Line 157  C         (1 + dt * K * d_zz) theta[n] = Line 168  C         (1 + dt * K * d_zz) theta[n] =
168  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
169  C---  C---
170    
171    #ifdef ALLOW_AUTODIFF_TAMC
172    C--   dummy statement to end declaration part
173          ikey = 1
174    #endif /* ALLOW_AUTODIFF_TAMC */
175    
176  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
177  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
178  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
# Line 168  C     uninitialised but inert locations. Line 184  C     uninitialised but inert locations.
184          yA(i,j)      = 0. _d 0          yA(i,j)      = 0. _d 0
185          uTrans(i,j)  = 0. _d 0          uTrans(i,j)  = 0. _d 0
186          vTrans(i,j)  = 0. _d 0          vTrans(i,j)  = 0. _d 0
187          aTerm(i,j)   = 0. _d 0          DO k=1,Nr
188          xTerm(i,j)   = 0. _d 0           phiHyd(i,j,k)  = 0. _d 0
189          cTerm(i,j)   = 0. _d 0           KappaRU(i,j,k) = 0. _d 0
190          mTerm(i,j)   = 0. _d 0           KappaRV(i,j,k) = 0. _d 0
191          pTerm(i,j)   = 0. _d 0           sigmaX(i,j,k) = 0. _d 0
192          fZon(i,j)    = 0. _d 0           sigmaY(i,j,k) = 0. _d 0
193          fMer(i,j)    = 0. _d 0           sigmaR(i,j,k) = 0. _d 0
         DO K=1,nZ  
          pH (i,j,k)  = 0. _d 0  
          K13(i,j,k) = 0. _d 0  
          K23(i,j,k) = 0. _d 0  
          K33(i,j,k) = 0. _d 0  
          KappaZT(i,j,k) = 0. _d 0  
194          ENDDO          ENDDO
195          rhokm1(i,j)  = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
196          rhok  (i,j)  = 0. _d 0          rhok   (i,j) = 0. _d 0
197          rhokp1(i,j)  = 0. _d 0          maskC  (i,j) = 0. _d 0
         rhotmp(i,j)  = 0. _d 0  
         buoyKM1(i,j) = 0. _d 0  
         buoyK  (i,j) = 0. _d 0  
         maskC (i,j)  = 0. _d 0  
198         ENDDO         ENDDO
199        ENDDO        ENDDO
200    
201    
202    #ifdef ALLOW_AUTODIFF_TAMC
203    C--   HPF directive to help TAMC
204    CHPF$ INDEPENDENT
205    #endif /* ALLOW_AUTODIFF_TAMC */
206    
207        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
208    
209    #ifdef ALLOW_AUTODIFF_TAMC
210    C--    HPF directive to help TAMC
211    CHPF$  INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
212    CHPF$&                  ,phiHyd,utrans,vtrans,maskc,xA,yA
213    CHPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV
214    CHPF$&                  )
215    #endif /* ALLOW_AUTODIFF_TAMC */
216    
217         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
218    
219    #ifdef ALLOW_AUTODIFF_TAMC
220              act1 = bi - myBxLo(myThid)
221              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
222    
223              act2 = bj - myByLo(myThid)
224              max2 = myByHi(myThid) - myByLo(myThid) + 1
225    
226              act3 = myThid - 1
227              max3 = nTx*nTy
228    
229              act4 = ikey_dynamics - 1
230    
231              ikey = (act1 + 1) + act2*max1
232         &                      + act3*max1*max2
233         &                      + act4*max1*max2*max3
234    #endif /* ALLOW_AUTODIFF_TAMC */
235    
236  C--     Set up work arrays that need valid initial values  C--     Set up work arrays that need valid initial values
237          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
238           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
239            wTrans(i,j)  = 0. _d 0            rTrans(i,j)   = 0. _d 0
240            wVel  (i,j,1) = 0. _d 0            rVel  (i,j,1) = 0. _d 0
241            wVel  (i,j,2) = 0. _d 0            rVel  (i,j,2) = 0. _d 0
242            fVerT(i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
243            fVerT(i,j,2) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
244            fVerS(i,j,1) = 0. _d 0            fVerS (i,j,1) = 0. _d 0
245            fVerS(i,j,2) = 0. _d 0            fVerS (i,j,2) = 0. _d 0
246            fVerU(i,j,1) = 0. _d 0            fVerU (i,j,1) = 0. _d 0
247            fVerU(i,j,2) = 0. _d 0            fVerU (i,j,2) = 0. _d 0
248            fVerV(i,j,1) = 0. _d 0            fVerV (i,j,1) = 0. _d 0
249            fVerV(i,j,2) = 0. _d 0            fVerV (i,j,2) = 0. _d 0
250            pH(i,j,1) = 0. _d 0           ENDDO
251            K13(i,j,1) = 0. _d 0          ENDDO
252            K23(i,j,1) = 0. _d 0  
253            K33(i,j,1) = 0. _d 0          DO k=1,Nr
254            KapGM(i,j) = GMkbackground           DO j=1-OLy,sNy+OLy
255              DO i=1-OLx,sNx+OLx
256    #ifdef INCLUDE_CONVECT_CALL
257               ConvectCount(i,j,k) = 0.
258    #endif
259               KappaRT(i,j,k) = 0. _d 0
260               KappaRS(i,j,k) = 0. _d 0
261              ENDDO
262           ENDDO           ENDDO
263          ENDDO          ENDDO
264    
# Line 222  C--     Set up work arrays that need val Line 267  C--     Set up work arrays that need val
267          jMin = 1-OLy+1          jMin = 1-OLy+1
268          jMax = sNy+OLy          jMax = sNy+OLy
269    
         K = 1  
         BOTTOM_LAYER = K .EQ. Nz  
270    
271  C--     Calculate gradient of surface pressure  C--     Start of diagnostic loop
272          CALL GRAD_PSURF(          DO k=Nr,1,-1
273       I       bi,bj,iMin,iMax,jMin,jMax,  
274       O       pSurfX,pSurfY,  #ifdef ALLOW_AUTODIFF_TAMC
275       I       myThid)  C? Patrick, is this formula correct now that we change the loop range?
276    C? Do we still need this?
277  C--     Update fields in top level according to tendency terms           kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
278          CALL CORRECTION_STEP(  #endif /* ALLOW_AUTODIFF_TAMC */
279       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid)  
280    C--       Integrate continuity vertically for vertical velocity
281          IF ( .NOT. BOTTOM_LAYER ) THEN            CALL INTEGRATE_FOR_W(
282  C--      Update fields in layer below according to tendency terms       I                         bi, bj, k, uVel, vVel,
283           CALL CORRECTION_STEP(       O                         wVel,
284       I        bi,bj,iMin,iMax,jMin,jMax,K+1,pSurfX,pSurfY,myTime,myThid)       I                         myThid )
285          ENDIF  
286  C--     Density of 1st level (below W(1)) reference to level 1  #ifdef    ALLOW_OBCS
287          CALL FIND_RHO(  #ifdef    ALLOW_NONHYDROSTATIC
288       I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  C--       Calculate future values on open boundaries
289       O     rhoKm1,            IF (useOBCS.AND.nonHydrostatic) THEN
290       I     myThid )              CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
291              ENDIF
292          IF ( .NOT. BOTTOM_LAYER ) THEN  #endif    /* ALLOW_NONHYDROSTATIC */
293    #endif    /* ALLOW_OBCS */
294  C--      Check static stability with layer below  
295  C        and mix as needed.  C--       Calculate gradients of potential density for isoneutral
296           CALL FIND_RHO(  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
297       I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,  c         IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
298       O      rhoKp1,            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
299       I      myThid )              CALL FIND_RHO(
300           CALL CONVECT(       I        bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
301       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,       I        theta, salt,
302       I       myTime,myIter,myThid)       O        rhoK,
303  C--      Recompute density after mixing       I        myThid )
304           CALL FIND_RHO(              IF (k.GT.1) CALL FIND_RHO(
305       I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,       I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
306       O      rhoKm1,       I        theta, salt,
307       I      myThid )       O        rhoKm1,
308         I        myThid )
309                CALL GRAD_SIGMA(
310         I             bi, bj, iMin, iMax, jMin, jMax, k,
311         I             rhoK, rhoKm1, rhoK,
312         O             sigmaX, sigmaY, sigmaR,
313         I             myThid )
314              ENDIF
315    
316    C--       Implicit Vertical Diffusion for Convection
317    c ==> should use sigmaR !!!
318              IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
319                CALL CALC_IVDC(
320         I        bi, bj, iMin, iMax, jMin, jMax, k,
321         I        rhoKm1, rhoK,
322         U        ConvectCount, KappaRT, KappaRS,
323         I        myTime, myIter, myThid)
324              END IF
325    
326    C--     end of diagnostic k loop (Nr:1)
327            ENDDO
328    
329    #ifdef  ALLOW_OBCS
330    C--     Calculate future values on open boundaries
331            IF (useOBCS) THEN
332              CALL OBCS_CALC( bi, bj, myTime+deltaT,
333         I            uVel, vVel, wVel, theta, salt,
334         I            myThid )
335          ENDIF          ENDIF
336    #endif  /* ALLOW_OBCS */
337    
338  C--     Calculate buoyancy  C--     Determines forcing terms based on external fields
339          CALL CALC_BUOY(  C       relaxation terms, etc.
340       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,          CALL EXTERNAL_FORCING_SURF(
341       O      buoyKm1,       I             bi, bj, iMin, iMax, jMin, jMax,
342       I      myThid )       I             myThid )
   
 C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
         CALL CALC_PHI_HYD(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,  
      U      phiHyd,  
      I      myThid )  
   
         DO K=2,Nz  
   
          BOTTOM_LAYER = K .EQ. Nz  
          IF ( .NOT. BOTTOM_LAYER ) THEN  
 C--       Update fields in layer below according to tendency terms  
           CALL CORRECTION_STEP(  
      I         bi,bj,iMin,iMax,jMin,jMax,K+1,pSurfX,pSurfY,myTime,myThid)  
          ENDIF  
 C--      Density of K level (below W(K)) reference to K level  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
      O      rhoK,  
      I      myThid )  
          IF ( .NOT. BOTTOM_LAYER ) THEN  
 C--       Check static stability with layer below  
 C         and mix as needed.  
 C--       Density of K+1 level (below W(K+1)) reference to K level  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,  
      O       rhoKp1,  
      I       myThid )  
           CALL CONVECT(  
      I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,  
      I        myTime,myIter,myThid)  
 C--       Recompute density after mixing  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,  
      O       rhoK,  
      I       myThid )  
          ENDIF  
343    
344  C--      Calculate buoyancy  #ifdef  ALLOW_GMREDI
345           CALL CALC_BUOY(  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
346       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,          IF (useGMRedi) THEN
347       O       buoyK,            DO k=1,Nr
348       I       myThid )              CALL GMREDI_CALC_TENSOR(
349         I             bi, bj, iMin, iMax, jMin, jMax, k,
350  C--      Integrate hydrostatic balance for pH with BC of pH(z=0)=0       I             sigmaX, sigmaY, sigmaR,
          CALL CALC_PHI_HYD(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,  
      U       phiHyd,  
      I       myThid )  
 C--      Calculate iso-neutral slopes for the GM/Redi parameterisation  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,  
      O      rhoTmp,  
      I      myThid )  
          CALL CALC_ISOSLOPES(  
      I             bi, bj, iMin, iMax, jMin, jMax, K,  
      I             rhoKm1, rhoK, rhotmp,  
      O             K13, K23, K33, KapGM,  
351       I             myThid )       I             myThid )
          DO J=jMin,jMax  
           DO I=iMin,iMax  
            rhoKm1(I,J) =rhoK(I,J)  
            buoyKm1(I,J)=buoyK(I,J)  
352            ENDDO            ENDDO
353           ENDDO  #ifdef ALLOW_AUTODIFF_TAMC
354            ELSE
355              DO k=1, Nr
356                CALL GMREDI_CALC_TENSOR_DUMMY(
357         I             bi, bj, iMin, iMax, jMin, jMax, k,
358         I             sigmaX, sigmaY, sigmaR,
359         I             myThid )
360              ENDDO
361    #endif /* ALLOW_AUTODIFF_TAMC */
362            ENDIF
363    #endif  /* ALLOW_GMREDI */
364    
365    #ifdef  ALLOW_KPP
366    C--     Compute KPP mixing coefficients
367            IF (useKPP) THEN
368              CALL KPP_CALC(
369         I                  bi, bj, myTime, myThid )
370            ENDIF
371    #endif  /* ALLOW_KPP */
372    
373          ENDDO ! K  #ifdef ALLOW_AUTODIFF_TAMC
374    CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
375    CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
376    CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
377    CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
378    CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
379    CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
380    #endif /* ALLOW_AUTODIFF_TAMC */
381    
382    #ifdef ALLOW_AIM
383    C       AIM - atmospheric intermediate model, physics package code.
384    C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
385            IF ( useAIM ) THEN
386             CALL TIMER_START('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)
387             CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
388             CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)
389            ENDIF
390    #endif /* ALLOW_AIM */
391    
392          DO K = Nz, 1, -1  
393           kM1  =max(1,k-1)   ! Points to level above k (=k-1)  C--     Start of thermodynamics loop
394           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above          DO k=Nr,1,-1
395           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer  
396           iMin = 1-OLx+2  C--       km1    Points to level above k (=k-1)
397           iMax = sNx+OLx-1  C--       kup    Cycles through 1,2 to point to layer above
398           jMin = 1-OLy+2  C--       kDown  Cycles through 2,1 to point to current layer
399           jMax = sNy+OLy-1  
400              km1  = MAX(1,k-1)
401              kup  = 1+MOD(k+1,2)
402              kDown= 1+MOD(k,2)
403    
404              iMin = 1-OLx+2
405              iMax = sNx+OLx-1
406              jMin = 1-OLy+2
407              jMax = sNy+OLy-1
408    
409    #ifdef ALLOW_AUTODIFF_TAMC
410    CPatrick Is this formula correct?
411             kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
412    CADJ STORE rvel  (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
413    CADJ STORE rTrans(:,:)       = comlev1_bibj_k, key = kkey, byte = isbyte
414    CADJ STORE KappaRT(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
415    CADJ STORE KappaRS(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
416    #endif /* ALLOW_AUTODIFF_TAMC */
417    
418  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
419           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
420       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
421       O        xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
422       I        myThid)       I        myThid)
423    
424    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
425  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
426           CALL CALC_DIFFUSIVITY(           CALL CALC_DIFFUSIVITY(
427       I        bi,bj,iMin,iMax,jMin,jMax,K,       I        bi,bj,iMin,iMax,jMin,jMax,k,
428       I        maskC,maskUp,KapGM,K33,       I        maskC,maskup,
429       O        KappaZT,KappaZS,       O        KappaRT,KappaRS,KappaRU,KappaRV,
430       I        myThid)       I        myThid)
431    #endif
432    
433  C--      Calculate accelerations in the momentum equations  C--      Calculate active tracer tendencies (gT,gS,...)
434           IF ( momStepping ) THEN  C        and step forward storing result in gTnm1, gSnm1, etc.
           CALL CALC_MOM_RHS(  
      I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
      I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,  
      I         phiHyd,  
      U         aTerm,xTerm,cTerm,mTerm,pTerm,  
      U         fZon, fMer, fVerU, fVerV,  
      I         myThid)  
          ENDIF  
   
 C--      Calculate active tracer tendencies  
435           IF ( tempStepping ) THEN           IF ( tempStepping ) THEN
436            CALL CALC_GT(             CALL CALC_GT(
437       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
438       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
439       I         K13,K23,KappaZT,KapGM,       I         KappaRT,
440       U         aTerm,xTerm,fZon,fMer,fVerT,       U         fVerT,
441       I         myThid)       I         myTime, myThid)
442               CALL TIMESTEP_TRACER(
443         I         bi,bj,iMin,iMax,jMin,jMax,k,
444         I         theta, gT,
445         U         gTnm1,
446         I         myIter, myThid)
447           ENDIF           ENDIF
448           IF ( saltStepping ) THEN           IF ( saltStepping ) THEN
449            CALL CALC_GS(             CALL CALC_GS(
450       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
451       I         xA,yA,uTrans,vTrans,wTrans,maskUp,maskC,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
452       I         K13,K23,KappaZS,KapGM,       I         KappaRS,
453       U         aTerm,xTerm,fZon,fMer,fVerS,       U         fVerS,
454       I         myThid)       I         myTime, myThid)
455               CALL TIMESTEP_TRACER(
456         I         bi,bj,iMin,iMax,jMin,jMax,k,
457         I         salt, gS,
458         U         gSnm1,
459         I         myIter, myThid)
460           ENDIF           ENDIF
461    
462  C--      Prediction step (step forward all model variables)  #ifdef   ALLOW_OBCS
463           CALL TIMESTEP(  C--      Apply open boundary conditions
464       I       bi,bj,iMin,iMax,jMin,jMax,K,           IF (useOBCS) THEN
465       I       myThid)             CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
466             END IF
467  C--      Diagnose barotropic divergence of predicted fields  #endif   /* ALLOW_OBCS */
468           CALL DIV_G(  
469       I       bi,bj,iMin,iMax,jMin,jMax,K,  C--      Freeze water
470       I       xA,yA,           IF (allowFreezing) THEN
471       I       myThid)  #ifdef ALLOW_AUTODIFF_TAMC
472    CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
473  C--      Cumulative diagnostic calculations (ie. time-averaging)  CADJ &   , key = kkey, byte = isbyte
474  #ifdef ALLOW_DIAGNOSTICS  #endif /* ALLOW_AUTODIFF_TAMC */
475           IF (taveFreq.GT.0.) THEN              CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
476            CALL DO_TIME_AVERAGES(           END IF
      I                           myTime, myIter, bi, bj, K, kUp, kDown,  
      I                           K13, K23, wVel, KapGM,  
      I                           myThid )  
          ENDIF  
 #endif  
477    
478          ENDDO ! K  C--     end of thermodynamic k loop (Nr:1)
479            ENDDO
480    
481    
482    #ifdef ALLOW_AUTODIFF_TAMC
483    CPatrick? What about this one?
484               maximpl = 6
485               iikey = (ikey-1)*maximpl
486    #endif /* ALLOW_AUTODIFF_TAMC */
487    
488  C--     Implicit diffusion  C--     Implicit diffusion
489          IF (implicitDiffusion) THEN          IF (implicitDiffusion) THEN
490           CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,  
491       I                  KappaZT,KappaZS,            IF (tempStepping) THEN
492       I                  myThid )  #ifdef ALLOW_AUTODIFF_TAMC
493                idkey = iikey + 1
494    #endif /* ALLOW_AUTODIFF_TAMC */
495                CALL IMPLDIFF(
496         I         bi, bj, iMin, iMax, jMin, jMax,
497         I         deltaTtracer, KappaRT, recip_HFacC,
498         U         gTNm1,
499         I         myThid )
500             ENDIF
501    
502             IF (saltStepping) THEN
503    #ifdef ALLOW_AUTODIFF_TAMC
504             idkey = iikey + 2
505    #endif /* ALLOW_AUTODIFF_TAMC */
506                CALL IMPLDIFF(
507         I         bi, bj, iMin, iMax, jMin, jMax,
508         I         deltaTtracer, KappaRS, recip_HFacC,
509         U         gSNm1,
510         I         myThid )
511             ENDIF
512    
513    #ifdef   ALLOW_OBCS
514    C--      Apply open boundary conditions
515             IF (useOBCS) THEN
516               DO K=1,Nr
517                 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
518               ENDDO
519             END IF
520    #endif   /* ALLOW_OBCS */
521    
522    C--     End If implicitDiffusion
523            ENDIF
524    
525    
526    
527    C--     Start of dynamics loop
528            DO k=1,Nr
529    
530    C--       km1    Points to level above k (=k-1)
531    C--       kup    Cycles through 1,2 to point to layer above
532    C--       kDown  Cycles through 2,1 to point to current layer
533    
534              km1  = MAX(1,k-1)
535              kup  = 1+MOD(k+1,2)
536              kDown= 1+MOD(k,2)
537    
538              iMin = 1-OLx+2
539              iMax = sNx+OLx-1
540              jMin = 1-OLy+2
541              jMax = sNy+OLy-1
542    
543    C--      Integrate hydrostatic balance for phiHyd with BC of
544    C        phiHyd(z=0)=0
545    C        distinguishe between Stagger and Non Stagger time stepping
546             IF (staggerTimeStep) THEN
547               CALL CALC_PHI_HYD(
548         I        bi,bj,iMin,iMax,jMin,jMax,k,
549         I        gTnm1, gSnm1,
550         U        phiHyd,
551         I        myThid )
552             ELSE
553               CALL CALC_PHI_HYD(
554         I        bi,bj,iMin,iMax,jMin,jMax,k,
555         I        theta, salt,
556         U        phiHyd,
557         I        myThid )
558             ENDIF
559    
560    C--      Calculate accelerations in the momentum equations (gU, gV, ...)
561    C        and step forward storing the result in gUnm1, gVnm1, etc...
562             IF ( momStepping ) THEN
563               CALL CALC_MOM_RHS(
564         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
565         I         phiHyd,KappaRU,KappaRV,
566         U         fVerU, fVerV,
567         I         myTime, myThid)
568               CALL TIMESTEP(
569         I         bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
570         I         myIter, myThid)
571    
572    #ifdef   ALLOW_OBCS
573    C--      Apply open boundary conditions
574             IF (useOBCS) THEN
575               CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
576             END IF
577    #endif   /* ALLOW_OBCS */
578    
579    #ifdef   ALLOW_AUTODIFF_TAMC
580    #ifdef   INCLUDE_CD_CODE
581             ELSE
582               DO j=1-OLy,sNy+OLy
583                 DO i=1-OLx,sNx+OLx
584                   guCD(i,j,k,bi,bj) = 0.0
585                   gvCD(i,j,k,bi,bj) = 0.0
586                 END DO
587               END DO
588    #endif   /* INCLUDE_CD_CODE */
589    #endif   /* ALLOW_AUTODIFF_TAMC */
590             ENDIF
591    
592    
593    C--     end of dynamics k loop (1:Nr)
594            ENDDO
595    
596    
597    
598    C--     Implicit viscosity
599            IF (implicitViscosity.AND.momStepping) THEN
600    #ifdef    ALLOW_AUTODIFF_TAMC
601              idkey = iikey + 3
602    #endif    /* ALLOW_AUTODIFF_TAMC */
603              CALL IMPLDIFF(
604         I         bi, bj, iMin, iMax, jMin, jMax,
605         I         deltaTmom, KappaRU,recip_HFacW,
606         U         gUNm1,
607         I         myThid )
608    #ifdef    ALLOW_AUTODIFF_TAMC
609              idkey = iikey + 4
610    #endif    /* ALLOW_AUTODIFF_TAMC */
611              CALL IMPLDIFF(
612         I         bi, bj, iMin, iMax, jMin, jMax,
613         I         deltaTmom, KappaRV,recip_HFacS,
614         U         gVNm1,
615         I         myThid )
616    
617    #ifdef   ALLOW_OBCS
618    C--      Apply open boundary conditions
619             IF (useOBCS) THEN
620               DO K=1,Nr
621                 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
622               ENDDO
623             END IF
624    #endif   /* ALLOW_OBCS */
625    
626    #ifdef    INCLUDE_CD_CODE
627    #ifdef    ALLOW_AUTODIFF_TAMC
628              idkey = iikey + 5
629    #endif    /* ALLOW_AUTODIFF_TAMC */
630              CALL IMPLDIFF(
631         I         bi, bj, iMin, iMax, jMin, jMax,
632         I         deltaTmom, KappaRU,recip_HFacW,
633         U         vVelD,
634         I         myThid )
635    #ifdef    ALLOW_AUTODIFF_TAMC
636              idkey = iikey + 6
637    #endif    /* ALLOW_AUTODIFF_TAMC */
638              CALL IMPLDIFF(
639         I         bi, bj, iMin, iMax, jMin, jMax,
640         I         deltaTmom, KappaRV,recip_HFacS,
641         U         uVelD,
642         I         myThid )
643    #endif    /* INCLUDE_CD_CODE */
644    C--     End If implicitViscosity.AND.momStepping
645          ENDIF          ENDIF
646    
647         ENDDO         ENDDO
648        ENDDO        ENDDO
649    
 C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),  
 C    &                           maxval(cg2d_x(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.),  
 C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)  
 C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),  
 C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)  
 C     write(0,*) 'dynamics: wVel(1) ',  
 C    &            minval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.),  
 C    &            maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.)  
 C     write(0,*) 'dynamics: wVel(2) ',  
 C    &            minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.),  
 C    &            maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.)  
 cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  
 cblk &                           maxval(K13(1:sNx,1:sNy,:))  
 cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  
 cblk &                           maxval(K23(1:sNx,1:sNy,:))  
 cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  
 cblk &                           maxval(K33(1:sNx,1:sNy,:))  
 C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),  
 C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))  
 C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),  
 C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))  
 C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),  
 C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))  
 C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),  
 C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))  
 C     write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.),  
 C    &                           maxval(pH/(Gravity*Rhonil))  
   
650        RETURN        RETURN
651        END        END

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

  ViewVC Help
Powered by ViewVC 1.1.22