/[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.40 - (hide annotations) (download)
Wed Dec 9 16:11:51 1998 UTC (25 years, 5 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint19, checkpoint20
Changes since 1.39: +2 -1 lines
Added IMPLICIT NONE in a lot of subroutines.
Also corrected the recip_Rhonil bug: we didn't set it in ini_parms.F

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

  ViewVC Help
Powered by ViewVC 1.1.22