/[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.37 - (hide annotations) (download)
Tue Nov 3 15:28:04 1998 UTC (25 years, 6 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint17
Changes since 1.36: +3 -3 lines
Partial changes to incoporate atmospheric configuration
Minor TAMC compliance changes
Included one-layer verification experiment exp0

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

  ViewVC Help
Powered by ViewVC 1.1.22