/[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.165 - (show annotations) (download)
Sat Aug 3 01:38:17 2013 UTC (10 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64m, checkpoint64l, checkpoint64n
Changes since 1.164: +32 -29 lines
skip the call to CALC_VISCOSITY if momViscosity=F

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

  ViewVC Help
Powered by ViewVC 1.1.22