/[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.175 - (show annotations) (download)
Tue Sep 9 22:32:09 2014 UTC (9 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65h, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e
Changes since 1.174: +2 -2 lines
Replace ALLOW_AUTODIFF_TAMC by ALLOW_AUTODIFF (except for tape/storage
 which are specific to TAF/TAMC).

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

  ViewVC Help
Powered by ViewVC 1.1.22