/[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.131 - (show annotations) (download)
Wed Mar 29 17:00:39 2006 UTC (18 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58d_post
Changes since 1.130: +15 -1 lines
Adding relevant headers for obcs+ptracers adjoint.

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

  ViewVC Help
Powered by ViewVC 1.1.22