/[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.10 - (show annotations) (download)
Wed Jan 17 15:05:59 2001 UTC (23 years, 4 months ago) by jmc
Branch: branch-atmos-merge
CVS Tags: branch-atmos-merge-shapiro, branch-atmos-merge-zonalfilt, branch-atmos-merge-phase5, branch-atmos-merge-phase7, branch-atmos-merge-phase6
Changes since 1.54.2.9: +13 -4 lines
enable a stagger time stepping of T,S and then U,V ;
   add phiHyd as argument of subroutine TIMESTEP

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

  ViewVC Help
Powered by ViewVC 1.1.22