/[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.127 - (show annotations) (download)
Thu Dec 15 21:09:00 2005 UTC (18 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58, checkpoint58a_post, checkpoint57z_post
Changes since 1.126: +2 -1 lines
isolate forward stepping of wVel in new S/R (previously part of calc_gw)

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

  ViewVC Help
Powered by ViewVC 1.1.22