/[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.4 by adcroft, Thu Apr 30 14:03:28 1998 UTC revision 1.80 by adcroft, Fri Aug 17 18:40:30 2001 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CPP_EEOPTIONS.h"  #include "CPP_OPTIONS.h"
5    
6        SUBROUTINE DYNAMICS(myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
7  C     /==========================================================\  C     /==========================================================\
8  C     | SUBROUTINE DYNAMICS                                      |  C     | SUBROUTINE DYNAMICS                                      |
9  C     | o Controlling routine for the explicit part of the model |  C     | o Controlling routine for the explicit part of the model |
# Line 20  C     | ===== Line 21  C     | =====
21  C     | C*P* comments indicating place holders for which code is |  C     | C*P* comments indicating place holders for which code is |
22  C     |      presently being developed.                          |  C     |      presently being developed.                          |
23  C     \==========================================================/  C     \==========================================================/
24          IMPLICIT NONE
25    
26  C     == Global variables ===  C     == Global variables ===
27  #include "SIZE.h"  #include "SIZE.h"
28  #include "EEPARAMS.h"  #include "EEPARAMS.h"
29  #include "CG2D.h"  #include "PARAMS.h"
30  #include "DYNVARS.h"  #include "DYNVARS.h"
31    #include "GRID.h"
32    #ifdef ALLOW_PASSIVE_TRACER
33    #include "TR1.h"
34    #endif
35    
36    #ifdef ALLOW_AUTODIFF_TAMC
37    # include "tamc.h"
38    # include "tamc_keys.h"
39    # include "FFIELDS.h"
40    # ifdef ALLOW_KPP
41    #  include "KPP.h"
42    # endif
43    # ifdef ALLOW_GMREDI
44    #  include "GMREDI.h"
45    # endif
46    #endif /* ALLOW_AUTODIFF_TAMC */
47    
48    #ifdef ALLOW_TIMEAVE
49    #include "TIMEAVE_STATV.h"
50    #endif
51    
52  C     == Routine arguments ==  C     == Routine arguments ==
53    C     myTime - Current time in simulation
54    C     myIter - Current iteration number in simulation
55  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
56          _RL myTime
57          INTEGER myIter
58        INTEGER myThid        INTEGER myThid
59    
60  C     == Local variables  C     == Local variables
61  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[STUV]               o fVer: Vertical flux term - note fVer
 C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  
 C                              o uTrans: Zonal transport  
 C                              o vTrans: Meridional transport  
 C                              o wTrans: Vertical transport  
 C     maskC,maskUp             o maskC: land/water mask for tracer cells  
 C                              o maskUp: land/water mask for W points  
 C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  
 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, and level above
66  C     jMin, jMax   are applied.  C     phiHyd         - Hydrostatic part of the potential phiHydi.
67    C                      In z coords phiHydiHyd is the hydrostatic
68    C                      Potential (=pressure/rho0) anomaly
69    C                      In p coords phiHydiHyd is the geopotential
70    C                      surface height anomaly.
71    C     phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
72    C     phiSurfY             or geopotentiel (atmos) in X and Y direction
73    C     iMin, iMax     - Ranges and sub-block indices on which calculations
74    C     jMin, jMax       are applied.
75  C     bi, bj  C     bi, bj
76  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
77  C                          are switched with layer to be the appropriate index  C     kDown, km1       are switched with layer to be the appropriate
78  C                          into fVerTerm  C                      index into fVerTerm.
79        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
80        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
81        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
87        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
88        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
89        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
92        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C This is currently used by IVDC and Diagnostics
93        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
94        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  
       _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
95        INTEGER iMin, iMax        INTEGER iMin, iMax
96        INTEGER jMin, jMax        INTEGER jMin, jMax
97        INTEGER bi, bj        INTEGER bi, bj
98        INTEGER i, j        INTEGER i, j
99        INTEGER k, kM1, kUp, kDown        INTEGER k, km1, kp1, kup, kDown
100    
101    Cjmc : add for phiHyd output <- but not working if multi tile per CPU
102    c     CHARACTER*(MAX_LEN_MBUF) suff
103    c     LOGICAL  DIFFERENT_MULTIPLE
104    c     EXTERNAL DIFFERENT_MULTIPLE
105    Cjmc(end)
106    
107    C---    The algorithm...
108    C
109    C       "Correction Step"
110    C       =================
111    C       Here we update the horizontal velocities with the surface
112    C       pressure such that the resulting flow is either consistent
113    C       with the free-surface evolution or the rigid-lid:
114    C         U[n] = U* + dt x d/dx P
115    C         V[n] = V* + dt x d/dy P
116    C
117    C       "Calculation of Gs"
118    C       ===================
119    C       This is where all the accelerations and tendencies (ie.
120    C       physics, parameterizations etc...) are calculated
121    C         rho = rho ( theta[n], salt[n] )
122    C         b   = b(rho, theta)
123    C         K31 = K31 ( rho )
124    C         Gu[n] = Gu( u[n], v[n], wVel, b, ... )
125    C         Gv[n] = Gv( u[n], v[n], wVel, b, ... )
126    C         Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
127    C         Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
128    C
129    C       "Time-stepping" or "Prediction"
130    C       ================================
131    C       The models variables are stepped forward with the appropriate
132    C       time-stepping scheme (currently we use Adams-Bashforth II)
133    C       - For momentum, the result is always *only* a "prediction"
134    C       in that the flow may be divergent and will be "corrected"
135    C       later with a surface pressure gradient.
136    C       - Normally for tracers the result is the new field at time
137    C       level [n+1} *BUT* in the case of implicit diffusion the result
138    C       is also *only* a prediction.
139    C       - We denote "predictors" with an asterisk (*).
140    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
141    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
142    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
143    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
144    C       With implicit diffusion:
145    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
146    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
147    C         (1 + dt * K * d_zz) theta[n] = theta*
148    C         (1 + dt * K * d_zz) salt[n] = salt*
149    C---
150    
151  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
152  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
# Line 94  C     point numbers. This prevents spuri Line 155  C     point numbers. This prevents spuri
155  C     uninitialised but inert locations.  C     uninitialised but inert locations.
156        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
157         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
158          xA(i,j)      = 0.*1. _d 37          DO k=1,Nr
159          yA(i,j)      = 0.*1. _d 37           phiHyd(i,j,k)  = 0. _d 0
160          uTrans(i,j)  = 0.*1. _d 37           KappaRU(i,j,k) = 0. _d 0
161          vTrans(i,j)  = 0.*1. _d 37           KappaRV(i,j,k) = 0. _d 0
162          aTerm(i,j)   = 0.*1. _d 37           sigmaX(i,j,k) = 0. _d 0
163          xTerm(i,j)   = 0.*1. _d 37           sigmaY(i,j,k) = 0. _d 0
164          cTerm(i,j)   = 0.*1. _d 37           sigmaR(i,j,k) = 0. _d 0
         mTerm(i,j)   = 0.*1. _d 37  
         pTerm(i,j)   = 0.*1. _d 37  
         fZon(i,j)    = 0.*1. _d 37  
         fMer(i,j)    = 0.*1. _d 37  
         DO K=1,nZ  
          pH (i,j,k)  = 0.*1. _d 37  
165          ENDDO          ENDDO
166          rhokm1(i,j)    = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
167          rhokp1(i,j)    = 0. _d 0          rhok   (i,j) = 0. _d 0
168         ENDDO          phiSurfX(i,j) = 0. _d 0
169        ENDDO          phiSurfY(i,j) = 0. _d 0
 C--   Set up work arrays that need valid initial values  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         wTrans(i,j)  = 0. _d 0  
         fVerT(i,j,1) = 0. _d 0  
         fVerT(i,j,2) = 0. _d 0  
         fVerS(i,j,1) = 0. _d 0  
         fVerS(i,j,2) = 0. _d 0  
         fVerU(i,j,1) = 0. _d 0  
         fVerU(i,j,2) = 0. _d 0  
         fVerV(i,j,1) = 0. _d 0  
         fVerV(i,j,2) = 0. _d 0  
170         ENDDO         ENDDO
171        ENDDO        ENDDO
172    
173    #ifdef ALLOW_AUTODIFF_TAMC
174    C--   HPF directive to help TAMC
175    CHPF$ INDEPENDENT
176    #endif /* ALLOW_AUTODIFF_TAMC */
177    
178        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
179    
180    #ifdef ALLOW_AUTODIFF_TAMC
181    C--    HPF directive to help TAMC
182    CHPF$  INDEPENDENT, NEW (fVerU,fVerV
183    CHPF$&                  ,phiHyd
184    CHPF$&                  ,KappaRU,KappaRV
185    CHPF$&                  )
186    #endif /* ALLOW_AUTODIFF_TAMC */
187    
188         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
189    
190  C--   Boundary condition on hydrostatic pressure is pH(z=0)=0  #ifdef ALLOW_AUTODIFF_TAMC
191              act1 = bi - myBxLo(myThid)
192              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
193    
194              act2 = bj - myByLo(myThid)
195              max2 = myByHi(myThid) - myByLo(myThid) + 1
196    
197              act3 = myThid - 1
198              max3 = nTx*nTy
199    
200              act4 = ikey_dynamics - 1
201    
202              ikey = (act1 + 1) + act2*max1
203         &                      + act3*max1*max2
204         &                      + act4*max1*max2*max3
205    #endif /* ALLOW_AUTODIFF_TAMC */
206    
207    C--     Set up work arrays that need valid initial values
208          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
209           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
210            pH(i,j,1) = 0. _d 0            fVerU  (i,j,1) = 0. _d 0
211              fVerU  (i,j,2) = 0. _d 0
212              fVerV  (i,j,1) = 0. _d 0
213              fVerV  (i,j,2) = 0. _d 0
214           ENDDO           ENDDO
215          ENDDO          ENDDO
216    
217          iMin = 1-OLx+1  C--     Start computation of dynamics
218          iMax = sNx+OLx          iMin = 1-OLx+2
219          jMin = 1-OLy+1          iMax = sNx+OLx-1
220          jMax = sNy+OLy          jMin = 1-OLy+2
221            jMax = sNy+OLy-1
222  C--     Calculate gradient of surface pressure  
223          CALL GRAD_PSURF(  #ifdef ALLOW_AUTODIFF_TAMC
224       I       bi,bj,iMin,iMax,jMin,jMax,  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
225       O       pSurfX,pSurfY,  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
226       I       myThid)  CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
227    #endif /* ALLOW_AUTODIFF_TAMC */
228  C--     Update fields in top level according to tendency terms  
229          CALL TIMESTEP(  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
230       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)  C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
231            IF (implicSurfPress.NE.1.) THEN
232          DO K=2,Nz            CALL CALC_GRAD_PHI_SURF(
233  C--     Update fields in Kth level according to tendency terms       I         bi,bj,iMin,iMax,jMin,jMax,
234          CALL TIMESTEP(       I         etaN,
235       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)       O         phiSurfX,phiSurfY,
236  C Density of K-1 level (above W(K)) reference to K level       I         myThid )                        
237           CALL FIND_RHO(          ENDIF
238       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',  
239       O      rhoKm1,  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
240       I      myThid )  C--      Calculate the total vertical diffusivity
241  C Density of K level (below W(K)) reference to K level          DO k=1,Nr
242           CALL FIND_RHO(           CALL CALC_VISCOSITY(
243       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',       I        bi,bj,iMin,iMax,jMin,jMax,k,
244       O      rhoKp1,       O        KappaRU,KappaRV,
      I      myThid )  
 C--     Calculate static stability and mix where convectively unstable  
          CALL CONVECT(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,myThid)  
 C Density of K-1 level (above W(K)) reference to K-1 level  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, 'LINEAR',  
      O      rhoKm1,  
      I      myThid )  
 C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
          CALL CALC_PH(  
      I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,  
      U       pH,  
      I       myThid )  
         ENDDO ! K  
   
 C Density of Nz level (bottom level) reference to Nz level  
          CALL FIND_RHO(  
      I      bi, bj, iMin, iMax, jMin, jMax,  Nz, Nz, 'LINEAR',  
      O      rhoKm1,  
      I      myThid )  
 C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  
          CALL CALC_PH(  
      I       bi,bj,iMin,iMax,jMin,jMax,Nz+1,rhoKm1,  
      U       pH,  
      I       myThid )  
   
         DO K = Nz, 1, -1  
          kM1  =max(1,k-1)   ! Points to level above k (=k-1)  
          kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above  
          kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer  
          iMin = 1-OLx+2  
          iMax = sNx+OLx-1  
          jMin = 1-OLy+2  
          jMax = sNy+OLy-1  
   
 C--      Get temporary terms used by tendency routines  
          CALL CALC_COMMON_FACTORS (  
      I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
      O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,  
245       I        myThid)       I        myThid)
246           ENDDO
247    #endif
248    
249  C--      Calculate accelerations in the momentum equations  C--     Start of dynamics loop
250           CALL CALC_MOM_RHS(          DO k=1,Nr
      I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,  
      I        xA,yA,uTrans,vTrans,wTrans,maskC,  
      I        pH,  
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
      I        myThid)  
251    
252  C--      Calculate active tracer tendencies  C--       km1    Points to level above k (=k-1)
253           CALL CALC_GT(  C--       kup    Cycles through 1,2 to point to layer above
254       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  C--       kDown  Cycles through 2,1 to point to current layer
255       I        xA,yA,uTrans,vTrans,wTrans,maskUp,  
256       U        aTerm,xTerm,fZon,fMer,fVerT,            km1  = MAX(1,k-1)
257       I        myThid)            kp1  = MIN(k+1,Nr)
258  Cdbg     CALL CALC_GS(            kup  = 1+MOD(k+1,2)
259  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            kDown= 1+MOD(k,2)
260  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  
261  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  #ifdef ALLOW_AUTODIFF_TAMC
262  Cdbg I        myThid)           kkey = (ikey-1)*Nr + k
263    #endif /* ALLOW_AUTODIFF_TAMC */
264    
265    C--      Integrate hydrostatic balance for phiHyd with BC of
266    C        phiHyd(z=0)=0
267    C        distinguishe between Stagger and Non Stagger time stepping
268             IF (staggerTimeStep) THEN
269               CALL CALC_PHI_HYD(
270         I        bi,bj,iMin,iMax,jMin,jMax,k,
271         I        gTnm1, gSnm1,
272         U        phiHyd,
273         I        myThid )
274             ELSE
275               CALL CALC_PHI_HYD(
276         I        bi,bj,iMin,iMax,jMin,jMax,k,
277         I        theta, salt,
278         U        phiHyd,
279         I        myThid )
280             ENDIF
281    
282    C--      Calculate accelerations in the momentum equations (gU, gV, ...)
283    C        and step forward storing the result in gUnm1, gVnm1, etc...
284             IF ( momStepping ) THEN
285    #ifndef DISABLE_MOM_FLUXFORM
286               IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
287         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
288         I         phiHyd,KappaRU,KappaRV,
289         U         fVerU, fVerV,
290         I         myTime, myIter, myThid)
291    #endif
292    #ifndef DISABLE_MOM_VECINV
293               IF (vectorInvariantMomentum) CALL MOM_VECINV(
294         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
295         I         phiHyd,KappaRU,KappaRV,
296         U         fVerU, fVerV,
297         I         myTime, myIter, myThid)
298    #endif
299               CALL TIMESTEP(
300         I         bi,bj,iMin,iMax,jMin,jMax,k,
301         I         phiHyd, phiSurfX, phiSurfY,
302         I         myIter, myThid)
303    
304    #ifdef   ALLOW_OBCS
305    C--      Apply open boundary conditions
306             IF (useOBCS) THEN
307               CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
308             END IF
309    #endif   /* ALLOW_OBCS */
310    
311    #ifdef   ALLOW_AUTODIFF_TAMC
312    #ifdef   INCLUDE_CD_CODE
313             ELSE
314               DO j=1-OLy,sNy+OLy
315                 DO i=1-OLx,sNx+OLx
316                   guCD(i,j,k,bi,bj) = 0.0
317                   gvCD(i,j,k,bi,bj) = 0.0
318                 END DO
319               END DO
320    #endif   /* INCLUDE_CD_CODE */
321    #endif   /* ALLOW_AUTODIFF_TAMC */
322             ENDIF
323    
324    
325    C--     end of dynamics k loop (1:Nr)
326          ENDDO          ENDDO
327    
328    
329    
330    C--     Implicit viscosity
331            IF (implicitViscosity.AND.momStepping) THEN
332    #ifdef    ALLOW_AUTODIFF_TAMC
333              idkey = iikey + 3
334    CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
335    #endif    /* ALLOW_AUTODIFF_TAMC */
336              CALL IMPLDIFF(
337         I         bi, bj, iMin, iMax, jMin, jMax,
338         I         deltaTmom, KappaRU,recip_HFacW,
339         U         gUNm1,
340         I         myThid )
341    #ifdef    ALLOW_AUTODIFF_TAMC
342              idkey = iikey + 4
343    CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
344    #endif    /* ALLOW_AUTODIFF_TAMC */
345              CALL IMPLDIFF(
346         I         bi, bj, iMin, iMax, jMin, jMax,
347         I         deltaTmom, KappaRV,recip_HFacS,
348         U         gVNm1,
349         I         myThid )
350    
351    #ifdef   ALLOW_OBCS
352    C--      Apply open boundary conditions
353             IF (useOBCS) THEN
354               DO K=1,Nr
355                 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
356               ENDDO
357             END IF
358    #endif   /* ALLOW_OBCS */
359    
360    #ifdef    INCLUDE_CD_CODE
361    #ifdef    ALLOW_AUTODIFF_TAMC
362              idkey = iikey + 5
363    CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
364    #endif    /* ALLOW_AUTODIFF_TAMC */
365              CALL IMPLDIFF(
366         I         bi, bj, iMin, iMax, jMin, jMax,
367         I         deltaTmom, KappaRU,recip_HFacW,
368         U         vVelD,
369         I         myThid )
370    #ifdef    ALLOW_AUTODIFF_TAMC
371              idkey = iikey + 6
372    CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
373    #endif    /* ALLOW_AUTODIFF_TAMC */
374              CALL IMPLDIFF(
375         I         bi, bj, iMin, iMax, jMin, jMax,
376         I         deltaTmom, KappaRV,recip_HFacS,
377         U         uVelD,
378         I         myThid )
379    #endif    /* INCLUDE_CD_CODE */
380    C--     End If implicitViscosity.AND.momStepping
381            ENDIF
382    
383    Cjmc : add for phiHyd output <- but not working if multi tile per CPU
384    c       IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
385    c    &  .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
386    c         WRITE(suff,'(I10.10)') myIter+1
387    c         CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
388    c       ENDIF
389    Cjmc(end)
390    
391    #ifdef ALLOW_TIMEAVE
392            IF (taveFreq.GT.0.) THEN
393              CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
394         I                              deltaTclock, bi, bj, myThid)
395              IF (ivdc_kappa.NE.0.) THEN
396                CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
397         I                              deltaTclock, bi, bj, myThid)
398              ENDIF
399            ENDIF
400    #endif /* ALLOW_TIMEAVE */
401    
402         ENDDO         ENDDO
403        ENDDO        ENDDO
404    
405    #ifndef DISABLE_DEBUGMODE
406          If (debugMode) THEN
407           CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
408           CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
409           CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
410           CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
411           CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
412           CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
413           CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
414           CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
415           CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
416           CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
417           CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
418           CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
419           CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
420           CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
421          ENDIF
422    #endif
423    
424        RETURN        RETURN
425        END        END

Legend:
Removed from v.1.4  
changed lines
  Added in v.1.80

  ViewVC Help
Powered by ViewVC 1.1.22