/[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.137 - (show annotations) (download)
Tue Jan 30 03:18:13 2007 UTC (17 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint58w_post, checkpoint59e, checkpoint59d, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint58v_post, checkpoint58x_post
Changes since 1.136: +1 -34 lines
Change storing for rstar adjoint.

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

  ViewVC Help
Powered by ViewVC 1.1.22