/[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.99 - (show annotations) (download)
Thu Oct 2 21:33:54 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51h_pre, checkpoint51i_pre, checkpoint51g_post
Changes since 1.98: +10 -1 lines
Bringing code up to date for AD
o remove some IF-statements which cause excessive dependencies
o provide interface for ADM*TLM

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.98 2003/07/08 15:00:26 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 # include "EOS.h"
84 # ifdef ALLOW_KPP
85 # include "KPP.h"
86 # endif
87 #endif /* ALLOW_AUTODIFF_TAMC */
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 phiHydC :: hydrostatic potential anomaly at cell center
129 C In z coords phiHyd is the hydrostatic potential
130 C (=pressure/rho0) anomaly
131 C In p coords phiHyd is the geopotential height anomaly.
132 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
133 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
134 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
135 C phiSurfY or geopotential (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 phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
147 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
148 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
151 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
152
153 INTEGER iMin, iMax
154 INTEGER jMin, jMax
155 INTEGER bi, bj
156 INTEGER i, j
157 INTEGER k, km1, kp1, kup, kDown
158
159 LOGICAL DIFFERENT_MULTIPLE
160 EXTERNAL DIFFERENT_MULTIPLE
161
162 C--- The algorithm...
163 C
164 C "Correction Step"
165 C =================
166 C Here we update the horizontal velocities with the surface
167 C pressure such that the resulting flow is either consistent
168 C with the free-surface evolution or the rigid-lid:
169 C U[n] = U* + dt x d/dx P
170 C V[n] = V* + dt x d/dy P
171 C
172 C "Calculation of Gs"
173 C ===================
174 C This is where all the accelerations and tendencies (ie.
175 C physics, parameterizations etc...) are calculated
176 C rho = rho ( theta[n], salt[n] )
177 C b = b(rho, theta)
178 C K31 = K31 ( rho )
179 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
180 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
181 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
182 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
183 C
184 C "Time-stepping" or "Prediction"
185 C ================================
186 C The models variables are stepped forward with the appropriate
187 C time-stepping scheme (currently we use Adams-Bashforth II)
188 C - For momentum, the result is always *only* a "prediction"
189 C in that the flow may be divergent and will be "corrected"
190 C later with a surface pressure gradient.
191 C - Normally for tracers the result is the new field at time
192 C level [n+1} *BUT* in the case of implicit diffusion the result
193 C is also *only* a prediction.
194 C - We denote "predictors" with an asterisk (*).
195 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
196 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
197 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
198 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
199 C With implicit diffusion:
200 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
201 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
202 C (1 + dt * K * d_zz) theta[n] = theta*
203 C (1 + dt * K * d_zz) salt[n] = salt*
204 C---
205 CEOP
206
207 C-- Call to routine for calculation of
208 C Eliassen-Palm-flux-forced U-tendency,
209 C if desired:
210 #ifdef INCLUDE_EP_FORCING_CODE
211 CALL CALC_EP_FORCING(myThid)
212 #endif
213
214 #ifdef ALLOW_AUTODIFF_TAMC
215 C-- HPF directive to help TAMC
216 CHPF$ INDEPENDENT
217 #endif /* ALLOW_AUTODIFF_TAMC */
218
219 DO bj=myByLo(myThid),myByHi(myThid)
220
221 #ifdef ALLOW_AUTODIFF_TAMC
222 C-- HPF directive to help TAMC
223 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
224 CHPF$& ,phiHydF
225 CHPF$& ,KappaRU,KappaRV
226 CHPF$& )
227 #endif /* ALLOW_AUTODIFF_TAMC */
228
229 DO bi=myBxLo(myThid),myBxHi(myThid)
230
231 #ifdef ALLOW_AUTODIFF_TAMC
232 act1 = bi - myBxLo(myThid)
233 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
234 act2 = bj - myByLo(myThid)
235 max2 = myByHi(myThid) - myByLo(myThid) + 1
236 act3 = myThid - 1
237 max3 = nTx*nTy
238 act4 = ikey_dynamics - 1
239 idynkey = (act1 + 1) + act2*max1
240 & + act3*max1*max2
241 & + act4*max1*max2*max3
242 #endif /* ALLOW_AUTODIFF_TAMC */
243
244 C-- Set up work arrays with valid (i.e. not NaN) values
245 C These inital values do not alter the numerical results. They
246 C just ensure that all memory references are to valid floating
247 C point numbers. This prevents spurious hardware signals due to
248 C uninitialised but inert locations.
249
250 DO k=1,Nr
251 DO j=1-OLy,sNy+OLy
252 DO i=1-OLx,sNx+OLx
253 KappaRU(i,j,k) = 0. _d 0
254 KappaRV(i,j,k) = 0. _d 0
255 #ifdef ALLOW_AUTODIFF_TAMC
256 cph(
257 c-- need some re-initialisation here to break dependencies
258 c-- totphihyd is assumed zero from ini_pressure, i.e.
259 c-- avoiding iterate pressure p = integral of (g*rho(p)*dz)
260 cph)
261 totPhiHyd(i,j,k,bi,bj) = 0. _d 0
262 gu(i,j,k,bi,bj) = 0. _d 0
263 gv(i,j,k,bi,bj) = 0. _d 0
264 #endif
265 ENDDO
266 ENDDO
267 ENDDO
268 DO j=1-OLy,sNy+OLy
269 DO i=1-OLx,sNx+OLx
270 fVerU (i,j,1) = 0. _d 0
271 fVerU (i,j,2) = 0. _d 0
272 fVerV (i,j,1) = 0. _d 0
273 fVerV (i,j,2) = 0. _d 0
274 phiHydF (i,j) = 0. _d 0
275 phiHydC (i,j) = 0. _d 0
276 dPhiHydX(i,j) = 0. _d 0
277 dPhiHydY(i,j) = 0. _d 0
278 phiSurfX(i,j) = 0. _d 0
279 phiSurfY(i,j) = 0. _d 0
280 ENDDO
281 ENDDO
282
283 C-- Start computation of dynamics
284 iMin = 0
285 iMax = sNx+1
286 jMin = 0
287 jMax = sNy+1
288
289 #ifdef ALLOW_AUTODIFF_TAMC
290 CADJ STORE wvel (:,:,:,bi,bj) =
291 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
292 #endif /* ALLOW_AUTODIFF_TAMC */
293
294 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
295 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
296 IF (implicSurfPress.NE.1.) THEN
297 CALL CALC_GRAD_PHI_SURF(
298 I bi,bj,iMin,iMax,jMin,jMax,
299 I etaN,
300 O phiSurfX,phiSurfY,
301 I myThid )
302 ENDIF
303
304 #ifdef ALLOW_AUTODIFF_TAMC
305 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
306 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
307 #ifdef ALLOW_KPP
308 CADJ STORE KPPviscAz (:,:,:,bi,bj)
309 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
310 #endif /* ALLOW_KPP */
311 #endif /* ALLOW_AUTODIFF_TAMC */
312
313 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
314 C-- Calculate the total vertical diffusivity
315 DO k=1,Nr
316 CALL CALC_VISCOSITY(
317 I bi,bj,iMin,iMax,jMin,jMax,k,
318 O KappaRU,KappaRV,
319 I myThid)
320 ENDDO
321 #endif
322
323 C-- Start of dynamics loop
324 DO k=1,Nr
325
326 C-- km1 Points to level above k (=k-1)
327 C-- kup Cycles through 1,2 to point to layer above
328 C-- kDown Cycles through 2,1 to point to current layer
329
330 km1 = MAX(1,k-1)
331 kp1 = MIN(k+1,Nr)
332 kup = 1+MOD(k+1,2)
333 kDown= 1+MOD(k,2)
334
335 #ifdef ALLOW_AUTODIFF_TAMC
336 kkey = (idynkey-1)*Nr + k
337 c
338 CADJ STORE totphihyd (:,:,k,bi,bj)
339 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
340 CADJ STORE gt (:,:,k,bi,bj)
341 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
342 CADJ STORE gs (:,:,k,bi,bj)
343 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
344 CADJ STORE theta (:,:,k,bi,bj)
345 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
346 CADJ STORE salt (:,:,k,bi,bj)
347 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
348 #endif /* ALLOW_AUTODIFF_TAMC */
349
350 C-- Integrate hydrostatic balance for phiHyd with BC of
351 C phiHyd(z=0)=0
352 C distinguishe between Stagger and Non Stagger time stepping
353 IF (staggerTimeStep) THEN
354 CALL CALC_PHI_HYD(
355 I bi,bj,iMin,iMax,jMin,jMax,k,
356 I gT, gS,
357 U phiHydF,
358 O phiHydC, dPhiHydX, dPhiHydY,
359 I myTime, myIter, myThid )
360 ELSE
361 CALL CALC_PHI_HYD(
362 I bi,bj,iMin,iMax,jMin,jMax,k,
363 I theta, salt,
364 U phiHydF,
365 O phiHydC, dPhiHydX, dPhiHydY,
366 I myTime, myIter, myThid )
367 ENDIF
368
369 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
370 C and step forward storing the result in gU, gV, etc...
371 IF ( momStepping ) THEN
372 #ifndef DISABLE_MOM_FLUXFORM
373 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
374 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
375 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
376 U fVerU, fVerV,
377 I myTime, myIter, myThid)
378 #endif
379 #ifndef DISABLE_MOM_VECINV
380 IF (vectorInvariantMomentum) CALL MOM_VECINV(
381 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
382 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
383 U fVerU, fVerV,
384 I myTime, myIter, myThid)
385 #endif
386 CALL TIMESTEP(
387 I bi,bj,iMin,iMax,jMin,jMax,k,
388 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
389 I myTime, myIter, myThid)
390
391 #ifdef ALLOW_OBCS
392 C-- Apply open boundary conditions
393 IF (useOBCS) THEN
394 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
395 ENDIF
396 #endif /* ALLOW_OBCS */
397
398 ENDIF
399
400
401 C-- end of dynamics k loop (1:Nr)
402 ENDDO
403
404 C-- Implicit viscosity
405 IF (implicitViscosity.AND.momStepping) THEN
406 #ifdef ALLOW_AUTODIFF_TAMC
407 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
408 #endif /* ALLOW_AUTODIFF_TAMC */
409 CALL IMPLDIFF(
410 I bi, bj, iMin, iMax, jMin, jMax,
411 I deltaTmom, KappaRU,recip_HFacW,
412 U gU,
413 I myThid )
414 #ifdef ALLOW_AUTODIFF_TAMC
415 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
416 #endif /* ALLOW_AUTODIFF_TAMC */
417 CALL IMPLDIFF(
418 I bi, bj, iMin, iMax, jMin, jMax,
419 I deltaTmom, KappaRV,recip_HFacS,
420 U gV,
421 I myThid )
422
423 #ifdef ALLOW_OBCS
424 C-- Apply open boundary conditions
425 IF (useOBCS) THEN
426 DO K=1,Nr
427 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
428 ENDDO
429 END IF
430 #endif /* ALLOW_OBCS */
431
432 #ifdef INCLUDE_CD_CODE
433 #ifdef ALLOW_AUTODIFF_TAMC
434 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
435 #endif /* ALLOW_AUTODIFF_TAMC */
436 CALL IMPLDIFF(
437 I bi, bj, iMin, iMax, jMin, jMax,
438 I deltaTmom, KappaRU,recip_HFacW,
439 U vVelD,
440 I myThid )
441 #ifdef ALLOW_AUTODIFF_TAMC
442 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, 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 ENDDO
454 ENDDO
455
456 Cml(
457 C In order to compare the variance of phiHydLow of a p/z-coordinate
458 C run with etaH of a z/p-coordinate run the drift of phiHydLow
459 C has to be removed by something like the following subroutine:
460 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
461 C & 'phiHydLow', myThid )
462 Cml)
463
464 #ifndef DISABLE_DEBUGMODE
465 If ( debugLevel .GE. debLevB ) 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