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

Legend:
Removed from v.1.2  
changed lines
  Added in v.1.85.2.1

  ViewVC Help
Powered by ViewVC 1.1.22