/[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.10 by adcroft, Thu May 28 16:19:50 1998 UTC revision 1.49 by heimbach, Fri Jun 9 02:45:04 2000 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    
3  #include "CPP_EEOPTIONS.h"  #include "CPP_OPTIONS.h"
4    
5        SUBROUTINE DYNAMICS(myTime, myIter, myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6  C     /==========================================================\  C     /==========================================================\
# Line 20  C     | ===== Line 20  C     | =====
20  C     | C*P* comments indicating place holders for which code is |  C     | C*P* comments indicating place holders for which code is |
21  C     |      presently being developed.                          |  C     |      presently being developed.                          |
22  C     \==========================================================/  C     \==========================================================/
23    c
24    c     changed: Patrick Heimbach heimbach@mit.edu 6-Jun-2000
25    c              - computation of ikey wrong for nTx,nTy > 1
26    c                and/or nsx,nsy > 1: act1 and act2 were
27    c                mixed up.
28    
29          IMPLICIT NONE
30    
31  C     == Global variables ===  C     == Global variables ===
32  #include "SIZE.h"  #include "SIZE.h"
# Line 27  C     == Global variables === Line 34  C     == Global variables ===
34  #include "CG2D.h"  #include "CG2D.h"
35  #include "PARAMS.h"  #include "PARAMS.h"
36  #include "DYNVARS.h"  #include "DYNVARS.h"
37    #include "GRID.h"
38    
39    #ifdef ALLOW_KPP
40    #include "KPPMIX.h"
41    #endif
42    
43    #ifdef ALLOW_AUTODIFF_TAMC
44    #include "tamc.h"
45    #include "tamc_keys.h"
46    #endif
47    
48  C     == Routine arguments ==  C     == Routine arguments ==
49  C     myTime - Current time in simulation  C     myTime - Current time in simulation
50  C     myIter - Current iteration number in simulation  C     myIter - Current iteration number in simulation
51  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
       INTEGER myThid  
52        _RL myTime        _RL myTime
53        INTEGER myIter        INTEGER myIter
54          INTEGER myThid
55    
56  C     == Local variables  C     == Local variables
57  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
58  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C     uTrans, vTrans, rTrans - Per block temporaries holding flow
59  C                              o uTrans: Zonal transport  C                              transport
60    C     rVel                     o uTrans: Zonal transport
61  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
62  C                              o wTrans: Vertical transport  C                              o rTrans: Vertical transport
63    C                              o rVel:   Vertical velocity at upper and
64    C                                        lower cell faces.
65  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
66  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
67  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in
# Line 57  C                              o fVer: V Line 77  C                              o fVer: V
77  C                                      is "pipelined" in the vertical  C                                      is "pipelined" in the vertical
78  C                                      so we need an fVer for each  C                                      so we need an fVer for each
79  C                                      variable.  C                                      variable.
80  C     iMin, iMax - Ranges and sub-block indices on which calculations  C     rhoK, rhoKM1   - Density at current level, level above and level
81  C     jMin, jMax   are applied.  C                      below.
82    C     rhoKP1                                                                  
83    C     buoyK, buoyKM1 - Buoyancy at current level and level above.
84    C     phiHyd         - Hydrostatic part of the potential phiHydi.
85    C                      In z coords phiHydiHyd is the hydrostatic
86    C                      pressure anomaly
87    C                      In p coords phiHydiHyd is the geopotential
88    C                      surface height
89    C                      anomaly.
90    C     etaSurfX,      - Holds surface elevation gradient in X and Y.
91    C     etaSurfY
92    C     K13, K23, K33  - Non-zero elements of small-angle approximation
93    C                      diffusion tensor.
94    C     KapGM          - Spatially varying Visbeck et. al mixing coeff.
95    C     KappaRT,       - Total diffusion in vertical for T and S.
96    C     KappaRS          (background + spatially varying, isopycnal term).
97    C     iMin, iMax     - Ranges and sub-block indices on which calculations
98    C     jMin, jMax       are applied.
99  C     bi, bj  C     bi, bj
100  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
101  C                          are switched with layer to be the appropriate index  C     kDown, kM1       are switched with layer to be the appropriate
102  C                          into fVerTerm  C                      index into fVerTerm.
103        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS xA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rVel    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
109        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111        _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112        _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL xTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113        _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL cTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114        _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL mTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pTerm   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL fZon    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fMer    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
119        _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
120        _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
121        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL fVerV   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
122        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
123        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124        _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL buoyK   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL rhotmp  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131          _RL K13     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
132          _RL K23     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
133          _RL K33     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
134          _RL KapGM   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135          _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
136          _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
137          _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
138          _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
139    
140    #ifdef INCLUDE_CONVECT_CALL
141          _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
142    #endif
143    
144        INTEGER iMin, iMax        INTEGER iMin, iMax
145        INTEGER jMin, jMax        INTEGER jMin, jMax
146        INTEGER bi, bj        INTEGER bi, bj
147        INTEGER i, j        INTEGER i, j
148        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
149          LOGICAL BOTTOM_LAYER
150    
151    #ifdef ALLOW_AUTODIFF_TAMC
152          INTEGER    isbyte
153          PARAMETER( isbyte = 4 )
154    
155          INTEGER act1, act2, act3, act4
156          INTEGER max1, max2, max3
157          INTEGER ikact, iikey,kkey
158          INTEGER maximpl
159    #endif
160    
161    C---    The algorithm...
162    C
163    C       "Correction Step"
164    C       =================
165    C       Here we update the horizontal velocities with the surface
166    C       pressure such that the resulting flow is either consistent
167    C       with the free-surface evolution or the rigid-lid:
168    C         U[n] = U* + dt x d/dx P
169    C         V[n] = V* + dt x d/dy P
170    C
171    C       "Calculation of Gs"
172    C       ===================
173    C       This is where all the accelerations and tendencies (ie.
174    C       phiHydysics, parameterizations etc...) are calculated
175    C         rVel = sum_r ( div. u[n] )
176    C         rho = rho ( theta[n], salt[n] )
177    C         b   = b(rho, theta)
178    C         K31 = K31 ( rho )
179    C         Gu[n] = Gu( u[n], v[n], rVel, b, ... )
180    C         Gv[n] = Gv( u[n], v[n], rVel, b, ... )
181    C         Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
182    C         Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
183    C
184    C       "Time-stepping" or "Prediction"
185    C       ================================
186    C       The models variables are stepped forward with the appropriate
187    C       time-stepping scheme (currently we use Adams-Bashforth II)
188    C       - For momentum, the result is always *only* a "prediction"
189    C       in that the flow may be divergent and will be "corrected"
190    C       later with a surface pressure gradient.
191    C       - Normally for tracers the result is the new field at time
192    C       level [n+1} *BUT* in the case of implicit diffusion the result
193    C       is also *only* a prediction.
194    C       - We denote "predictors" with an asterisk (*).
195    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
196    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
197    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
198    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
199    C       With implicit diffusion:
200    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
201    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
202    C         (1 + dt * K * d_zz) theta[n] = theta*
203    C         (1 + dt * K * d_zz) salt[n] = salt*
204    C---
205    
206    #ifdef ALLOW_AUTODIFF_TAMC
207    C--   dummy statement to end declaration part
208          ikey = 1
209    #endif
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 115  C     uninitialised but inert locations. Line 226  C     uninitialised but inert locations.
226          pTerm(i,j)   = 0. _d 0          pTerm(i,j)   = 0. _d 0
227          fZon(i,j)    = 0. _d 0          fZon(i,j)    = 0. _d 0
228          fMer(i,j)    = 0. _d 0          fMer(i,j)    = 0. _d 0
229          DO K=1,nZ          DO K=1,Nr
230           pH (i,j,k)  = 0. _d 0           phiHyd (i,j,k)  = 0. _d 0
231           K13(i,j,k) = 0. _d 0           K13(i,j,k)  = 0. _d 0
232           K23(i,j,k) = 0. _d 0           K23(i,j,k)  = 0. _d 0
233           K33(i,j,k) = 0. _d 0           K33(i,j,k)  = 0. _d 0
234             KappaRU(i,j,k) = 0. _d 0
235             KappaRV(i,j,k) = 0. _d 0
236          ENDDO          ENDDO
237          rhokm1(i,j)  = 0. _d 0          rhoKM1 (i,j) = 0. _d 0
238          rhokp1(i,j)  = 0. _d 0          rhok   (i,j) = 0. _d 0
239          rhotmp(i,j)  = 0. _d 0          rhoKP1 (i,j) = 0. _d 0
240            rhoTMP (i,j) = 0. _d 0
241            buoyKM1(i,j) = 0. _d 0
242            buoyK  (i,j) = 0. _d 0
243            maskC  (i,j) = 0. _d 0
244         ENDDO         ENDDO
245        ENDDO        ENDDO
246    
247    
248    #ifdef ALLOW_AUTODIFF_TAMC
249    C--   HPF directive to help TAMC
250    !HPF$ INDEPENDENT
251    #endif
252    
253        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
254    
255    #ifdef ALLOW_AUTODIFF_TAMC
256    C--    HPF directive to help TAMC
257    !HPF$  INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
258    !HPF$&                  ,phiHyd,K13,K23,K33,KapGM
259    !HPF$&                  ,utrans,vtrans,maskc,xA,yA
260    !HPF$&                  ,KappaRT,KappaRS,KappaRU,KappaRV
261    !HPF$&                  )
262    #endif
263    
264         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
265    
266  C--     Boundary condition on hydrostatic pressure is pH(z=0)=0  #ifdef ALLOW_AUTODIFF_TAMC
267              act1 = bi - myBxLo(myThid)
268              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
269    
270              act2 = bj - myByLo(myThid)
271              max2 = myByHi(myThid) - myByLo(myThid) + 1
272    
273              act3 = myThid - 1
274              max3 = nTx*nTy
275    
276              act4 = ikey_dynamics - 1
277    
278              ikey = (act1 + 1) + act2*max1
279         &                      + act3*max1*max2
280         &                      + act4*max1*max2*max3
281    #endif
282    
283    C--     Set up work arrays that need valid initial values
284          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
285           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
286            pH(i,j,1) = 0. _d 0            rTrans(i,j)   = 0. _d 0
287            K13(i,j,1) = 0. _d 0            rVel  (i,j,1) = 0. _d 0
288            K23(i,j,1) = 0. _d 0            rVel  (i,j,2) = 0. _d 0
289            K33(i,j,1) = 0. _d 0            fVerT (i,j,1) = 0. _d 0
290            KapGM(i,j) = 0. _d 0            fVerT (i,j,2) = 0. _d 0
291              fVerS (i,j,1) = 0. _d 0
292              fVerS (i,j,2) = 0. _d 0
293              fVerU (i,j,1) = 0. _d 0
294              fVerU (i,j,2) = 0. _d 0
295              fVerV (i,j,1) = 0. _d 0
296              fVerV (i,j,2) = 0. _d 0
297              phiHyd(i,j,1) = 0. _d 0
298              K13   (i,j,1) = 0. _d 0
299              K23   (i,j,1) = 0. _d 0
300              K33   (i,j,1) = 0. _d 0
301              KapGM (i,j)   = GMkbackground
302           ENDDO           ENDDO
303          ENDDO          ENDDO
304    
305  C--     Set up work arrays that need valid initial values          DO k=1,Nr
306          DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
307           DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
308            wTrans(i,j)  = 0. _d 0  #ifdef INCLUDE_CONVECT_CALL
309            fVerT(i,j,1) = 0. _d 0             ConvectCount(i,j,k) = 0.
310            fVerT(i,j,2) = 0. _d 0  #endif
311            fVerS(i,j,1) = 0. _d 0             KappaRT(i,j,k) = 0. _d 0
312            fVerS(i,j,2) = 0. _d 0             KappaRS(i,j,k) = 0. _d 0
313            fVerU(i,j,1) = 0. _d 0            ENDDO
           fVerU(i,j,2) = 0. _d 0  
           fVerV(i,j,1) = 0. _d 0  
           fVerV(i,j,2) = 0. _d 0  
314           ENDDO           ENDDO
315          ENDDO          ENDDO
316    
# Line 161  C--     Set up work arrays that need val Line 319  C--     Set up work arrays that need val
319          jMin = 1-OLy+1          jMin = 1-OLy+1
320          jMax = sNy+OLy          jMax = sNy+OLy
321    
322    
323            K = 1
324            BOTTOM_LAYER = K .EQ. Nr
325    
326    #ifdef DO_PIPELINED_CORRECTION_STEP
327  C--     Calculate gradient of surface pressure  C--     Calculate gradient of surface pressure
328          CALL GRAD_PSURF(          CALL CALC_GRAD_ETA_SURF(
329       I       bi,bj,iMin,iMax,jMin,jMax,       I       bi,bj,iMin,iMax,jMin,jMax,
330       O       pSurfX,pSurfY,       O       etaSurfX,etaSurfY,
331       I       myThid)       I       myThid)
   
332  C--     Update fields in top level according to tendency terms  C--     Update fields in top level according to tendency terms
333          CALL TIMESTEP(          CALL CORRECTION_STEP(
334       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,K,
335         I       etaSurfX,etaSurfY,myTime,myThid)
336    
337    #ifdef ALLOW_OBCS
338            IF (openBoundaries) THEN
339    #ifdef ALLOW_AUTODIFF_TAMC
340    CADJ STORE uvel (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
341    CADJ STORE vvel (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
342    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
343    CADJ STORE salt(:,:,k,bi,bj)   = comlev1_2d, key = ikey, byte = isbyte
344    #endif
345               CALL APPLY_OBCS1( bi, bj, K, myThid )
346            END IF
347    #endif
348    
349            IF ( .NOT. BOTTOM_LAYER ) THEN
350    C--      Update fields in layer below according to tendency terms
351             CALL CORRECTION_STEP(
352         I        bi,bj,iMin,iMax,jMin,jMax,K+1,
353         I        etaSurfX,etaSurfY,myTime,myThid)
354    #ifdef ALLOW_OBCS
355             IF (openBoundaries) THEN
356    #ifdef ALLOW_AUTODIFF_TAMC
357    CADJ STORE uvel (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
358    CADJ STORE vvel (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
359    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
360    CADJ STORE salt(:,:,k,bi,bj)   = comlev1_2d, key = ikey, byte = isbyte
361    #endif
362                CALL APPLY_OBCS1( bi, bj, K+1, myThid )
363             END IF
364    #endif
365            ENDIF
366    #endif
367  C--     Density of 1st level (below W(1)) reference to level 1  C--     Density of 1st level (below W(1)) reference to level 1
368    #ifdef  INCLUDE_FIND_RHO_CALL
369    #ifdef ALLOW_AUTODIFF_TAMC
370    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
371    CADJ STORE salt (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
372    #endif
373          CALL FIND_RHO(          CALL FIND_RHO(
374       I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,       I     bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
375       O     rhoKm1,       O     rhoKm1,
376       I     myThid )       I     myThid )
377  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  #endif
378          CALL CALC_PH(  
379       I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,          IF (       (.NOT. BOTTOM_LAYER)
380       U      pH,  #ifdef ALLOW_KPP
381         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
382    #endif
383         &     ) THEN
384    C--      Check static stability with layer below
385    C--      and mix as needed.
386    #ifdef  INCLUDE_FIND_RHO_CALL
387    #ifdef ALLOW_AUTODIFF_TAMC
388    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
389    CADJ STORE salt (:,:,k,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
390    #endif
391             CALL FIND_RHO(
392         I      bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
393         O      rhoKp1,
394       I      myThid )       I      myThid )
395          DO J=1-Oly,sNy+Oly  #endif
          DO I=1-Olx,sNx+Olx  
           rhoKp1(I,J)=rhoKm1(I,J)  
          ENDDO  
         ENDDO  
396    
397          DO K=2,Nz  #ifdef  INCLUDE_CONVECT_CALL
398  C--     Update fields in Kth level according to tendency terms  
399          CALL TIMESTEP(  #ifdef ALLOW_AUTODIFF_TAMC
400       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)  CADJ STORE rhoKm1(:,:)  = comlev1_2d, key = ikey, byte = isbyte
401  C--     Density of K-1 level (above W(K)) reference to K-1 level  CADJ STORE rhoKp1(:,:)  = comlev1_2d, key = ikey, byte = isbyte
402  copt    CALL FIND_RHO(  #endif
403  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,           CALL CONVECT(
404  copt O     rhoKm1,       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
405  copt I     myThid )       U       ConvectCount,
406  C       rhoKm1=rhoKp1       I       myTime,myIter,myThid)
407          DO J=1-Oly,sNy+Oly  #ifdef ALLOW_AUTODIFF_TAMC
408           DO I=1-Olx,sNx+Olx  CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)
409            rhoKm1(I,J)=rhoKp1(I,J)  CADJ &     = comlev1_2d, key = ikey, byte = isbyte
410           ENDDO  CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)
411          ENDDO  CADJ &     = comlev1_2d, key = ikey, byte = isbyte
412  C--     Density of K level (below W(K)) reference to K level  #endif
413          CALL FIND_RHO(  
414       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  #endif
415       O     rhoKp1,  
416       I     myThid )  C--      Implicit Vertical Diffusion for Convection
417  C--     Density of K-1 level (above W(K)) reference to K level           IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
418          CALL FIND_RHO(       I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
419       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,       U       ConvectCount, KappaRT, KappaRS,
420       O     rhotmp,       I       myTime,myIter,myThid)
421       I     myThid )  CRG: do we need do store STORE KappaRT, KappaRS ?
422  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
423          CALL CALC_ISOSLOPES(  C--      Recompute density after mixing
424       I            bi, bj, iMin, iMax, jMin, jMax, K,  #ifdef  INCLUDE_FIND_RHO_CALL
425       I            rhoKm1, rhoKp1, rhotmp,           CALL FIND_RHO(
426       O            K13, K23, K33, KapGM,       I      bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
427       I            myThid )       O      rhoKm1,
 C--     Calculate static stability and mix where convectively unstable  
         CALL CONVECT(  
      I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,  
      I      myTime,myIter,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, eosType,  
      O     rhoKm1,  
      I     myThid )  
 C--     Density of K level (below W(K)) referenced to K level  
         CALL FIND_RHO(  
      I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,  
      O     rhoKp1,  
      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,rhoKp1,  
      U      pH,  
428       I      myThid )       I      myThid )
429    #endif
430            ENDIF
431    C--     Calculate buoyancy
432            CALL CALC_BUOYANCY(
433         I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
434         O      buoyKm1,
435         I      myThid )
436    C--     Integrate hydrostatic balance for phiHyd with BC of
437    C--     phiHyd(z=0)=0
438            CALL CALC_PHI_HYD(
439         I      bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
440         U      phiHyd,
441         I      myThid )
442    
443    C----------------------------------------------
444    C--     start of downward loop
445    C----------------------------------------------
446            DO K=2,Nr
447    
448    #ifdef ALLOW_AUTODIFF_TAMC
449             kkey = ikact*(Nr-2+1) + (k-2) + 1
450    #endif
451    
452             BOTTOM_LAYER = K .EQ. Nr
453    
454    #ifdef DO_PIPELINED_CORRECTION_STEP
455             IF ( .NOT. BOTTOM_LAYER ) THEN
456    C--       Update fields in layer below according to tendency terms
457              CALL CORRECTION_STEP(
458         I         bi,bj,iMin,iMax,jMin,jMax,K+1,
459         I         etaSurfX,etaSurfY,myTime,myThid)
460    #ifdef ALLOW_OBCS
461              IF (openBoundaries) THEN
462    #ifdef ALLOW_AUTODIFF_TAMC
463    CADJ STORE uvel (:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
464    CADJ STORE vvel (:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
465    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
466    CADJ STORE salt(:,:,k,bi,bj)   = comlev1_2d, key = ikey, byte = isbyte
467    #endif
468                 CALL APPLY_OBCS1( bi, bj, K+1, myThid )
469              END IF
470    #endif
471             ENDIF
472    #endif
473    
474    C--      Density of K level (below W(K)) reference to K level
475    #ifdef  INCLUDE_FIND_RHO_CALL
476    #ifdef ALLOW_AUTODIFF_TAMC
477    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
478    CADJ STORE salt (:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
479    #endif
480             CALL FIND_RHO(
481         I      bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
482         O      rhoK,
483         I      myThid )
484    #endif
485             IF (       (.NOT. BOTTOM_LAYER)
486    #ifdef ALLOW_KPP
487         &       .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
488    #endif
489         &      ) THEN
490    C--       Check static stability with layer below and mix as needed.
491    C--       Density of K+1 level (below W(K+1)) reference to K level.
492    #ifdef  INCLUDE_FIND_RHO_CALL
493    #ifdef ALLOW_AUTODIFF_TAMC
494    CADJ STORE theta(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
495    CADJ STORE salt (:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
496    #endif
497              CALL FIND_RHO(
498         I       bi, bj, iMin, iMax, jMin, jMax,  K+1, K, eosType,
499         O       rhoKp1,
500         I       myThid )
501    #endif
502    
503    #ifdef ALLOW_AUTODIFF_TAMC
504    CADJ STORE rhok  (:,:)   = comlev1_3d, key = kkey, byte = isbyte
505    CADJ STORE rhoKm1(:,:)   = comlev1_3d, key = kkey, byte = isbyte
506    CADJ STORE rhoKp1(:,:)   = comlev1_3d, key = kkey, byte = isbyte
507    #endif
508    
509    #ifdef  INCLUDE_CONVECT_CALL
510              CALL CONVECT(
511         I        bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
512         U        ConvectCount,
513         I        myTime,myIter,myThid)
514    #ifdef ALLOW_AUTODIFF_TAMC
515    CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)
516    CADJ &     = comlev1_3d, key = kkey, byte = isbyte
517    CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)
518    CADJ &     = comlev1_3d, key = kkey, byte = isbyte
519    #endif
520    #endif
521    
522    C--      Implicit Vertical Diffusion for Convection
523             IF (ivdc_kappa.NE.0.) THEN
524                CALL CALC_IVDC(
525         I       bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
526         U       ConvectCount, KappaRT, KappaRS,
527         I       myTime,myIter,myThid)
528    CRG: do we need do store STORE KappaRT, KappaRS ?
529             END IF
530    
531    C--       Recompute density after mixing
532    #ifdef  INCLUDE_FIND_RHO_CALL
533              CALL FIND_RHO(
534         I       bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
535         O       rhoK,
536         I       myThid )
537    #endif
538             ENDIF
539    C--      Calculate buoyancy
540             CALL CALC_BUOYANCY(
541         I       bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
542         O       buoyK,
543         I       myThid )
544    C--      Integrate hydrostatic balance for phiHyd with BC of
545    C--      phiHyd(z=0)=0
546             CALL CALC_PHI_HYD(
547         I        bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
548         U        phiHyd,
549         I        myThid )
550    C--      Calculate iso-neutral slopes for the GM/Redi parameterisation
551    #ifdef  INCLUDE_FIND_RHO_CALL
552             CALL FIND_RHO(
553         I        bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
554         O        rhoTmp,
555         I        myThid )
556    #endif
557    
558    #ifdef  INCLUDE_CALC_ISOSLOPES_CALL
559    #ifdef ALLOW_AUTODIFF_TAMC
560    CADJ STORE rhoTmp(:,:)  = comlev1_3d, key = kkey, byte = isbyte
561    CADJ STORE rhok  (:,:)  = comlev1_3d, key = kkey, byte = isbyte
562    CADJ STORE rhoKm1(:,:)  = comlev1_3d, key = kkey, byte = isbyte
563    CADJ STORE kapgm (:,:)  = comlev1_3d, key = kkey, byte = isbyte
564    #endif
565             CALL CALC_ISOSLOPES(
566         I        bi, bj, iMin, iMax, jMin, jMax, K,
567         I        rhoKm1, rhoK, rhotmp,
568         O        K13, K23, K33, KapGM,
569         I        myThid )
570    #endif
571    
572             DO J=jMin,jMax
573              DO I=iMin,iMax
574    #ifdef  INCLUDE_FIND_RHO_CALL
575               rhoKm1 (I,J) = rhoK(I,J)
576    #endif
577               buoyKm1(I,J) = buoyK(I,J)
578              ENDDO
579             ENDDO
580          ENDDO          ENDDO
581    C--     end of k loop
582    
583    #ifdef ALLOW_AUTODIFF_TAMC
584    CADJ STORE theta(:,:,:,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
585    CADJ STORE salt (:,:,:,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
586    CADJ STORE uvel (:,:,:,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
587    CADJ STORE vvel (:,:,:,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
588    #endif
589    
590    #ifdef ALLOW_KPP
591    C----------------------------------------------
592    C--     Compute KPP mixing coefficients
593    C----------------------------------------------
594            IF (usingKPPmixing) THEN
595    #ifdef ALLOW_AUTODIFF_TAMC
596    CADJ STORE fu  (:,:  ,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
597    CADJ STORE fv  (:,:  ,bi,bj)  = comlev1_2d, key = ikey, byte = isbyte
598    #endif
599             CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
600         I          , myThid)
601             CALL KVMIX(
602         I               bi, bj, myTime, myThid )
603             CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
604         I        , myThid)
605            ENDIF
606    #endif
607    
608    C----------------------------------------------
609    C--     start of upward loop
610    C----------------------------------------------
611            DO K = Nr, 1, -1
612    
         DO K = Nz, 1, -1  
613           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
614           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
615           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer           kDown=1+MOD(k,2)   ! Cycles through 2,1 to point to current layer
616    
617           iMin = 1-OLx+2           iMin = 1-OLx+2
618           iMax = sNx+OLx-1           iMax = sNx+OLx-1
619           jMin = 1-OLy+2           jMin = 1-OLy+2
620           jMax = sNy+OLy-1           jMax = sNy+OLy-1
621    
622    #ifdef ALLOW_AUTODIFF_TAMC
623             kkey = ikact*(Nr-1+1) + (k-1) + 1
624    #endif
625    
626    #ifdef ALLOW_AUTODIFF_TAMC
627    CADJ STORE rvel  (:,:,kDown)  = comlev1_3d, key = kkey, byte = isbyte
628    CADJ STORE rTrans(:,:)        = comlev1_3d, key = kkey, byte = isbyte
629    CADJ STORE KappaRT(:,:,:)     = comlev1_3d, key = kkey, byte = isbyte
630    CADJ STORE KappaRS(:,:,:)     = comlev1_3d, key = kkey, byte = isbyte
631    #endif
632    
633  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
634           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
635       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
636       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
637       I        myThid)       I        myThid)
638    
639    #ifdef ALLOW_OBCS
640            IF (openBoundaries) THEN
641             CALL APPLY_OBCS3( bi, bj, K, Kup, rTrans, rVel, myThid )
642            ENDIF
643    #endif
644    
645    #ifdef  INCLUDE_CALC_DIFFUSIVITY_CALL
646    C--      Calculate the total vertical diffusivity
647             CALL CALC_DIFFUSIVITY(
648         I        bi,bj,iMin,iMax,jMin,jMax,K,
649         I        maskC,maskUp,KapGM,K33,
650         O        KappaRT,KappaRS,KappaRU,KappaRV,
651         I        myThid)
652    #endif
653  C--      Calculate accelerations in the momentum equations  C--      Calculate accelerations in the momentum equations
654           IF ( momStepping ) THEN           IF ( momStepping ) THEN
655            CALL CALC_MOM_RHS(            CALL CALC_MOM_RHS(
656       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
657       I         xA,yA,uTrans,vTrans,wTrans,maskC,       I         xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
658       I         pH,       I         phiHyd,KappaRU,KappaRV,
659       U         aTerm,xTerm,cTerm,mTerm,pTerm,       U         aTerm,xTerm,cTerm,mTerm,pTerm,
660       U         fZon, fMer, fVerU, fVerV,       U         fZon, fMer, fVerU, fVerV,
661       I         myThid)       I         myTime, myThid)
662    #ifdef ALLOW_AUTODIFF_TAMC
663    #ifdef INCLUDE_CD_CODE
664             ELSE
665                DO j=1-OLy,sNy+OLy
666                   DO i=1-OLx,sNx+OLx
667                      guCD(i,j,k,bi,bj) = 0.0
668                      gvCD(i,j,k,bi,bj) = 0.0
669                   END DO
670                END DO
671    #endif
672    #endif
673           ENDIF           ENDIF
   
674  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
675           IF ( tempStepping ) THEN           IF ( tempStepping ) THEN
676            CALL CALC_GT(            CALL CALC_GT(
677       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
678       I         xA,yA,uTrans,vTrans,wTrans,maskUp,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
679       I         K13,K23,K33,KapGM,       I         K13,K23,KappaRT,KapGM,
680       U         aTerm,xTerm,fZon,fMer,fVerT,       U         aTerm,xTerm,fZon,fMer,fVerT,
681       I         myThid)       I         myTime, myThid)
682             ENDIF
683             IF ( saltStepping ) THEN
684              CALL CALC_GS(
685         I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
686         I         xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
687         I         K13,K23,KappaRS,KapGM,
688         U         aTerm,xTerm,fZon,fMer,fVerS,
689         I         myTime, myThid)
690             ENDIF
691    #ifdef ALLOW_OBCS
692    C--      Calculate future values on open boundaries
693             IF (openBoundaries) THEN
694    Caja      CALL CYCLE_OBCS( K, bi, bj, myThid )
695              CALL SET_OBCS( K, bi, bj, myTime+deltaTclock, myThid )
696           ENDIF           ENDIF
697  Cdbg     CALL CALC_GS(  #endif
698  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  C--      Prediction step (step forward all model variables)
699  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,           CALL TIMESTEP(
700  Cdbg I        K13,K23,K33,KapGM,       I       bi,bj,iMin,iMax,jMin,jMax,K,
701  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,       I       myIter, myThid)
702  Cdbg I        myThid)  #ifdef ALLOW_OBCS
703    C--      Apply open boundary conditions
704             IF (openBoundaries) THEN
705    #ifdef ALLOW_AUTODIFF_TAMC
706    CADJ STORE gunm1(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
707    CADJ STORE gvnm1(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
708    CADJ STORE gwnm1(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
709    #endif
710                CALL APPLY_OBCS2( bi, bj, K, myThid )
711             END IF
712    #endif
713    C--      Freeze water
714             IF (allowFreezing) THEN
715    #ifdef ALLOW_AUTODIFF_TAMC
716    CADJ STORE gTNm1(:,:,k,bi,bj)  = comlev1_3d, key = kkey, byte = isbyte
717    #endif
718                CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
719             END IF
720    
721    #ifdef DIVG_IN_DYNAMICS
722    C--      Diagnose barotropic divergence of predicted fields
723             CALL CALC_DIV_GHAT(
724         I       bi,bj,iMin,iMax,jMin,jMax,K,
725         I       xA,yA,
726         I       myThid)
727    #endif /* DIVG_IN_DYNAMICS */
728    
729          ENDDO  C--      Cumulative diagnostic calculations (ie. time-averaging)
730    #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
731             IF (taveFreq.GT.0.) THEN
732              CALL DO_TIME_AVERAGES(
733         I                           myTime, myIter, bi, bj, K, kUp, kDown,
734         I                           K13, K23, rVel, KapGM, ConvectCount,
735         I                           myThid )
736             ENDIF
737    #endif
738    
739    
740            ENDDO ! K
741    
742    C--     Implicit diffusion
743            IF (implicitDiffusion) THEN
744    
745    #ifdef ALLOW_AUTODIFF_TAMC
746               maximpl = 6
747               iikey = ikact*maximpl
748    #endif
749    
750             IF (tempStepping) THEN
751    #ifdef ALLOW_AUTODIFF_TAMC
752                idkey = iikey + 1
753    #endif
754                CALL IMPLDIFF(
755         I         bi, bj, iMin, iMax, jMin, jMax,
756         I         deltaTtracer, KappaRT,recip_HFacC,
757         U         gTNm1,
758         I         myThid )
759             END IF
760    
761             IF (saltStepping) THEN
762    #ifdef ALLOW_AUTODIFF_TAMC
763             idkey = iikey + 2
764    #endif
765                CALL IMPLDIFF(
766         I         bi, bj, iMin, iMax, jMin, jMax,
767         I         deltaTtracer, KappaRS,recip_HFacC,
768         U         gSNm1,
769         I         myThid )
770             END IF
771    
772            ENDIF ! implicitDiffusion
773    
774    C--     Implicit viscosity
775            IF (implicitViscosity) THEN
776    
777             IF (momStepping) THEN
778    #ifdef ALLOW_AUTODIFF_TAMC
779             idkey = iikey + 3
780    #endif
781              CALL IMPLDIFF(
782         I         bi, bj, iMin, iMax, jMin, jMax,
783         I         deltaTmom, KappaRU,recip_HFacW,
784         U         gUNm1,
785         I         myThid )
786    #ifdef ALLOW_AUTODIFF_TAMC
787             idkey = iikey + 4
788    #endif
789              CALL IMPLDIFF(
790         I         bi, bj, iMin, iMax, jMin, jMax,
791         I         deltaTmom, KappaRV,recip_HFacS,
792         U         gVNm1,
793         I         myThid )
794    
795    #ifdef INCLUDE_CD_CODE
796    
797    #ifdef ALLOW_AUTODIFF_TAMC
798             idkey = iikey + 5
799    #endif
800              CALL IMPLDIFF(
801         I         bi, bj, iMin, iMax, jMin, jMax,
802         I         deltaTmom, KappaRU,recip_HFacW,
803         U         vVelD,
804         I         myThid )
805    #ifdef ALLOW_AUTODIFF_TAMC
806            idkey = iikey + 6
807    #endif
808              CALL IMPLDIFF(
809         I         bi, bj, iMin, iMax, jMin, jMax,
810         I         deltaTmom, KappaRV,recip_HFacS,
811         U         uVelD,
812         I         myThid )
813    
814    #endif
815    
816             ENDIF ! momStepping
817            ENDIF ! implicitViscosity
818    
819         ENDDO         ENDDO
820        ENDDO        ENDDO
821    
822  !dbg  write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)  C     write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
823  !dbg  write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
824  !dbg &                         maxval(uVel(1:sNx,1:sNy,:,:,:))  C     write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
825  !dbg  write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),  C    &                           maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
826  !dbg &                         maxval(vVel(1:sNx,1:sNy,:,:,:))  C     write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
827  !dbg  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),  C    &                           maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
828  !dbg &                         maxval(K13(1:sNx,1:sNy,:))  C     write(0,*) 'dynamics: rVel(1) ',
829  !dbg  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),  C    &            minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
830  !dbg &                         maxval(K23(1:sNx,1:sNy,:))  C    &            maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
831  !dbg  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),  C     write(0,*) 'dynamics: rVel(2) ',
832  !dbg &                         maxval(K33(1:sNx,1:sNy,:))  C    &            minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
833  !dbg  write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),  C    &            maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
834  !dbg &                         maxval(gT(1:sNx,1:sNy,:,:,:))  cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
835  !dbg  write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),  cblk &                           maxval(K13(1:sNx,1:sNy,:))
836  !dbg &                         maxval(Theta(1:sNx,1:sNy,:,:,:))  cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
837  !dbg  write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),  cblk &                           maxval(K23(1:sNx,1:sNy,:))
838  !dbg &                          maxval(pH/(Gravity*Rhonil))  cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
839    cblk &                           maxval(K33(1:sNx,1:sNy,:))
840    C     write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
841    C    &                           maxval(gT(1:sNx,1:sNy,:,:,:))
842    C     write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
843    C    &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
844    C     write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
845    C    &                           maxval(gS(1:sNx,1:sNy,:,:,:))
846    C     write(0,*) 'dynamics: S  ',minval(salt(1:sNx,1:sNy,:,:,:)),
847    C    &                           maxval(salt(1:sNx,1:sNy,:,:,:))
848    C     write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
849    C    &                           maxval(phiHyd/(Gravity*Rhonil))
850    C     CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
851    C    &Nr, 1, myThid )
852    C     CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
853    C    &Nr, 1, myThid )
854    C     CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
855    C    &Nr, 1, myThid )
856    C     CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
857    C    &Nr, 1, myThid )
858    C     CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
859    C    &Nr, 1, myThid )
860    
861    
862        RETURN        RETURN
863        END        END

Legend:
Removed from v.1.10  
changed lines
  Added in v.1.49

  ViewVC Help
Powered by ViewVC 1.1.22