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

Legend:
Removed from v.1.21  
changed lines
  Added in v.1.54.2.7

  ViewVC Help
Powered by ViewVC 1.1.22