/[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.146 - (show annotations) (download)
Fri May 14 23:21:02 2010 UTC (14 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62g, checkpoint62i, checkpoint62h
Changes since 1.145: +3 -1 lines
move initialisation of dPhiHydX,dPhiHydY (= output of S/R CALC_GRAD_PHI_HYD)
 from dynamics.F to calc_grad_phi_hyd.F

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

  ViewVC Help
Powered by ViewVC 1.1.22