/[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.83.2.5 - (show annotations) (download)
Thu Jul 11 14:24:26 2002 UTC (21 years, 10 months ago) by heimbach
Branch: release1
CVS Tags: release1_p13_pre, release1_p13, release1_p8, release1_p9, release1_p5, release1_p6, release1_p7, release1_p11, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, release1_p12_pre
Branch point for: release1_50yr
Changes since 1.83.2.4: +11 -7 lines
o added Eliassen-Palm flux hook (dynamics)
o removed unused TAF keys iikey, idkey (dynamics,thermodynamics)

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.83.2.4 2002/04/17 01:38:00 heimbach Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: DYNAMICS
8 C !INTERFACE:
9 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
10 C !DESCRIPTION: \bv
11 C *==========================================================*
12 C | SUBROUTINE DYNAMICS
13 C | o Controlling routine for the explicit part of the model
14 C | dynamics.
15 C *==========================================================*
16 C | This routine evaluates the "dynamics" terms for each
17 C | block of ocean in turn. Because the blocks of ocean have
18 C | overlap regions they are independent of one another.
19 C | If terms involving lateral integrals are needed in this
20 C | routine care will be needed. Similarly finite-difference
21 C | operations with stencils wider than the overlap region
22 C | require special consideration.
23 C | The algorithm...
24 C |
25 C | "Correction Step"
26 C | =================
27 C | Here we update the horizontal velocities with the surface
28 C | pressure such that the resulting flow is either consistent
29 C | with the free-surface evolution or the rigid-lid:
30 C | U[n] = U* + dt x d/dx P
31 C | V[n] = V* + dt x d/dy P
32 C |
33 C | "Calculation of Gs"
34 C | ===================
35 C | This is where all the accelerations and tendencies (ie.
36 C | physics, parameterizations etc...) are calculated
37 C | rho = rho ( theta[n], salt[n] )
38 C | b = b(rho, theta)
39 C | K31 = K31 ( rho )
40 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
41 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
42 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
43 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
44 C |
45 C | "Time-stepping" or "Prediction"
46 C | ================================
47 C | The models variables are stepped forward with the appropriate
48 C | time-stepping scheme (currently we use Adams-Bashforth II)
49 C | - For momentum, the result is always *only* a "prediction"
50 C | in that the flow may be divergent and will be "corrected"
51 C | later with a surface pressure gradient.
52 C | - Normally for tracers the result is the new field at time
53 C | level [n+1} *BUT* in the case of implicit diffusion the result
54 C | is also *only* a prediction.
55 C | - We denote "predictors" with an asterisk (*).
56 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
57 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
58 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
59 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
60 C | With implicit diffusion:
61 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63 C | (1 + dt * K * d_zz) theta[n] = theta*
64 C | (1 + dt * K * d_zz) salt[n] = salt*
65 C |
66 C *==========================================================*
67 C \ev
68 C !USES:
69 IMPLICIT NONE
70 C == Global variables ===
71 #include "SIZE.h"
72 #include "EEPARAMS.h"
73 #include "PARAMS.h"
74 #include "DYNVARS.h"
75 #include "GRID.h"
76 #ifdef ALLOW_PASSIVE_TRACER
77 #include "TR1.h"
78 #endif
79 #ifdef ALLOW_AUTODIFF_TAMC
80 # include "tamc.h"
81 # include "tamc_keys.h"
82 # include "FFIELDS.h"
83 # ifdef ALLOW_KPP
84 # include "KPP.h"
85 # endif
86 # ifdef ALLOW_GMREDI
87 # include "GMREDI.h"
88 # endif
89 #endif /* ALLOW_AUTODIFF_TAMC */
90 #ifdef ALLOW_TIMEAVE
91 #include "TIMEAVE_STATV.h"
92 #endif
93
94 C !CALLING SEQUENCE:
95 C DYNAMICS()
96 C |
97 C |-- CALC_GRAD_PHI_SURF
98 C |
99 C |-- CALC_VISCOSITY
100 C |
101 C |-- CALC_PHI_HYD
102 C |
103 C |-- MOM_FLUXFORM
104 C |
105 C |-- MOM_VECINV
106 C |
107 C |-- TIMESTEP
108 C |
109 C |-- OBCS_APPLY_UV
110 C |
111 C |-- IMPLDIFF
112 C |
113 C |-- OBCS_APPLY_UV
114 C |
115 C |-- CALL TIMEAVE_CUMUL_1T
116 C |-- CALL DEBUG_STATS_RL
117
118 C !INPUT/OUTPUT PARAMETERS:
119 C == Routine arguments ==
120 C myTime - Current time in simulation
121 C myIter - Current iteration number in simulation
122 C myThid - Thread number for this instance of the routine.
123 _RL myTime
124 INTEGER myIter
125 INTEGER myThid
126
127 C !LOCAL VARIABLES:
128 C == Local variables
129 C fVer[STUV] o fVer: Vertical flux term - note fVer
130 C is "pipelined" in the vertical
131 C so we need an fVer for each
132 C variable.
133 C rhoK, rhoKM1 - Density at current level, and level above
134 C phiHyd - Hydrostatic part of the potential phiHydi.
135 C In z coords phiHydiHyd is the hydrostatic
136 C Potential (=pressure/rho0) anomaly
137 C In p coords phiHydiHyd is the geopotential
138 C surface height anomaly.
139 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
140 C phiSurfY or geopotentiel (atmos) in X and Y direction
141 C iMin, iMax - Ranges and sub-block indices on which calculations
142 C jMin, jMax are applied.
143 C bi, bj
144 C k, kup, - Index for layer above and below. kup and kDown
145 C kDown, km1 are switched with layer to be the appropriate
146 C index into fVerTerm.
147 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
148 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
150 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
155 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156
157 INTEGER iMin, iMax
158 INTEGER jMin, jMax
159 INTEGER bi, bj
160 INTEGER i, j
161 INTEGER k, km1, kp1, kup, kDown
162
163 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
164 c CHARACTER*(MAX_LEN_MBUF) suff
165 c LOGICAL DIFFERENT_MULTIPLE
166 c EXTERNAL DIFFERENT_MULTIPLE
167 Cjmc(end)
168
169 C--- The algorithm...
170 C
171 C "Correction Step"
172 C =================
173 C Here we update the horizontal velocities with the surface
174 C pressure such that the resulting flow is either consistent
175 C with the free-surface evolution or the rigid-lid:
176 C U[n] = U* + dt x d/dx P
177 C V[n] = V* + dt x d/dy P
178 C
179 C "Calculation of Gs"
180 C ===================
181 C This is where all the accelerations and tendencies (ie.
182 C physics, parameterizations etc...) are calculated
183 C rho = rho ( theta[n], salt[n] )
184 C b = b(rho, theta)
185 C K31 = K31 ( rho )
186 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
187 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
188 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
189 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
190 C
191 C "Time-stepping" or "Prediction"
192 C ================================
193 C The models variables are stepped forward with the appropriate
194 C time-stepping scheme (currently we use Adams-Bashforth II)
195 C - For momentum, the result is always *only* a "prediction"
196 C in that the flow may be divergent and will be "corrected"
197 C later with a surface pressure gradient.
198 C - Normally for tracers the result is the new field at time
199 C level [n+1} *BUT* in the case of implicit diffusion the result
200 C is also *only* a prediction.
201 C - We denote "predictors" with an asterisk (*).
202 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
203 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
204 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
205 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
206 C With implicit diffusion:
207 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
208 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209 C (1 + dt * K * d_zz) theta[n] = theta*
210 C (1 + dt * K * d_zz) salt[n] = salt*
211 C---
212 CEOP
213
214 C-- Set up work arrays with valid (i.e. not NaN) values
215 C These inital values do not alter the numerical results. They
216 C just ensure that all memory references are to valid floating
217 C point numbers. This prevents spurious hardware signals due to
218 C uninitialised but inert locations.
219 DO j=1-OLy,sNy+OLy
220 DO i=1-OLx,sNx+OLx
221 rhoKM1 (i,j) = 0. _d 0
222 rhok (i,j) = 0. _d 0
223 phiSurfX(i,j) = 0. _d 0
224 phiSurfY(i,j) = 0. _d 0
225 ENDDO
226 ENDDO
227
228 C-- Call to routine for calculation of
229 C Eliassen-Palm-flux-forced U-tendency,
230 C if desired:
231 #ifdef INCLUDE_EP_FORCING_CODE
232 CALL CALC_EP_FORCING(myThid)
233 #endif
234
235 #ifdef ALLOW_AUTODIFF_TAMC
236 C-- HPF directive to help TAMC
237 CHPF$ INDEPENDENT
238 #endif /* ALLOW_AUTODIFF_TAMC */
239
240 DO bj=myByLo(myThid),myByHi(myThid)
241
242 #ifdef ALLOW_AUTODIFF_TAMC
243 C-- HPF directive to help TAMC
244 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
245 CHPF$& ,phiHyd
246 CHPF$& ,KappaRU,KappaRV
247 CHPF$& )
248 #endif /* ALLOW_AUTODIFF_TAMC */
249
250 DO bi=myBxLo(myThid),myBxHi(myThid)
251
252 #ifdef ALLOW_AUTODIFF_TAMC
253 act1 = bi - myBxLo(myThid)
254 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
255 act2 = bj - myByLo(myThid)
256 max2 = myByHi(myThid) - myByLo(myThid) + 1
257 act3 = myThid - 1
258 max3 = nTx*nTy
259 act4 = ikey_dynamics - 1
260 ikey = (act1 + 1) + act2*max1
261 & + act3*max1*max2
262 & + act4*max1*max2*max3
263 #endif /* ALLOW_AUTODIFF_TAMC */
264
265 C-- Set up work arrays that need valid initial values
266 DO j=1-OLy,sNy+OLy
267 DO i=1-OLx,sNx+OLx
268 DO k=1,Nr
269 phiHyd(i,j,k) = 0. _d 0
270 KappaRU(i,j,k) = 0. _d 0
271 KappaRV(i,j,k) = 0. _d 0
272 ENDDO
273 fVerU (i,j,1) = 0. _d 0
274 fVerU (i,j,2) = 0. _d 0
275 fVerV (i,j,1) = 0. _d 0
276 fVerV (i,j,2) = 0. _d 0
277 ENDDO
278 ENDDO
279
280 C-- Start computation of dynamics
281 iMin = 1-OLx+2
282 iMax = sNx+OLx-1
283 jMin = 1-OLy+2
284 jMax = sNy+OLy-1
285
286 #ifdef ALLOW_AUTODIFF_TAMC
287 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
288 #endif /* ALLOW_AUTODIFF_TAMC */
289
290 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
291 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
292 IF (implicSurfPress.NE.1.) THEN
293 CALL CALC_GRAD_PHI_SURF(
294 I bi,bj,iMin,iMax,jMin,jMax,
295 I etaN,
296 O phiSurfX,phiSurfY,
297 I myThid )
298 ENDIF
299
300 #ifdef ALLOW_AUTODIFF_TAMC
301 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
302 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
303 #ifdef ALLOW_KPP
304 CADJ STORE KPPviscAz (:,:,:,bi,bj)
305 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
306 #endif /* ALLOW_KPP */
307 #endif /* ALLOW_AUTODIFF_TAMC */
308
309 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
310 C-- Calculate the total vertical diffusivity
311 DO k=1,Nr
312 CALL CALC_VISCOSITY(
313 I bi,bj,iMin,iMax,jMin,jMax,k,
314 O KappaRU,KappaRV,
315 I myThid)
316 ENDDO
317 #endif
318
319 C-- Start of dynamics loop
320 DO k=1,Nr
321
322 C-- km1 Points to level above k (=k-1)
323 C-- kup Cycles through 1,2 to point to layer above
324 C-- kDown Cycles through 2,1 to point to current layer
325
326 km1 = MAX(1,k-1)
327 kp1 = MIN(k+1,Nr)
328 kup = 1+MOD(k+1,2)
329 kDown= 1+MOD(k,2)
330
331 #ifdef ALLOW_AUTODIFF_TAMC
332 kkey = (ikey-1)*Nr + k
333 #endif /* ALLOW_AUTODIFF_TAMC */
334
335 C-- Integrate hydrostatic balance for phiHyd with BC of
336 C phiHyd(z=0)=0
337 C distinguishe between Stagger and Non Stagger time stepping
338 IF (staggerTimeStep) THEN
339 CALL CALC_PHI_HYD(
340 I bi,bj,iMin,iMax,jMin,jMax,k,
341 I gT, gS,
342 U phiHyd,
343 I myThid )
344 ELSE
345 CALL CALC_PHI_HYD(
346 I bi,bj,iMin,iMax,jMin,jMax,k,
347 I theta, salt,
348 U phiHyd,
349 I myThid )
350 ENDIF
351
352 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
353 C and step forward storing the result in gUnm1, gVnm1, etc...
354 IF ( momStepping ) THEN
355 #ifndef DISABLE_MOM_FLUXFORM
356 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
357 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
358 I phiHyd,KappaRU,KappaRV,
359 U fVerU, fVerV,
360 I myTime, myIter, myThid)
361 #endif
362 #ifndef DISABLE_MOM_VECINV
363 IF (vectorInvariantMomentum) CALL MOM_VECINV(
364 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
365 I phiHyd,KappaRU,KappaRV,
366 U fVerU, fVerV,
367 I myTime, myIter, myThid)
368 #endif
369 CALL TIMESTEP(
370 I bi,bj,iMin,iMax,jMin,jMax,k,
371 I phiHyd, phiSurfX, phiSurfY,
372 I myIter, myThid)
373
374 #ifdef ALLOW_OBCS
375 C-- Apply open boundary conditions
376 IF (useOBCS) THEN
377 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
378 END IF
379 #endif /* ALLOW_OBCS */
380
381 #ifdef ALLOW_AUTODIFF_TAMC
382 #ifdef INCLUDE_CD_CODE
383 ELSE
384 DO j=1-OLy,sNy+OLy
385 DO i=1-OLx,sNx+OLx
386 guCD(i,j,k,bi,bj) = 0.0
387 gvCD(i,j,k,bi,bj) = 0.0
388 END DO
389 END DO
390 #endif /* INCLUDE_CD_CODE */
391 #endif /* ALLOW_AUTODIFF_TAMC */
392 ENDIF
393
394
395 C-- end of dynamics k loop (1:Nr)
396 ENDDO
397
398 C-- Implicit viscosity
399 IF (implicitViscosity.AND.momStepping) THEN
400 #ifdef ALLOW_AUTODIFF_TAMC
401 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
402 #endif /* ALLOW_AUTODIFF_TAMC */
403 CALL IMPLDIFF(
404 I bi, bj, iMin, iMax, jMin, jMax,
405 I deltaTmom, KappaRU,recip_HFacW,
406 U gUNm1,
407 I myThid )
408 #ifdef ALLOW_AUTODIFF_TAMC
409 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
410 #endif /* ALLOW_AUTODIFF_TAMC */
411 CALL IMPLDIFF(
412 I bi, bj, iMin, iMax, jMin, jMax,
413 I deltaTmom, KappaRV,recip_HFacS,
414 U gVNm1,
415 I myThid )
416
417 #ifdef ALLOW_OBCS
418 C-- Apply open boundary conditions
419 IF (useOBCS) THEN
420 DO K=1,Nr
421 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
422 ENDDO
423 END IF
424 #endif /* ALLOW_OBCS */
425
426 #ifdef INCLUDE_CD_CODE
427 #ifdef ALLOW_AUTODIFF_TAMC
428 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
429 #endif /* ALLOW_AUTODIFF_TAMC */
430 CALL IMPLDIFF(
431 I bi, bj, iMin, iMax, jMin, jMax,
432 I deltaTmom, KappaRU,recip_HFacW,
433 U vVelD,
434 I myThid )
435 #ifdef ALLOW_AUTODIFF_TAMC
436 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
437 #endif /* ALLOW_AUTODIFF_TAMC */
438 CALL IMPLDIFF(
439 I bi, bj, iMin, iMax, jMin, jMax,
440 I deltaTmom, KappaRV,recip_HFacS,
441 U uVelD,
442 I myThid )
443 #endif /* INCLUDE_CD_CODE */
444 C-- End If implicitViscosity.AND.momStepping
445 ENDIF
446
447 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
448 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
449 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
450 c WRITE(suff,'(I10.10)') myIter+1
451 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
452 c ENDIF
453 Cjmc(end)
454
455 #ifdef ALLOW_TIMEAVE
456 IF (taveFreq.GT.0.) THEN
457 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
458 I deltaTclock, bi, bj, myThid)
459 ENDIF
460 #endif /* ALLOW_TIMEAVE */
461
462 ENDDO
463 ENDDO
464
465 #ifndef DISABLE_DEBUGMODE
466 If (debugMode) THEN
467 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
468 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
469 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
470 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
471 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
472 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
473 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
474 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
475 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
476 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
477 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
478 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
479 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
480 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
481 ENDIF
482 #endif
483
484 RETURN
485 END

  ViewVC Help
Powered by ViewVC 1.1.22