/[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.3 - (show annotations) (download)
Sun Feb 4 14:38:51 2001 UTC (23 years, 3 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint36, checkpoint35
Changes since 1.2: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

1 C $Header: /u/gcmpack/models/MITgcmUV/verification/aim.5l_Equatorial_Channel/code/dynamics.F,v 1.2 2001/02/02 21:36:30 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
7 C /==========================================================\
8 C | SUBROUTINE DYNAMICS |
9 C | o Controlling routine for the explicit part of the model |
10 C | dynamics. |
11 C |==========================================================|
12 C | This routine evaluates the "dynamics" terms for each |
13 C | block of ocean in turn. Because the blocks of ocean have |
14 C | overlap regions they are independent of one another. |
15 C | If terms involving lateral integrals are needed in this |
16 C | routine care will be needed. Similarly finite-difference |
17 C | operations with stencils wider than the overlap region |
18 C | require special consideration. |
19 C | Notes |
20 C | ===== |
21 C | C*P* comments indicating place holders for which code is |
22 C | presently being developed. |
23 C \==========================================================/
24 IMPLICIT NONE
25
26 C == Global variables ===
27 #include "SIZE.h"
28 #include "EEPARAMS.h"
29 #include "CG2D.h"
30 #include "PARAMS.h"
31 #include "DYNVARS.h"
32 #include "GRID.h"
33
34 #ifdef ALLOW_AUTODIFF_TAMC
35 # include "tamc.h"
36 # include "tamc_keys.h"
37 #endif /* ALLOW_AUTODIFF_TAMC */
38
39 #ifdef ALLOW_KPP
40 # include "KPP.h"
41 #endif
42
43 C == Routine arguments ==
44 C myTime - Current time in simulation
45 C myIter - Current iteration number in simulation
46 C myThid - Thread number for this instance of the routine.
47 _RL myTime
48 INTEGER myIter
49 INTEGER myThid
50
51 C == Local variables
52 C xA, yA - Per block temporaries holding face areas
53 C uTrans, vTrans, rTrans - Per block temporaries holding flow
54 C transport
55 C rVel o uTrans: Zonal transport
56 C o vTrans: Meridional transport
57 C o rTrans: Vertical transport
58 C o rVel: Vertical velocity at upper and
59 C lower cell faces.
60 C maskC,maskUp o maskC: land/water mask for tracer cells
61 C o maskUp: land/water mask for W points
62 C fVer[STUV] o fVer: Vertical flux term - note fVer
63 C is "pipelined" in the vertical
64 C so we need an fVer for each
65 C variable.
66 C rhoK, rhoKM1 - Density at current level, and level above
67 C phiHyd - Hydrostatic part of the potential phiHydi.
68 C In z coords phiHydiHyd is the hydrostatic
69 C pressure anomaly
70 C In p coords phiHydiHyd is the geopotential
71 C surface height
72 C anomaly.
73 C etaSurfX, - Holds surface elevation gradient in X and Y.
74 C etaSurfY
75 C KappaRT, - Total diffusion in vertical for T and S.
76 C KappaRS (background + spatially varying, isopycnal term).
77 C iMin, iMax - Ranges and sub-block indices on which calculations
78 C jMin, jMax are applied.
79 C bi, bj
80 C k, kup, - Index for layer above and below. kup and kDown
81 C kDown, km1 are switched with layer to be the appropriate
82 C index into fVerTerm.
83 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
89 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
92 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
93 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
94 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
95 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
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 #ifdef ALLOW_AIM
368 C AIM - atmospheric intermediate model, physics package code.
369 IF ( useAIM ) THEN
370 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
371 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
372 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
373 ENDIF
374 #endif /* ALLOW_AIM */
375
376
377
378 C-- Start of thermodynamics loop
379 DO k=Nr,1,-1
380
381 C-- km1 Points to level above k (=k-1)
382 C-- kup Cycles through 1,2 to point to layer above
383 C-- kDown Cycles through 2,1 to point to current layer
384
385 km1 = MAX(1,k-1)
386 kup = 1+MOD(k+1,2)
387 kDown= 1+MOD(k,2)
388
389 iMin = 1-OLx+2
390 iMax = sNx+OLx-1
391 jMin = 1-OLy+2
392 jMax = sNy+OLy-1
393
394 #ifdef ALLOW_AUTODIFF_TAMC
395 CPatrick Is this formula correct?
396 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
397 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
398 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
399 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
400 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
401 #endif /* ALLOW_AUTODIFF_TAMC */
402
403 C-- Get temporary terms used by tendency routines
404 CALL CALC_COMMON_FACTORS (
405 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
406 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
407 I myThid)
408
409 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
410 C-- Calculate the total vertical diffusivity
411 CALL CALC_DIFFUSIVITY(
412 I bi,bj,iMin,iMax,jMin,jMax,k,
413 I maskC,maskup,
414 O KappaRT,KappaRS,KappaRU,KappaRV,
415 I myThid)
416 #endif
417
418 C-- Calculate active tracer tendencies (gT,gS,...)
419 C and step forward storing result in gTnm1, gSnm1, etc.
420 IF ( tempStepping ) THEN
421 CALL CALC_GT(
422 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
423 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
424 I KappaRT,
425 U fVerT,
426 I myTime, myThid)
427 CALL TIMESTEP_TRACER(
428 I bi,bj,iMin,iMax,jMin,jMax,k,
429 I theta, gT,
430 U gTnm1,
431 I myIter, myThid)
432 ENDIF
433 IF ( saltStepping ) THEN
434 CALL CALC_GS(
435 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
436 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
437 I KappaRS,
438 U fVerS,
439 I myTime, myThid)
440 CALL TIMESTEP_TRACER(
441 I bi,bj,iMin,iMax,jMin,jMax,k,
442 I salt, gS,
443 U gSnm1,
444 I myIter, myThid)
445 ENDIF
446
447 #ifdef ALLOW_OBCS
448 C-- Apply open boundary conditions
449 IF (openBoundaries) THEN
450 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
451 END IF
452 #endif /* ALLOW_OBCS */
453
454 C-- Freeze water
455 IF (allowFreezing) THEN
456 #ifdef ALLOW_AUTODIFF_TAMC
457 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
458 CADJ & , key = kkey, byte = isbyte
459 #endif /* ALLOW_AUTODIFF_TAMC */
460 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
461 END IF
462
463 C-- end of thermodynamic k loop (Nr:1)
464 ENDDO
465
466
467 #ifdef ALLOW_AUTODIFF_TAMC
468 CPatrick? What about this one?
469 maximpl = 6
470 iikey = (ikey-1)*maximpl
471 #endif /* ALLOW_AUTODIFF_TAMC */
472
473 C-- Implicit diffusion
474 IF (implicitDiffusion) THEN
475
476 IF (tempStepping) THEN
477 #ifdef ALLOW_AUTODIFF_TAMC
478 idkey = iikey + 1
479 #endif /* ALLOW_AUTODIFF_TAMC */
480 CALL IMPLDIFF(
481 I bi, bj, iMin, iMax, jMin, jMax,
482 I deltaTtracer, KappaRT, recip_HFacC,
483 U gTNm1,
484 I myThid )
485 ENDIF
486
487 IF (saltStepping) THEN
488 #ifdef ALLOW_AUTODIFF_TAMC
489 idkey = iikey + 2
490 #endif /* ALLOW_AUTODIFF_TAMC */
491 CALL IMPLDIFF(
492 I bi, bj, iMin, iMax, jMin, jMax,
493 I deltaTtracer, KappaRS, recip_HFacC,
494 U gSNm1,
495 I myThid )
496 ENDIF
497
498 #ifdef ALLOW_OBCS
499 C-- Apply open boundary conditions
500 IF (openBoundaries) THEN
501 DO K=1,Nr
502 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
503 ENDDO
504 END IF
505 #endif /* ALLOW_OBCS */
506
507 C-- End If implicitDiffusion
508 ENDIF
509
510
511
512 C-- Start of dynamics loop
513 DO k=1,Nr
514
515 C-- km1 Points to level above k (=k-1)
516 C-- kup Cycles through 1,2 to point to layer above
517 C-- kDown Cycles through 2,1 to point to current layer
518
519 km1 = MAX(1,k-1)
520 kup = 1+MOD(k+1,2)
521 kDown= 1+MOD(k,2)
522
523 iMin = 1-OLx+2
524 iMax = sNx+OLx-1
525 jMin = 1-OLy+2
526 jMax = sNy+OLy-1
527
528 C-- Integrate hydrostatic balance for phiHyd with BC of
529 C phiHyd(z=0)=0
530 C distinguishe between Stagger and Non Stagger time stepping
531 IF (staggerTimeStep) THEN
532 CALL CALC_PHI_HYD(
533 I bi,bj,iMin,iMax,jMin,jMax,k,
534 I gTnm1, gSnm1,
535 U phiHyd,
536 I myThid )
537 ELSE
538 CALL CALC_PHI_HYD(
539 I bi,bj,iMin,iMax,jMin,jMax,k,
540 I theta, salt,
541 U phiHyd,
542 I myThid )
543 ENDIF
544
545 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
546 C and step forward storing the result in gUnm1, gVnm1, etc...
547 IF ( momStepping ) THEN
548 CALL CALC_MOM_RHS(
549 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
550 I phiHyd,KappaRU,KappaRV,
551 U fVerU, fVerV,
552 I myTime, myThid)
553 CALL TIMESTEP(
554 I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
555 I myIter, myThid)
556
557 #ifdef ALLOW_OBCS
558 C-- Apply open boundary conditions
559 IF (openBoundaries) THEN
560 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
561 END IF
562 #endif /* ALLOW_OBCS */
563
564 #ifdef ALLOW_AUTODIFF_TAMC
565 #ifdef INCLUDE_CD_CODE
566 ELSE
567 DO j=1-OLy,sNy+OLy
568 DO i=1-OLx,sNx+OLx
569 guCD(i,j,k,bi,bj) = 0.0
570 gvCD(i,j,k,bi,bj) = 0.0
571 END DO
572 END DO
573 #endif /* INCLUDE_CD_CODE */
574 #endif /* ALLOW_AUTODIFF_TAMC */
575 ENDIF
576
577
578 C-- end of dynamics k loop (1:Nr)
579 ENDDO
580
581
582
583 C-- Implicit viscosity
584 IF (implicitViscosity.AND.momStepping) THEN
585 #ifdef ALLOW_AUTODIFF_TAMC
586 idkey = iikey + 3
587 #endif /* ALLOW_AUTODIFF_TAMC */
588 CALL IMPLDIFF(
589 I bi, bj, iMin, iMax, jMin, jMax,
590 I deltaTmom, KappaRU,recip_HFacW,
591 U gUNm1,
592 I myThid )
593 #ifdef ALLOW_AUTODIFF_TAMC
594 idkey = iikey + 4
595 #endif /* ALLOW_AUTODIFF_TAMC */
596 CALL IMPLDIFF(
597 I bi, bj, iMin, iMax, jMin, jMax,
598 I deltaTmom, KappaRV,recip_HFacS,
599 U gVNm1,
600 I myThid )
601
602 #ifdef ALLOW_OBCS
603 C-- Apply open boundary conditions
604 IF (openBoundaries) THEN
605 DO K=1,Nr
606 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
607 ENDDO
608 END IF
609 #endif /* ALLOW_OBCS */
610
611 #ifdef INCLUDE_CD_CODE
612 #ifdef ALLOW_AUTODIFF_TAMC
613 idkey = iikey + 5
614 #endif /* ALLOW_AUTODIFF_TAMC */
615 CALL IMPLDIFF(
616 I bi, bj, iMin, iMax, jMin, jMax,
617 I deltaTmom, KappaRU,recip_HFacW,
618 U vVelD,
619 I myThid )
620 #ifdef ALLOW_AUTODIFF_TAMC
621 idkey = iikey + 6
622 #endif /* ALLOW_AUTODIFF_TAMC */
623 CALL IMPLDIFF(
624 I bi, bj, iMin, iMax, jMin, jMax,
625 I deltaTmom, KappaRV,recip_HFacS,
626 U uVelD,
627 I myThid )
628 #endif /* INCLUDE_CD_CODE */
629 C-- End If implicitViscosity.AND.momStepping
630 ENDIF
631
632 ENDDO
633 ENDDO
634
635 RETURN
636 END
637
638
639 C-- Cumulative diagnostic calculations (ie. time-averaging)
640 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
641 c IF (taveFreq.GT.0.) THEN
642 c CALL DO_TIME_AVERAGES(
643 c I myTime, myIter, bi, bj, k, kup, kDown,
644 c I ConvectCount,
645 c I myThid )
646 c ENDIF
647 #endif
648

  ViewVC Help
Powered by ViewVC 1.1.22