/[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.84 - (show annotations) (download)
Fri Nov 16 03:25:40 2001 UTC (22 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint43a-release1mods, release1-branch_tutorials, chkpt44a_post, chkpt44c_pre, release1-branch-end, checkpoint44b_post, chkpt44a_pre, checkpoint44b_pre, checkpoint44, chkpt44c_post, release1-branch_branchpoint
Branch point for: release1-branch
Changes since 1.83: +1 -9 lines
fix diagnostic of convective adjustment (IVDC)
 (broken since thermo-/dynamics split)

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.83 2001/09/27 20:12:10 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 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
157 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
158 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
159
160 INTEGER iMin, iMax
161 INTEGER jMin, jMax
162 INTEGER bi, bj
163 INTEGER i, j
164 INTEGER k, km1, kp1, kup, kDown
165
166 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
167 c CHARACTER*(MAX_LEN_MBUF) suff
168 c LOGICAL DIFFERENT_MULTIPLE
169 c EXTERNAL DIFFERENT_MULTIPLE
170 Cjmc(end)
171
172 C--- The algorithm...
173 C
174 C "Correction Step"
175 C =================
176 C Here we update the horizontal velocities with the surface
177 C pressure such that the resulting flow is either consistent
178 C with the free-surface evolution or the rigid-lid:
179 C U[n] = U* + dt x d/dx P
180 C V[n] = V* + dt x d/dy P
181 C
182 C "Calculation of Gs"
183 C ===================
184 C This is where all the accelerations and tendencies (ie.
185 C physics, parameterizations etc...) are calculated
186 C rho = rho ( theta[n], salt[n] )
187 C b = b(rho, theta)
188 C K31 = K31 ( rho )
189 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
190 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
191 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
192 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
193 C
194 C "Time-stepping" or "Prediction"
195 C ================================
196 C The models variables are stepped forward with the appropriate
197 C time-stepping scheme (currently we use Adams-Bashforth II)
198 C - For momentum, the result is always *only* a "prediction"
199 C in that the flow may be divergent and will be "corrected"
200 C later with a surface pressure gradient.
201 C - Normally for tracers the result is the new field at time
202 C level [n+1} *BUT* in the case of implicit diffusion the result
203 C is also *only* a prediction.
204 C - We denote "predictors" with an asterisk (*).
205 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
206 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
207 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
208 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209 C With implicit diffusion:
210 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
211 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
212 C (1 + dt * K * d_zz) theta[n] = theta*
213 C (1 + dt * K * d_zz) salt[n] = salt*
214 C---
215 CEOP
216
217 C-- Set up work arrays with valid (i.e. not NaN) values
218 C These inital values do not alter the numerical results. They
219 C just ensure that all memory references are to valid floating
220 C point numbers. This prevents spurious hardware signals due to
221 C uninitialised but inert locations.
222 DO j=1-OLy,sNy+OLy
223 DO i=1-OLx,sNx+OLx
224 DO k=1,Nr
225 phiHyd(i,j,k) = 0. _d 0
226 KappaRU(i,j,k) = 0. _d 0
227 KappaRV(i,j,k) = 0. _d 0
228 sigmaX(i,j,k) = 0. _d 0
229 sigmaY(i,j,k) = 0. _d 0
230 sigmaR(i,j,k) = 0. _d 0
231 ENDDO
232 rhoKM1 (i,j) = 0. _d 0
233 rhok (i,j) = 0. _d 0
234 phiSurfX(i,j) = 0. _d 0
235 phiSurfY(i,j) = 0. _d 0
236 ENDDO
237 ENDDO
238
239 #ifdef ALLOW_AUTODIFF_TAMC
240 C-- HPF directive to help TAMC
241 CHPF$ INDEPENDENT
242 #endif /* ALLOW_AUTODIFF_TAMC */
243
244 DO bj=myByLo(myThid),myByHi(myThid)
245
246 #ifdef ALLOW_AUTODIFF_TAMC
247 C-- HPF directive to help TAMC
248 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
249 CHPF$& ,phiHyd
250 CHPF$& ,KappaRU,KappaRV
251 CHPF$& )
252 #endif /* ALLOW_AUTODIFF_TAMC */
253
254 DO bi=myBxLo(myThid),myBxHi(myThid)
255
256 #ifdef ALLOW_AUTODIFF_TAMC
257 act1 = bi - myBxLo(myThid)
258 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
259 act2 = bj - myByLo(myThid)
260 max2 = myByHi(myThid) - myByLo(myThid) + 1
261 act3 = myThid - 1
262 max3 = nTx*nTy
263 act4 = ikey_dynamics - 1
264 ikey = (act1 + 1) + act2*max1
265 & + act3*max1*max2
266 & + act4*max1*max2*max3
267 #endif /* ALLOW_AUTODIFF_TAMC */
268
269 C-- Set up work arrays that need valid initial values
270 DO j=1-OLy,sNy+OLy
271 DO i=1-OLx,sNx+OLx
272 fVerU (i,j,1) = 0. _d 0
273 fVerU (i,j,2) = 0. _d 0
274 fVerV (i,j,1) = 0. _d 0
275 fVerV (i,j,2) = 0. _d 0
276 ENDDO
277 ENDDO
278
279 C-- Start computation of dynamics
280 iMin = 1-OLx+2
281 iMax = sNx+OLx-1
282 jMin = 1-OLy+2
283 jMax = sNy+OLy-1
284
285 #ifdef ALLOW_AUTODIFF_TAMC
286 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
287 #endif /* ALLOW_AUTODIFF_TAMC */
288
289 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
290 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
291 IF (implicSurfPress.NE.1.) THEN
292 CALL CALC_GRAD_PHI_SURF(
293 I bi,bj,iMin,iMax,jMin,jMax,
294 I etaN,
295 O phiSurfX,phiSurfY,
296 I myThid )
297 ENDIF
298
299 #ifdef ALLOW_AUTODIFF_TAMC
300 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
301 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
302 #ifdef ALLOW_KPP
303 CADJ STORE KPPviscAz (:,:,:,bi,bj)
304 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
305 #endif /* ALLOW_KPP */
306 #endif /* ALLOW_AUTODIFF_TAMC */
307
308 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
309 C-- Calculate the total vertical diffusivity
310 DO k=1,Nr
311 CALL CALC_VISCOSITY(
312 I bi,bj,iMin,iMax,jMin,jMax,k,
313 O KappaRU,KappaRV,
314 I myThid)
315 ENDDO
316 #endif
317
318 C-- Start of dynamics loop
319 DO k=1,Nr
320
321 C-- km1 Points to level above k (=k-1)
322 C-- kup Cycles through 1,2 to point to layer above
323 C-- kDown Cycles through 2,1 to point to current layer
324
325 km1 = MAX(1,k-1)
326 kp1 = MIN(k+1,Nr)
327 kup = 1+MOD(k+1,2)
328 kDown= 1+MOD(k,2)
329
330 #ifdef ALLOW_AUTODIFF_TAMC
331 kkey = (ikey-1)*Nr + k
332 #endif /* ALLOW_AUTODIFF_TAMC */
333
334 C-- Integrate hydrostatic balance for phiHyd with BC of
335 C phiHyd(z=0)=0
336 C distinguishe between Stagger and Non Stagger time stepping
337 IF (staggerTimeStep) THEN
338 CALL CALC_PHI_HYD(
339 I bi,bj,iMin,iMax,jMin,jMax,k,
340 I gT, gS,
341 U phiHyd,
342 I myThid )
343 ELSE
344 CALL CALC_PHI_HYD(
345 I bi,bj,iMin,iMax,jMin,jMax,k,
346 I theta, salt,
347 U phiHyd,
348 I myThid )
349 ENDIF
350
351 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
352 C and step forward storing the result in gUnm1, gVnm1, etc...
353 IF ( momStepping ) THEN
354 #ifndef DISABLE_MOM_FLUXFORM
355 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
356 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
357 I phiHyd,KappaRU,KappaRV,
358 U fVerU, fVerV,
359 I myTime, myIter, myThid)
360 #endif
361 #ifndef DISABLE_MOM_VECINV
362 IF (vectorInvariantMomentum) CALL MOM_VECINV(
363 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
364 I phiHyd,KappaRU,KappaRV,
365 U fVerU, fVerV,
366 I myTime, myIter, myThid)
367 #endif
368 CALL TIMESTEP(
369 I bi,bj,iMin,iMax,jMin,jMax,k,
370 I phiHyd, phiSurfX, phiSurfY,
371 I myIter, myThid)
372
373 #ifdef ALLOW_OBCS
374 C-- Apply open boundary conditions
375 IF (useOBCS) THEN
376 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
377 END IF
378 #endif /* ALLOW_OBCS */
379
380 #ifdef ALLOW_AUTODIFF_TAMC
381 #ifdef INCLUDE_CD_CODE
382 ELSE
383 DO j=1-OLy,sNy+OLy
384 DO i=1-OLx,sNx+OLx
385 guCD(i,j,k,bi,bj) = 0.0
386 gvCD(i,j,k,bi,bj) = 0.0
387 END DO
388 END DO
389 #endif /* INCLUDE_CD_CODE */
390 #endif /* ALLOW_AUTODIFF_TAMC */
391 ENDIF
392
393
394 C-- end of dynamics k loop (1:Nr)
395 ENDDO
396
397
398
399 C-- Implicit viscosity
400 IF (implicitViscosity.AND.momStepping) THEN
401 #ifdef ALLOW_AUTODIFF_TAMC
402 idkey = iikey + 3
403 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
404 #endif /* ALLOW_AUTODIFF_TAMC */
405 CALL IMPLDIFF(
406 I bi, bj, iMin, iMax, jMin, jMax,
407 I deltaTmom, KappaRU,recip_HFacW,
408 U gUNm1,
409 I myThid )
410 #ifdef ALLOW_AUTODIFF_TAMC
411 idkey = iikey + 4
412 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
413 #endif /* ALLOW_AUTODIFF_TAMC */
414 CALL IMPLDIFF(
415 I bi, bj, iMin, iMax, jMin, jMax,
416 I deltaTmom, KappaRV,recip_HFacS,
417 U gVNm1,
418 I myThid )
419
420 #ifdef ALLOW_OBCS
421 C-- Apply open boundary conditions
422 IF (useOBCS) THEN
423 DO K=1,Nr
424 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
425 ENDDO
426 END IF
427 #endif /* ALLOW_OBCS */
428
429 #ifdef INCLUDE_CD_CODE
430 #ifdef ALLOW_AUTODIFF_TAMC
431 idkey = iikey + 5
432 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
433 #endif /* ALLOW_AUTODIFF_TAMC */
434 CALL IMPLDIFF(
435 I bi, bj, iMin, iMax, jMin, jMax,
436 I deltaTmom, KappaRU,recip_HFacW,
437 U vVelD,
438 I myThid )
439 #ifdef ALLOW_AUTODIFF_TAMC
440 idkey = iikey + 6
441 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
442 #endif /* ALLOW_AUTODIFF_TAMC */
443 CALL IMPLDIFF(
444 I bi, bj, iMin, iMax, jMin, jMax,
445 I deltaTmom, KappaRV,recip_HFacS,
446 U uVelD,
447 I myThid )
448 #endif /* INCLUDE_CD_CODE */
449 C-- End If implicitViscosity.AND.momStepping
450 ENDIF
451
452 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
453 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
454 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
455 c WRITE(suff,'(I10.10)') myIter+1
456 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
457 c ENDIF
458 Cjmc(end)
459
460 #ifdef ALLOW_TIMEAVE
461 IF (taveFreq.GT.0.) THEN
462 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
463 I deltaTclock, bi, bj, myThid)
464 ENDIF
465 #endif /* ALLOW_TIMEAVE */
466
467 ENDDO
468 ENDDO
469
470 #ifndef DISABLE_DEBUGMODE
471 If (debugMode) THEN
472 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
473 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
474 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
475 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
476 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
477 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
478 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
479 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
480 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
481 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
482 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
483 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
484 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
485 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
486 ENDIF
487 #endif
488
489 RETURN
490 END

  ViewVC Help
Powered by ViewVC 1.1.22