/[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.29 - (hide annotations) (download)
Thu Aug 20 19:26:40 1998 UTC (25 years, 9 months ago) by cnh
Branch: MAIN
Changes since 1.28: +8 -7 lines
Isomorphism consistency changes

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

  ViewVC Help
Powered by ViewVC 1.1.22