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