/[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.82 - (show annotations) (download)
Wed Sep 26 18:09:14 2001 UTC (22 years, 8 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint41
Changes since 1.81: +91 -21 lines
Bringing comments up to data and formatting for document extraction.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.81 2001/09/19 02:43:27 adcroft 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 TIMEAVE_CUMULATE
117 C |-- CALL DEBUG_STATS_RL
118
119 C !INPUT/OUTPUT PARAMETERS:
120 C == Routine arguments ==
121 C myTime - Current time in simulation
122 C myIter - Current iteration number in simulation
123 C myThid - Thread number for this instance of the routine.
124 _RL myTime
125 INTEGER myIter
126 INTEGER myThid
127
128 C !LOCAL VARIABLES:
129 C == Local variables
130 C fVer[STUV] o fVer: Vertical flux term - note fVer
131 C is "pipelined" in the vertical
132 C so we need an fVer for each
133 C variable.
134 C rhoK, rhoKM1 - Density at current level, and level above
135 C phiHyd - Hydrostatic part of the potential phiHydi.
136 C In z coords phiHydiHyd is the hydrostatic
137 C Potential (=pressure/rho0) anomaly
138 C In p coords phiHydiHyd is the geopotential
139 C surface height anomaly.
140 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
141 C phiSurfY or geopotentiel (atmos) in X and Y direction
142 C iMin, iMax - Ranges and sub-block indices on which calculations
143 C jMin, jMax are applied.
144 C bi, bj
145 C k, kup, - Index for layer above and below. kup and kDown
146 C kDown, km1 are switched with layer to be the appropriate
147 C index into fVerTerm.
148 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
150 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
151 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
157 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
158 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
159 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
160
161 C This is currently used by IVDC and Diagnostics
162 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
163
164 INTEGER iMin, iMax
165 INTEGER jMin, jMax
166 INTEGER bi, bj
167 INTEGER i, j
168 INTEGER k, km1, kp1, kup, kDown
169
170 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
171 c CHARACTER*(MAX_LEN_MBUF) suff
172 c LOGICAL DIFFERENT_MULTIPLE
173 c EXTERNAL DIFFERENT_MULTIPLE
174 Cjmc(end)
175
176 C--- The algorithm...
177 C
178 C "Correction Step"
179 C =================
180 C Here we update the horizontal velocities with the surface
181 C pressure such that the resulting flow is either consistent
182 C with the free-surface evolution or the rigid-lid:
183 C U[n] = U* + dt x d/dx P
184 C V[n] = V* + dt x d/dy P
185 C
186 C "Calculation of Gs"
187 C ===================
188 C This is where all the accelerations and tendencies (ie.
189 C physics, parameterizations etc...) are calculated
190 C rho = rho ( theta[n], salt[n] )
191 C b = b(rho, theta)
192 C K31 = K31 ( rho )
193 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
194 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
195 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
196 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
197 C
198 C "Time-stepping" or "Prediction"
199 C ================================
200 C The models variables are stepped forward with the appropriate
201 C time-stepping scheme (currently we use Adams-Bashforth II)
202 C - For momentum, the result is always *only* a "prediction"
203 C in that the flow may be divergent and will be "corrected"
204 C later with a surface pressure gradient.
205 C - Normally for tracers the result is the new field at time
206 C level [n+1} *BUT* in the case of implicit diffusion the result
207 C is also *only* a prediction.
208 C - We denote "predictors" with an asterisk (*).
209 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
210 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
211 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
212 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213 C With implicit diffusion:
214 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
215 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
216 C (1 + dt * K * d_zz) theta[n] = theta*
217 C (1 + dt * K * d_zz) salt[n] = salt*
218 C---
219 CEOP
220
221 C-- Set up work arrays with valid (i.e. not NaN) values
222 C These inital values do not alter the numerical results. They
223 C just ensure that all memory references are to valid floating
224 C point numbers. This prevents spurious hardware signals due to
225 C uninitialised but inert locations.
226 DO j=1-OLy,sNy+OLy
227 DO i=1-OLx,sNx+OLx
228 DO k=1,Nr
229 phiHyd(i,j,k) = 0. _d 0
230 KappaRU(i,j,k) = 0. _d 0
231 KappaRV(i,j,k) = 0. _d 0
232 sigmaX(i,j,k) = 0. _d 0
233 sigmaY(i,j,k) = 0. _d 0
234 sigmaR(i,j,k) = 0. _d 0
235 ENDDO
236 rhoKM1 (i,j) = 0. _d 0
237 rhok (i,j) = 0. _d 0
238 phiSurfX(i,j) = 0. _d 0
239 phiSurfY(i,j) = 0. _d 0
240 ENDDO
241 ENDDO
242
243 #ifdef ALLOW_AUTODIFF_TAMC
244 C-- HPF directive to help TAMC
245 CHPF$ INDEPENDENT
246 #endif /* ALLOW_AUTODIFF_TAMC */
247
248 DO bj=myByLo(myThid),myByHi(myThid)
249
250 #ifdef ALLOW_AUTODIFF_TAMC
251 C-- HPF directive to help TAMC
252 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
253 CHPF$& ,phiHyd
254 CHPF$& ,KappaRU,KappaRV
255 CHPF$& )
256 #endif /* ALLOW_AUTODIFF_TAMC */
257
258 DO bi=myBxLo(myThid),myBxHi(myThid)
259
260 #ifdef ALLOW_AUTODIFF_TAMC
261 act1 = bi - myBxLo(myThid)
262 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
263
264 act2 = bj - myByLo(myThid)
265 max2 = myByHi(myThid) - myByLo(myThid) + 1
266
267 act3 = myThid - 1
268 max3 = nTx*nTy
269
270 act4 = ikey_dynamics - 1
271
272 ikey = (act1 + 1) + act2*max1
273 & + act3*max1*max2
274 & + act4*max1*max2*max3
275 #endif /* ALLOW_AUTODIFF_TAMC */
276
277 C-- Set up work arrays that need valid initial values
278 DO j=1-OLy,sNy+OLy
279 DO i=1-OLx,sNx+OLx
280 fVerU (i,j,1) = 0. _d 0
281 fVerU (i,j,2) = 0. _d 0
282 fVerV (i,j,1) = 0. _d 0
283 fVerV (i,j,2) = 0. _d 0
284 ENDDO
285 ENDDO
286
287 C-- Start computation of dynamics
288 iMin = 1-OLx+2
289 iMax = sNx+OLx-1
290 jMin = 1-OLy+2
291 jMax = sNy+OLy-1
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 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
297 #endif /* ALLOW_AUTODIFF_TAMC */
298
299 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
300 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
301 IF (implicSurfPress.NE.1.) THEN
302 CALL CALC_GRAD_PHI_SURF(
303 I bi,bj,iMin,iMax,jMin,jMax,
304 I etaN,
305 O phiSurfX,phiSurfY,
306 I myThid )
307 ENDIF
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
399
400 C-- Implicit viscosity
401 IF (implicitViscosity.AND.momStepping) THEN
402 #ifdef ALLOW_AUTODIFF_TAMC
403 idkey = iikey + 3
404 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
405 #endif /* ALLOW_AUTODIFF_TAMC */
406 CALL IMPLDIFF(
407 I bi, bj, iMin, iMax, jMin, jMax,
408 I deltaTmom, KappaRU,recip_HFacW,
409 U gUNm1,
410 I myThid )
411 #ifdef ALLOW_AUTODIFF_TAMC
412 idkey = iikey + 4
413 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
414 #endif /* ALLOW_AUTODIFF_TAMC */
415 CALL IMPLDIFF(
416 I bi, bj, iMin, iMax, jMin, jMax,
417 I deltaTmom, KappaRV,recip_HFacS,
418 U gVNm1,
419 I myThid )
420
421 #ifdef ALLOW_OBCS
422 C-- Apply open boundary conditions
423 IF (useOBCS) THEN
424 DO K=1,Nr
425 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
426 ENDDO
427 END IF
428 #endif /* ALLOW_OBCS */
429
430 #ifdef INCLUDE_CD_CODE
431 #ifdef ALLOW_AUTODIFF_TAMC
432 idkey = iikey + 5
433 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
434 #endif /* ALLOW_AUTODIFF_TAMC */
435 CALL IMPLDIFF(
436 I bi, bj, iMin, iMax, jMin, jMax,
437 I deltaTmom, KappaRU,recip_HFacW,
438 U vVelD,
439 I myThid )
440 #ifdef ALLOW_AUTODIFF_TAMC
441 idkey = iikey + 6
442 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
443 #endif /* ALLOW_AUTODIFF_TAMC */
444 CALL IMPLDIFF(
445 I bi, bj, iMin, iMax, jMin, jMax,
446 I deltaTmom, KappaRV,recip_HFacS,
447 U uVelD,
448 I myThid )
449 #endif /* INCLUDE_CD_CODE */
450 C-- End If implicitViscosity.AND.momStepping
451 ENDIF
452
453 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
454 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
455 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
456 c WRITE(suff,'(I10.10)') myIter+1
457 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
458 c ENDIF
459 Cjmc(end)
460
461 #ifdef ALLOW_TIMEAVE
462 IF (taveFreq.GT.0.) THEN
463 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
464 I deltaTclock, bi, bj, myThid)
465 IF (ivdc_kappa.NE.0.) THEN
466 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
467 I deltaTclock, bi, bj, myThid)
468 ENDIF
469 ENDIF
470 #endif /* ALLOW_TIMEAVE */
471
472 ENDDO
473 ENDDO
474
475 #ifndef DISABLE_DEBUGMODE
476 If (debugMode) THEN
477 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
478 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
479 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
480 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
481 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
482 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
483 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
484 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
485 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
486 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
487 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
488 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
489 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
490 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
491 ENDIF
492 #endif
493
494 RETURN
495 END

  ViewVC Help
Powered by ViewVC 1.1.22