/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.30 - (show annotations) (download)
Thu Aug 20 20:05:01 1998 UTC (25 years, 8 months ago) by cnh
Branch: MAIN
Changes since 1.29: +103 -110 lines
Isomorphism consistency changes

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.29 1998/08/20 19:26:40 cnh Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6 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 #include "PARAMS.h"
29 #include "DYNVARS.h"
30
31 C == Routine arguments ==
32 C myTime - Current time in simulation
33 C myIter - Current iteration number in simulation
34 C myThid - Thread number for this instance of the routine.
35 INTEGER myThid
36 _RL myTime
37 INTEGER myIter
38
39 C == Local variables
40 C xA, yA - Per block temporaries holding face areas
41 C uTrans, vTrans, rTrans - Per block temporaries holding flow transport
42 C rVel o uTrans: Zonal transport
43 C o vTrans: Meridional transport
44 C o rTrans: Vertical transport
45 C o rVel: Vertical velocity at upper and lower
46 C cell faces.
47 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 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
68 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 C bi, bj
79 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 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
102 _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 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
111 _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
112 _RL K33 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
113 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114 _RL KappaZT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
115 _RL KappaZS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
116
117 INTEGER iMin, iMax
118 INTEGER jMin, jMax
119 INTEGER bi, bj
120 INTEGER i, j
121 INTEGER k, kM1, kUp, kDown
122 LOGICAL BOTTOM_LAYER
123
124 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 C physics, parameterizations etc...) are calculated
138 C rVel = sum_r ( div. u[n] )
139 C rho = rho ( theta[n], salt[n] )
140 C b = b(rho, theta)
141 C K31 = K31 ( rho )
142 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 C
147 C "Time-stepping" or "Prediction"
148 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 C With implicit diffusion:
163 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 C (1 + dt * K * d_zz) theta[n] = theta*
166 C (1 + dt * K * d_zz) salt[n] = salt*
167 C---
168
169 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 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 DO K=1,nZ
188 pH (i,j,k) = 0. _d 0
189 K13(i,j,k) = 0. _d 0
190 K23(i,j,k) = 0. _d 0
191 K33(i,j,k) = 0. _d 0
192 KappaZT(i,j,k) = 0. _d 0
193 ENDDO
194 rhoKM1 (i,j) = 0. _d 0
195 rhok (i,j) = 0. _d 0
196 rhoKP1 (i,j) = 0. _d 0
197 rhoTMP (i,j) = 0. _d 0
198 buoyKM1(i,j) = 0. _d 0
199 buoyK (i,j) = 0. _d 0
200 maskC (i,j) = 0. _d 0
201 ENDDO
202 ENDDO
203
204 DO bj=myByLo(myThid),myByHi(myThid)
205 DO bi=myBxLo(myThid),myBxHi(myThid)
206
207 C-- Set up work arrays that need valid initial values
208 DO j=1-OLy,sNy+OLy
209 DO i=1-OLx,sNx+OLx
210 rTrans(i,j) = 0. _d 0
211 rVel (i,j,1) = 0. _d 0
212 rVel (i,j,2) = 0. _d 0
213 fVerT (i,j,1) = 0. _d 0
214 fVerT (i,j,2) = 0. _d 0
215 fVerS (i,j,1) = 0. _d 0
216 fVerS (i,j,2) = 0. _d 0
217 fVerU (i,j,1) = 0. _d 0
218 fVerU (i,j,2) = 0. _d 0
219 fVerV (i,j,1) = 0. _d 0
220 fVerV (i,j,2) = 0. _d 0
221 phiHyd(i,j,1) = 0. _d 0
222 K13 (i,j,1) = 0. _d 0
223 K23 (i,j,1) = 0. _d 0
224 K33 (i,j,1) = 0. _d 0
225 KapGM (i,j) = GMkbackground
226 ENDDO
227 ENDDO
228
229 iMin = 1-OLx+1
230 iMax = sNx+OLx
231 jMin = 1-OLy+1
232 jMax = sNy+OLy
233
234 K = 1
235 BOTTOM_LAYER = K .EQ. Nz
236
237 C-- Calculate gradient of surface pressure
238 CALL CALC_GRAD_ETA_SURF(
239 I bi,bj,iMin,iMax,jMin,jMax,
240 O etaSurfX,etaSurfY,
241 I myThid)
242 C-- Update fields in top level according to tendency terms
243 CALL CORRECTION_STEP(
244 I bi,bj,iMin,iMax,jMin,jMax,K,
245 I etaSurfX,etaSurfY,myTime,myThid)
246 IF ( .NOT. BOTTOM_LAYER ) THEN
247 C-- Update fields in layer below according to tendency terms
248 CALL CORRECTION_STEP(
249 I bi,bj,iMin,iMax,jMin,jMax,K+1,
250 I etaSurfX,etaSurfY,myTime,myThid)
251 ENDIF
252 C-- Density of 1st level (below W(1)) reference to level 1
253 CALL FIND_RHO(
254 I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
255 O rhoKm1,
256 I myThid )
257
258 IF ( .NOT. BOTTOM_LAYER ) THEN
259 C-- Check static stability with layer below
260 C-- and mix as needed.
261 CALL FIND_RHO(
262 I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
263 O rhoKp1,
264 I myThid )
265 CALL CONVECT(
266 I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
267 I myTime,myIter,myThid)
268 C-- Recompute density after mixing
269 CALL FIND_RHO(
270 I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
271 O rhoKm1,
272 I myThid )
273 ENDIF
274 C-- Calculate buoyancy
275 CALL CALC_BUOY(
276 I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
277 O buoyKm1,
278 I myThid )
279 C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0
280 CALL CALC_PHI_HYD(
281 I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
282 U phiHyd,
283 I myThid )
284
285 DO K=2,Nz
286 BOTTOM_LAYER = K .EQ. Nz
287 IF ( .NOT. BOTTOM_LAYER ) THEN
288 C-- Update fields in layer below according to tendency terms
289 CALL CORRECTION_STEP(
290 I bi,bj,iMin,iMax,jMin,jMax,K+1,
291 I etaSurfX,etaSurfY,myTime,myThid)
292 ENDIF
293 C-- Density of K level (below W(K)) reference to K level
294 CALL FIND_RHO(
295 I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
296 O rhoK,
297 I myThid )
298 IF ( .NOT. BOTTOM_LAYER ) THEN
299 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 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 C-- Calculate buoyancy
315 CALL CALC_BUOY(
316 I bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
317 O buoyK,
318 I myThid )
319 C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0
320 CALL CALC_PHI_HYD(
321 I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
322 U phiHyd,
323 I myThid )
324 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
325 CALL FIND_RHO(
326 I bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
327 O rhoTmp,
328 I myThid )
329 CALL CALC_ISOSLOPES(
330 I bi, bj, iMin, iMax, jMin, jMax, K,
331 I rhoKm1, rhoK, rhotmp,
332 O K13, K23, K33, KapGM,
333 I myThid )
334 DO J=jMin,jMax
335 DO I=iMin,iMax
336 rhoKm1 (I,J) = rhoK(I,J)
337 buoyKm1(I,J) = buoyK(I,J)
338 ENDDO
339 ENDDO
340 ENDDO ! K
341
342 DO K = Nz, 1, -1
343
344 kM1 =max(1,k-1) ! Points to level above k (=k-1)
345 kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
346 kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer
347 iMin = 1-OLx+2
348 iMax = sNx+OLx-1
349 jMin = 1-OLy+2
350 jMax = sNy+OLy-1
351
352 C-- Get temporary terms used by tendency routines
353 CALL CALC_COMMON_FACTORS (
354 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
355 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
356 I myThid)
357 C-- Calculate the total vertical diffusivity
358 CALL CALC_DIFFUSIVITY(
359 I bi,bj,iMin,iMax,jMin,jMax,K,
360 I maskC,maskUp,KapGM,K33,
361 O KappaZT,KappaZS,
362 I myThid)
363 C-- Calculate accelerations in the momentum equations
364 IF ( momStepping ) THEN
365 CALL CALC_MOM_RHS(
366 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
367 I xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
368 I phiHyd,
369 U aTerm,xTerm,cTerm,mTerm,pTerm,
370 U fZon, fMer, fVerU, fVerV,
371 I myThid)
372 ENDIF
373 C-- Calculate active tracer tendencies
374 IF ( tempStepping ) THEN
375 CALL CALC_GT(
376 I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
377 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
378 I K13,K23,KappaRT,KapGM,
379 U aTerm,xTerm,fZon,fMer,fVerT,
380 I myThid)
381 ENDIF
382 IF ( saltStepping ) THEN
383 CALL CALC_GS(
384 I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
385 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
386 I K13,K23,KappaRS,KapGM,
387 U aTerm,xTerm,fZon,fMer,fVerS,
388 I myThid)
389 ENDIF
390 C-- Prediction step (step forward all model variables)
391 CALL TIMESTEP(
392 I bi,bj,iMin,iMax,jMin,jMax,K,
393 I myThid)
394 C-- Diagnose barotropic divergence of predicted fields
395 CALL CALC_DIV_G(
396 I bi,bj,iMin,iMax,jMin,jMax,K,
397 I xA,yA,
398 I myThid)
399
400 C-- Cumulative diagnostic calculations (ie. time-averaging)
401 #ifdef ALLOW_DIAGNOSTICS
402 IF (taveFreq.GT.0.) THEN
403 CALL DO_TIME_AVERAGES(
404 I myTime, myIter, bi, bj, K, kUp, kDown,
405 I K13, K23, rVel, KapGM,
406 I myThid )
407 ENDIF
408 #endif
409
410 ENDDO ! K
411
412 C-- Implicit diffusion
413 IF (implicitDiffusion) THEN
414 CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
415 I KappaZT,KappaZS,
416 I myThid )
417 ENDIF
418
419 ENDDO
420 ENDDO
421
422 C write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
423 C & maxval(cg2d_x(1:sNx,1:sNy,:,:))
424 C write(0,*) 'dynamics: U ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
425 C & maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
426 C write(0,*) 'dynamics: V ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
427 C & maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
428 C write(0,*) 'dynamics: rVel(1) ',
429 C & minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
430 C & maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
431 C write(0,*) 'dynamics: rVel(2) ',
432 C & minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
433 C & maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
434 cblk write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
435 cblk & maxval(K13(1:sNx,1:sNy,:))
436 cblk write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
437 cblk & maxval(K23(1:sNx,1:sNy,:))
438 cblk write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
439 cblk & maxval(K33(1:sNx,1:sNy,:))
440 C write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
441 C & maxval(gT(1:sNx,1:sNy,:,:,:))
442 C write(0,*) 'dynamics: T ',minval(Theta(1:sNx,1:sNy,:,:,:)),
443 C & maxval(Theta(1:sNx,1:sNy,:,:,:))
444 C write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
445 C & maxval(gS(1:sNx,1:sNy,:,:,:))
446 C write(0,*) 'dynamics: S ',minval(salt(1:sNx,1:sNy,:,:,:)),
447 C & maxval(salt(1:sNx,1:sNy,:,:,:))
448 C write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.),
449 C & maxval(pH/(Gravity*Rhonil))
450
451 RETURN
452 END

  ViewVC Help
Powered by ViewVC 1.1.22