/[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.4.4 - (show annotations) (download)
Wed Feb 19 22:43:51 2003 UTC (21 years, 3 months ago) by dimitri
Branch: ecco-branch
CVS Tags: icebear5, icebear4, icebear3, ecco_c44_e27
Branch point for: icebear
Changes since 1.83.4.3: +0 -4 lines
o Removed spurious iikey and idkey

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

  ViewVC Help
Powered by ViewVC 1.1.22