/[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.132 - (show annotations) (download)
Wed May 3 23:34:41 2006 UTC (18 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58e_post
Changes since 1.131: +43 -6 lines
o Now rstar adjoint.

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

  ViewVC Help
Powered by ViewVC 1.1.22