/[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.68 by adcroft, Tue May 29 14:01:37 2001 UTC revision 1.153 by jmc, Thu Nov 18 00:57:24 2010 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2  C $Name$  C $Name$
3    
4    #include "PACKAGES_CONFIG.h"
5  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
6    #ifdef ALLOW_OBCS
7    # include "OBCS_OPTIONS.h"
8    #endif
9    
10    #undef DYNAMICS_GUGV_EXCH_CHECK
11    
12    CBOP
13    C     !ROUTINE: DYNAMICS
14    C     !INTERFACE:
15        SUBROUTINE DYNAMICS(myTime, myIter, myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
16  C     /==========================================================\  C     !DESCRIPTION: \bv
17  C     | SUBROUTINE DYNAMICS                                      |  C     *==========================================================*
18  C     | o Controlling routine for the explicit part of the model |  C     | SUBROUTINE DYNAMICS
19  C     |   dynamics.                                              |  C     | o Controlling routine for the explicit part of the model
20  C     |==========================================================|  C     |   dynamics.
21  C     | This routine evaluates the "dynamics" terms for each     |  C     *==========================================================*
22  C     | block of ocean in turn. Because the blocks of ocean have |  C     | This routine evaluates the "dynamics" terms for each
23  C     | overlap regions they are independent of one another.     |  C     | block of ocean in turn. Because the blocks of ocean have
24  C     | If terms involving lateral integrals are needed in this  |  C     | overlap regions they are independent of one another.
25  C     | routine care will be needed. Similarly finite-difference |  C     | If terms involving lateral integrals are needed in this
26  C     | operations with stencils wider than the overlap region   |  C     | routine care will be needed. Similarly finite-difference
27  C     | require special consideration.                           |  C     | operations with stencils wider than the overlap region
28  C     | Notes                                                    |  C     | require special consideration.
29  C     | =====                                                    |  C     | The algorithm...
30  C     | C*P* comments indicating place holders for which code is |  C     |
31  C     |      presently being developed.                          |  C     | "Correction Step"
32  C     \==========================================================/  C     | =================
33    C     | Here we update the horizontal velocities with the surface
34    C     | pressure such that the resulting flow is either consistent
35    C     | with the free-surface evolution or the rigid-lid:
36    C     |   U[n] = U* + dt x d/dx P
37    C     |   V[n] = V* + dt x d/dy P
38    C     |   W[n] = W* + dt x d/dz P  (NH mode)
39    C     |
40    C     | "Calculation of Gs"
41    C     | ===================
42    C     | This is where all the accelerations and tendencies (ie.
43    C     | physics, parameterizations etc...) are calculated
44    C     |   rho = rho ( theta[n], salt[n] )
45    C     |   b   = b(rho, theta)
46    C     |   K31 = K31 ( rho )
47    C     |   Gu[n] = Gu( u[n], v[n], wVel, b, ... )
48    C     |   Gv[n] = Gv( u[n], v[n], wVel, b, ... )
49    C     |   Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
50    C     |   Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
51    C     |
52    C     | "Time-stepping" or "Prediction"
53    C     | ================================
54    C     | The models variables are stepped forward with the appropriate
55    C     | time-stepping scheme (currently we use Adams-Bashforth II)
56    C     | - For momentum, the result is always *only* a "prediction"
57    C     | in that the flow may be divergent and will be "corrected"
58    C     | later with a surface pressure gradient.
59    C     | - Normally for tracers the result is the new field at time
60    C     | level [n+1} *BUT* in the case of implicit diffusion the result
61    C     | is also *only* a prediction.
62    C     | - We denote "predictors" with an asterisk (*).
63    C     |   U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
64    C     |   V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
65    C     |   theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66    C     |   salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
67    C     | With implicit diffusion:
68    C     |   theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69    C     |   salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
70    C     |   (1 + dt * K * d_zz) theta[n] = theta*
71    C     |   (1 + dt * K * d_zz) salt[n] = salt*
72    C     |
73    C     *==========================================================*
74    C     \ev
75    C     !USES:
76        IMPLICIT NONE        IMPLICIT NONE
   
77  C     == Global variables ===  C     == Global variables ===
78  #include "SIZE.h"  #include "SIZE.h"
79  #include "EEPARAMS.h"  #include "EEPARAMS.h"
80  #include "PARAMS.h"  #include "PARAMS.h"
81  #include "DYNVARS.h"  #include "DYNVARS.h"
82    #ifdef ALLOW_CD_CODE
83    #include "CD_CODE_VARS.h"
84    #endif
85  #include "GRID.h"  #include "GRID.h"
   
86  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
87  # include "tamc.h"  # include "tamc.h"
88  # include "tamc_keys.h"  # include "tamc_keys.h"
89  # include "FFIELDS.h"  # include "FFIELDS.h"
90    # include "EOS.h"
91  # ifdef ALLOW_KPP  # ifdef ALLOW_KPP
92  #  include "KPP.h"  #  include "KPP.h"
93  # endif  # endif
94  # ifdef ALLOW_GMREDI  # ifdef ALLOW_PTRACERS
95  #  include "GMREDI.h"  #  include "PTRACERS_SIZE.h"
96    #  include "PTRACERS_FIELDS.h"
97    # endif
98    # ifdef ALLOW_OBCS
99    #  include "OBCS.h"
100    #  ifdef ALLOW_PTRACERS
101    #   include "OBCS_PTRACERS.h"
102    #  endif
103    # endif
104    # ifdef ALLOW_MOM_FLUXFORM
105    #  include "MOM_FLUXFORM.h"
106  # endif  # endif
107  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
108    
109  #ifdef ALLOW_TIMEAVE  C     !CALLING SEQUENCE:
110  #include "TIMEAVE_STATV.h"  C     DYNAMICS()
111  #endif  C      |
112    C      |-- CALC_EP_FORCING
113    C      |
114    C      |-- CALC_GRAD_PHI_SURF
115    C      |
116    C      |-- CALC_VISCOSITY
117    C      |
118    C      |-- CALC_PHI_HYD
119    C      |
120    C      |-- MOM_FLUXFORM
121    C      |
122    C      |-- MOM_VECINV
123    C      |
124    C      |-- TIMESTEP
125    C      |
126    C      |-- OBCS_APPLY_UV
127    C      |
128    C      |-- MOM_U_IMPLICIT_R
129    C      |-- MOM_V_IMPLICIT_R
130    C      |
131    C      |-- IMPLDIFF
132    C      |
133    C      |-- OBCS_APPLY_UV
134    C      |
135    C      |-- CALC_GW
136    C      |
137    C      |-- DIAGNOSTICS_FILL
138    C      |-- DEBUG_STATS_RL
139    
140    C     !INPUT/OUTPUT PARAMETERS:
141  C     == Routine arguments ==  C     == Routine arguments ==
142  C     myTime - Current time in simulation  C     myTime :: Current time in simulation
143  C     myIter - Current iteration number in simulation  C     myIter :: Current iteration number in simulation
144  C     myThid - Thread number for this instance of the routine.  C     myThid :: Thread number for this instance of the routine.
145        _RL myTime        _RL myTime
146        INTEGER myIter        INTEGER myIter
147        INTEGER myThid        INTEGER myThid
148    
149    C     !FUNCTIONS:
150    #ifdef ALLOW_DIAGNOSTICS
151          LOGICAL  DIAGNOSTICS_IS_ON
152          EXTERNAL DIAGNOSTICS_IS_ON
153    #endif
154    
155    C     !LOCAL VARIABLES:
156  C     == Local variables  C     == Local variables
157  C     xA, yA                 - Per block temporaries holding face areas  C     fVer[UV]               o fVer: Vertical flux term - note fVer
158  C     uTrans, vTrans, rTrans - Per block temporaries holding flow  C                                    is "pipelined" in the vertical
159  C                              transport  C                                    so we need an fVer for each
160  C                              o uTrans: Zonal transport  C                                    variable.
161  C                              o vTrans: Meridional transport  C     phiHydC    :: hydrostatic potential anomaly at cell center
162  C                              o rTrans: Vertical transport  C                   In z coords phiHyd is the hydrostatic potential
163  C     maskUp                   o maskUp: land/water mask for W points  C                      (=pressure/rho0) anomaly
164  C     fVer[STUV]               o fVer: Vertical flux term - note fVer  C                   In p coords phiHyd is the geopotential height anomaly.
165  C                                      is "pipelined" in the vertical  C     phiHydF    :: hydrostatic potential anomaly at middle between 2 centers
166  C                                      so we need an fVer for each  C     dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
167  C                                      variable.  C     phiSurfX,  ::  gradient of Surface potential (Pressure/rho, ocean)
168  C     rhoK, rhoKM1   - Density at current level, and level above  C     phiSurfY             or geopotential (atmos) in X and Y direction
169  C     phiHyd         - Hydrostatic part of the potential phiHydi.  C     guDissip   :: dissipation tendency (all explicit terms), u component
170  C                      In z coords phiHydiHyd is the hydrostatic  C     gvDissip   :: dissipation tendency (all explicit terms), v component
171  C                      Potential (=pressure/rho0) anomaly  C     KappaRU    :: vertical viscosity
172  C                      In p coords phiHydiHyd is the geopotential  C     KappaRV    :: vertical viscosity
 C                      surface height anomaly.  
 C     phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)  
 C     phiSurfY             or geopotentiel (atmos) in X and Y direction  
 C     KappaRT,       - Total diffusion in vertical for T and S.  
 C     KappaRS          (background + spatially varying, isopycnal term).  
173  C     iMin, iMax     - Ranges and sub-block indices on which calculations  C     iMin, iMax     - Ranges and sub-block indices on which calculations
174  C     jMin, jMax       are applied.  C     jMin, jMax       are applied.
175  C     bi, bj  C     bi, bj
176  C     k, kup,        - Index for layer above and below. kup and kDown  C     k, kup,        - Index for layer above and below. kup and kDown
177  C     kDown, km1       are switched with layer to be the appropriate  C     kDown, km1       are switched with layer to be the appropriate
178  C                      index into fVerTerm.  C                      index into fVerTerm.
 C     tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf.  
       _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS maskUp  (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)  
179        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
180        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
181        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
182        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
183        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
184          _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
185        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
186        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
187        _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
188        _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
189        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
190        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
       _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL tauAB  
   
 C This is currently used by IVDC and Diagnostics  
       _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
191    
192        INTEGER iMin, iMax        INTEGER iMin, iMax
193        INTEGER jMin, jMax        INTEGER jMin, jMax
194        INTEGER bi, bj        INTEGER bi, bj
195        INTEGER i, j        INTEGER i, j
196        INTEGER k, km1, kup, kDown        INTEGER k, km1, kp1, kup, kDown
197    
198    #ifdef ALLOW_DIAGNOSTICS
199          LOGICAL dPhiHydDiagIsOn
200          _RL tmpFac
201    #endif /* ALLOW_DIAGNOSTICS */
202    
203    
 Cjmc : add for phiHyd output <- but not working if multi tile per CPU  
 c     CHARACTER*(MAX_LEN_MBUF) suff  
 c     LOGICAL  DIFFERENT_MULTIPLE  
 c     EXTERNAL DIFFERENT_MULTIPLE  
 Cjmc(end)  
   
204  C---    The algorithm...  C---    The algorithm...
205  C  C
206  C       "Correction Step"  C       "Correction Step"
# Line 165  C         salt* = salt[n] + dt x ( 3/2 G Line 244  C         salt* = salt[n] + dt x ( 3/2 G
244  C         (1 + dt * K * d_zz) theta[n] = theta*  C         (1 + dt * K * d_zz) theta[n] = theta*
245  C         (1 + dt * K * d_zz) salt[n] = salt*  C         (1 + dt * K * d_zz) salt[n] = salt*
246  C---  C---
247    CEOP
248    
249  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_DEBUG
250  C--   dummy statement to end declaration part        IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
251        ikey = 1  #endif
 #endif /* ALLOW_AUTODIFF_TAMC */  
252    
253  C--   Set up work arrays with valid (i.e. not NaN) values  #ifdef ALLOW_DIAGNOSTICS
254  C     These inital values do not alter the numerical results. They        dPhiHydDiagIsOn = .FALSE.
255  C     just ensure that all memory references are to valid floating        IF ( useDiagnostics )
256  C     point numbers. This prevents spurious hardware signals due to       &  dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
257  C     uninitialised but inert locations.       &               .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
258        DO j=1-OLy,sNy+OLy  #endif
        DO i=1-OLx,sNx+OLx  
         xA(i,j)      = 0. _d 0  
         yA(i,j)      = 0. _d 0  
         uTrans(i,j)  = 0. _d 0  
         vTrans(i,j)  = 0. _d 0  
         DO k=1,Nr  
          phiHyd(i,j,k)  = 0. _d 0  
          KappaRU(i,j,k) = 0. _d 0  
          KappaRV(i,j,k) = 0. _d 0  
          sigmaX(i,j,k) = 0. _d 0  
          sigmaY(i,j,k) = 0. _d 0  
          sigmaR(i,j,k) = 0. _d 0  
         ENDDO  
         rhoKM1 (i,j) = 0. _d 0  
         rhok   (i,j) = 0. _d 0  
         phiSurfX(i,j) = 0. _d 0  
         phiSurfY(i,j) = 0. _d 0  
        ENDDO  
       ENDDO  
259    
260    C-- Call to routine for calculation of
261    C   Eliassen-Palm-flux-forced U-tendency,
262    C   if desired:
263    #ifdef INCLUDE_EP_FORCING_CODE
264          CALL CALC_EP_FORCING(myThid)
265    #endif
266    
267  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
268  C--   HPF directive to help TAMC  C--   HPF directive to help TAMC
# Line 207  CHPF$ INDEPENDENT Line 273  CHPF$ INDEPENDENT
273    
274  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
275  C--    HPF directive to help TAMC  C--    HPF directive to help TAMC
276  CHPF$  INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV  CHPF$  INDEPENDENT, NEW (fVerU,fVerV
277  CHPF$&                  ,phiHyd,utrans,vtrans,xA,yA  CHPF$&                  ,phiHydF
278  CHPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV  CHPF$&                  ,KappaRU,KappaRV
279  CHPF$&                  )  CHPF$&                  )
280  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
281    
# Line 218  CHPF$&                  ) Line 284  CHPF$&                  )
284  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
285            act1 = bi - myBxLo(myThid)            act1 = bi - myBxLo(myThid)
286            max1 = myBxHi(myThid) - myBxLo(myThid) + 1            max1 = myBxHi(myThid) - myBxLo(myThid) + 1
   
287            act2 = bj - myByLo(myThid)            act2 = bj - myByLo(myThid)
288            max2 = myByHi(myThid) - myByLo(myThid) + 1            max2 = myByHi(myThid) - myByLo(myThid) + 1
   
289            act3 = myThid - 1            act3 = myThid - 1
290            max3 = nTx*nTy            max3 = nTx*nTy
   
291            act4 = ikey_dynamics - 1            act4 = ikey_dynamics - 1
292              idynkey = (act1 + 1) + act2*max1
           ikey = (act1 + 1) + act2*max1  
293       &                      + act3*max1*max2       &                      + act3*max1*max2
294       &                      + act4*max1*max2*max3       &                      + act4*max1*max2*max3
295  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
296    
297  C--     Set up work arrays that need valid initial values  C--   Set up work arrays with valid (i.e. not NaN) values
298          DO j=1-OLy,sNy+OLy  C     These inital values do not alter the numerical results. They
299           DO i=1-OLx,sNx+OLx  C     just ensure that all memory references are to valid floating
300            rTrans(i,j)   = 0. _d 0  C     point numbers. This prevents spurious hardware signals due to
301            fVerT (i,j,1) = 0. _d 0  C     uninitialised but inert locations.
           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  
302    
303    #ifdef ALLOW_AUTODIFF_TAMC
304          DO k=1,Nr          DO k=1,Nr
305           DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
306            DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
307  C This is currently also used by IVDC and Diagnostics             KappaRU(i,j,k) = 0. _d 0
308             ConvectCount(i,j,k) = 0.             KappaRV(i,j,k) = 0. _d 0
309             KappaRT(i,j,k) = 0. _d 0  cph(
310             KappaRS(i,j,k) = 0. _d 0  c--   need some re-initialisation here to break dependencies
311    cph)
312               gU(i,j,k,bi,bj) = 0. _d 0
313               gV(i,j,k,bi,bj) = 0. _d 0
314            ENDDO            ENDDO
315           ENDDO           ENDDO
316          ENDDO          ENDDO
   
         iMin = 1-OLx+1  
         iMax = sNx+OLx  
         jMin = 1-OLy+1  
         jMax = sNy+OLy  
   
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 C--     Start of diagnostic loop  
         DO k=Nr,1,-1  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 C? Patrick, is this formula correct now that we change the loop range?  
 C? Do we still need this?  
 cph kkey formula corrected.  
 cph Needed for rhok, rhokm1, in the case useGMREDI.  
          kkey = (ikey-1)*Nr + k  
 CADJ STORE rhokm1(:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
 CADJ STORE rhok  (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
317  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
318            DO j=1-OLy,sNy+OLy
319  C--       Integrate continuity vertically for vertical velocity           DO i=1-OLx,sNx+OLx
320            CALL INTEGRATE_FOR_W(            fVerU  (i,j,1) = 0. _d 0
321       I                         bi, bj, k, uVel, vVel,            fVerU  (i,j,2) = 0. _d 0
322       O                         wVel,            fVerV  (i,j,1) = 0. _d 0
323       I                         myThid )            fVerV  (i,j,2) = 0. _d 0
324              phiHydF (i,j)  = 0. _d 0
325  #ifdef    ALLOW_OBCS            phiHydC (i,j)  = 0. _d 0
326  #ifdef    ALLOW_NONHYDROSTATIC  #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
327  C--       Apply OBC to W if in N-H mode            dPhiHydX(i,j)  = 0. _d 0
328            IF (useOBCS.AND.nonHydrostatic) THEN            dPhiHydY(i,j)  = 0. _d 0
329              CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )  #endif
330            ENDIF            phiSurfX(i,j)  = 0. _d 0
331  #endif    /* ALLOW_NONHYDROSTATIC */            phiSurfY(i,j)  = 0. _d 0
332  #endif    /* ALLOW_OBCS */            guDissip(i,j)  = 0. _d 0
333              gvDissip(i,j)  = 0. _d 0
334  C--       Calculate gradients of potential density for isoneutral  #ifdef ALLOW_AUTODIFF_TAMC
335  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)            phiHydLow(i,j,bi,bj) = 0. _d 0
336  c         IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN  # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
337            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN  #  ifndef DISABLE_RSTAR_CODE
338  #ifdef ALLOW_AUTODIFF_TAMC            dWtransC(i,j,bi,bj) = 0. _d 0
339  CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte            dWtransU(i,j,bi,bj) = 0. _d 0
340  CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte            dWtransV(i,j,bi,bj) = 0. _d 0
341  #endif /* ALLOW_AUTODIFF_TAMC */  #  endif
342              CALL FIND_RHO(  # endif
343       I        bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,  #endif
344       I        theta, salt,           ENDDO
      O        rhoK,  
      I        myThid )  
             IF (k.GT.1) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte  
 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
              CALL FIND_RHO(  
      I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,  
      I        theta, salt,  
      O        rhoKm1,  
      I        myThid )  
             ENDIF  
             CALL GRAD_SIGMA(  
      I             bi, bj, iMin, iMax, jMin, jMax, k,  
      I             rhoK, rhoKm1, rhoK,  
      O             sigmaX, sigmaY, sigmaR,  
      I             myThid )  
           ENDIF  
   
 C--       Implicit Vertical Diffusion for Convection  
 c ==> should use sigmaR !!!  
           IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN  
             CALL CALC_IVDC(  
      I        bi, bj, iMin, iMax, jMin, jMax, k,  
      I        rhoKm1, rhoK,  
      U        ConvectCount, KappaRT, KappaRS,  
      I        myTime, myIter, myThid)  
           ENDIF  
   
 C--     end of diagnostic k loop (Nr:1)  
345          ENDDO          ENDDO
346    
347  #ifdef ALLOW_AUTODIFF_TAMC  C--     Start computation of dynamics
348  cph avoids recomputation of integrate_for_w          iMin = 0
349  CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte          iMax = sNx+1
350  #endif /* ALLOW_AUTODIFF_TAMC */          jMin = 0
351            jMax = sNy+1
 #ifdef  ALLOW_OBCS  
 C--     Calculate future values on open boundaries  
         IF (useOBCS) THEN  
           CALL OBCS_CALC( bi, bj, myTime+deltaT,  
      I            uVel, vVel, wVel, theta, salt,  
      I            myThid )  
         ENDIF  
 #endif  /* ALLOW_OBCS */  
   
 C--     Determines forcing terms based on external fields  
 C       relaxation terms, etc.  
         CALL EXTERNAL_FORCING_SURF(  
      I             bi, bj, iMin, iMax, jMin, jMax,  
      I             myThid )  
 #ifdef ALLOW_AUTODIFF_TAMC  
 cph needed for KPP  
 CADJ STORE surfacetendencyU(:,:,bi,bj)  
 CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  
 CADJ STORE surfacetendencyV(:,:,bi,bj)  
 CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  
 CADJ STORE surfacetendencyS(:,:,bi,bj)  
 CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  
 CADJ STORE surfacetendencyT(:,:,bi,bj)  
 CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 #ifdef  ALLOW_GMREDI  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte  
 CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte  
 CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
 C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
         IF (useGMRedi) THEN  
           DO k=1,Nr  
             CALL GMREDI_CALC_TENSOR(  
      I             bi, bj, iMin, iMax, jMin, jMax, k,  
      I             sigmaX, sigmaY, sigmaR,  
      I             myThid )  
           ENDDO  
 #ifdef ALLOW_AUTODIFF_TAMC  
         ELSE  
           DO k=1, Nr  
             CALL GMREDI_CALC_TENSOR_DUMMY(  
      I             bi, bj, iMin, iMax, jMin, jMax, k,  
      I             sigmaX, sigmaY, sigmaR,  
      I             myThid )  
           ENDDO  
 #endif /* ALLOW_AUTODIFF_TAMC */  
         ENDIF  
352    
353  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
354  CADJ STORE Kwx(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE wvel (:,:,:,bi,bj) =
355  CADJ STORE Kwy(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  CADJ &     comlev1_bibj, key=idynkey, byte=isbyte
 CADJ STORE Kwz(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  
356  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
357    
358  #endif  /* ALLOW_GMREDI */  C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
359    C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)
360  #ifdef  ALLOW_KPP          IF (implicSurfPress.NE.1.) THEN
361  C--     Compute KPP mixing coefficients            CALL CALC_GRAD_PHI_SURF(
362          IF (useKPP) THEN       I         bi,bj,iMin,iMax,jMin,jMax,
363            CALL KPP_CALC(       I         etaN,
364       I                  bi, bj, myTime, myThid )       O         phiSurfX,phiSurfY,
365  #ifdef ALLOW_AUTODIFF_TAMC       I         myThid )
         ELSE  
           CALL KPP_CALC_DUMMY(  
      I                  bi, bj, myTime, myThid )  
 #endif /* ALLOW_AUTODIFF_TAMC */  
366          ENDIF          ENDIF
367    
368  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
369  CADJ STORE KPPghat   (:,:,:,bi,bj)  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
370  CADJ &   , KPPviscAz (:,:,:,bi,bj)  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
371  CADJ &   , KPPdiffKzT(:,:,:,bi,bj)  #ifdef ALLOW_KPP
372  CADJ &   , KPPdiffKzS(:,:,:,bi,bj)  CADJ STORE KPPviscAz (:,:,:,bi,bj)
373  CADJ &   , KPPfrac   (:,:  ,bi,bj)  CADJ &                 = comlev1_bibj, key=idynkey, byte=isbyte
374  CADJ &                 = comlev1_bibj, key=ikey, byte=isbyte  #endif /* ALLOW_KPP */
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 #endif  /* ALLOW_KPP */  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 #ifdef ALLOW_AIM  
 C       AIM - atmospheric intermediate model, physics package code.  
 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics  
         IF ( useAIM ) THEN  
          CALL TIMER_START('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)  
          CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )  
          CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)  
         ENDIF  
 #endif /* ALLOW_AIM */  
   
   
 C--     Start of thermodynamics loop  
         DO k=Nr,1,-1  
 #ifdef ALLOW_AUTODIFF_TAMC  
 C? Patrick Is this formula correct?  
 cph Yes, but I rewrote it.  
 cph Also, the KappaR? need the index and subscript k!  
          kkey = (ikey-1)*Nr + k  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 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  
   
 C--      Get temporary terms used by tendency routines  
          CALL CALC_COMMON_FACTORS (  
      I        bi,bj,iMin,iMax,jMin,jMax,k,  
      O        xA,yA,uTrans,vTrans,rTrans,maskUp,  
      I        myThid)  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE KappaRT(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte  
 CADJ STORE KappaRS(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte  
375  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
376    
377  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL  #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
378  C--      Calculate the total vertical diffusivity  C--     Calculate the total vertical viscosity
379           CALL CALC_DIFFUSIVITY(          CALL CALC_VISCOSITY(
380       I        bi,bj,iMin,iMax,jMin,jMax,k,       I            bi,bj, iMin,iMax,jMin,jMax,
381       I        maskUp,       O            KappaRU, KappaRV,
382       O        KappaRT,KappaRS,KappaRU,KappaRV,       I            myThid )
383       I        myThid)  #else
384  #endif          DO k=1,Nr
385             DO j=1-OLy,sNy+OLy
386  C--      Calculate active tracer tendencies (gT,gS,...)            DO i=1-OLx,sNx+OLx
387  C        and step forward storing result in gTnm1, gSnm1, etc.             KappaRU(i,j,k) = 0. _d 0
388           IF ( tempStepping ) THEN             KappaRV(i,j,k) = 0. _d 0
389             CALL CALC_GT(            ENDDO
390       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,           ENDDO
      I         xA,yA,uTrans,vTrans,rTrans,maskUp,  
      I         KappaRT,  
      U         fVerT,  
      I         myTime, myThid)  
            tauAB = 0.5d0 + abEps  
            CALL TIMESTEP_TRACER(  
      I         bi,bj,iMin,iMax,jMin,jMax,k,tauAB,  
      I         theta, gT,  
      U         gTnm1,  
      I         myIter, myThid)  
          ENDIF  
          IF ( saltStepping ) THEN  
            CALL CALC_GS(  
      I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,  
      I         xA,yA,uTrans,vTrans,rTrans,maskUp,  
      I         KappaRS,  
      U         fVerS,  
      I         myTime, myThid)  
            tauAB = 0.5d0 + abEps  
            CALL TIMESTEP_TRACER(  
      I         bi,bj,iMin,iMax,jMin,jMax,k,tauAB,  
      I         salt, gS,  
      U         gSnm1,  
      I         myIter, myThid)  
          ENDIF  
   
 #ifdef   ALLOW_OBCS  
 C--      Apply open boundary conditions  
          IF (useOBCS) THEN  
            CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )  
          END IF  
 #endif   /* ALLOW_OBCS */  
   
 C--      Freeze water  
          IF (allowFreezing) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k  
 CADJ &   , key = kkey, byte = isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
             CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )  
          END IF  
   
 C--     end of thermodynamic k loop (Nr:1)  
391          ENDDO          ENDDO
392    #endif
393    
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 C? Patrick? What about this one?  
 cph Keys iikey and idkey don't seem to be needed  
 cph since storing occurs on different tape for each  
 cph impldiff call anyways.  
 cph Thus, common block comlev1_impl isn't needed either.  
 cph Storing below needed in the case useGMREDI.  
         iikey = (ikey-1)*maximpl  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 C--     Implicit diffusion  
         IF (implicitDiffusion) THEN  
   
          IF (tempStepping) THEN  
 #ifdef ALLOW_AUTODIFF_TAMC  
             idkey = iikey + 1  
 CADJ STORE gTNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
             CALL IMPLDIFF(  
      I         bi, bj, iMin, iMax, jMin, jMax,  
      I         deltaTtracer, KappaRT, recip_HFacC,  
      U         gTNm1,  
      I         myThid )  
          ENDIF  
   
          IF (saltStepping) THEN  
394  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
395           idkey = iikey + 2  CADJ STORE KappaRU(:,:,:)
396  CADJ STORE gSNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  CADJ &     = comlev1_bibj, key=idynkey, byte=isbyte
397    CADJ STORE KappaRV(:,:,:)
398    CADJ &     = comlev1_bibj, key=idynkey, byte=isbyte
399  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
             CALL IMPLDIFF(  
      I         bi, bj, iMin, iMax, jMin, jMax,  
      I         deltaTtracer, KappaRS, recip_HFacC,  
      U         gSNm1,  
      I         myThid )  
          ENDIF  
   
 #ifdef   ALLOW_OBCS  
 C--      Apply open boundary conditions  
          IF (useOBCS) THEN  
            DO K=1,Nr  
              CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )  
            ENDDO  
          END IF  
 #endif   /* ALLOW_OBCS */  
   
 C--     End If implicitDiffusion  
         ENDIF  
   
 C--     Start computation of dynamics  
         iMin = 1-OLx+2  
         iMax = sNx+OLx-1  
         jMin = 1-OLy+2  
         jMax = sNy+OLy-1  
   
 C--     Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)  
 C       (note: this loop will be replaced by CALL CALC_GRAD_ETA)  
         IF (implicSurfPress.NE.1.) THEN  
           CALL CALC_GRAD_PHI_SURF(  
      I         bi,bj,iMin,iMax,jMin,jMax,  
      I         etaN,  
      O         phiSurfX,phiSurfY,  
      I         myThid )                          
         ENDIF  
400    
401  C--     Start of dynamics loop  C--     Start of dynamics loop
402          DO k=1,Nr          DO k=1,Nr
# Line 617  C--       kup    Cycles through 1,2 to p Line 406  C--       kup    Cycles through 1,2 to p
406  C--       kDown  Cycles through 2,1 to point to current layer  C--       kDown  Cycles through 2,1 to point to current layer
407    
408            km1  = MAX(1,k-1)            km1  = MAX(1,k-1)
409              kp1  = MIN(k+1,Nr)
410            kup  = 1+MOD(k+1,2)            kup  = 1+MOD(k+1,2)
411            kDown= 1+MOD(k,2)            kDown= 1+MOD(k,2)
412    
413  C--      Integrate hydrostatic balance for phiHyd with BC of  #ifdef ALLOW_AUTODIFF_TAMC
414             kkey = (idynkey-1)*Nr + k
415    c
416    CADJ STORE totphihyd (:,:,k,bi,bj)
417    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
418    CADJ STORE phihydlow (:,:,bi,bj)
419    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
420    CADJ STORE theta (:,:,k,bi,bj)
421    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
422    CADJ STORE salt  (:,:,k,bi,bj)
423    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
424    CADJ STORE gt(:,:,k,bi,bj)
425    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
426    CADJ STORE gs(:,:,k,bi,bj)
427    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
428    # ifdef NONLIN_FRSURF
429    cph-test
430    CADJ STORE  phiHydC (:,:)
431    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
432    CADJ STORE  phiHydF (:,:)
433    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
434    CADJ STORE  gudissip (:,:)
435    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
436    CADJ STORE  gvdissip (:,:)
437    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
438    CADJ STORE  fVerU (:,:,:)
439    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
440    CADJ STORE  fVerV (:,:,:)
441    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
442    CADJ STORE gu(:,:,k,bi,bj)
443    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
444    CADJ STORE gv(:,:,k,bi,bj)
445    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
446    #  ifndef ALLOW_ADAMSBASHFORTH_3
447    CADJ STORE gunm1(:,:,k,bi,bj)
448    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
449    CADJ STORE gvnm1(:,:,k,bi,bj)
450    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
451    #  else
452    CADJ STORE gunm(:,:,k,bi,bj,1)
453    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
454    CADJ STORE gunm(:,:,k,bi,bj,2)
455    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
456    CADJ STORE gvnm(:,:,k,bi,bj,1)
457    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
458    CADJ STORE gvnm(:,:,k,bi,bj,2)
459    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
460    #  endif
461    #  ifdef ALLOW_CD_CODE
462    CADJ STORE unm1(:,:,k,bi,bj)
463    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
464    CADJ STORE vnm1(:,:,k,bi,bj)
465    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
466    CADJ STORE uVelD(:,:,k,bi,bj)
467    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
468    CADJ STORE vVelD(:,:,k,bi,bj)
469    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
470    #  endif
471    # endif
472    # ifdef ALLOW_DEPTH_CONTROL
473    CADJ STORE  fVerU (:,:,:)
474    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
475    CADJ STORE  fVerV (:,:,:)
476    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
477    # endif
478    #endif /* ALLOW_AUTODIFF_TAMC */
479    
480    C--      Integrate hydrostatic balance for phiHyd with BC of
481  C        phiHyd(z=0)=0  C        phiHyd(z=0)=0
482  C        distinguishe between Stagger and Non Stagger time stepping           IF ( implicitIntGravWave ) THEN
          IF (staggerTimeStep) THEN  
483             CALL CALC_PHI_HYD(             CALL CALC_PHI_HYD(
484       I        bi,bj,iMin,iMax,jMin,jMax,k,       I        bi,bj,iMin,iMax,jMin,jMax,k,
485       I        gTnm1, gSnm1,       I        gT, gS,
486       U        phiHyd,       U        phiHydF,
487       I        myThid )       O        phiHydC, dPhiHydX, dPhiHydY,
488         I        myTime, myIter, myThid )
489           ELSE           ELSE
490             CALL CALC_PHI_HYD(             CALL CALC_PHI_HYD(
491       I        bi,bj,iMin,iMax,jMin,jMax,k,       I        bi,bj,iMin,iMax,jMin,jMax,k,
492       I        theta, salt,       I        theta, salt,
493       U        phiHyd,       U        phiHydF,
494       I        myThid )       O        phiHydC, dPhiHydX, dPhiHydY,
495         I        myTime, myIter, myThid )
496             ENDIF
497    #ifdef ALLOW_DIAGNOSTICS
498             IF ( dPhiHydDiagIsOn ) THEN
499               tmpFac = -1. _d 0
500               CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
501         &                           'Um_dPHdx', k, 1, 2, bi, bj, myThid )
502               CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
503         &                           'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
504           ENDIF           ENDIF
505    #endif /* ALLOW_DIAGNOSTICS */
506    
507  C--      Calculate accelerations in the momentum equations (gU, gV, ...)  C--      Calculate accelerations in the momentum equations (gU, gV, ...)
508  C        and step forward storing the result in gUnm1, gVnm1, etc...  C        and step forward storing the result in gU, gV, etc...
509           IF ( momStepping ) THEN           IF ( momStepping ) THEN
510             CALL CALC_MOM_RHS(  #ifdef ALLOW_AUTODIFF_TAMC
511    # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
512    #  ifndef DISABLE_RSTAR_CODE
513    CADJ STORE dWtransC(:,:,bi,bj)
514    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
515    CADJ STORE dWtransU(:,:,bi,bj)
516    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
517    CADJ STORE dWtransV(:,:,bi,bj)
518    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
519    #  endif
520    # endif
521    #endif
522               IF (.NOT. vectorInvariantMomentum) THEN
523    #ifdef ALLOW_MOM_FLUXFORM
524    C
525                  CALL MOM_FLUXFORM(
526       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,       I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
527       I         phiHyd,KappaRU,KappaRV,       I         KappaRU, KappaRV,
528       U         fVerU, fVerV,       U         fVerU, fVerV,
529       I         myTime, myThid)       O         guDissip, gvDissip,
530         I         myTime, myIter, myThid)
531    #endif
532               ELSE
533    #ifdef ALLOW_MOM_VECINV
534    C
535    # ifdef ALLOW_AUTODIFF_TAMC
536    #  ifdef NONLIN_FRSURF
537    CADJ STORE fVerU(:,:,:)
538    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
539    CADJ STORE fVerV(:,:,:)
540    CADJ &     = comlev1_bibj_k, key=kkey, byte=isbyte
541    #  endif
542    # endif /* ALLOW_AUTODIFF_TAMC */
543    C
544                 CALL MOM_VECINV(
545         I         bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
546         I         KappaRU, KappaRV,
547         U         fVerU, fVerV,
548         O         guDissip, gvDissip,
549         I         myTime, myIter, myThid)
550    #endif
551               ENDIF
552    C
553             CALL TIMESTEP(             CALL TIMESTEP(
554       I         bi,bj,iMin,iMax,jMin,jMax,k,       I         bi,bj,iMin,iMax,jMin,jMax,k,
555       I         phiHyd, phiSurfX, phiSurfY,       I         dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
556       I         myIter, myThid)       I         guDissip, gvDissip,
557         I         myTime, myIter, myThid)
558    
559  #ifdef   ALLOW_OBCS  #ifdef   ALLOW_OBCS
560  C--      Apply open boundary conditions  C--      Apply open boundary conditions
561           IF (useOBCS) THEN  c          IF (useOBCS) THEN
562             CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )  c            CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
563           END IF  c          ENDIF
564  #endif   /* ALLOW_OBCS */  #endif   /* ALLOW_OBCS */
565    
 #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   /* INCLUDE_CD_CODE */  
 #endif   /* ALLOW_AUTODIFF_TAMC */  
566           ENDIF           ENDIF
567    
   
568  C--     end of dynamics k loop (1:Nr)  C--     end of dynamics k loop (1:Nr)
569          ENDDO          ENDDO
570    
571    C--     Implicit Vertical advection & viscosity
572    #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
573  C--     Implicit viscosity       defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
574          IF (implicitViscosity.AND.momStepping) THEN          IF ( momImplVertAdv ) THEN
575              CALL MOM_U_IMPLICIT_R( kappaRU,
576         I                           bi, bj, myTime, myIter, myThid )
577              CALL MOM_V_IMPLICIT_R( kappaRV,
578         I                           bi, bj, myTime, myIter, myThid )
579            ELSEIF ( implicitViscosity ) THEN
580    #else /* INCLUDE_IMPLVERTADV_CODE */
581            IF     ( implicitViscosity ) THEN
582    #endif /* INCLUDE_IMPLVERTADV_CODE */
583  #ifdef    ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
584            idkey = iikey + 3  CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
585  CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
586  #endif    /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
587            CALL IMPLDIFF(            CALL IMPLDIFF(
588       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
589       I         deltaTmom, KappaRU,recip_HFacW,       I         -1, KappaRU,recip_HFacW,
590       U         gUNm1,       U         gU,
591       I         myThid )       I         myThid )
592  #ifdef    ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
593            idkey = iikey + 4  CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
594  CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
595  #endif    /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
596            CALL IMPLDIFF(            CALL IMPLDIFF(
597       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
598       I         deltaTmom, KappaRV,recip_HFacS,       I         -2, KappaRV,recip_HFacS,
599       U         gVNm1,       U         gV,
600       I         myThid )       I         myThid )
601            ENDIF
602    
603  #ifdef   ALLOW_OBCS  #ifdef   ALLOW_OBCS
604  C--      Apply open boundary conditions  C--      Apply open boundary conditions
605           IF (useOBCS) THEN  c       IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
606             DO K=1,Nr  c          DO K=1,Nr
607               CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )  c            CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
608             ENDDO  c          ENDDO
609           END IF          IF ( useOBCS ) THEN
610              CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
611            ENDIF
612  #endif   /* ALLOW_OBCS */  #endif   /* ALLOW_OBCS */
613    
614  #ifdef    INCLUDE_CD_CODE  #ifdef    ALLOW_CD_CODE
615            IF (implicitViscosity.AND.useCDscheme) THEN
616  #ifdef    ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
617            idkey = iikey + 5  CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  
618  #endif    /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
619            CALL IMPLDIFF(            CALL IMPLDIFF(
620       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
621       I         deltaTmom, KappaRU,recip_HFacW,       I         0, KappaRU,recip_HFacW,
622       U         vVelD,       U         vVelD,
623       I         myThid )       I         myThid )
624  #ifdef    ALLOW_AUTODIFF_TAMC  #ifdef    ALLOW_AUTODIFF_TAMC
625            idkey = iikey + 6  CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  
626  #endif    /* ALLOW_AUTODIFF_TAMC */  #endif    /* ALLOW_AUTODIFF_TAMC */
627            CALL IMPLDIFF(            CALL IMPLDIFF(
628       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
629       I         deltaTmom, KappaRV,recip_HFacS,       I         0, KappaRV,recip_HFacS,
630       U         uVelD,       U         uVelD,
631       I         myThid )       I         myThid )
 #endif    /* INCLUDE_CD_CODE */  
 C--     End If implicitViscosity.AND.momStepping  
632          ENDIF          ENDIF
633    #endif    /* ALLOW_CD_CODE */
634  Cjmc : add for phiHyd output <- but not working if multi tile per CPU  C--     End implicit Vertical advection & viscosity
635  c       IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)  
636  c    &  .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
637  c         WRITE(suff,'(I10.10)') myIter+1  
638  c         CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)  #ifdef ALLOW_NONHYDROSTATIC
639  c       ENDIF  C--   Step forward W field in N-H algorithm
640  Cjmc(end)          IF ( nonHydrostatic ) THEN
641    #ifdef ALLOW_DEBUG
642  #ifdef ALLOW_TIMEAVE           IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
643          IF (taveFreq.GT.0.) THEN  #endif
644            CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,           CALL TIMER_START('CALC_GW          [DYNAMICS]',myThid)
645       I                              deltaTclock, bi, bj, myThid)           CALL CALC_GW(
646            IF (ivdc_kappa.NE.0.) THEN       I                 bi,bj, KappaRU, KappaRV,
647              CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,       I                 myTime, myIter, myThid )
      I                              deltaTclock, bi, bj, myThid)  
           ENDIF  
648          ENDIF          ENDIF
649  #endif /* ALLOW_TIMEAVE */          IF ( nonHydrostatic.OR.implicitIntGravWave )
650         &   CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
651            IF ( nonHydrostatic )
652         &   CALL TIMER_STOP ('CALC_GW          [DYNAMICS]',myThid)
653    #endif
654    
655    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
656    
657    C-    end of bi,bj loops
658         ENDDO         ENDDO
659        ENDDO        ENDDO
660    
661    #ifdef ALLOW_OBCS
662          IF (useOBCS) THEN
663           CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
664          ENDIF
665    #endif
666    
667    Cml(
668    C     In order to compare the variance of phiHydLow of a p/z-coordinate
669    C     run with etaH of a z/p-coordinate run the drift of phiHydLow
670    C     has to be removed by something like the following subroutine:
671    C      CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
672    C     &                     'phiHydLow', myTime, myThid )
673    Cml)
674    
675    #ifdef ALLOW_DIAGNOSTICS
676          IF ( useDiagnostics ) THEN
677    
678           CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD  ',0,Nr,0,1,1,myThid)
679           CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT  ',0, 1,0,1,1,myThid)
680    
681           tmpFac = 1. _d 0
682           CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
683         &                                 'PHIHYDSQ',0,Nr,0,1,1,myThid)
684    
685           CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
686         &                                 'PHIBOTSQ',0, 1,0,1,1,myThid)
687    
688          ENDIF
689    #endif /* ALLOW_DIAGNOSTICS */
690    
691    #ifdef ALLOW_DEBUG
692          IF ( debugLevel .GE. debLevB ) THEN
693           CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
694           CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
695           CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
696           CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
697           CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
698           CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
699           CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
700           CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
701           CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
702           CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
703    #ifndef ALLOW_ADAMSBASHFORTH_3
704           CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
705           CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
706           CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
707           CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
708    #endif
709          ENDIF
710    #endif
711    
712    #ifdef DYNAMICS_GUGV_EXCH_CHECK
713    C- jmc: For safety checking only: This Exchange here should not change
714    C       the solution. If solution changes, it means something is wrong,
715    C       but it does not mean that it is less wrong with this exchange.
716          IF ( debugLevel .GT. debLevB ) THEN
717           CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
718          ENDIF
719    #endif
720    
721    #ifdef ALLOW_DEBUG
722          IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
723    #endif
724    
725        RETURN        RETURN
726        END        END

Legend:
Removed from v.1.68  
changed lines
  Added in v.1.153

  ViewVC Help
Powered by ViewVC 1.1.22