/[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.26 - (hide annotations) (download)
Wed Aug 19 16:20:49 1998 UTC (25 years, 9 months ago) by cnh
Branch: MAIN
Changes since 1.25: +35 -16 lines
Changes to support r as vertical coordinate

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

  ViewVC Help
Powered by ViewVC 1.1.22