/[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.173 - (show annotations) (download)
Fri Aug 15 19:22:06 2014 UTC (9 years, 10 months ago) by jmc
Branch: MAIN
Changes since 1.172: +1 -7 lines
remove gT,gS

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

  ViewVC Help
Powered by ViewVC 1.1.22