/[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.39 - (hide annotations) (download)
Tue Dec 8 19:44:28 1998 UTC (25 years, 5 months ago) by adcroft
Branch: MAIN
Changes since 1.38: +5 -1 lines
Implementation of Open Boundaries:
 o new source code: ini_obcs.F set_obcs.F apply_obcs1.F apply_obcs2.F
                    OBCS.h
 o modified code at a few points, key changes are in
    dynamcis.F the_model_main.F and ini_cg2d.F
 o documentation in OBCS.h and doc/OpenBound.*

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

  ViewVC Help
Powered by ViewVC 1.1.22