/[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.166 - (show annotations) (download)
Sun Sep 15 14:28:31 2013 UTC (10 years, 8 months ago) by m_bates
Branch: MAIN
CVS Tags: checkpoint64p, checkpoint64o
Changes since 1.165: +7 -1 lines
Changes to implement a residual model.  Also, calculation of the mean velocity from the residual and bolus.

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.165 2013/08/03 01:38:17 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6 #ifdef ALLOW_OBCS
7 # include "OBCS_OPTIONS.h"
8 #endif
9
10 #undef DYNAMICS_GUGV_EXCH_CHECK
11
12 CBOP
13 C !ROUTINE: DYNAMICS
14 C !INTERFACE:
15 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE DYNAMICS
19 C | o Controlling routine for the explicit part of the model
20 C | dynamics.
21 C *==========================================================*
22 C | This routine evaluates the "dynamics" terms for each
23 C | block of ocean in turn. Because the blocks of ocean have
24 C | overlap regions they are independent of one another.
25 C | If terms involving lateral integrals are needed in this
26 C | routine care will be needed. Similarly finite-difference
27 C | operations with stencils wider than the overlap region
28 C | require special consideration.
29 C | The algorithm...
30 C |
31 C | "Correction Step"
32 C | =================
33 C | Here we update the horizontal velocities with the surface
34 C | pressure such that the resulting flow is either consistent
35 C | with the free-surface evolution or the rigid-lid:
36 C | U[n] = U* + dt x d/dx P
37 C | V[n] = V* + dt x d/dy P
38 C | W[n] = W* + dt x d/dz P (NH mode)
39 C |
40 C | "Calculation of Gs"
41 C | ===================
42 C | This is where all the accelerations and tendencies (ie.
43 C | physics, parameterizations etc...) are calculated
44 C | rho = rho ( theta[n], salt[n] )
45 C | b = b(rho, theta)
46 C | K31 = K31 ( rho )
47 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
48 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
49 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
50 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
51 C |
52 C | "Time-stepping" or "Prediction"
53 C | ================================
54 C | The models variables are stepped forward with the appropriate
55 C | time-stepping scheme (currently we use Adams-Bashforth II)
56 C | - For momentum, the result is always *only* a "prediction"
57 C | in that the flow may be divergent and will be "corrected"
58 C | later with a surface pressure gradient.
59 C | - Normally for tracers the result is the new field at time
60 C | level [n+1} *BUT* in the case of implicit diffusion the result
61 C | is also *only* a prediction.
62 C | - We denote "predictors" with an asterisk (*).
63 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
64 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
65 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
67 C | With implicit diffusion:
68 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
70 C | (1 + dt * K * d_zz) theta[n] = theta*
71 C | (1 + dt * K * d_zz) salt[n] = salt*
72 C |
73 C *==========================================================*
74 C \ev
75 C !USES:
76 IMPLICIT NONE
77 C == Global variables ===
78 #include "SIZE.h"
79 #include "EEPARAMS.h"
80 #include "PARAMS.h"
81 #include "GRID.h"
82 #include "DYNVARS.h"
83 #ifdef ALLOW_CD_CODE
84 # include "CD_CODE_VARS.h"
85 #endif
86 #ifdef ALLOW_AUTODIFF_TAMC
87 # include "tamc.h"
88 # include "tamc_keys.h"
89 # include "FFIELDS.h"
90 # include "EOS.h"
91 # ifdef ALLOW_KPP
92 # include "KPP.h"
93 # endif
94 # ifdef ALLOW_PTRACERS
95 # include "PTRACERS_SIZE.h"
96 # include "PTRACERS_FIELDS.h"
97 # endif
98 # ifdef ALLOW_OBCS
99 # include "OBCS_PARAMS.h"
100 # include "OBCS_FIELDS.h"
101 # ifdef ALLOW_PTRACERS
102 # include "OBCS_PTRACERS.h"
103 # endif
104 # endif
105 # ifdef ALLOW_MOM_FLUXFORM
106 # include "MOM_FLUXFORM.h"
107 # endif
108 #endif /* ALLOW_AUTODIFF_TAMC */
109
110 C !CALLING SEQUENCE:
111 C DYNAMICS()
112 C |
113 C |-- CALC_EP_FORCING
114 C |
115 C |-- CALC_GRAD_PHI_SURF
116 C |
117 C |-- CALC_VISCOSITY
118 C |
119 C |-- CALC_EDDY_STRESS
120 C |
121 C |-- CALC_PHI_HYD
122 C |
123 C |-- MOM_FLUXFORM
124 C |
125 C |-- MOM_VECINV
126 C |
127 C |-- TIMESTEP
128 C |
129 C |-- MOM_U_IMPLICIT_R
130 C |-- MOM_V_IMPLICIT_R
131 C |
132 C |-- IMPLDIFF
133 C |
134 C |-- OBCS_APPLY_UV
135 C |
136 C |-- CALC_GW
137 C |
138 C |-- DIAGNOSTICS_FILL
139 C |-- DEBUG_STATS_RL
140
141 C !INPUT/OUTPUT PARAMETERS:
142 C == Routine arguments ==
143 C myTime :: Current time in simulation
144 C myIter :: Current iteration number in simulation
145 C myThid :: Thread number for this instance of the routine.
146 _RL myTime
147 INTEGER myIter
148 INTEGER myThid
149
150 C !FUNCTIONS:
151 #ifdef ALLOW_DIAGNOSTICS
152 LOGICAL DIAGNOSTICS_IS_ON
153 EXTERNAL DIAGNOSTICS_IS_ON
154 #endif
155
156 C !LOCAL VARIABLES:
157 C == Local variables
158 C fVer[UV] o fVer: Vertical flux term - note fVer
159 C is "pipelined" in the vertical
160 C so we need an fVer for each
161 C variable.
162 C phiHydC :: hydrostatic potential anomaly at cell center
163 C In z coords phiHyd is the hydrostatic potential
164 C (=pressure/rho0) anomaly
165 C In p coords phiHyd is the geopotential height anomaly.
166 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
167 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
168 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
169 C phiSurfY or geopotential (atmos) in X and Y direction
170 C guDissip :: dissipation tendency (all explicit terms), u component
171 C gvDissip :: dissipation tendency (all explicit terms), v component
172 C KappaRU :: vertical viscosity for velocity U-component
173 C KappaRV :: vertical viscosity for velocity V-component
174 C iMin, iMax :: Ranges and sub-block indices on which calculations
175 C jMin, jMax are applied.
176 C bi, bj :: tile indices
177 C k :: current level index
178 C km1, kp1 :: index of level above (k-1) and below (k+1)
179 C kUp, kDown :: Index for interface above and below. kUp and kDown are
180 C are switched with k to be the appropriate index into fVerU,V
181 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
182 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
183 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
184 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
185 _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
186 _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
187 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
188 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
189 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
190 _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
191 _RL KappaRU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
192 _RL KappaRV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
193
194 INTEGER iMin, iMax
195 INTEGER jMin, jMax
196 INTEGER bi, bj
197 INTEGER i, j
198 INTEGER k, km1, kp1, kUp, kDown
199
200 #ifdef ALLOW_DIAGNOSTICS
201 LOGICAL dPhiHydDiagIsOn
202 _RL tmpFac
203 #endif /* ALLOW_DIAGNOSTICS */
204
205 C--- The algorithm...
206 C
207 C "Correction Step"
208 C =================
209 C Here we update the horizontal velocities with the surface
210 C pressure such that the resulting flow is either consistent
211 C with the free-surface evolution or the rigid-lid:
212 C U[n] = U* + dt x d/dx P
213 C V[n] = V* + dt x d/dy P
214 C
215 C "Calculation of Gs"
216 C ===================
217 C This is where all the accelerations and tendencies (ie.
218 C physics, parameterizations etc...) are calculated
219 C rho = rho ( theta[n], salt[n] )
220 C b = b(rho, theta)
221 C K31 = K31 ( rho )
222 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
223 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
224 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
225 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
226 C
227 C "Time-stepping" or "Prediction"
228 C ================================
229 C The models variables are stepped forward with the appropriate
230 C time-stepping scheme (currently we use Adams-Bashforth II)
231 C - For momentum, the result is always *only* a "prediction"
232 C in that the flow may be divergent and will be "corrected"
233 C later with a surface pressure gradient.
234 C - Normally for tracers the result is the new field at time
235 C level [n+1} *BUT* in the case of implicit diffusion the result
236 C is also *only* a prediction.
237 C - We denote "predictors" with an asterisk (*).
238 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
239 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
240 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
241 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
242 C With implicit diffusion:
243 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
244 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
245 C (1 + dt * K * d_zz) theta[n] = theta*
246 C (1 + dt * K * d_zz) salt[n] = salt*
247 C---
248 CEOP
249
250 #ifdef ALLOW_DEBUG
251 IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
252 #endif
253
254 #ifdef ALLOW_DIAGNOSTICS
255 dPhiHydDiagIsOn = .FALSE.
256 IF ( useDiagnostics )
257 & dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
258 & .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
259 #endif
260
261 C-- Call to routine for calculation of Eliassen-Palm-flux-forced
262 C U-tendency, if desired:
263 #ifdef INCLUDE_EP_FORCING_CODE
264 CALL CALC_EP_FORCING(myThid)
265 #endif
266
267 #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
268 CALL DUMMY_IN_DYNAMICS( myTime, myIter, myThid )
269 #endif
270
271 #ifdef ALLOW_AUTODIFF_TAMC
272 C-- HPF directive to help TAMC
273 CHPF$ INDEPENDENT
274 #endif /* ALLOW_AUTODIFF_TAMC */
275
276 DO bj=myByLo(myThid),myByHi(myThid)
277
278 #ifdef ALLOW_AUTODIFF_TAMC
279 C-- HPF directive to help TAMC
280 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
281 CHPF$& ,phiHydF
282 CHPF$& ,KappaRU,KappaRV
283 CHPF$& )
284 #endif /* ALLOW_AUTODIFF_TAMC */
285
286 DO bi=myBxLo(myThid),myBxHi(myThid)
287
288 #ifdef ALLOW_AUTODIFF_TAMC
289 act1 = bi - myBxLo(myThid)
290 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
291 act2 = bj - myByLo(myThid)
292 max2 = myByHi(myThid) - myByLo(myThid) + 1
293 act3 = myThid - 1
294 max3 = nTx*nTy
295 act4 = ikey_dynamics - 1
296 idynkey = (act1 + 1) + act2*max1
297 & + act3*max1*max2
298 & + act4*max1*max2*max3
299 #endif /* ALLOW_AUTODIFF_TAMC */
300
301 C-- Set up work arrays with valid (i.e. not NaN) values
302 C These initial values do not alter the numerical results. They
303 C just ensure that all memory references are to valid floating
304 C point numbers. This prevents spurious hardware signals due to
305 C uninitialised but inert locations.
306
307 #ifdef ALLOW_AUTODIFF_TAMC
308 DO k=1,Nr
309 DO j=1-OLy,sNy+OLy
310 DO i=1-OLx,sNx+OLx
311 cph(
312 c-- need some re-initialisation here to break dependencies
313 cph)
314 gU(i,j,k,bi,bj) = 0. _d 0
315 gV(i,j,k,bi,bj) = 0. _d 0
316 ENDDO
317 ENDDO
318 ENDDO
319 #endif /* ALLOW_AUTODIFF_TAMC */
320 DO j=1-OLy,sNy+OLy
321 DO i=1-OLx,sNx+OLx
322 fVerU (i,j,1) = 0. _d 0
323 fVerU (i,j,2) = 0. _d 0
324 fVerV (i,j,1) = 0. _d 0
325 fVerV (i,j,2) = 0. _d 0
326 phiHydF (i,j) = 0. _d 0
327 phiHydC (i,j) = 0. _d 0
328 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
329 dPhiHydX(i,j) = 0. _d 0
330 dPhiHydY(i,j) = 0. _d 0
331 #endif
332 phiSurfX(i,j) = 0. _d 0
333 phiSurfY(i,j) = 0. _d 0
334 guDissip(i,j) = 0. _d 0
335 gvDissip(i,j) = 0. _d 0
336 #ifdef ALLOW_AUTODIFF_TAMC
337 phiHydLow(i,j,bi,bj) = 0. _d 0
338 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
339 # ifndef DISABLE_RSTAR_CODE
340 # ifndef ALLOW_AUTODIFF_OPENAD
341 dWtransC(i,j,bi,bj) = 0. _d 0
342 dWtransU(i,j,bi,bj) = 0. _d 0
343 dWtransV(i,j,bi,bj) = 0. _d 0
344 # endif
345 # endif
346 # endif
347 #endif
348 ENDDO
349 ENDDO
350
351 C-- Start computation of dynamics
352 iMin = 0
353 iMax = sNx+1
354 jMin = 0
355 jMax = sNy+1
356
357 #ifdef ALLOW_AUTODIFF_TAMC
358 CADJ STORE wVel (:,:,:,bi,bj) =
359 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
360 #endif /* ALLOW_AUTODIFF_TAMC */
361
362 C-- Explicit part of the Surface Potential Gradient (add in TIMESTEP)
363 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
364 IF (implicSurfPress.NE.1.) THEN
365 CALL CALC_GRAD_PHI_SURF(
366 I bi,bj,iMin,iMax,jMin,jMax,
367 I etaN,
368 O phiSurfX,phiSurfY,
369 I myThid )
370 ENDIF
371
372 #ifdef ALLOW_AUTODIFF_TAMC
373 CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
374 CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
375 #ifdef ALLOW_KPP
376 CADJ STORE KPPviscAz (:,:,:,bi,bj)
377 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
378 #endif /* ALLOW_KPP */
379 #endif /* ALLOW_AUTODIFF_TAMC */
380
381 #if (defined INCLUDE_CALC_DIFFUSIVITY_CALL) && !(defined ALLOW_AUTODIFF)
382 IF ( .NOT.momViscosity ) THEN
383 #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL and not ALLOW_AUTODIFF */
384 DO k=1,Nr
385 DO j=1-OLy,sNy+OLy
386 DO i=1-OLx,sNx+OLx
387 KappaRU(i,j,k) = 0. _d 0
388 KappaRV(i,j,k) = 0. _d 0
389 ENDDO
390 ENDDO
391 ENDDO
392 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
393 C-- Calculate the total vertical viscosity
394 #ifdef ALLOW_AUTODIFF
395 IF ( momViscosity ) THEN
396 #else
397 ELSE
398 #endif
399 CALL CALC_VISCOSITY(
400 I bi,bj, iMin,iMax,jMin,jMax,
401 O KappaRU, KappaRV,
402 I myThid )
403 ENDIF
404 #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
405
406 #ifdef ALLOW_AUTODIFF_TAMC
407 CADJ STORE KappaRU(:,:,:)
408 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
409 CADJ STORE KappaRV(:,:,:)
410 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
411 #endif /* ALLOW_AUTODIFF_TAMC */
412
413 #ifdef ALLOW_OBCS
414 C-- For Stevens boundary conditions velocities need to be extrapolated
415 C (copied) to a narrow strip outside the domain
416 IF ( useOBCS ) THEN
417 CALL OBCS_COPY_UV_N(
418 U uVel(1-OLx,1-OLy,1,bi,bj),
419 U vVel(1-OLx,1-OLy,1,bi,bj),
420 I Nr, bi, bj, myThid )
421 ENDIF
422 #endif /* ALLOW_OBCS */
423
424 #ifdef ALLOW_EDDYPSI
425 CALL CALC_EDDY_STRESS(bi,bj,myThid)
426 #endif
427
428 C-- Start of dynamics loop
429 DO k=1,Nr
430
431 C-- km1 Points to level above k (=k-1)
432 C-- kup Cycles through 1,2 to point to layer above
433 C-- kDown Cycles through 2,1 to point to current layer
434
435 km1 = MAX(1,k-1)
436 kp1 = MIN(k+1,Nr)
437 kup = 1+MOD(k+1,2)
438 kDown= 1+MOD(k,2)
439
440 #ifdef ALLOW_AUTODIFF_TAMC
441 kkey = (idynkey-1)*Nr + k
442 c
443 CADJ STORE totPhiHyd (:,:,k,bi,bj)
444 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
445 CADJ STORE phiHydLow (:,:,bi,bj)
446 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
447 CADJ STORE theta (:,:,k,bi,bj)
448 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
449 CADJ STORE salt (:,:,k,bi,bj)
450 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
451 CADJ STORE gT(:,:,k,bi,bj)
452 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
453 CADJ STORE gS(:,:,k,bi,bj)
454 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
455 # ifdef NONLIN_FRSURF
456 cph-test
457 CADJ STORE phiHydC (:,:)
458 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
459 CADJ STORE phiHydF (:,:)
460 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
461 CADJ STORE guDissip (:,:)
462 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
463 CADJ STORE gvDissip (:,:)
464 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
465 CADJ STORE fVerU (:,:,:)
466 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
467 CADJ STORE fVerV (:,:,:)
468 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
469 CADJ STORE gU(:,:,k,bi,bj)
470 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
471 CADJ STORE gV(:,:,k,bi,bj)
472 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
473 # ifndef ALLOW_ADAMSBASHFORTH_3
474 CADJ STORE guNm1(:,:,k,bi,bj)
475 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
476 CADJ STORE gvNm1(:,:,k,bi,bj)
477 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
478 # else
479 CADJ STORE guNm(:,:,k,bi,bj,1)
480 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
481 CADJ STORE guNm(:,:,k,bi,bj,2)
482 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
483 CADJ STORE gvNm(:,:,k,bi,bj,1)
484 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
485 CADJ STORE gvNm(:,:,k,bi,bj,2)
486 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
487 # endif
488 # ifdef ALLOW_CD_CODE
489 CADJ STORE uNM1(:,:,k,bi,bj)
490 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
491 CADJ STORE vNM1(:,:,k,bi,bj)
492 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
493 CADJ STORE uVelD(:,:,k,bi,bj)
494 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
495 CADJ STORE vVelD(:,:,k,bi,bj)
496 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
497 # endif
498 # endif
499 # ifdef ALLOW_DEPTH_CONTROL
500 CADJ STORE fVerU (:,:,:)
501 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
502 CADJ STORE fVerV (:,:,:)
503 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
504 # endif
505 #endif /* ALLOW_AUTODIFF_TAMC */
506
507 C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
508 IF ( implicitIntGravWave ) THEN
509 CALL CALC_PHI_HYD(
510 I bi,bj,iMin,iMax,jMin,jMax,k,
511 I gT, gS,
512 U phiHydF,
513 O phiHydC, dPhiHydX, dPhiHydY,
514 I myTime, myIter, myThid )
515 ELSE
516 CALL CALC_PHI_HYD(
517 I bi,bj,iMin,iMax,jMin,jMax,k,
518 I theta, salt,
519 U phiHydF,
520 O phiHydC, dPhiHydX, dPhiHydY,
521 I myTime, myIter, myThid )
522 ENDIF
523 #ifdef ALLOW_DIAGNOSTICS
524 IF ( dPhiHydDiagIsOn ) THEN
525 tmpFac = -1. _d 0
526 CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
527 & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
528 CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
529 & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
530 ENDIF
531 #endif /* ALLOW_DIAGNOSTICS */
532
533 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
534 C and step forward storing the result in gU, gV, etc...
535 IF ( momStepping ) THEN
536 #ifdef ALLOW_AUTODIFF_TAMC
537 # ifdef NONLIN_FRSURF
538 # if (defined ALLOW_MOM_FLUXFORM) && !(defined DISABLE_RSTAR_CODE)
539 CADJ STORE dWtransC(:,:,bi,bj)
540 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
541 CADJ STORE dWtransU(:,:,bi,bj)
542 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
543 CADJ STORE dWtransV(:,:,bi,bj)
544 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
545 # endif
546 CADJ STORE fVerU(:,:,:)
547 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
548 CADJ STORE fVerV(:,:,:)
549 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
550 # endif /* NONLIN_FRSURF */
551 #endif /* ALLOW_AUTODIFF_TAMC */
552 IF (.NOT. vectorInvariantMomentum) THEN
553 #ifdef ALLOW_MOM_FLUXFORM
554 CALL MOM_FLUXFORM(
555 I bi,bj,k,iMin,iMax,jMin,jMax,
556 I KappaRU, KappaRV,
557 U fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
558 O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
559 O guDissip, gvDissip,
560 I myTime, myIter, myThid)
561 #endif
562 ELSE
563 #ifdef ALLOW_MOM_VECINV
564 CALL MOM_VECINV(
565 I bi,bj,k,iMin,iMax,jMin,jMax,
566 I KappaRU, KappaRV,
567 I fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
568 O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
569 O guDissip, gvDissip,
570 I myTime, myIter, myThid)
571 #endif
572 ENDIF
573
574 CALL TIMESTEP(
575 I bi,bj,iMin,iMax,jMin,jMax,k,
576 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
577 I guDissip, gvDissip,
578 I myTime, myIter, myThid)
579
580 ENDIF
581
582 C-- end of dynamics k loop (1:Nr)
583 ENDDO
584
585 C-- Implicit Vertical advection & viscosity
586 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
587 defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
588 IF ( momImplVertAdv ) THEN
589 CALL MOM_U_IMPLICIT_R( kappaRU,
590 I bi, bj, myTime, myIter, myThid )
591 CALL MOM_V_IMPLICIT_R( kappaRV,
592 I bi, bj, myTime, myIter, myThid )
593 ELSEIF ( implicitViscosity ) THEN
594 #else /* INCLUDE_IMPLVERTADV_CODE */
595 IF ( implicitViscosity ) THEN
596 #endif /* INCLUDE_IMPLVERTADV_CODE */
597 #ifdef ALLOW_AUTODIFF_TAMC
598 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
599 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
600 #endif /* ALLOW_AUTODIFF_TAMC */
601 CALL IMPLDIFF(
602 I bi, bj, iMin, iMax, jMin, jMax,
603 I -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
604 U gU,
605 I myThid )
606 #ifdef ALLOW_AUTODIFF_TAMC
607 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
608 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
609 #endif /* ALLOW_AUTODIFF_TAMC */
610 CALL IMPLDIFF(
611 I bi, bj, iMin, iMax, jMin, jMax,
612 I -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
613 U gV,
614 I myThid )
615 ENDIF
616
617 #ifdef ALLOW_OBCS
618 C-- Apply open boundary conditions
619 IF ( useOBCS ) THEN
620 C-- but first save intermediate velocities to be used in the
621 C next time step for the Stevens boundary conditions
622 CALL OBCS_SAVE_UV_N(
623 I bi, bj, iMin, iMax, jMin, jMax, 0,
624 I gU, gV, myThid )
625 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
626 ENDIF
627 #endif /* ALLOW_OBCS */
628
629 #ifdef ALLOW_CD_CODE
630 IF (implicitViscosity.AND.useCDscheme) THEN
631 #ifdef ALLOW_AUTODIFF_TAMC
632 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
633 #endif /* ALLOW_AUTODIFF_TAMC */
634 CALL IMPLDIFF(
635 I bi, bj, iMin, iMax, jMin, jMax,
636 I 0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
637 U vVelD,
638 I myThid )
639 #ifdef ALLOW_AUTODIFF_TAMC
640 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
641 #endif /* ALLOW_AUTODIFF_TAMC */
642 CALL IMPLDIFF(
643 I bi, bj, iMin, iMax, jMin, jMax,
644 I 0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
645 U uVelD,
646 I myThid )
647 ENDIF
648 #endif /* ALLOW_CD_CODE */
649 C-- End implicit Vertical advection & viscosity
650
651 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
652
653 #ifdef ALLOW_NONHYDROSTATIC
654 C-- Step forward W field in N-H algorithm
655 IF ( nonHydrostatic ) THEN
656 #ifdef ALLOW_DEBUG
657 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
658 #endif
659 CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
660 CALL CALC_GW(
661 I bi,bj, KappaRU, KappaRV,
662 I myTime, myIter, myThid )
663 ENDIF
664 IF ( nonHydrostatic.OR.implicitIntGravWave )
665 & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
666 IF ( nonHydrostatic )
667 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
668 #endif
669
670 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
671
672 C- end of bi,bj loops
673 ENDDO
674 ENDDO
675
676 #ifdef ALLOW_OBCS
677 IF (useOBCS) THEN
678 CALL OBCS_EXCHANGES( myThid )
679 ENDIF
680 #endif
681
682 Cml(
683 C In order to compare the variance of phiHydLow of a p/z-coordinate
684 C run with etaH of a z/p-coordinate run the drift of phiHydLow
685 C has to be removed by something like the following subroutine:
686 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
687 C & 'phiHydLow', myTime, myThid )
688 Cml)
689
690 #ifdef ALLOW_DIAGNOSTICS
691 IF ( useDiagnostics ) THEN
692
693 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
694 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
695
696 tmpFac = 1. _d 0
697 CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
698 & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
699
700 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
701 & 'PHIBOTSQ',0, 1,0,1,1,myThid)
702
703 ENDIF
704 #endif /* ALLOW_DIAGNOSTICS */
705
706 #ifdef ALLOW_DEBUG
707 IF ( debugLevel .GE. debLevD ) THEN
708 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
709 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
710 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
711 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
712 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
713 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
714 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
715 CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
716 CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
717 CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
718 #ifndef ALLOW_ADAMSBASHFORTH_3
719 CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
720 CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
721 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
722 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
723 #endif
724 ENDIF
725 #endif
726
727 #ifdef DYNAMICS_GUGV_EXCH_CHECK
728 C- jmc: For safety checking only: This Exchange here should not change
729 C the solution. If solution changes, it means something is wrong,
730 C but it does not mean that it is less wrong with this exchange.
731 IF ( debugLevel .GE. debLevE ) THEN
732 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
733 ENDIF
734 #endif
735
736 #ifdef ALLOW_DEBUG
737 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
738 #endif
739
740 RETURN
741 END

  ViewVC Help
Powered by ViewVC 1.1.22