/[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.14 - (show annotations) (download)
Mon Jun 8 21:43:01 1998 UTC (25 years, 11 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint6
Changes since 1.13: +9 -5 lines
Merge of GM Redi and spherical polar and inplicit diffusion
and CD. Everything for a global run is now included, however,
still some discrepancies with GM Redi.

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

  ViewVC Help
Powered by ViewVC 1.1.22