/[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.28 - (show annotations) (download)
Thu Aug 20 19:25:05 1998 UTC (25 years, 9 months ago) by cnh
Branch: MAIN
Changes since 1.27: +2 -2 lines
Isomorphism consistency changes

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

  ViewVC Help
Powered by ViewVC 1.1.22