/[MITgcm]/MITgcm/verification/aim.5l_Equatorial_Channel/code/dynamics.F
ViewVC logotype

Contents of /MITgcm/verification/aim.5l_Equatorial_Channel/code/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (show annotations) (download)
Fri Feb 2 21:36:30 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
Changes since 1.1: +647 -0 lines
Merged changes from branch "branch-atmos-merge" into MAIN (checkpoint34)
 - substantial modifications to algorithm sequence (dynamics.F)
 - packaged OBCS, Shapiro filter, Zonal filter, Atmospheric Physics

1 C $Header: /u/gcmpack/models/MITgcmUV/verification/aim.5l_Equatorial_Channel/code/Attic/dynamics.F,v 1.1.2.1 2001/01/31 03:16:11 cnh 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 #ifdef ALLOW_AIM
367 C AIM - atmospheric intermediate model, physics package code.
368 IF ( useAIM ) THEN
369 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
370 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
371 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
372 ENDIF
373 #endif /* ALLOW_AIM */
374
375
376
377 C-- Start of thermodynamics loop
378 DO k=Nr,1,-1
379
380 C-- km1 Points to level above k (=k-1)
381 C-- kup Cycles through 1,2 to point to layer above
382 C-- kDown Cycles through 2,1 to point to current layer
383
384 km1 = MAX(1,k-1)
385 kup = 1+MOD(k+1,2)
386 kDown= 1+MOD(k,2)
387
388 iMin = 1-OLx+2
389 iMax = sNx+OLx-1
390 jMin = 1-OLy+2
391 jMax = sNy+OLy-1
392
393 #ifdef ALLOW_AUTODIFF_TAMC
394 CPatrick Is this formula correct?
395 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
396 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
397 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
398 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
399 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
400 #endif /* ALLOW_AUTODIFF_TAMC */
401
402 C-- Get temporary terms used by tendency routines
403 CALL CALC_COMMON_FACTORS (
404 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
405 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
406 I myThid)
407
408 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
409 C-- Calculate the total vertical diffusivity
410 CALL CALC_DIFFUSIVITY(
411 I bi,bj,iMin,iMax,jMin,jMax,k,
412 I maskC,maskup,
413 O KappaRT,KappaRS,KappaRU,KappaRV,
414 I myThid)
415 #endif
416
417 C-- Calculate active tracer tendencies (gT,gS,...)
418 C and step forward storing result in gTnm1, gSnm1, etc.
419 IF ( tempStepping ) THEN
420 CALL CALC_GT(
421 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
422 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
423 I KappaRT,
424 U fVerT,
425 I myTime, myThid)
426 CALL TIMESTEP_TRACER(
427 I bi,bj,iMin,iMax,jMin,jMax,k,
428 I theta, gT,
429 U gTnm1,
430 I myIter, myThid)
431 ENDIF
432 IF ( saltStepping ) THEN
433 CALL CALC_GS(
434 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
435 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
436 I KappaRS,
437 U fVerS,
438 I myTime, myThid)
439 CALL TIMESTEP_TRACER(
440 I bi,bj,iMin,iMax,jMin,jMax,k,
441 I salt, gS,
442 U gSnm1,
443 I myIter, myThid)
444 ENDIF
445
446 #ifdef ALLOW_OBCS
447 C-- Apply open boundary conditions
448 IF (openBoundaries) THEN
449 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
450 END IF
451 #endif /* ALLOW_OBCS */
452
453 C-- Freeze water
454 IF (allowFreezing) THEN
455 #ifdef ALLOW_AUTODIFF_TAMC
456 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
457 CADJ & , key = kkey, byte = isbyte
458 #endif /* ALLOW_AUTODIFF_TAMC */
459 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
460 END IF
461
462 C-- end of thermodynamic k loop (Nr:1)
463 ENDDO
464
465
466 #ifdef ALLOW_AUTODIFF_TAMC
467 CPatrick? What about this one?
468 maximpl = 6
469 iikey = (ikey-1)*maximpl
470 #endif /* ALLOW_AUTODIFF_TAMC */
471
472 C-- Implicit diffusion
473 IF (implicitDiffusion) THEN
474
475 IF (tempStepping) THEN
476 #ifdef ALLOW_AUTODIFF_TAMC
477 idkey = iikey + 1
478 #endif /* ALLOW_AUTODIFF_TAMC */
479 CALL IMPLDIFF(
480 I bi, bj, iMin, iMax, jMin, jMax,
481 I deltaTtracer, KappaRT, recip_HFacC,
482 U gTNm1,
483 I myThid )
484 ENDIF
485
486 IF (saltStepping) THEN
487 #ifdef ALLOW_AUTODIFF_TAMC
488 idkey = iikey + 2
489 #endif /* ALLOW_AUTODIFF_TAMC */
490 CALL IMPLDIFF(
491 I bi, bj, iMin, iMax, jMin, jMax,
492 I deltaTtracer, KappaRS, recip_HFacC,
493 U gSNm1,
494 I myThid )
495 ENDIF
496
497 #ifdef ALLOW_OBCS
498 C-- Apply open boundary conditions
499 IF (openBoundaries) THEN
500 DO K=1,Nr
501 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
502 ENDDO
503 END IF
504 #endif /* ALLOW_OBCS */
505
506 C-- End If implicitDiffusion
507 ENDIF
508
509
510
511 C-- Start of dynamics loop
512 DO k=1,Nr
513
514 C-- km1 Points to level above k (=k-1)
515 C-- kup Cycles through 1,2 to point to layer above
516 C-- kDown Cycles through 2,1 to point to current layer
517
518 km1 = MAX(1,k-1)
519 kup = 1+MOD(k+1,2)
520 kDown= 1+MOD(k,2)
521
522 iMin = 1-OLx+2
523 iMax = sNx+OLx-1
524 jMin = 1-OLy+2
525 jMax = sNy+OLy-1
526
527 C-- Integrate hydrostatic balance for phiHyd with BC of
528 C phiHyd(z=0)=0
529 C distinguishe between Stagger and Non Stagger time stepping
530 IF (staggerTimeStep) THEN
531 CALL CALC_PHI_HYD(
532 I bi,bj,iMin,iMax,jMin,jMax,k,
533 I gTnm1, gSnm1,
534 U phiHyd,
535 I myThid )
536 ELSE
537 CALL CALC_PHI_HYD(
538 I bi,bj,iMin,iMax,jMin,jMax,k,
539 I theta, salt,
540 U phiHyd,
541 I myThid )
542 ENDIF
543
544 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
545 C and step forward storing the result in gUnm1, gVnm1, etc...
546 IF ( momStepping ) THEN
547 CALL CALC_MOM_RHS(
548 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
549 I phiHyd,KappaRU,KappaRV,
550 U fVerU, fVerV,
551 I myTime, myThid)
552 CALL TIMESTEP(
553 I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
554 I myIter, myThid)
555
556 #ifdef ALLOW_OBCS
557 C-- Apply open boundary conditions
558 IF (openBoundaries) THEN
559 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
560 END IF
561 #endif /* ALLOW_OBCS */
562
563 #ifdef ALLOW_AUTODIFF_TAMC
564 #ifdef INCLUDE_CD_CODE
565 ELSE
566 DO j=1-OLy,sNy+OLy
567 DO i=1-OLx,sNx+OLx
568 guCD(i,j,k,bi,bj) = 0.0
569 gvCD(i,j,k,bi,bj) = 0.0
570 END DO
571 END DO
572 #endif /* INCLUDE_CD_CODE */
573 #endif /* ALLOW_AUTODIFF_TAMC */
574 ENDIF
575
576
577 C-- end of dynamics k loop (1:Nr)
578 ENDDO
579
580
581
582 C-- Implicit viscosity
583 IF (implicitViscosity.AND.momStepping) THEN
584 #ifdef ALLOW_AUTODIFF_TAMC
585 idkey = iikey + 3
586 #endif /* ALLOW_AUTODIFF_TAMC */
587 CALL IMPLDIFF(
588 I bi, bj, iMin, iMax, jMin, jMax,
589 I deltaTmom, KappaRU,recip_HFacW,
590 U gUNm1,
591 I myThid )
592 #ifdef ALLOW_AUTODIFF_TAMC
593 idkey = iikey + 4
594 #endif /* ALLOW_AUTODIFF_TAMC */
595 CALL IMPLDIFF(
596 I bi, bj, iMin, iMax, jMin, jMax,
597 I deltaTmom, KappaRV,recip_HFacS,
598 U gVNm1,
599 I myThid )
600
601 #ifdef ALLOW_OBCS
602 C-- Apply open boundary conditions
603 IF (openBoundaries) THEN
604 DO K=1,Nr
605 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
606 ENDDO
607 END IF
608 #endif /* ALLOW_OBCS */
609
610 #ifdef INCLUDE_CD_CODE
611 #ifdef ALLOW_AUTODIFF_TAMC
612 idkey = iikey + 5
613 #endif /* ALLOW_AUTODIFF_TAMC */
614 CALL IMPLDIFF(
615 I bi, bj, iMin, iMax, jMin, jMax,
616 I deltaTmom, KappaRU,recip_HFacW,
617 U vVelD,
618 I myThid )
619 #ifdef ALLOW_AUTODIFF_TAMC
620 idkey = iikey + 6
621 #endif /* ALLOW_AUTODIFF_TAMC */
622 CALL IMPLDIFF(
623 I bi, bj, iMin, iMax, jMin, jMax,
624 I deltaTmom, KappaRV,recip_HFacS,
625 U uVelD,
626 I myThid )
627 #endif /* INCLUDE_CD_CODE */
628 C-- End If implicitViscosity.AND.momStepping
629 ENDIF
630
631 ENDDO
632 ENDDO
633
634 RETURN
635 END
636
637
638 C-- Cumulative diagnostic calculations (ie. time-averaging)
639 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
640 c IF (taveFreq.GT.0.) THEN
641 c CALL DO_TIME_AVERAGES(
642 c I myTime, myIter, bi, bj, k, kup, kDown,
643 c I ConvectCount,
644 c I myThid )
645 c ENDIF
646 #endif
647

  ViewVC Help
Powered by ViewVC 1.1.22