/[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.133 - (show annotations) (download)
Wed May 31 19:53:28 2006 UTC (17 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58f_post, checkpoint58g_post
Changes since 1.132: +4 -2 lines
Enable variable grid visc. for adjoint, but still exclude
calcLeith and calcSmag

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

  ViewVC Help
Powered by ViewVC 1.1.22