/[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.54.2.11 - (show annotations) (download)
Thu Jan 25 19:43:32 2001 UTC (23 years, 3 months ago) by adcroft
Branch: branch-atmos-merge
Changes since 1.54.2.10: +3 -4 lines
o Removed array phiHydInterface. Why have more arguments than are necessary.
It made the "finite volume" integration easier but wasn't used for
the default energy conserving method. The fv is still available in
comments but has been coded without the phiHydInterface array.
o Put a safety "IF" in front of k+1 references.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.54.2.10 2001/01/17 15:05:59 jmc Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6 C /==========================================================\
7 C | SUBROUTINE DYNAMICS |
8 C | o Controlling routine for the explicit part of the model |
9 C | dynamics. |
10 C |==========================================================|
11 C | This routine evaluates the "dynamics" terms for each |
12 C | block of ocean in turn. Because the blocks of ocean have |
13 C | overlap regions they are independent of one another. |
14 C | If terms involving lateral integrals are needed in this |
15 C | routine care will be needed. Similarly finite-difference |
16 C | operations with stencils wider than the overlap region |
17 C | require special consideration. |
18 C | Notes |
19 C | ===== |
20 C | C*P* comments indicating place holders for which code is |
21 C | presently being developed. |
22 C \==========================================================/
23 IMPLICIT NONE
24
25 C == Global variables ===
26 #include "SIZE.h"
27 #include "EEPARAMS.h"
28 #include "CG2D.h"
29 #include "PARAMS.h"
30 #include "DYNVARS.h"
31 #include "GRID.h"
32
33 #ifdef ALLOW_AUTODIFF_TAMC
34 # include "tamc.h"
35 # include "tamc_keys.h"
36 #endif /* ALLOW_AUTODIFF_TAMC */
37
38 #ifdef ALLOW_KPP
39 # include "KPP.h"
40 #endif
41
42 C == Routine arguments ==
43 C myTime - Current time in simulation
44 C myIter - Current iteration number in simulation
45 C myThid - Thread number for this instance of the routine.
46 _RL myTime
47 INTEGER myIter
48 INTEGER myThid
49
50 C == Local variables
51 C xA, yA - Per block temporaries holding face areas
52 C uTrans, vTrans, rTrans - Per block temporaries holding flow
53 C transport
54 C rVel o uTrans: Zonal transport
55 C o vTrans: Meridional transport
56 C o rTrans: Vertical transport
57 C o rVel: Vertical velocity at upper and
58 C lower cell faces.
59 C maskC,maskUp o maskC: land/water mask for tracer cells
60 C o maskUp: land/water mask for W points
61 C fVer[STUV] o fVer: Vertical flux term - note fVer
62 C is "pipelined" in the vertical
63 C so we need an fVer for each
64 C variable.
65 C rhoK, rhoKM1 - Density at current level, and level above
66 C phiHyd - Hydrostatic part of the potential phiHydi.
67 C In z coords phiHydiHyd is the hydrostatic
68 C pressure anomaly
69 C In p coords phiHydiHyd is the geopotential
70 C surface height
71 C anomaly.
72 C etaSurfX, - Holds surface elevation gradient in X and Y.
73 C etaSurfY
74 C KappaRT, - Total diffusion in vertical for T and S.
75 C KappaRS (background + spatially varying, isopycnal term).
76 C iMin, iMax - Ranges and sub-block indices on which calculations
77 C jMin, jMax are applied.
78 C bi, bj
79 C k, kup, - Index for layer above and below. kup and kDown
80 C kDown, km1 are switched with layer to be the appropriate
81 C index into fVerTerm.
82 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
88 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
91 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
92 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
93 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
94 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
98 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
99 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
100 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
101 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
103 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
104
105 C This is currently also used by IVDC and Diagnostics
106 C #ifdef INCLUDE_CONVECT_CALL
107 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
108 C #endif
109
110 INTEGER iMin, iMax
111 INTEGER jMin, jMax
112 INTEGER bi, bj
113 INTEGER i, j
114 INTEGER k, km1, kup, kDown
115
116 #ifdef ALLOW_AUTODIFF_TAMC
117 INTEGER isbyte
118 PARAMETER( isbyte = 4 )
119
120 INTEGER act1, act2, act3, act4
121 INTEGER max1, max2, max3
122 INTEGER iikey, kkey
123 INTEGER maximpl
124 #endif /* ALLOW_AUTODIFF_TAMC */
125
126 C--- The algorithm...
127 C
128 C "Correction Step"
129 C =================
130 C Here we update the horizontal velocities with the surface
131 C pressure such that the resulting flow is either consistent
132 C with the free-surface evolution or the rigid-lid:
133 C U[n] = U* + dt x d/dx P
134 C V[n] = V* + dt x d/dy P
135 C
136 C "Calculation of Gs"
137 C ===================
138 C This is where all the accelerations and tendencies (ie.
139 C physics, parameterizations etc...) are calculated
140 C rVel = sum_r ( div. u[n] )
141 C rho = rho ( theta[n], salt[n] )
142 C b = b(rho, theta)
143 C K31 = K31 ( rho )
144 C Gu[n] = Gu( u[n], v[n], rVel, b, ... )
145 C Gv[n] = Gv( u[n], v[n], rVel, b, ... )
146 C Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
147 C Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
148 C
149 C "Time-stepping" or "Prediction"
150 C ================================
151 C The models variables are stepped forward with the appropriate
152 C time-stepping scheme (currently we use Adams-Bashforth II)
153 C - For momentum, the result is always *only* a "prediction"
154 C in that the flow may be divergent and will be "corrected"
155 C later with a surface pressure gradient.
156 C - Normally for tracers the result is the new field at time
157 C level [n+1} *BUT* in the case of implicit diffusion the result
158 C is also *only* a prediction.
159 C - We denote "predictors" with an asterisk (*).
160 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
161 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
162 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
163 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
164 C With implicit diffusion:
165 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
166 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
167 C (1 + dt * K * d_zz) theta[n] = theta*
168 C (1 + dt * K * d_zz) salt[n] = salt*
169 C---
170
171 #ifdef ALLOW_AUTODIFF_TAMC
172 C-- dummy statement to end declaration part
173 ikey = 1
174 #endif /* ALLOW_AUTODIFF_TAMC */
175
176 C-- Set up work arrays with valid (i.e. not NaN) values
177 C These inital values do not alter the numerical results. They
178 C just ensure that all memory references are to valid floating
179 C point numbers. This prevents spurious hardware signals due to
180 C uninitialised but inert locations.
181 DO j=1-OLy,sNy+OLy
182 DO i=1-OLx,sNx+OLx
183 xA(i,j) = 0. _d 0
184 yA(i,j) = 0. _d 0
185 uTrans(i,j) = 0. _d 0
186 vTrans(i,j) = 0. _d 0
187 DO k=1,Nr
188 KappaRU(i,j,k) = 0. _d 0
189 KappaRV(i,j,k) = 0. _d 0
190 sigmaX(i,j,k) = 0. _d 0
191 sigmaY(i,j,k) = 0. _d 0
192 sigmaR(i,j,k) = 0. _d 0
193 ENDDO
194 rhoKM1 (i,j) = 0. _d 0
195 rhok (i,j) = 0. _d 0
196 maskC (i,j) = 0. _d 0
197 ENDDO
198 ENDDO
199
200
201 #ifdef ALLOW_AUTODIFF_TAMC
202 C-- HPF directive to help TAMC
203 CHPF$ INDEPENDENT
204 #endif /* ALLOW_AUTODIFF_TAMC */
205
206 DO bj=myByLo(myThid),myByHi(myThid)
207
208 #ifdef ALLOW_AUTODIFF_TAMC
209 C-- HPF directive to help TAMC
210 CHPF$ INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
211 CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA
212 CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
213 CHPF$& )
214 #endif /* ALLOW_AUTODIFF_TAMC */
215
216 DO bi=myBxLo(myThid),myBxHi(myThid)
217
218 #ifdef ALLOW_AUTODIFF_TAMC
219 act1 = bi - myBxLo(myThid)
220 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
221
222 act2 = bj - myByLo(myThid)
223 max2 = myByHi(myThid) - myByLo(myThid) + 1
224
225 act3 = myThid - 1
226 max3 = nTx*nTy
227
228 act4 = ikey_dynamics - 1
229
230 ikey = (act1 + 1) + act2*max1
231 & + act3*max1*max2
232 & + act4*max1*max2*max3
233 #endif /* ALLOW_AUTODIFF_TAMC */
234
235 C-- Set up work arrays that need valid initial values
236 DO j=1-OLy,sNy+OLy
237 DO i=1-OLx,sNx+OLx
238 rTrans(i,j) = 0. _d 0
239 rVel (i,j,1) = 0. _d 0
240 rVel (i,j,2) = 0. _d 0
241 fVerT (i,j,1) = 0. _d 0
242 fVerT (i,j,2) = 0. _d 0
243 fVerS (i,j,1) = 0. _d 0
244 fVerS (i,j,2) = 0. _d 0
245 fVerU (i,j,1) = 0. _d 0
246 fVerU (i,j,2) = 0. _d 0
247 fVerV (i,j,1) = 0. _d 0
248 fVerV (i,j,2) = 0. _d 0
249 ENDDO
250 ENDDO
251
252 DO k=1,Nr
253 DO j=1-OLy,sNy+OLy
254 DO i=1-OLx,sNx+OLx
255 #ifdef INCLUDE_CONVECT_CALL
256 ConvectCount(i,j,k) = 0.
257 #endif
258 KappaRT(i,j,k) = 0. _d 0
259 KappaRS(i,j,k) = 0. _d 0
260 ENDDO
261 ENDDO
262 ENDDO
263
264 iMin = 1-OLx+1
265 iMax = sNx+OLx
266 jMin = 1-OLy+1
267 jMax = sNy+OLy
268
269
270 C-- Start of diagnostic loop
271 DO k=Nr,1,-1
272
273 #ifdef ALLOW_AUTODIFF_TAMC
274 C? Patrick, is this formula correct now that we change the loop range?
275 C? Do we still need this?
276 kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
277 #endif /* ALLOW_AUTODIFF_TAMC */
278
279 C-- Integrate continuity vertically for vertical velocity
280 CALL INTEGRATE_FOR_W(
281 I bi, bj, k, uVel, vVel,
282 O wVel,
283 I myThid )
284
285 #ifdef ALLOW_OBCS
286 C-- Calculate future values on open boundaries
287 IF (openBoundaries) THEN
288 #ifdef ALLOW_NONHYDROSTATIC
289 IF (nonHydrostatic) THEN
290 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
291 ENDIF
292 #endif /* ALLOW_NONHYDROSTATIC */
293 CALL OBCS_CALC( bi, bj, k, myTime+deltaT, myThid )
294 ENDIF
295 #endif /* ALLOW_OBCS */
296
297 C-- Calculate gradients of potential density for isoneutral
298 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
299 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
300 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
301 CALL FIND_RHO(
302 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
303 I theta, salt,
304 O rhoK,
305 I myThid )
306 IF (k.GT.1) CALL FIND_RHO(
307 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
308 I theta, salt,
309 O rhoKm1,
310 I myThid )
311 CALL GRAD_SIGMA(
312 I bi, bj, iMin, iMax, jMin, jMax, k,
313 I rhoK, rhoKm1, rhoK,
314 O sigmaX, sigmaY, sigmaR,
315 I myThid )
316 ENDIF
317
318 C-- Implicit Vertical Diffusion for Convection
319 c ==> should use sigmaR !!!
320 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
321 CALL CALC_IVDC(
322 I bi, bj, iMin, iMax, jMin, jMax, k,
323 I rhoKm1, rhoK,
324 U ConvectCount, KappaRT, KappaRS,
325 I myTime, myIter, myThid)
326 END IF
327
328 C-- end of diagnostic k loop (Nr:1)
329 ENDDO
330
331 C-- Determines forcing terms based on external fields
332 C relaxation terms, etc.
333 CALL EXTERNAL_FORCING_SURF(
334 I bi, bj, iMin, iMax, jMin, jMax,
335 I myThid )
336
337 #ifdef ALLOW_GMREDI
338 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
339 IF (useGMRedi) THEN
340 DO k=1,Nr
341 CALL GMREDI_CALC_TENSOR(
342 I bi, bj, iMin, iMax, jMin, jMax, k,
343 I sigmaX, sigmaY, sigmaR,
344 I myThid )
345 ENDDO
346 ENDIF
347 #endif /* ALLOW_GMREDI */
348
349 #ifdef ALLOW_KPP
350 C-- Compute KPP mixing coefficients
351 IF (useKPP) THEN
352 CALL KPP_CALC(
353 I bi, bj, myTime, myThid )
354 ENDIF
355 #endif /* ALLOW_KPP */
356
357 #ifdef ALLOW_AUTODIFF_TAMC
358 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
359 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
360 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
361 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
362 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
363 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
364 #endif /* ALLOW_AUTODIFF_TAMC */
365
366
367
368 C-- Start of thermodynamics loop
369 DO k=Nr,1,-1
370
371 C-- km1 Points to level above k (=k-1)
372 C-- kup Cycles through 1,2 to point to layer above
373 C-- kDown Cycles through 2,1 to point to current layer
374
375 km1 = MAX(1,k-1)
376 kup = 1+MOD(k+1,2)
377 kDown= 1+MOD(k,2)
378
379 iMin = 1-OLx+2
380 iMax = sNx+OLx-1
381 jMin = 1-OLy+2
382 jMax = sNy+OLy-1
383
384 #ifdef ALLOW_AUTODIFF_TAMC
385 CPatrick Is this formula correct?
386 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
387 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
388 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
389 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
390 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
391 #endif /* ALLOW_AUTODIFF_TAMC */
392
393 C-- Get temporary terms used by tendency routines
394 CALL CALC_COMMON_FACTORS (
395 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
396 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
397 I myThid)
398
399 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
400 C-- Calculate the total vertical diffusivity
401 CALL CALC_DIFFUSIVITY(
402 I bi,bj,iMin,iMax,jMin,jMax,k,
403 I maskC,maskup,
404 O KappaRT,KappaRS,KappaRU,KappaRV,
405 I myThid)
406 #endif
407
408 C-- Calculate active tracer tendencies (gT,gS,...)
409 C and step forward storing result in gTnm1, gSnm1, etc.
410 IF ( tempStepping ) THEN
411 CALL CALC_GT(
412 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
413 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
414 I KappaRT,
415 U fVerT,
416 I myTime, myThid)
417 CALL TIMESTEP_TRACER(
418 I bi,bj,iMin,iMax,jMin,jMax,k,
419 I theta, gT,
420 U gTnm1,
421 I myIter, myThid)
422 ENDIF
423 IF ( saltStepping ) THEN
424 CALL CALC_GS(
425 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
426 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
427 I KappaRS,
428 U fVerS,
429 I myTime, myThid)
430 CALL TIMESTEP_TRACER(
431 I bi,bj,iMin,iMax,jMin,jMax,k,
432 I salt, gS,
433 U gSnm1,
434 I myIter, myThid)
435 ENDIF
436
437 #ifdef ALLOW_OBCS
438 C-- Apply open boundary conditions
439 IF (openBoundaries) THEN
440 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
441 END IF
442 #endif /* ALLOW_OBCS */
443
444 C-- Freeze water
445 IF (allowFreezing) THEN
446 #ifdef ALLOW_AUTODIFF_TAMC
447 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
448 CADJ & , key = kkey, byte = isbyte
449 #endif /* ALLOW_AUTODIFF_TAMC */
450 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
451 END IF
452
453 C-- end of thermodynamic k loop (Nr:1)
454 ENDDO
455
456
457 #ifdef ALLOW_AUTODIFF_TAMC
458 CPatrick? What about this one?
459 maximpl = 6
460 iikey = (ikey-1)*maximpl
461 #endif /* ALLOW_AUTODIFF_TAMC */
462
463 C-- Implicit diffusion
464 IF (implicitDiffusion) THEN
465
466 IF (tempStepping) THEN
467 #ifdef ALLOW_AUTODIFF_TAMC
468 idkey = iikey + 1
469 #endif /* ALLOW_AUTODIFF_TAMC */
470 CALL IMPLDIFF(
471 I bi, bj, iMin, iMax, jMin, jMax,
472 I deltaTtracer, KappaRT, recip_HFacC,
473 U gTNm1,
474 I myThid )
475 ENDIF
476
477 IF (saltStepping) THEN
478 #ifdef ALLOW_AUTODIFF_TAMC
479 idkey = iikey + 2
480 #endif /* ALLOW_AUTODIFF_TAMC */
481 CALL IMPLDIFF(
482 I bi, bj, iMin, iMax, jMin, jMax,
483 I deltaTtracer, KappaRS, recip_HFacC,
484 U gSNm1,
485 I myThid )
486 ENDIF
487
488 #ifdef ALLOW_OBCS
489 C-- Apply open boundary conditions
490 IF (openBoundaries) THEN
491 DO K=1,Nr
492 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
493 ENDDO
494 END IF
495 #endif /* ALLOW_OBCS */
496
497 C-- End If implicitDiffusion
498 ENDIF
499
500
501
502 C-- Start of dynamics loop
503 DO k=1,Nr
504
505 C-- km1 Points to level above k (=k-1)
506 C-- kup Cycles through 1,2 to point to layer above
507 C-- kDown Cycles through 2,1 to point to current layer
508
509 km1 = MAX(1,k-1)
510 kup = 1+MOD(k+1,2)
511 kDown= 1+MOD(k,2)
512
513 iMin = 1-OLx+2
514 iMax = sNx+OLx-1
515 jMin = 1-OLy+2
516 jMax = sNy+OLy-1
517
518 C-- Integrate hydrostatic balance for phiHyd with BC of
519 C phiHyd(z=0)=0
520 C distinguishe between Stagger and Non Stagger time stepping
521 IF (staggerTimeStep) THEN
522 CALL CALC_PHI_HYD(
523 I bi,bj,iMin,iMax,jMin,jMax,k,
524 I gTnm1, gSnm1,
525 U phiHyd,
526 I myThid )
527 ELSE
528 CALL CALC_PHI_HYD(
529 I bi,bj,iMin,iMax,jMin,jMax,k,
530 I theta, salt,
531 U phiHyd,
532 I myThid )
533 ENDIF
534
535 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
536 C and step forward storing the result in gUnm1, gVnm1, etc...
537 IF ( momStepping ) THEN
538 CALL CALC_MOM_RHS(
539 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
540 I phiHyd,KappaRU,KappaRV,
541 U fVerU, fVerV,
542 I myTime, myThid)
543 CALL TIMESTEP(
544 I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
545 I myIter, myThid)
546
547 #ifdef ALLOW_OBCS
548 C-- Apply open boundary conditions
549 IF (openBoundaries) THEN
550 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
551 END IF
552 #endif /* ALLOW_OBCS */
553
554 #ifdef ALLOW_AUTODIFF_TAMC
555 #ifdef INCLUDE_CD_CODE
556 ELSE
557 DO j=1-OLy,sNy+OLy
558 DO i=1-OLx,sNx+OLx
559 guCD(i,j,k,bi,bj) = 0.0
560 gvCD(i,j,k,bi,bj) = 0.0
561 END DO
562 END DO
563 #endif /* INCLUDE_CD_CODE */
564 #endif /* ALLOW_AUTODIFF_TAMC */
565 ENDIF
566
567
568 C-- end of dynamics k loop (1:Nr)
569 ENDDO
570
571
572
573 C-- Implicit viscosity
574 IF (implicitViscosity.AND.momStepping) THEN
575 #ifdef ALLOW_AUTODIFF_TAMC
576 idkey = iikey + 3
577 #endif /* ALLOW_AUTODIFF_TAMC */
578 CALL IMPLDIFF(
579 I bi, bj, iMin, iMax, jMin, jMax,
580 I deltaTmom, KappaRU,recip_HFacW,
581 U gUNm1,
582 I myThid )
583 #ifdef ALLOW_AUTODIFF_TAMC
584 idkey = iikey + 4
585 #endif /* ALLOW_AUTODIFF_TAMC */
586 CALL IMPLDIFF(
587 I bi, bj, iMin, iMax, jMin, jMax,
588 I deltaTmom, KappaRV,recip_HFacS,
589 U gVNm1,
590 I myThid )
591
592 #ifdef ALLOW_OBCS
593 C-- Apply open boundary conditions
594 IF (openBoundaries) THEN
595 DO K=1,Nr
596 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
597 ENDDO
598 END IF
599 #endif /* ALLOW_OBCS */
600
601 #ifdef INCLUDE_CD_CODE
602 #ifdef ALLOW_AUTODIFF_TAMC
603 idkey = iikey + 5
604 #endif /* ALLOW_AUTODIFF_TAMC */
605 CALL IMPLDIFF(
606 I bi, bj, iMin, iMax, jMin, jMax,
607 I deltaTmom, KappaRU,recip_HFacW,
608 U vVelD,
609 I myThid )
610 #ifdef ALLOW_AUTODIFF_TAMC
611 idkey = iikey + 6
612 #endif /* ALLOW_AUTODIFF_TAMC */
613 CALL IMPLDIFF(
614 I bi, bj, iMin, iMax, jMin, jMax,
615 I deltaTmom, KappaRV,recip_HFacS,
616 U uVelD,
617 I myThid )
618 #endif /* INCLUDE_CD_CODE */
619 C-- End If implicitViscosity.AND.momStepping
620 ENDIF
621
622 ENDDO
623 ENDDO
624
625 RETURN
626 END
627
628
629 C-- Cumulative diagnostic calculations (ie. time-averaging)
630 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
631 c IF (taveFreq.GT.0.) THEN
632 c CALL DO_TIME_AVERAGES(
633 c I myTime, myIter, bi, bj, k, kup, kDown,
634 c I ConvectCount,
635 c I myThid )
636 c ENDIF
637 #endif
638

  ViewVC Help
Powered by ViewVC 1.1.22