/[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.157 - (show annotations) (download)
Tue May 24 20:25:33 2011 UTC (13 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62y
Changes since 1.156: +2 -2 lines
split "OBCS.h" into 4 separated header files (OBCS_PARAMS,GRID,FIELDS,SEAICE)

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.156 2011/05/23 00:41:09 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_FIELDS.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 (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
251 #endif
252
253 #ifdef ALLOW_DIAGNOSTICS
254 dPhiHydDiagIsOn = .FALSE.
255 IF ( useDiagnostics )
256 & dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
257 & .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
258 #endif
259
260 C-- Call to routine for calculation of
261 C Eliassen-Palm-flux-forced U-tendency,
262 C if desired:
263 #ifdef INCLUDE_EP_FORCING_CODE
264 CALL CALC_EP_FORCING(myThid)
265 #endif
266
267 #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
268 CALL DUMMY_IN_DYNAMICS( mytime, myiter, myThid )
269 #endif
270
271 #ifdef ALLOW_AUTODIFF_TAMC
272 C-- HPF directive to help TAMC
273 CHPF$ INDEPENDENT
274 #endif /* ALLOW_AUTODIFF_TAMC */
275
276 DO bj=myByLo(myThid),myByHi(myThid)
277
278 #ifdef ALLOW_AUTODIFF_TAMC
279 C-- HPF directive to help TAMC
280 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
281 CHPF$& ,phiHydF
282 CHPF$& ,KappaRU,KappaRV
283 CHPF$& )
284 #endif /* ALLOW_AUTODIFF_TAMC */
285
286 DO bi=myBxLo(myThid),myBxHi(myThid)
287
288 #ifdef ALLOW_AUTODIFF_TAMC
289 act1 = bi - myBxLo(myThid)
290 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
291 act2 = bj - myByLo(myThid)
292 max2 = myByHi(myThid) - myByLo(myThid) + 1
293 act3 = myThid - 1
294 max3 = nTx*nTy
295 act4 = ikey_dynamics - 1
296 idynkey = (act1 + 1) + act2*max1
297 & + act3*max1*max2
298 & + act4*max1*max2*max3
299 #endif /* ALLOW_AUTODIFF_TAMC */
300
301 C-- Set up work arrays with valid (i.e. not NaN) values
302 C These inital values do not alter the numerical results. They
303 C just ensure that all memory references are to valid floating
304 C point numbers. This prevents spurious hardware signals due to
305 C uninitialised but inert locations.
306
307 #ifdef ALLOW_AUTODIFF_TAMC
308 DO k=1,Nr
309 DO j=1-OLy,sNy+OLy
310 DO i=1-OLx,sNx+OLx
311 KappaRU(i,j,k) = 0. _d 0
312 KappaRV(i,j,k) = 0. _d 0
313 cph(
314 c-- need some re-initialisation here to break dependencies
315 cph)
316 gU(i,j,k,bi,bj) = 0. _d 0
317 gV(i,j,k,bi,bj) = 0. _d 0
318 ENDDO
319 ENDDO
320 ENDDO
321 #endif /* ALLOW_AUTODIFF_TAMC */
322 DO j=1-OLy,sNy+OLy
323 DO i=1-OLx,sNx+OLx
324 fVerU (i,j,1) = 0. _d 0
325 fVerU (i,j,2) = 0. _d 0
326 fVerV (i,j,1) = 0. _d 0
327 fVerV (i,j,2) = 0. _d 0
328 phiHydF (i,j) = 0. _d 0
329 phiHydC (i,j) = 0. _d 0
330 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
331 dPhiHydX(i,j) = 0. _d 0
332 dPhiHydY(i,j) = 0. _d 0
333 #endif
334 phiSurfX(i,j) = 0. _d 0
335 phiSurfY(i,j) = 0. _d 0
336 guDissip(i,j) = 0. _d 0
337 gvDissip(i,j) = 0. _d 0
338 #ifdef ALLOW_AUTODIFF_TAMC
339 phiHydLow(i,j,bi,bj) = 0. _d 0
340 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
341 # ifndef DISABLE_RSTAR_CODE
342 dWtransC(i,j,bi,bj) = 0. _d 0
343 dWtransU(i,j,bi,bj) = 0. _d 0
344 dWtransV(i,j,bi,bj) = 0. _d 0
345 # endif
346 # endif
347 #endif
348 ENDDO
349 ENDDO
350
351 C-- Start computation of dynamics
352 iMin = 0
353 iMax = sNx+1
354 jMin = 0
355 jMax = sNy+1
356
357 #ifdef ALLOW_AUTODIFF_TAMC
358 CADJ STORE wvel (:,:,:,bi,bj) =
359 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
360 #endif /* ALLOW_AUTODIFF_TAMC */
361
362 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
363 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
364 IF (implicSurfPress.NE.1.) THEN
365 CALL CALC_GRAD_PHI_SURF(
366 I bi,bj,iMin,iMax,jMin,jMax,
367 I etaN,
368 O phiSurfX,phiSurfY,
369 I myThid )
370 ENDIF
371
372 #ifdef ALLOW_AUTODIFF_TAMC
373 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
374 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
375 #ifdef ALLOW_KPP
376 CADJ STORE KPPviscAz (:,:,:,bi,bj)
377 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
378 #endif /* ALLOW_KPP */
379 #endif /* ALLOW_AUTODIFF_TAMC */
380
381 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
382 C-- Calculate the total vertical viscosity
383 CALL CALC_VISCOSITY(
384 I bi,bj, iMin,iMax,jMin,jMax,
385 O KappaRU, KappaRV,
386 I myThid )
387 #else
388 DO k=1,Nr
389 DO j=1-OLy,sNy+OLy
390 DO i=1-OLx,sNx+OLx
391 KappaRU(i,j,k) = 0. _d 0
392 KappaRV(i,j,k) = 0. _d 0
393 ENDDO
394 ENDDO
395 ENDDO
396 #endif
397
398 #ifdef ALLOW_AUTODIFF_TAMC
399 CADJ STORE KappaRU(:,:,:)
400 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
401 CADJ STORE KappaRV(:,:,:)
402 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
403 #endif /* ALLOW_AUTODIFF_TAMC */
404
405 C-- Start of dynamics loop
406 DO k=1,Nr
407
408 C-- km1 Points to level above k (=k-1)
409 C-- kup Cycles through 1,2 to point to layer above
410 C-- kDown Cycles through 2,1 to point to current layer
411
412 km1 = MAX(1,k-1)
413 kp1 = MIN(k+1,Nr)
414 kup = 1+MOD(k+1,2)
415 kDown= 1+MOD(k,2)
416
417 #ifdef ALLOW_AUTODIFF_TAMC
418 kkey = (idynkey-1)*Nr + k
419 c
420 CADJ STORE totphihyd (:,:,k,bi,bj)
421 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
422 CADJ STORE phihydlow (:,:,bi,bj)
423 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
424 CADJ STORE theta (:,:,k,bi,bj)
425 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
426 CADJ STORE salt (:,:,k,bi,bj)
427 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
428 CADJ STORE gt(:,:,k,bi,bj)
429 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
430 CADJ STORE gs(:,:,k,bi,bj)
431 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
432 # ifdef NONLIN_FRSURF
433 cph-test
434 CADJ STORE phiHydC (:,:)
435 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
436 CADJ STORE phiHydF (:,:)
437 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
438 CADJ STORE gudissip (:,:)
439 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
440 CADJ STORE gvdissip (:,:)
441 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
442 CADJ STORE fVerU (:,:,:)
443 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
444 CADJ STORE fVerV (:,:,:)
445 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
446 CADJ STORE gu(:,:,k,bi,bj)
447 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
448 CADJ STORE gv(:,:,k,bi,bj)
449 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
450 # ifndef ALLOW_ADAMSBASHFORTH_3
451 CADJ STORE gunm1(:,:,k,bi,bj)
452 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
453 CADJ STORE gvnm1(:,:,k,bi,bj)
454 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
455 # else
456 CADJ STORE gunm(:,:,k,bi,bj,1)
457 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
458 CADJ STORE gunm(:,:,k,bi,bj,2)
459 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
460 CADJ STORE gvnm(:,:,k,bi,bj,1)
461 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
462 CADJ STORE gvnm(:,:,k,bi,bj,2)
463 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
464 # endif
465 # ifdef ALLOW_CD_CODE
466 CADJ STORE unm1(:,:,k,bi,bj)
467 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
468 CADJ STORE vnm1(:,:,k,bi,bj)
469 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
470 CADJ STORE uVelD(:,:,k,bi,bj)
471 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
472 CADJ STORE vVelD(:,:,k,bi,bj)
473 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
474 # endif
475 # endif
476 # ifdef ALLOW_DEPTH_CONTROL
477 CADJ STORE fVerU (:,:,:)
478 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
479 CADJ STORE fVerV (:,:,:)
480 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
481 # endif
482 #endif /* ALLOW_AUTODIFF_TAMC */
483
484 C-- Integrate hydrostatic balance for phiHyd with BC of
485 C phiHyd(z=0)=0
486 IF ( implicitIntGravWave ) THEN
487 CALL CALC_PHI_HYD(
488 I bi,bj,iMin,iMax,jMin,jMax,k,
489 I gT, gS,
490 U phiHydF,
491 O phiHydC, dPhiHydX, dPhiHydY,
492 I myTime, myIter, myThid )
493 ELSE
494 CALL CALC_PHI_HYD(
495 I bi,bj,iMin,iMax,jMin,jMax,k,
496 I theta, salt,
497 U phiHydF,
498 O phiHydC, dPhiHydX, dPhiHydY,
499 I myTime, myIter, myThid )
500 ENDIF
501 #ifdef ALLOW_DIAGNOSTICS
502 IF ( dPhiHydDiagIsOn ) THEN
503 tmpFac = -1. _d 0
504 CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
505 & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
506 CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
507 & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
508 ENDIF
509 #endif /* ALLOW_DIAGNOSTICS */
510
511 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
512 C and step forward storing the result in gU, gV, etc...
513 IF ( momStepping ) THEN
514 #ifdef ALLOW_AUTODIFF_TAMC
515 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
516 # ifndef DISABLE_RSTAR_CODE
517 CADJ STORE dWtransC(:,:,bi,bj)
518 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
519 CADJ STORE dWtransU(:,:,bi,bj)
520 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
521 CADJ STORE dWtransV(:,:,bi,bj)
522 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
523 # endif
524 # endif
525 #endif
526 IF (.NOT. vectorInvariantMomentum) THEN
527 #ifdef ALLOW_MOM_FLUXFORM
528 C
529 CALL MOM_FLUXFORM(
530 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
531 I KappaRU, KappaRV,
532 U fVerU, fVerV,
533 O guDissip, gvDissip,
534 I myTime, myIter, myThid)
535 #endif
536 ELSE
537 #ifdef ALLOW_MOM_VECINV
538 C
539 # ifdef ALLOW_AUTODIFF_TAMC
540 # ifdef NONLIN_FRSURF
541 CADJ STORE fVerU(:,:,:)
542 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
543 CADJ STORE fVerV(:,:,:)
544 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
545 # endif
546 # endif /* ALLOW_AUTODIFF_TAMC */
547 C
548 CALL MOM_VECINV(
549 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
550 I KappaRU, KappaRV,
551 U fVerU, fVerV,
552 O guDissip, gvDissip,
553 I myTime, myIter, myThid)
554 #endif
555 ENDIF
556 C
557 CALL TIMESTEP(
558 I bi,bj,iMin,iMax,jMin,jMax,k,
559 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
560 I guDissip, gvDissip,
561 I myTime, myIter, myThid)
562
563 ENDIF
564
565 C-- end of dynamics k loop (1:Nr)
566 ENDDO
567
568 C-- Implicit Vertical advection & viscosity
569 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
570 defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
571 IF ( momImplVertAdv ) THEN
572 CALL MOM_U_IMPLICIT_R( kappaRU,
573 I bi, bj, myTime, myIter, myThid )
574 CALL MOM_V_IMPLICIT_R( kappaRV,
575 I bi, bj, myTime, myIter, myThid )
576 ELSEIF ( implicitViscosity ) THEN
577 #else /* INCLUDE_IMPLVERTADV_CODE */
578 IF ( implicitViscosity ) THEN
579 #endif /* INCLUDE_IMPLVERTADV_CODE */
580 #ifdef ALLOW_AUTODIFF_TAMC
581 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
582 CADJ STORE gU(:,:,:,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 -1, KappaRU,recip_HFacW,
587 U gU,
588 I myThid )
589 #ifdef ALLOW_AUTODIFF_TAMC
590 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
591 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
592 #endif /* ALLOW_AUTODIFF_TAMC */
593 CALL IMPLDIFF(
594 I bi, bj, iMin, iMax, jMin, jMax,
595 I -2, KappaRV,recip_HFacS,
596 U gV,
597 I myThid )
598 ENDIF
599
600 #ifdef ALLOW_OBCS
601 C-- Apply open boundary conditions
602 IF ( useOBCS ) THEN
603 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
604 ENDIF
605 #endif /* ALLOW_OBCS */
606
607 #ifdef ALLOW_CD_CODE
608 IF (implicitViscosity.AND.useCDscheme) THEN
609 #ifdef ALLOW_AUTODIFF_TAMC
610 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
611 #endif /* ALLOW_AUTODIFF_TAMC */
612 CALL IMPLDIFF(
613 I bi, bj, iMin, iMax, jMin, jMax,
614 I 0, KappaRU,recip_HFacW,
615 U vVelD,
616 I myThid )
617 #ifdef ALLOW_AUTODIFF_TAMC
618 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
619 #endif /* ALLOW_AUTODIFF_TAMC */
620 CALL IMPLDIFF(
621 I bi, bj, iMin, iMax, jMin, jMax,
622 I 0, KappaRV,recip_HFacS,
623 U uVelD,
624 I myThid )
625 ENDIF
626 #endif /* ALLOW_CD_CODE */
627 C-- End implicit Vertical advection & viscosity
628
629 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
630
631 #ifdef ALLOW_NONHYDROSTATIC
632 C-- Step forward W field in N-H algorithm
633 IF ( nonHydrostatic ) THEN
634 #ifdef ALLOW_DEBUG
635 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
636 #endif
637 CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
638 CALL CALC_GW(
639 I bi,bj, KappaRU, KappaRV,
640 I myTime, myIter, myThid )
641 ENDIF
642 IF ( nonHydrostatic.OR.implicitIntGravWave )
643 & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
644 IF ( nonHydrostatic )
645 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
646 #endif
647
648 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
649
650 C- end of bi,bj loops
651 ENDDO
652 ENDDO
653
654 #ifdef ALLOW_OBCS
655 IF (useOBCS) THEN
656 CALL OBCS_EXCHANGES( myThid )
657 ENDIF
658 #endif
659
660 Cml(
661 C In order to compare the variance of phiHydLow of a p/z-coordinate
662 C run with etaH of a z/p-coordinate run the drift of phiHydLow
663 C has to be removed by something like the following subroutine:
664 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
665 C & 'phiHydLow', myTime, myThid )
666 Cml)
667
668 #ifdef ALLOW_DIAGNOSTICS
669 IF ( useDiagnostics ) THEN
670
671 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
672 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
673
674 tmpFac = 1. _d 0
675 CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
676 & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
677
678 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
679 & 'PHIBOTSQ',0, 1,0,1,1,myThid)
680
681 ENDIF
682 #endif /* ALLOW_DIAGNOSTICS */
683
684 #ifdef ALLOW_DEBUG
685 IF ( debugLevel .GE. debLevB ) THEN
686 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
687 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
688 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
689 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
690 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
691 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
692 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
693 CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
694 CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
695 CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
696 #ifndef ALLOW_ADAMSBASHFORTH_3
697 CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
698 CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
699 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
700 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
701 #endif
702 ENDIF
703 #endif
704
705 #ifdef DYNAMICS_GUGV_EXCH_CHECK
706 C- jmc: For safety checking only: This Exchange here should not change
707 C the solution. If solution changes, it means something is wrong,
708 C but it does not mean that it is less wrong with this exchange.
709 IF ( debugLevel .GT. debLevB ) THEN
710 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
711 ENDIF
712 #endif
713
714 #ifdef ALLOW_DEBUG
715 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
716 #endif
717
718 RETURN
719 END

  ViewVC Help
Powered by ViewVC 1.1.22