/[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.169 - (show annotations) (download)
Thu Feb 6 23:16:46 2014 UTC (10 years, 3 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint64u
Changes since 1.168: +3 -1 lines
- avoid recomputation.

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

  ViewVC Help
Powered by ViewVC 1.1.22