/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Annotation of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.45 - (hide annotations) (download)
Thu Aug 26 17:47:37 1999 UTC (24 years, 9 months ago) by adcroft
Branch: MAIN
Changes since 1.44: +33 -4 lines
Added IVDC (Implicit Vertical Diffusion Convection).
Also facilitated a "convection counter" that is output through "diags".

1 adcroft 1.45 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.44 1999/07/28 16:32:12 adcroft Exp $
2 cnh 1.1
3 adcroft 1.24 #include "CPP_OPTIONS.h"
4 cnh 1.1
5 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6 cnh 1.1 C /==========================================================\
7     C | SUBROUTINE DYNAMICS |
8     C | o Controlling routine for the explicit part of the model |
9     C | dynamics. |
10     C |==========================================================|
11     C | This routine evaluates the "dynamics" terms for each |
12     C | block of ocean in turn. Because the blocks of ocean have |
13     C | overlap regions they are independent of one another. |
14     C | If terms involving lateral integrals are needed in this |
15     C | routine care will be needed. Similarly finite-difference |
16     C | operations with stencils wider than the overlap region |
17     C | require special consideration. |
18     C | Notes |
19     C | ===== |
20     C | C*P* comments indicating place holders for which code is |
21     C | presently being developed. |
22     C \==========================================================/
23 adcroft 1.40 IMPLICIT NONE
24 cnh 1.1
25     C == Global variables ===
26     #include "SIZE.h"
27     #include "EEPARAMS.h"
28     #include "CG2D.h"
29 adcroft 1.6 #include "PARAMS.h"
30 adcroft 1.3 #include "DYNVARS.h"
31 adcroft 1.42 #include "GRID.h"
32     #ifdef ALLOW_KPP
33     #include "KPPMIX.h"
34     #endif
35 cnh 1.1
36     C == Routine arguments ==
37 cnh 1.8 C myTime - Current time in simulation
38     C myIter - Current iteration number in simulation
39 cnh 1.1 C myThid - Thread number for this instance of the routine.
40     INTEGER myThid
41 cnh 1.8 _RL myTime
42     INTEGER myIter
43 cnh 1.1
44     C == Local variables
45     C xA, yA - Per block temporaries holding face areas
46 cnh 1.38 C uTrans, vTrans, rTrans - Per block temporaries holding flow
47     C transport
48 cnh 1.30 C rVel o uTrans: Zonal transport
49 cnh 1.1 C o vTrans: Meridional transport
50 cnh 1.30 C o rTrans: Vertical transport
51 cnh 1.38 C o rVel: Vertical velocity at upper and
52     C lower cell faces.
53 cnh 1.1 C maskC,maskUp o maskC: land/water mask for tracer cells
54     C o maskUp: land/water mask for W points
55     C aTerm, xTerm, cTerm - Work arrays for holding separate terms in
56     C mTerm, pTerm, tendency equations.
57     C fZon, fMer, fVer[STUV] o aTerm: Advection term
58     C o xTerm: Mixing term
59     C o cTerm: Coriolis term
60     C o mTerm: Metric term
61     C o pTerm: Pressure term
62     C o fZon: Zonal flux term
63     C o fMer: Meridional flux term
64     C o fVer: Vertical flux term - note fVer
65     C is "pipelined" in the vertical
66     C so we need an fVer for each
67     C variable.
68 cnh 1.38 C rhoK, rhoKM1 - Density at current level, level above and level
69     C below.
70 cnh 1.26 C rhoKP1
71     C buoyK, buoyKM1 - Buoyancy at current level and level above.
72 cnh 1.31 C phiHyd - Hydrostatic part of the potential phiHydi.
73 cnh 1.38 C In z coords phiHydiHyd is the hydrostatic
74     C pressure anomaly
75     C In p coords phiHydiHyd is the geopotential
76     C surface height
77 cnh 1.30 C anomaly.
78     C etaSurfX, - Holds surface elevation gradient in X and Y.
79     C etaSurfY
80     C K13, K23, K33 - Non-zero elements of small-angle approximation
81     C diffusion tensor.
82     C KapGM - Spatially varying Visbeck et. al mixing coeff.
83     C KappaRT, - Total diffusion in vertical for T and S.
84 cnh 1.38 C KappaRS (background + spatially varying, isopycnal term).
85 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
86     C jMin, jMax are applied.
87 cnh 1.1 C bi, bj
88 cnh 1.30 C k, kUp, - Index for layer above and below. kUp and kDown
89 cnh 1.38 C kDown, kM1 are switched with layer to be the appropriate
90     C index into fVerTerm.
91 cnh 1.30 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96     _RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
97     _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99     _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100     _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101     _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102     _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103     _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
107     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
108     _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
109     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
110 cnh 1.31 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111 cnh 1.30 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112     _RL rhokp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114     _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115     _RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116     _RL rhotmp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 cnh 1.29 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118     _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 cnh 1.31 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
120     _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
121     _RL K33 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
122 cnh 1.30 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123 cnh 1.31 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
124     _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
125 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
126     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
127 adcroft 1.12
128 adcroft 1.45 #ifdef INCLUDE_CONVECT_CALL
129     _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
130     #endif
131    
132 cnh 1.1 INTEGER iMin, iMax
133     INTEGER jMin, jMax
134     INTEGER bi, bj
135     INTEGER i, j
136     INTEGER k, kM1, kUp, kDown
137 cnh 1.19 LOGICAL BOTTOM_LAYER
138 cnh 1.1
139 adcroft 1.11 C--- The algorithm...
140     C
141     C "Correction Step"
142     C =================
143     C Here we update the horizontal velocities with the surface
144     C pressure such that the resulting flow is either consistent
145     C with the free-surface evolution or the rigid-lid:
146     C U[n] = U* + dt x d/dx P
147     C V[n] = V* + dt x d/dy P
148     C
149     C "Calculation of Gs"
150     C ===================
151     C This is where all the accelerations and tendencies (ie.
152 cnh 1.31 C phiHydysics, parameterizations etc...) are calculated
153 cnh 1.27 C rVel = sum_r ( div. u[n] )
154 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
155 cnh 1.27 C b = b(rho, theta)
156 adcroft 1.11 C K31 = K31 ( rho )
157 cnh 1.27 C Gu[n] = Gu( u[n], v[n], rVel, b, ... )
158     C Gv[n] = Gv( u[n], v[n], rVel, b, ... )
159     C Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
160     C Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
161 adcroft 1.11 C
162 adcroft 1.12 C "Time-stepping" or "Prediction"
163 adcroft 1.11 C ================================
164     C The models variables are stepped forward with the appropriate
165     C time-stepping scheme (currently we use Adams-Bashforth II)
166     C - For momentum, the result is always *only* a "prediction"
167     C in that the flow may be divergent and will be "corrected"
168     C later with a surface pressure gradient.
169     C - Normally for tracers the result is the new field at time
170     C level [n+1} *BUT* in the case of implicit diffusion the result
171     C is also *only* a prediction.
172     C - We denote "predictors" with an asterisk (*).
173     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
174     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
175     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
176     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
177 adcroft 1.12 C With implicit diffusion:
178 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
179     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
180 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
181     C (1 + dt * K * d_zz) salt[n] = salt*
182 adcroft 1.11 C---
183    
184 cnh 1.1 C-- Set up work arrays with valid (i.e. not NaN) values
185     C These inital values do not alter the numerical results. They
186     C just ensure that all memory references are to valid floating
187     C point numbers. This prevents spurious hardware signals due to
188     C uninitialised but inert locations.
189     DO j=1-OLy,sNy+OLy
190     DO i=1-OLx,sNx+OLx
191 adcroft 1.5 xA(i,j) = 0. _d 0
192     yA(i,j) = 0. _d 0
193     uTrans(i,j) = 0. _d 0
194     vTrans(i,j) = 0. _d 0
195     aTerm(i,j) = 0. _d 0
196     xTerm(i,j) = 0. _d 0
197     cTerm(i,j) = 0. _d 0
198     mTerm(i,j) = 0. _d 0
199     pTerm(i,j) = 0. _d 0
200     fZon(i,j) = 0. _d 0
201     fMer(i,j) = 0. _d 0
202 cnh 1.31 DO K=1,Nr
203     phiHyd (i,j,k) = 0. _d 0
204 cnh 1.30 K13(i,j,k) = 0. _d 0
205     K23(i,j,k) = 0. _d 0
206     K33(i,j,k) = 0. _d 0
207 adcroft 1.45 KappaRU(i,j,k) = 0. _d 0
208     KappaRV(i,j,k) = 0. _d 0
209 cnh 1.1 ENDDO
210 cnh 1.30 rhoKM1 (i,j) = 0. _d 0
211     rhok (i,j) = 0. _d 0
212     rhoKP1 (i,j) = 0. _d 0
213     rhoTMP (i,j) = 0. _d 0
214 cnh 1.26 buoyKM1(i,j) = 0. _d 0
215     buoyK (i,j) = 0. _d 0
216 cnh 1.30 maskC (i,j) = 0. _d 0
217 cnh 1.1 ENDDO
218     ENDDO
219    
220 cnh 1.35
221 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
222     DO bi=myBxLo(myThid),myBxHi(myThid)
223    
224 cnh 1.7 C-- Set up work arrays that need valid initial values
225     DO j=1-OLy,sNy+OLy
226     DO i=1-OLx,sNx+OLx
227 cnh 1.27 rTrans(i,j) = 0. _d 0
228     rVel (i,j,1) = 0. _d 0
229     rVel (i,j,2) = 0. _d 0
230 cnh 1.30 fVerT (i,j,1) = 0. _d 0
231     fVerT (i,j,2) = 0. _d 0
232     fVerS (i,j,1) = 0. _d 0
233     fVerS (i,j,2) = 0. _d 0
234     fVerU (i,j,1) = 0. _d 0
235     fVerU (i,j,2) = 0. _d 0
236     fVerV (i,j,1) = 0. _d 0
237     fVerV (i,j,2) = 0. _d 0
238 cnh 1.27 phiHyd(i,j,1) = 0. _d 0
239 cnh 1.30 K13 (i,j,1) = 0. _d 0
240     K23 (i,j,1) = 0. _d 0
241     K33 (i,j,1) = 0. _d 0
242     KapGM (i,j) = GMkbackground
243 cnh 1.7 ENDDO
244     ENDDO
245    
246 adcroft 1.45 DO k=1,Nr
247     DO j=1-OLy,sNy+OLy
248     DO i=1-OLx,sNx+OLx
249     #ifdef INCLUDE_CONVECT_CALL
250     ConvectCount(i,j,k) = 0.
251     #endif
252     KappaRT(i,j,k) = 0. _d 0
253     KappaRS(i,j,k) = 0. _d 0
254     ENDDO
255     ENDDO
256     ENDDO
257    
258 cnh 1.1 iMin = 1-OLx+1
259     iMax = sNx+OLx
260     jMin = 1-OLy+1
261     jMax = sNy+OLy
262 cnh 1.35
263 cnh 1.1
264 cnh 1.19 K = 1
265 cnh 1.31 BOTTOM_LAYER = K .EQ. Nr
266 cnh 1.19
267 cnh 1.38 #ifdef DO_PIPELINED_CORRECTION_STEP
268 adcroft 1.4 C-- Calculate gradient of surface pressure
269 cnh 1.28 CALL CALC_GRAD_ETA_SURF(
270 adcroft 1.4 I bi,bj,iMin,iMax,jMin,jMax,
271 cnh 1.29 O etaSurfX,etaSurfY,
272 adcroft 1.4 I myThid)
273     C-- Update fields in top level according to tendency terms
274 adcroft 1.11 CALL CORRECTION_STEP(
275 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K,
276     I etaSurfX,etaSurfY,myTime,myThid)
277 adcroft 1.43 #ifdef ALLOW_OBCS
278 adcroft 1.39 IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K, myThid )
279 adcroft 1.43 #endif
280 adcroft 1.21 IF ( .NOT. BOTTOM_LAYER ) THEN
281     C-- Update fields in layer below according to tendency terms
282     CALL CORRECTION_STEP(
283 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K+1,
284     I etaSurfX,etaSurfY,myTime,myThid)
285 adcroft 1.43 #ifdef ALLOW_OBCS
286 adcroft 1.39 IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
287 adcroft 1.43 #endif
288 adcroft 1.21 ENDIF
289 cnh 1.38 #endif
290 cnh 1.7 C-- Density of 1st level (below W(1)) reference to level 1
291 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
292 cnh 1.7 CALL FIND_RHO(
293 cnh 1.19 I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
294 cnh 1.7 O rhoKm1,
295     I myThid )
296 cnh 1.38 #endif
297 cnh 1.19
298 adcroft 1.42 IF ( (.NOT. BOTTOM_LAYER)
299     #ifdef ALLOW_KPP
300     & .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
301     #endif
302     & ) THEN
303 cnh 1.19 C-- Check static stability with layer below
304 cnh 1.30 C-- and mix as needed.
305 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
306 cnh 1.19 CALL FIND_RHO(
307     I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
308     O rhoKp1,
309     I myThid )
310 cnh 1.38 #endif
311     #ifdef INCLUDE_CONVECT_CALL
312 cnh 1.19 CALL CONVECT(
313     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
314 adcroft 1.45 U ConvectCount,
315 cnh 1.19 I myTime,myIter,myThid)
316 cnh 1.38 #endif
317 adcroft 1.45 C-- Implicit Vertical Diffusion for Convection
318     IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
319     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
320     U ConvectCount, KappaRT, KappaRS,
321     I myTime,myIter,myThid)
322 cnh 1.19 C-- Recompute density after mixing
323 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
324 cnh 1.19 CALL FIND_RHO(
325     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
326     O rhoKm1,
327     I myThid )
328 cnh 1.38 #endif
329 cnh 1.19 ENDIF
330 cnh 1.26 C-- Calculate buoyancy
331 cnh 1.32 CALL CALC_BUOYANCY(
332 cnh 1.26 I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
333     O buoyKm1,
334     I myThid )
335 cnh 1.38 C-- Integrate hydrostatic balance for phiHyd with BC of
336     C-- phiHyd(z=0)=0
337 cnh 1.26 CALL CALC_PHI_HYD(
338     I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
339     U phiHyd,
340 adcroft 1.5 I myThid )
341    
342 cnh 1.31 DO K=2,Nr
343     BOTTOM_LAYER = K .EQ. Nr
344 cnh 1.38 #ifdef DO_PIPELINED_CORRECTION_STEP
345 adcroft 1.21 IF ( .NOT. BOTTOM_LAYER ) THEN
346     C-- Update fields in layer below according to tendency terms
347     CALL CORRECTION_STEP(
348 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K+1,
349     I etaSurfX,etaSurfY,myTime,myThid)
350 adcroft 1.43 #ifdef ALLOW_OBCS
351 adcroft 1.39 IF (openBoundaries) CALL APPLY_OBCS1( bi, bj, K+1, myThid )
352 adcroft 1.43 #endif
353 adcroft 1.21 ENDIF
354 cnh 1.38 #endif
355 cnh 1.19 C-- Density of K level (below W(K)) reference to K level
356 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
357 cnh 1.19 CALL FIND_RHO(
358     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
359     O rhoK,
360     I myThid )
361 cnh 1.38 #endif
362 adcroft 1.42 IF ( (.NOT. BOTTOM_LAYER)
363     #ifdef ALLOW_KPP
364     & .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
365     #endif
366     & ) THEN
367 cnh 1.27 C-- Check static stability with layer below and mix as needed.
368     C-- Density of K+1 level (below W(K+1)) reference to K level.
369 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
370 cnh 1.19 CALL FIND_RHO(
371     I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
372     O rhoKp1,
373     I myThid )
374 cnh 1.38 #endif
375     #ifdef INCLUDE_CONVECT_CALL
376 cnh 1.19 CALL CONVECT(
377     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
378 adcroft 1.45 U ConvectCount,
379 cnh 1.19 I myTime,myIter,myThid)
380 cnh 1.38 #endif
381 adcroft 1.45 C-- Implicit Vertical Diffusion for Convection
382     IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
383     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
384     U ConvectCount, KappaRT, KappaRS,
385     I myTime,myIter,myThid)
386 cnh 1.19 C-- Recompute density after mixing
387 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
388 cnh 1.19 CALL FIND_RHO(
389     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
390     O rhoK,
391     I myThid )
392 cnh 1.38 #endif
393 cnh 1.19 ENDIF
394 cnh 1.26 C-- Calculate buoyancy
395 cnh 1.32 CALL CALC_BUOYANCY(
396 cnh 1.26 I bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
397     O buoyK,
398     I myThid )
399 cnh 1.38 C-- Integrate hydrostatic balance for phiHyd with BC of
400     C-- phiHyd(z=0)=0
401 cnh 1.26 CALL CALC_PHI_HYD(
402 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
403     U phiHyd,
404     I myThid )
405 cnh 1.19 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
406 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
407 cnh 1.19 CALL FIND_RHO(
408 cnh 1.30 I bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
409     O rhoTmp,
410     I myThid )
411 cnh 1.38 #endif
412     #ifdef INCLUDE_CALC_ISOSLOPES_CALL
413 cnh 1.19 CALL CALC_ISOSLOPES(
414 cnh 1.30 I bi, bj, iMin, iMax, jMin, jMax, K,
415     I rhoKm1, rhoK, rhotmp,
416     O K13, K23, K33, KapGM,
417     I myThid )
418 cnh 1.38 #endif
419 cnh 1.19 DO J=jMin,jMax
420     DO I=iMin,iMax
421 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
422 cnh 1.27 rhoKm1 (I,J) = rhoK(I,J)
423 cnh 1.38 #endif
424 cnh 1.27 buoyKm1(I,J) = buoyK(I,J)
425 cnh 1.19 ENDDO
426 adcroft 1.10 ENDDO
427 adcroft 1.11 ENDDO ! K
428    
429 adcroft 1.42 #ifdef ALLOW_KPP
430     C-- Compute KPP mixing coefficients
431     IF (usingKPPmixing) THEN
432     CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
433     I , myThid)
434     CALL KVMIX(
435     I bi, bj, myTime, myThid )
436     CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
437     I , myThid)
438     ENDIF
439     #endif
440    
441 cnh 1.31 DO K = Nr, 1, -1
442 cnh 1.30
443 cnh 1.1 kM1 =max(1,k-1) ! Points to level above k (=k-1)
444     kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
445     kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer
446     iMin = 1-OLx+2
447     iMax = sNx+OLx-1
448     jMin = 1-OLy+2
449     jMax = sNy+OLy-1
450    
451     C-- Get temporary terms used by tendency routines
452     CALL CALC_COMMON_FACTORS (
453     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
454 cnh 1.30 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
455 cnh 1.1 I myThid)
456 cnh 1.38 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
457 adcroft 1.12 C-- Calculate the total vertical diffusivity
458     CALL CALC_DIFFUSIVITY(
459     I bi,bj,iMin,iMax,jMin,jMax,K,
460     I maskC,maskUp,KapGM,K33,
461 adcroft 1.42 O KappaRT,KappaRS,KappaRU,KappaRV,
462 adcroft 1.12 I myThid)
463 cnh 1.38 #endif
464 cnh 1.1 C-- Calculate accelerations in the momentum equations
465 cnh 1.9 IF ( momStepping ) THEN
466     CALL CALC_MOM_RHS(
467     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
468 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
469 adcroft 1.42 I phiHyd,KappaRU,KappaRV,
470 cnh 1.9 U aTerm,xTerm,cTerm,mTerm,pTerm,
471     U fZon, fMer, fVerU, fVerV,
472 cnh 1.38 I myTime, myThid)
473 cnh 1.9 ENDIF
474 cnh 1.1 C-- Calculate active tracer tendencies
475 cnh 1.9 IF ( tempStepping ) THEN
476     CALL CALC_GT(
477     I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
478 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
479     I K13,K23,KappaRT,KapGM,
480 cnh 1.9 U aTerm,xTerm,fZon,fMer,fVerT,
481 cnh 1.37 I myTime, myThid)
482 cnh 1.9 ENDIF
483 adcroft 1.18 IF ( saltStepping ) THEN
484     CALL CALC_GS(
485     I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
486 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
487     I K13,K23,KappaRS,KapGM,
488 adcroft 1.18 U aTerm,xTerm,fZon,fMer,fVerS,
489 cnh 1.37 I myTime, myThid)
490 adcroft 1.18 ENDIF
491 adcroft 1.11 C-- Prediction step (step forward all model variables)
492     CALL TIMESTEP(
493     I bi,bj,iMin,iMax,jMin,jMax,K,
494     I myThid)
495 adcroft 1.43 #ifdef ALLOW_OBCS
496 adcroft 1.41 C-- Apply open boundary conditions
497 adcroft 1.39 IF (openBoundaries) CALL APPLY_OBCS2( bi, bj, K, myThid )
498 adcroft 1.43 #endif
499 adcroft 1.41 C-- Freeze water
500     IF (allowFreezing)
501     & CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
502 adcroft 1.11 C-- Diagnose barotropic divergence of predicted fields
503 cnh 1.31 CALL CALC_DIV_GHAT(
504 adcroft 1.11 I bi,bj,iMin,iMax,jMin,jMax,K,
505     I xA,yA,
506     I myThid)
507 adcroft 1.23
508     C-- Cumulative diagnostic calculations (ie. time-averaging)
509 cnh 1.38 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
510 adcroft 1.23 IF (taveFreq.GT.0.) THEN
511     CALL DO_TIME_AVERAGES(
512     I myTime, myIter, bi, bj, K, kUp, kDown,
513 adcroft 1.45 I K13, K23, rVel, KapGM, ConvectCount,
514 adcroft 1.23 I myThid )
515     ENDIF
516     #endif
517 adcroft 1.45
518 adcroft 1.11
519     ENDDO ! K
520 adcroft 1.12
521     C-- Implicit diffusion
522     IF (implicitDiffusion) THEN
523 adcroft 1.42 IF (tempStepping) CALL IMPLDIFF(
524     I bi, bj, iMin, iMax, jMin, jMax,
525     I deltaTtracer, KappaRT,recip_HFacC,
526     U gTNm1,
527     I myThid )
528     IF (saltStepping) CALL IMPLDIFF(
529     I bi, bj, iMin, iMax, jMin, jMax,
530     I deltaTtracer, KappaRS,recip_HFacC,
531     U gSNm1,
532     I myThid )
533 adcroft 1.44 ENDIF ! implicitDiffusion
534     C-- Implicit viscosity
535     IF (implicitViscosity) THEN
536 adcroft 1.42 IF (momStepping) THEN
537     CALL IMPLDIFF(
538     I bi, bj, iMin, iMax, jMin, jMax,
539     I deltaTmom, KappaRU,recip_HFacW,
540     U gUNm1,
541     I myThid )
542     CALL IMPLDIFF(
543     I bi, bj, iMin, iMax, jMin, jMax,
544     I deltaTmom, KappaRV,recip_HFacS,
545     U gVNm1,
546     I myThid )
547     #ifdef INCLUDE_CD_CODE
548     CALL IMPLDIFF(
549     I bi, bj, iMin, iMax, jMin, jMax,
550     I deltaTmom, KappaRU,recip_HFacW,
551     U vVelD,
552     I myThid )
553     CALL IMPLDIFF(
554     I bi, bj, iMin, iMax, jMin, jMax,
555     I deltaTmom, KappaRV,recip_HFacS,
556     U uVelD,
557     I myThid )
558     #endif
559     ENDIF ! momStepping
560 adcroft 1.44 ENDIF ! implicitViscosity
561 cnh 1.1
562     ENDDO
563     ENDDO
564 adcroft 1.6
565 cnh 1.19 C write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
566     C & maxval(cg2d_x(1:sNx,1:sNy,:,:))
567 cnh 1.20 C write(0,*) 'dynamics: U ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
568 adcroft 1.21 C & maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
569 cnh 1.20 C write(0,*) 'dynamics: V ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
570 adcroft 1.21 C & maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
571 cnh 1.30 C write(0,*) 'dynamics: rVel(1) ',
572     C & minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
573     C & maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
574     C write(0,*) 'dynamics: rVel(2) ',
575     C & minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
576     C & maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
577 adcroft 1.15 cblk write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
578     cblk & maxval(K13(1:sNx,1:sNy,:))
579     cblk write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
580     cblk & maxval(K23(1:sNx,1:sNy,:))
581     cblk write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
582     cblk & maxval(K33(1:sNx,1:sNy,:))
583 cnh 1.19 C write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
584     C & maxval(gT(1:sNx,1:sNy,:,:,:))
585     C write(0,*) 'dynamics: T ',minval(Theta(1:sNx,1:sNy,:,:,:)),
586     C & maxval(Theta(1:sNx,1:sNy,:,:,:))
587     C write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
588     C & maxval(gS(1:sNx,1:sNy,:,:,:))
589     C write(0,*) 'dynamics: S ',minval(salt(1:sNx,1:sNy,:,:,:)),
590     C & maxval(salt(1:sNx,1:sNy,:,:,:))
591 cnh 1.31 C write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
592     C & maxval(phiHyd/(Gravity*Rhonil))
593 cnh 1.36 C CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
594     C &Nr, 1, myThid )
595     C CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
596     C &Nr, 1, myThid )
597     C CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
598     C &Nr, 1, myThid )
599     C CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
600 cnh 1.38 C &Nr, 1, myThid )
601     C CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
602 cnh 1.36 C &Nr, 1, myThid )
603    
604 cnh 1.1
605     RETURN
606     END

  ViewVC Help
Powered by ViewVC 1.1.22