/[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.53 by heimbach, Mon Sep 11 23:07:29 2000 UTC revision 1.62 by jmc, Wed Feb 14 22:51:27 2001 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
5    
# Line 33  C     == Global variables === Line 34  C     == Global variables ===
34  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
35  # include "tamc.h"  # include "tamc.h"
36  # include "tamc_keys.h"  # include "tamc_keys.h"
 # ifdef ALLOW_KPP  
 #  include "KPP.h"  
 # endif  
37  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
38    
39    #ifdef ALLOW_KPP
40    # include "KPP.h"
41    #endif
42    
43    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
44    #include "AVER.h"
45    #endif
46    
47  C     == Routine arguments ==  C     == Routine arguments ==
48  C     myTime - Current time in simulation  C     myTime - Current time in simulation
49  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
# Line 50  C     == Local variables Line 56  C     == Local variables
56  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
57  C     uTrans, vTrans, rTrans - Per block temporaries holding flow  C     uTrans, vTrans, rTrans - Per block temporaries holding flow
58  C                              transport  C                              transport
59  C     rVel                     o uTrans: Zonal transport  C                              o uTrans: Zonal transport
60  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
61  C                              o rTrans: Vertical transport  C                              o rTrans: Vertical transport
 C                              o rVel:   Vertical velocity at upper and  
 C                                        lower cell faces.  
62  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
63  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
64  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  
65  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
66  C                                      so we need an fVer for each  C                                      so we need an fVer for each
67  C                                      variable.  C                                      variable.
68  C     rhoK, rhoKM1   - Density at current level, level above and level  C     rhoK, rhoKM1   - Density at current level, and level above
 C                      below.  
 C     rhoKP1                                                                    
 C     buoyK, buoyKM1 - Buoyancy at current level and level above.  
69  C     phiHyd         - Hydrostatic part of the potential phiHydi.  C     phiHyd         - Hydrostatic part of the potential phiHydi.
70  C                      In z coords phiHydiHyd is the hydrostatic  C                      In z coords phiHydiHyd is the hydrostatic
71  C                      pressure anomaly  C                      pressure anomaly
# Line 95  C                      index into fVerTe Line 87  C                      index into fVerTe
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 rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
90        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
92        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
93        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
94        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
95        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
96        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
98        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
99        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
100        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
101        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
# Line 126  C                      index into fVerTe Line 104  C                      index into fVerTe
104        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
105        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
106    
107  C This is currently also used by IVDC and Diagnostics  C This is currently used by IVDC and Diagnostics
 C #ifdef INCLUDE_CONVECT_CALL  
108        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
 C #endif  
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
       LOGICAL BOTTOM_LAYER  
115    
116    Cjmc : add for phiHyd output <- but not working if multi tile per CPU
117    c     CHARACTER*(MAX_LEN_MBUF) suff
118    c     LOGICAL  DIFFERENT_MULTIPLE
119    c     EXTERNAL DIFFERENT_MULTIPLE
120    Cjmc(end)
121    
122  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
123        INTEGER    isbyte        INTEGER    isbyte
124        PARAMETER( isbyte = 4 )        PARAMETER( isbyte = 4 )
# Line 162  C       "Calculation of Gs" Line 143  C       "Calculation of Gs"
143  C       ===================  C       ===================
144  C       This is where all the accelerations and tendencies (ie.  C       This is where all the accelerations and tendencies (ie.
145  C       physics, parameterizations etc...) are calculated  C       physics, parameterizations etc...) are calculated
 C         rVel = sum_r ( div. u[n] )  
146  C         rho = rho ( theta[n], salt[n] )  C         rho = rho ( theta[n], salt[n] )
147  C         b   = b(rho, theta)  C         b   = b(rho, theta)
148  C         K31 = K31 ( rho )  C         K31 = K31 ( rho )
149  C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )  C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
150  C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )  C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
151  C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )  C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
152  C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
153  C  C
154  C       "Time-stepping" or "Prediction"  C       "Time-stepping" or "Prediction"
155  C       ================================  C       ================================
# Line 209  C     uninitialised but inert locations. Line 189  C     uninitialised but inert locations.
189          yA(i,j)      = 0. _d 0          yA(i,j)      = 0. _d 0
190          uTrans(i,j)  = 0. _d 0          uTrans(i,j)  = 0. _d 0
191          vTrans(i,j)  = 0. _d 0          vTrans(i,j)  = 0. _d 0
         aTerm(i,j)   = 0. _d 0  
         xTerm(i,j)   = 0. _d 0  
         cTerm(i,j)   = 0. _d 0  
         mTerm(i,j)   = 0. _d 0  
         pTerm(i,j)   = 0. _d 0  
         fZon(i,j)    = 0. _d 0  
         fMer(i,j)    = 0. _d 0  
192          DO k=1,Nr          DO k=1,Nr
193           phiHyd (i,j,k)  = 0. _d 0           phiHyd(i,j,k)  = 0. _d 0
194           KappaRU(i,j,k) = 0. _d 0           KappaRU(i,j,k) = 0. _d 0
195           KappaRV(i,j,k) = 0. _d 0           KappaRV(i,j,k) = 0. _d 0
196           sigmaX(i,j,k) = 0. _d 0           sigmaX(i,j,k) = 0. _d 0
# Line 226  C     uninitialised but inert locations. Line 199  C     uninitialised but inert locations.
199          ENDDO          ENDDO
200          rhoKM1 (i,j) = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
201          rhok   (i,j) = 0. _d 0          rhok   (i,j) = 0. _d 0
         rhoKP1 (i,j) = 0. _d 0  
         rhoTMP (i,j) = 0. _d 0  
         buoyKM1(i,j) = 0. _d 0  
         buoyK  (i,j) = 0. _d 0  
202          maskC  (i,j) = 0. _d 0          maskC  (i,j) = 0. _d 0
203         ENDDO         ENDDO
204        ENDDO        ENDDO
# Line 244  CHPF$ INDEPENDENT Line 213  CHPF$ INDEPENDENT
213    
214  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
215  C--    HPF directive to help TAMC  C--    HPF directive to help TAMC
216  CHPF$  INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV  CHPF$  INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV
217  CHPF$&                  ,phiHyd,utrans,vtrans,maskc,xA,yA  CHPF$&                  ,phiHyd,utrans,vtrans,maskc,xA,yA
218  CHPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV  CHPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV
219  CHPF$&                  )  CHPF$&                  )
# Line 273  C--     Set up work arrays that need val Line 242  C--     Set up work arrays that need val
242          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
243           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
244            rTrans(i,j)   = 0. _d 0            rTrans(i,j)   = 0. _d 0
           rVel  (i,j,1) = 0. _d 0  
           rVel  (i,j,2) = 0. _d 0  
245            fVerT (i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
246            fVerT (i,j,2) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
247            fVerS (i,j,1) = 0. _d 0            fVerS (i,j,1) = 0. _d 0
# Line 283  C--     Set up work arrays that need val Line 250  C--     Set up work arrays that need val
250            fVerU (i,j,2) = 0. _d 0            fVerU (i,j,2) = 0. _d 0
251            fVerV (i,j,1) = 0. _d 0            fVerV (i,j,1) = 0. _d 0
252            fVerV (i,j,2) = 0. _d 0            fVerV (i,j,2) = 0. _d 0
           phiHyd(i,j,1) = 0. _d 0  
253           ENDDO           ENDDO
254          ENDDO          ENDDO
255    
256          DO k=1,Nr          DO k=1,Nr
257           DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
258            DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
259  #ifdef INCLUDE_CONVECT_CALL  C This is currently also used by IVDC and Diagnostics
260             ConvectCount(i,j,k) = 0.             ConvectCount(i,j,k) = 0.
 #endif  
261             KappaRT(i,j,k) = 0. _d 0             KappaRT(i,j,k) = 0. _d 0
262             KappaRS(i,j,k) = 0. _d 0             KappaRS(i,j,k) = 0. _d 0
263            ENDDO            ENDDO
# Line 305  C--     Set up work arrays that need val Line 270  C--     Set up work arrays that need val
270          jMax = sNy+OLy          jMax = sNy+OLy
271    
272    
273          k = 1  C--     Start of diagnostic loop
274          BOTTOM_LAYER = k .EQ. Nr          DO k=Nr,1,-1
   
 #ifdef DO_PIPELINED_CORRECTION_STEP  
 C--     Calculate gradient of surface pressure  
         CALL CALC_GRAD_ETA_SURF(  
      I       bi,bj,iMin,iMax,jMin,jMax,  
      O       etaSurfX,etaSurfY,  
      I       myThid)  
 C--     Update fields in top level according to tendency terms  
         CALL CORRECTION_STEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,k,  
      I       etaSurfX,etaSurfY,myTime,myThid)  
   
 #ifdef ALLOW_OBCS  
         IF (openBoundaries) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
            CALL APPLY_OBCS1( bi, bj, k, myThid )  
         END IF  
 #endif  
   
         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,  
      I        etaSurfX,etaSurfY,myTime,myThid)  
 #ifdef ALLOW_OBCS  
          IF (openBoundaries) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
             CALL APPLY_OBCS1( bi, bj, k+1, myThid )  
          END IF  
 #endif  
         ENDIF  
 #endif  
 C--     Density of 1st level (below W(1)) reference to level 1  
 #ifdef  INCLUDE_FIND_RHO_CALL  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,  
      O     rhoKm1,  
      I     myThid )  
 #endif  
   
         IF (.NOT. BOTTOM_LAYER) THEN  
   
 C--      Check static stability with layer below  
 C--      and mix as needed.  
 #ifdef  INCLUDE_FIND_RHO_CALL  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj  
 CADJ &   , key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj  
 CADJ &   , key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax, k+1, k, eosType,  
      O      rhoKp1,  
      I      myThid )  
 #endif  
   
 #ifdef  INCLUDE_CONVECT_CALL  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE rhoKm1(:,:)  = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE rhoKp1(:,:)  = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
          CALL CONVECT(  
      I       bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,  
      U       ConvectCount,  
      I       myTime,myIter,myThid)  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)  
 CADJ &     = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)  
 CADJ &     = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 #endif  
   
 C--      Implicit Vertical Diffusion for Convection  
          IF (ivdc_kappa.NE.0.) THEN  
             CALL CALC_IVDC(  
      I       bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,  
      U       ConvectCount, KappaRT, KappaRS,  
      I       myTime,myIter,myThid)  
          ENDIF  
   
 C--      Recompute density after mixing  
 #ifdef  INCLUDE_FIND_RHO_CALL  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,  
      O      rhoKm1,  
      I      myThid )  
 #endif  
         ENDIF  
   
 C--     Calculate buoyancy  
         CALL CALC_BUOYANCY(  
      I      bi,bj,iMin,iMax,jMin,jMax,k,rhoKm1,  
      O      buoyKm1,  
      I      myThid )  
   
 C--     Integrate hydrostatic balance for phiHyd with BC of  
 C--     phiHyd(z=0)=0  
         CALL CALC_PHI_HYD(  
      I      bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyKm1,  
      U      phiHyd,  
      I      myThid )  
   
 #ifdef ALLOW_GMREDI  
         CALL GRAD_SIGMA(  
      I            bi, bj, iMin, iMax, jMin, jMax, k,  
      I            rhoKm1, rhoKm1, rhoKm1,  
      O            sigmaX, sigmaY, sigmaR,  
      I            myThid )  
 #endif  
   
 C--     Start of downward loop  
         DO k=2,Nr  
275    
276  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
277    C? Patrick, is this formula correct now that we change the loop range?
278    C? Do we still need this?
279           kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1           kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
280  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
281    
282           BOTTOM_LAYER = k .EQ. Nr  C--       Integrate continuity vertically for vertical velocity
283              CALL INTEGRATE_FOR_W(
284  #ifdef DO_PIPELINED_CORRECTION_STEP       I                         bi, bj, k, uVel, vVel,
285           IF ( .NOT. BOTTOM_LAYER ) THEN       O                         wVel,
286  C--       Update fields in layer below according to tendency terms       I                         myThid )
287            CALL CORRECTION_STEP(  
288       I         bi,bj,iMin,iMax,jMin,jMax,k+1,  #ifdef    ALLOW_OBCS
289       I         etaSurfX,etaSurfY,myTime,myThid)  #ifdef    ALLOW_NONHYDROSTATIC
290  #ifdef ALLOW_OBCS  C--       Apply OBC to W if in N-H mode
291            IF (openBoundaries) THEN            IF (useOBCS.AND.nonHydrostatic) THEN
292  #ifdef ALLOW_AUTODIFF_TAMC              CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
293  CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj_k            ENDIF
294  CADJ &   , key = kkey, byte = isbyte  #endif    /* ALLOW_NONHYDROSTATIC */
295  CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj_k  #endif    /* ALLOW_OBCS */
296  CADJ &   , key = kkey, byte = isbyte  
297  CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k  C--       Calculate gradients of potential density for isoneutral
298  CADJ &   , key = kkey, byte = isbyte  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
299  CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k  c         IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
300  CADJ &   , key = kkey, byte = isbyte            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
301  #endif /* ALLOW_AUTODIFF_TAMC */              CALL FIND_RHO(
302               CALL APPLY_OBCS1( bi, bj, k+1, myThid )       I        bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
303            END IF       I        theta, salt,
304  #endif       O        rhoK,
          ENDIF  
 #endif /* DO_PIPELINED_CORRECTION_STEP */  
   
 C--      Density of k level (below W(k)) reference to k level  
 #ifdef  INCLUDE_FIND_RHO_CALL  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  k, k, eosType,  
      O      rhoK,  
      I      myThid )  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE rhoK(:,:)  = comlev1_bibj_k, key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
 #endif  
   
          IF (.NOT. BOTTOM_LAYER) THEN  
   
 C--       Check static stability with layer below and mix as needed.  
 C--       Density of k+1 level (below W(k+1)) reference to k level.  
 #ifdef  INCLUDE_FIND_RHO_CALL  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax,  k+1, k, eosType,  
      O       rhoKp1,  
      I       myThid )  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE rhoKp1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
 #endif  
   
 #ifdef  INCLUDE_CONVECT_CALL  
           CALL CONVECT(  
      I        bi,bj,iMin,iMax,jMin,jMax,k+1,rhoK,rhoKp1,  
      U        ConvectCount,  
      I        myTime,myIter,myThid)  
   
 #endif  
   
 C--      Implicit Vertical Diffusion for Convection  
          IF (ivdc_kappa.NE.0.) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE rhoKm1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
             CALL CALC_IVDC(  
      I       bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,  
      U       ConvectCount, KappaRT, KappaRS,  
      I       myTime,myIter,myThid)  
          END IF  
   
 C--       Recompute density after mixing  
 #ifdef  INCLUDE_FIND_RHO_CALL  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
           CALL FIND_RHO(  
      I       bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,  
      O       rhoK,  
      I       myThid )  
 #endif  
   
 C--            IF (.NOT. BOTTOM_LAYER) ends here  
          ENDIF  
   
 C--      Calculate buoyancy  
          CALL CALC_BUOYANCY(  
      I       bi,bj,iMin,iMax,jMin,jMax,k,rhoK,  
      O       buoyK,  
      I       myThid )  
   
 C--      Integrate hydrostatic balance for phiHyd with BC of  
 C--      phiHyd(z=0)=0  
          CALL CALC_PHI_HYD(  
      I        bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyK,  
      U        phiHyd,  
305       I        myThid )       I        myThid )
306                IF (k.GT.1) CALL FIND_RHO(
 #ifdef  INCLUDE_FIND_RHO_CALL  
 C--      Calculate iso-neutral slopes for the GM/Redi parameterisation  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
          CALL FIND_RHO(  
307       I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,       I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
308       O        rhoTmp,       I        theta, salt,
309         O        rhoKm1,
310       I        myThid )       I        myThid )
311  #endif              CALL GRAD_SIGMA(
   
   
 #ifdef ALLOW_GMREDI  
          CALL GRAD_SIGMA(  
312       I             bi, bj, iMin, iMax, jMin, jMax, k,       I             bi, bj, iMin, iMax, jMin, jMax, k,
313       I             rhoK, rhotmp, rhoK,       I             rhoK, rhoKm1, rhoK,
314       O             sigmaX, sigmaY, sigmaR,       O             sigmaX, sigmaY, sigmaR,
315       I             myThid )       I             myThid )
316  #endif            ENDIF
317    
318           DO J=jMin,jMax  C--       Implicit Vertical Diffusion for Convection
319            DO I=iMin,iMax  c ==> should use sigmaR !!!
320  #ifdef  INCLUDE_FIND_RHO_CALL            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
321             rhoKm1 (I,J) = rhoK(I,J)              CALL CALC_IVDC(
322  #endif       I        bi, bj, iMin, iMax, jMin, jMax, k,
323             buoyKm1(I,J) = buoyK(I,J)       I        rhoKm1, rhoK,
324            ENDDO       U        ConvectCount, KappaRT, KappaRS,
325           ENDDO       I        myTime, myIter, myThid)
326              ENDIF
327    
328  C--     end of k loop  C--     end of diagnostic k loop (Nr:1)
329          ENDDO          ENDDO
330    
331  #ifdef ALLOW_GMREDI  #ifdef  ALLOW_OBCS
332    C--     Calculate future values on open boundaries
333            IF (useOBCS) THEN
334              CALL OBCS_CALC( bi, bj, myTime+deltaT,
335         I            uVel, vVel, wVel, theta, salt,
336         I            myThid )
337            ENDIF
338    #endif  /* ALLOW_OBCS */
339    
340    C--     Determines forcing terms based on external fields
341    C       relaxation terms, etc.
342            CALL EXTERNAL_FORCING_SURF(
343         I             bi, bj, iMin, iMax, jMin, jMax,
344         I             myThid )
345    
346    #ifdef  ALLOW_GMREDI
347    C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
348          IF (useGMRedi) THEN          IF (useGMRedi) THEN
349            DO k=1, Nr            DO k=1,Nr
350              CALL GMREDI_CALC_TENSOR(              CALL GMREDI_CALC_TENSOR(
351       I             bi, bj, iMin, iMax, jMin, jMax, k,       I             bi, bj, iMin, iMax, jMin, jMax, k,
352       I             sigmaX, sigmaY, sigmaR,       I             sigmaX, sigmaY, sigmaR,
353       I             myThid )       I             myThid )
354            ENDDO            ENDDO
355    #ifdef ALLOW_AUTODIFF_TAMC
356            ELSE
357              DO k=1, Nr
358                CALL GMREDI_CALC_TENSOR_DUMMY(
359         I             bi, bj, iMin, iMax, jMin, jMax, k,
360         I             sigmaX, sigmaY, sigmaR,
361         I             myThid )
362              ENDDO
363    #endif /* ALLOW_AUTODIFF_TAMC */
364          ENDIF          ENDIF
365  #endif  #endif  /* ALLOW_GMREDI */
366    
367    #ifdef  ALLOW_KPP
368    C--     Compute KPP mixing coefficients
369            IF (useKPP) THEN
370              CALL KPP_CALC(
371         I                  bi, bj, myTime, myThid )
372            ENDIF
373    #endif  /* ALLOW_KPP */
374    
375  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
376  CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
377  CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte
   
378  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
379  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
380  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
381  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
   
 C--     dummy initialization to break data flow because  
 C--     calc_div_ghat has a condition for initialization  
         DO J=jMin,jMax  
            DO I=iMin,iMax  
               cg2d_b(i,j,bi,bj) = 0.0  
            ENDDO  
         ENDDO  
382  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
383    
384    #ifdef ALLOW_AIM
385    C       AIM - atmospheric intermediate model, physics package code.
386    C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
387            IF ( useAIM ) THEN
388             CALL TIMER_START('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)
389             CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
390             CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)
391            ENDIF
392    #endif /* ALLOW_AIM */
393    
 #ifdef ALLOW_KPP  
 C--   Compute KPP mixing coefficients  
         IF (useKPP) THEN  
394    
395            CALL TIMER_START('KPP_CALC               [DYNAMICS]', myThid)  C--     Start of thermodynamics loop
396            CALL KPP_CALC(          DO k=Nr,1,-1
      I                 bi, bj, myTime, myThid )  
           CALL TIMER_STOP ('KPP_CALC               [DYNAMICS]', myThid)  
397    
398  #ifdef ALLOW_AUTODIFF_TAMC  C--       km1    Points to level above k (=k-1)
399          ELSE  C--       kup    Cycles through 1,2 to point to layer above
400            DO j=1-OLy,sNy+OLy  C--       kDown  Cycles through 2,1 to point to current layer
             DO i=1-OLx,sNx+OLx  
               KPPhbl (i,j,bi,bj) = 1.0  
               KPPfrac(i,j,bi,bj) = 0.0  
               DO k = 1,Nr  
                  KPPghat   (i,j,k,bi,bj) = 0.0  
                  KPPviscAz (i,j,k,bi,bj) = viscAz  
                  KPPdiffKzT(i,j,k,bi,bj) = diffKzT  
                  KPPdiffKzS(i,j,k,bi,bj) = diffKzS  
               ENDDO  
             ENDDO  
           ENDDO  
 #endif /* ALLOW_AUTODIFF_TAMC */  
         ENDIF  
401    
402  #ifdef ALLOW_AUTODIFF_TAMC            km1  = MAX(1,k-1)
403  CADJ STORE KPPghat   (:,:,:,bi,bj)            kup  = 1+MOD(k+1,2)
404  CADJ &   , KPPviscAz (:,:,:,bi,bj)            kDown= 1+MOD(k,2)
405  CADJ &   , KPPdiffKzT(:,:,:,bi,bj)  
406  CADJ &   , KPPdiffKzS(:,:,:,bi,bj)            iMin = 1-OLx+2
407  CADJ &   , KPPfrac   (:,:  ,bi,bj)            iMax = sNx+OLx-1
408  CADJ &                 = comlev1_bibj, key = ikey, byte = isbyte            jMin = 1-OLy+2
409  #endif /* ALLOW_AUTODIFF_TAMC */            jMax = sNy+OLy-1
   
 #endif /* ALLOW_KPP */  
   
 C--     Start of upward loop  
         DO k = Nr, 1, -1  
   
 C--      km1    Points to level above k (=k-1)  
 C--      kup    Cycles through 1,2 to point to layer above  
 C--      kDown  Cycles through 2,1 to point to current layer  
   
          km1  =max(1,k-1)  
          kup  =1+MOD(k+1,2)  
          kDown=1+MOD(k,2)  
   
          iMin = 1-OLx+2  
          iMax = sNx+OLx-1  
          jMin = 1-OLy+2  
          jMax = sNy+OLy-1  
410    
411  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
412    CPatrick Is this formula correct?
413           kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1           kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
   
 CADJ STORE rvel  (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte  
414  CADJ STORE rTrans(:,:)       = comlev1_bibj_k, key = kkey, byte = isbyte  CADJ STORE rTrans(:,:)       = comlev1_bibj_k, key = kkey, byte = isbyte
415  CADJ STORE KappaRT(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte  CADJ STORE KappaRT(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
416  CADJ STORE KappaRS(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte  CADJ STORE KappaRS(:,:,:)    = comlev1_bibj_k, key = kkey, byte = isbyte
# Line 688  CADJ STORE KappaRS(:,:,:)    = comlev1_b Line 419  CADJ STORE KappaRS(:,:,:)    = comlev1_b
419  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
420           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
421       I        bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
422       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,maskC,maskUp,
423       I        myThid)       I        myThid)
424    
 #ifdef ALLOW_OBCS  
         IF (openBoundaries) THEN  
          CALL APPLY_OBCS3( bi, bj, k, kup, rTrans, rVel, myThid )  
         ENDIF  
 #endif  
   
425  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
426  C--      Calculate the total vertical diffusivity  C--      Calculate the total vertical diffusivity
427           CALL CALC_DIFFUSIVITY(           CALL CALC_DIFFUSIVITY(
428       I        bi,bj,iMin,iMax,jMin,jMax,k,       I        bi,bj,iMin,iMax,jMin,jMax,k,
429       I        maskC,maskUp,       I        maskC,maskup,
430       O        KappaRT,KappaRS,KappaRU,KappaRV,       O        KappaRT,KappaRS,KappaRU,KappaRV,
431       I        myThid)       I        myThid)
432  #endif  #endif
433  C--      Calculate accelerations in the momentum equations  
434           IF ( momStepping ) THEN  C--      Calculate active tracer tendencies (gT,gS,...)
435            CALL CALC_MOM_RHS(  C        and step forward storing result in gTnm1, gSnm1, etc.
      I         bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,  
      I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,  
      I         phiHyd,KappaRU,KappaRV,  
      U         aTerm,xTerm,cTerm,mTerm,pTerm,  
      U         fZon, fMer, fVerU, fVerV,  
      I         myTime, myThid)  
 #ifdef ALLOW_AUTODIFF_TAMC  
 #ifdef INCLUDE_CD_CODE  
          ELSE  
             DO j=1-OLy,sNy+OLy  
                DO i=1-OLx,sNx+OLx  
                   guCD(i,j,k,bi,bj) = 0.0  
                   gvCD(i,j,k,bi,bj) = 0.0  
                END DO  
             END DO  
 #endif  
 #endif /* ALLOW_AUTODIFF_TAMC */  
          ENDIF  
 C--      Calculate active tracer tendencies  
436           IF ( tempStepping ) THEN           IF ( tempStepping ) THEN
437            CALL CALC_GT(             CALL CALC_GT(
438       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
439       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
440       I         KappaRT,       I         KappaRT,
441       U         aTerm,xTerm,fZon,fMer,fVerT,       U         fVerT,
442       I         myTime, myThid)       I         myTime, myThid)
443               CALL TIMESTEP_TRACER(
444         I         bi,bj,iMin,iMax,jMin,jMax,k,
445         I         theta, gT,
446         U         gTnm1,
447         I         myIter, myThid)
448           ENDIF           ENDIF
449           IF ( saltStepping ) THEN           IF ( saltStepping ) THEN
450            CALL CALC_GS(             CALL CALC_GS(
451       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
452       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
453       I         KappaRS,       I         KappaRS,
454       U         aTerm,xTerm,fZon,fMer,fVerS,       U         fVerS,
455       I         myTime, myThid)       I         myTime, myThid)
456               CALL TIMESTEP_TRACER(
457         I         bi,bj,iMin,iMax,jMin,jMax,k,
458         I         salt, gS,
459         U         gSnm1,
460         I         myIter, myThid)
461           ENDIF           ENDIF
462  #ifdef ALLOW_OBCS  
463  C--      Calculate future values on open boundaries  #ifdef   ALLOW_OBCS
          IF (openBoundaries) THEN  
 Caja      CALL CYCLE_OBCS( k, bi, bj, myThid )  
           CALL SET_OBCS( k, bi, bj, myTime+deltaTclock, myThid )  
          ENDIF  
 #endif  
 C--      Prediction step (step forward all model variables)  
          CALL TIMESTEP(  
      I       bi,bj,iMin,iMax,jMin,jMax,k,  
      I       myIter, myThid)  
 #ifdef ALLOW_OBCS  
464  C--      Apply open boundary conditions  C--      Apply open boundary conditions
465           IF (openBoundaries) THEN           IF (useOBCS) THEN
466  #ifdef ALLOW_AUTODIFF_TAMC             CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
 CADJ STORE gunm1(:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE gvnm1(:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 CADJ STORE gwnm1(:,:,k,bi,bj) = comlev1_bibj_k, key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
 CADJ &    
             CALL APPLY_OBCS2( bi, bj, k, myThid )  
467           END IF           END IF
468  #endif  #endif   /* ALLOW_OBCS */
469    
470  C--      Freeze water  C--      Freeze water
471           IF (allowFreezing) THEN           IF (allowFreezing) THEN
472  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 777  CADJ &   , key = kkey, byte = isbyte Line 476  CADJ &   , key = kkey, byte = isbyte
476              CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )              CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
477           END IF           END IF
478    
479  #ifdef DIVG_IN_DYNAMICS  C--     end of thermodynamic k loop (Nr:1)
 C--      Diagnose barotropic divergence of predicted fields  
          CALL CALC_DIV_GHAT(  
      I       bi,bj,iMin,iMax,jMin,jMax,k,  
      I       xA,yA,  
      I       myThid)  
 #endif /* DIVG_IN_DYNAMICS */  
   
 C--      Cumulative diagnostic calculations (ie. time-averaging)  
 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE  
          IF (taveFreq.GT.0.) THEN  
           CALL DO_TIME_AVERAGES(  
      I                           myTime, myIter, bi, bj, k, kup, kDown,  
      I                           rVel, ConvectCount,  
      I                           myThid )  
          ENDIF  
 #endif  
   
   
 C--     k loop  
480          ENDDO          ENDDO
481    
482    
483  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
484    CPatrick? What about this one?
485             maximpl = 6             maximpl = 6
486             iikey = (ikey-1)*maximpl             iikey = (ikey-1)*maximpl
487  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
# Line 813  C--     Implicit diffusion Line 495  C--     Implicit diffusion
495  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
496              CALL IMPLDIFF(              CALL IMPLDIFF(
497       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
498       I         deltaTtracer, KappaRT,recip_HFacC,       I         deltaTtracer, KappaRT, recip_HFacC,
499       U         gTNm1,       U         gTNm1,
500       I         myThid )       I         myThid )
501           END IF           ENDIF
502    
503           IF (saltStepping) THEN           IF (saltStepping) THEN
504  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 824  C--     Implicit diffusion Line 506  C--     Implicit diffusion
506  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
507              CALL IMPLDIFF(              CALL IMPLDIFF(
508       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
509       I         deltaTtracer, KappaRS,recip_HFacC,       I         deltaTtracer, KappaRS, recip_HFacC,
510       U         gSNm1,       U         gSNm1,
511       I         myThid )       I         myThid )
512             ENDIF
513    
514    #ifdef   ALLOW_OBCS
515    C--      Apply open boundary conditions
516             IF (useOBCS) THEN
517               DO K=1,Nr
518                 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
519               ENDDO
520           END IF           END IF
521    #endif   /* ALLOW_OBCS */
522    
523  C--     implicitDiffusion  C--     End If implicitDiffusion
524          ENDIF          ENDIF
525    
 C--     Implicit viscosity  
         IF (implicitViscosity) THEN  
526    
527           IF (momStepping) THEN  
528  #ifdef ALLOW_AUTODIFF_TAMC  C--     Start of dynamics loop
529           idkey = iikey + 3          DO k=1,Nr
530  #endif /* ALLOW_AUTODIFF_TAMC */  
531    C--       km1    Points to level above k (=k-1)
532    C--       kup    Cycles through 1,2 to point to layer above
533    C--       kDown  Cycles through 2,1 to point to current layer
534    
535              km1  = MAX(1,k-1)
536              kup  = 1+MOD(k+1,2)
537              kDown= 1+MOD(k,2)
538    
539              iMin = 1-OLx+2
540              iMax = sNx+OLx-1
541              jMin = 1-OLy+2
542              jMax = sNy+OLy-1
543    
544    C--      Integrate hydrostatic balance for phiHyd with BC of
545    C        phiHyd(z=0)=0
546    C        distinguishe between Stagger and Non Stagger time stepping
547             IF (staggerTimeStep) THEN
548               CALL CALC_PHI_HYD(
549         I        bi,bj,iMin,iMax,jMin,jMax,k,
550         I        gTnm1, gSnm1,
551         U        phiHyd,
552         I        myThid )
553             ELSE
554               CALL CALC_PHI_HYD(
555         I        bi,bj,iMin,iMax,jMin,jMax,k,
556         I        theta, salt,
557         U        phiHyd,
558         I        myThid )
559             ENDIF
560    
561    C--      Calculate accelerations in the momentum equations (gU, gV, ...)
562    C        and step forward storing the result in gUnm1, gVnm1, etc...
563             IF ( momStepping ) THEN
564               CALL CALC_MOM_RHS(
565         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
566         I         phiHyd,KappaRU,KappaRV,
567         U         fVerU, fVerV,
568         I         myTime, myThid)
569               CALL TIMESTEP(
570         I         bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
571         I         myIter, myThid)
572    
573    #ifdef   ALLOW_OBCS
574    C--      Apply open boundary conditions
575             IF (useOBCS) THEN
576               CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
577             END IF
578    #endif   /* ALLOW_OBCS */
579    
580    #ifdef   ALLOW_AUTODIFF_TAMC
581    #ifdef   INCLUDE_CD_CODE
582             ELSE
583               DO j=1-OLy,sNy+OLy
584                 DO i=1-OLx,sNx+OLx
585                   guCD(i,j,k,bi,bj) = 0.0
586                   gvCD(i,j,k,bi,bj) = 0.0
587                 END DO
588               END DO
589    #endif   /* INCLUDE_CD_CODE */
590    #endif   /* ALLOW_AUTODIFF_TAMC */
591             ENDIF
592    
593    
594    C--     end of dynamics k loop (1:Nr)
595            ENDDO
596    
597    
598    
599    C--     Implicit viscosity
600            IF (implicitViscosity.AND.momStepping) THEN
601    #ifdef    ALLOW_AUTODIFF_TAMC
602              idkey = iikey + 3
603    #endif    /* ALLOW_AUTODIFF_TAMC */
604            CALL IMPLDIFF(            CALL IMPLDIFF(
605       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
606       I         deltaTmom, KappaRU,recip_HFacW,       I         deltaTmom, KappaRU,recip_HFacW,
607       U         gUNm1,       U         gUNm1,
608       I         myThid )       I         myThid )
609  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
610           idkey = iikey + 4            idkey = iikey + 4
611  #endif /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
612            CALL IMPLDIFF(            CALL IMPLDIFF(
613       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
614       I         deltaTmom, KappaRV,recip_HFacS,       I         deltaTmom, KappaRV,recip_HFacS,
615       U         gVNm1,       U         gVNm1,
616       I         myThid )       I         myThid )
617    
618  #ifdef INCLUDE_CD_CODE  #ifdef   ALLOW_OBCS
619    C--      Apply open boundary conditions
620             IF (useOBCS) THEN
621               DO K=1,Nr
622                 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
623               ENDDO
624             END IF
625    #endif   /* ALLOW_OBCS */
626    
627  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef    INCLUDE_CD_CODE
628           idkey = iikey + 5  #ifdef    ALLOW_AUTODIFF_TAMC
629  #endif /* ALLOW_AUTODIFF_TAMC */            idkey = iikey + 5
630    #endif    /* ALLOW_AUTODIFF_TAMC */
631            CALL IMPLDIFF(            CALL IMPLDIFF(
632       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
633       I         deltaTmom, KappaRU,recip_HFacW,       I         deltaTmom, KappaRU,recip_HFacW,
634       U         vVelD,       U         vVelD,
635       I         myThid )       I         myThid )
636  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
637          idkey = iikey + 6            idkey = iikey + 6
638  #endif /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
639            CALL IMPLDIFF(            CALL IMPLDIFF(
640       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
641       I         deltaTmom, KappaRV,recip_HFacS,       I         deltaTmom, KappaRV,recip_HFacS,
642       U         uVelD,       U         uVelD,
643       I         myThid )       I         myThid )
644    #endif    /* INCLUDE_CD_CODE */
645  #endif  C--     End If implicitViscosity.AND.momStepping
   
 C--      momStepping  
          ENDIF  
   
 C--     implicitViscosity  
646          ENDIF          ENDIF
647    
648    Cjmc : add for phiHyd output <- but not working if multi tile per CPU
649    c       IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
650    c    &  .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
651    c         WRITE(suff,'(I10.10)') myIter+1
652    c         CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
653    c       ENDIF
654    Cjmc(end)
655    
656    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
657            IF (taveFreq.GT.0.) THEN
658             DO K=1,Nr
659              CALL TIMEAVER_1FLD_XYZ(phiHyd, phiHydtave,
660         I                              deltaTclock, bi, bj, K, myThid)
661              IF (ivdc_kappa.NE.0.) THEN
662                CALL TIMEAVER_1FLD_XYZ(ConvectCount, ConvectCountTave,
663         I                              deltaTclock, bi, bj, K, myThid)
664              ENDIF
665             ENDDO
666            ENDIF
667    #endif /* INCLUDE_DIAGNOSTICS_INTERFACE_CODE */
668    
669         ENDDO         ENDDO
670        ENDDO        ENDDO
671    
   
672        RETURN        RETURN
673        END        END

Legend:
Removed from v.1.53  
changed lines
  Added in v.1.62

  ViewVC Help
Powered by ViewVC 1.1.22