/[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.85 - (show annotations) (download)
Fri Feb 15 21:29:04 2002 UTC (22 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint44e_post, checkpoint44f_post, chkpt44d_post, checkpoint44e_pre, checkpoint44h_pre, ecco_c44_e18, ecco_c44_e17, checkpoint44g_post, checkpoint44f_pre
Branch point for: release1_final
Changes since 1.84: +1 -7 lines
Removed arrays which are no longer needed.

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.84 2001/11/16 03:25:40 jmc 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 DO k=1,Nr
222 phiHyd(i,j,k) = 0. _d 0
223 KappaRU(i,j,k) = 0. _d 0
224 KappaRV(i,j,k) = 0. _d 0
225 ENDDO
226 rhoKM1 (i,j) = 0. _d 0
227 rhok (i,j) = 0. _d 0
228 phiSurfX(i,j) = 0. _d 0
229 phiSurfY(i,j) = 0. _d 0
230 ENDDO
231 ENDDO
232
233 #ifdef ALLOW_AUTODIFF_TAMC
234 C-- HPF directive to help TAMC
235 CHPF$ INDEPENDENT
236 #endif /* ALLOW_AUTODIFF_TAMC */
237
238 DO bj=myByLo(myThid),myByHi(myThid)
239
240 #ifdef ALLOW_AUTODIFF_TAMC
241 C-- HPF directive to help TAMC
242 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
243 CHPF$& ,phiHyd
244 CHPF$& ,KappaRU,KappaRV
245 CHPF$& )
246 #endif /* ALLOW_AUTODIFF_TAMC */
247
248 DO bi=myBxLo(myThid),myBxHi(myThid)
249
250 #ifdef ALLOW_AUTODIFF_TAMC
251 act1 = bi - myBxLo(myThid)
252 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
253 act2 = bj - myByLo(myThid)
254 max2 = myByHi(myThid) - myByLo(myThid) + 1
255 act3 = myThid - 1
256 max3 = nTx*nTy
257 act4 = ikey_dynamics - 1
258 ikey = (act1 + 1) + act2*max1
259 & + act3*max1*max2
260 & + act4*max1*max2*max3
261 #endif /* ALLOW_AUTODIFF_TAMC */
262
263 C-- Set up work arrays that need valid initial values
264 DO j=1-OLy,sNy+OLy
265 DO i=1-OLx,sNx+OLx
266 fVerU (i,j,1) = 0. _d 0
267 fVerU (i,j,2) = 0. _d 0
268 fVerV (i,j,1) = 0. _d 0
269 fVerV (i,j,2) = 0. _d 0
270 ENDDO
271 ENDDO
272
273 C-- Start computation of dynamics
274 iMin = 1-OLx+2
275 iMax = sNx+OLx-1
276 jMin = 1-OLy+2
277 jMax = sNy+OLy-1
278
279 #ifdef ALLOW_AUTODIFF_TAMC
280 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
281 #endif /* ALLOW_AUTODIFF_TAMC */
282
283 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
284 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
285 IF (implicSurfPress.NE.1.) THEN
286 CALL CALC_GRAD_PHI_SURF(
287 I bi,bj,iMin,iMax,jMin,jMax,
288 I etaN,
289 O phiSurfX,phiSurfY,
290 I myThid )
291 ENDIF
292
293 #ifdef ALLOW_AUTODIFF_TAMC
294 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
295 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
296 #ifdef ALLOW_KPP
297 CADJ STORE KPPviscAz (:,:,:,bi,bj)
298 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
299 #endif /* ALLOW_KPP */
300 #endif /* ALLOW_AUTODIFF_TAMC */
301
302 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
303 C-- Calculate the total vertical diffusivity
304 DO k=1,Nr
305 CALL CALC_VISCOSITY(
306 I bi,bj,iMin,iMax,jMin,jMax,k,
307 O KappaRU,KappaRV,
308 I myThid)
309 ENDDO
310 #endif
311
312 C-- Start of dynamics loop
313 DO k=1,Nr
314
315 C-- km1 Points to level above k (=k-1)
316 C-- kup Cycles through 1,2 to point to layer above
317 C-- kDown Cycles through 2,1 to point to current layer
318
319 km1 = MAX(1,k-1)
320 kp1 = MIN(k+1,Nr)
321 kup = 1+MOD(k+1,2)
322 kDown= 1+MOD(k,2)
323
324 #ifdef ALLOW_AUTODIFF_TAMC
325 kkey = (ikey-1)*Nr + k
326 #endif /* ALLOW_AUTODIFF_TAMC */
327
328 C-- Integrate hydrostatic balance for phiHyd with BC of
329 C phiHyd(z=0)=0
330 C distinguishe between Stagger and Non Stagger time stepping
331 IF (staggerTimeStep) THEN
332 CALL CALC_PHI_HYD(
333 I bi,bj,iMin,iMax,jMin,jMax,k,
334 I gT, gS,
335 U phiHyd,
336 I myThid )
337 ELSE
338 CALL CALC_PHI_HYD(
339 I bi,bj,iMin,iMax,jMin,jMax,k,
340 I theta, salt,
341 U phiHyd,
342 I myThid )
343 ENDIF
344
345 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
346 C and step forward storing the result in gUnm1, gVnm1, etc...
347 IF ( momStepping ) THEN
348 #ifndef DISABLE_MOM_FLUXFORM
349 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
350 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
351 I phiHyd,KappaRU,KappaRV,
352 U fVerU, fVerV,
353 I myTime, myIter, myThid)
354 #endif
355 #ifndef DISABLE_MOM_VECINV
356 IF (vectorInvariantMomentum) CALL MOM_VECINV(
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 CALL TIMESTEP(
363 I bi,bj,iMin,iMax,jMin,jMax,k,
364 I phiHyd, phiSurfX, phiSurfY,
365 I myIter, myThid)
366
367 #ifdef ALLOW_OBCS
368 C-- Apply open boundary conditions
369 IF (useOBCS) THEN
370 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
371 END IF
372 #endif /* ALLOW_OBCS */
373
374 #ifdef ALLOW_AUTODIFF_TAMC
375 #ifdef INCLUDE_CD_CODE
376 ELSE
377 DO j=1-OLy,sNy+OLy
378 DO i=1-OLx,sNx+OLx
379 guCD(i,j,k,bi,bj) = 0.0
380 gvCD(i,j,k,bi,bj) = 0.0
381 END DO
382 END DO
383 #endif /* INCLUDE_CD_CODE */
384 #endif /* ALLOW_AUTODIFF_TAMC */
385 ENDIF
386
387
388 C-- end of dynamics k loop (1:Nr)
389 ENDDO
390
391
392
393 C-- Implicit viscosity
394 IF (implicitViscosity.AND.momStepping) THEN
395 #ifdef ALLOW_AUTODIFF_TAMC
396 idkey = iikey + 3
397 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
398 #endif /* ALLOW_AUTODIFF_TAMC */
399 CALL IMPLDIFF(
400 I bi, bj, iMin, iMax, jMin, jMax,
401 I deltaTmom, KappaRU,recip_HFacW,
402 U gUNm1,
403 I myThid )
404 #ifdef ALLOW_AUTODIFF_TAMC
405 idkey = iikey + 4
406 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
407 #endif /* ALLOW_AUTODIFF_TAMC */
408 CALL IMPLDIFF(
409 I bi, bj, iMin, iMax, jMin, jMax,
410 I deltaTmom, KappaRV,recip_HFacS,
411 U gVNm1,
412 I myThid )
413
414 #ifdef ALLOW_OBCS
415 C-- Apply open boundary conditions
416 IF (useOBCS) THEN
417 DO K=1,Nr
418 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
419 ENDDO
420 END IF
421 #endif /* ALLOW_OBCS */
422
423 #ifdef INCLUDE_CD_CODE
424 #ifdef ALLOW_AUTODIFF_TAMC
425 idkey = iikey + 5
426 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
427 #endif /* ALLOW_AUTODIFF_TAMC */
428 CALL IMPLDIFF(
429 I bi, bj, iMin, iMax, jMin, jMax,
430 I deltaTmom, KappaRU,recip_HFacW,
431 U vVelD,
432 I myThid )
433 #ifdef ALLOW_AUTODIFF_TAMC
434 idkey = iikey + 6
435 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
436 #endif /* ALLOW_AUTODIFF_TAMC */
437 CALL IMPLDIFF(
438 I bi, bj, iMin, iMax, jMin, jMax,
439 I deltaTmom, KappaRV,recip_HFacS,
440 U uVelD,
441 I myThid )
442 #endif /* INCLUDE_CD_CODE */
443 C-- End If implicitViscosity.AND.momStepping
444 ENDIF
445
446 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
447 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
448 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
449 c WRITE(suff,'(I10.10)') myIter+1
450 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
451 c ENDIF
452 Cjmc(end)
453
454 #ifdef ALLOW_TIMEAVE
455 IF (taveFreq.GT.0.) THEN
456 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
457 I deltaTclock, bi, bj, myThid)
458 ENDIF
459 #endif /* ALLOW_TIMEAVE */
460
461 ENDDO
462 ENDDO
463
464 #ifndef DISABLE_DEBUGMODE
465 If (debugMode) THEN
466 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
467 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
468 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
469 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
470 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
471 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
472 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
473 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
474 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
475 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
476 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
477 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
478 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
479 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
480 ENDIF
481 #endif
482
483 RETURN
484 END

  ViewVC Help
Powered by ViewVC 1.1.22