/[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.5 - (show annotations) (download)
Tue Jan 9 16:21:18 2001 UTC (23 years, 4 months ago) by adcroft
Branch: branch-atmos-merge
Changes since 1.54.2.4: +4 -27 lines
Moved af,df,fZon and fMer temporary arrays into calc_gs() and calc_gt()
 - this uses slightly more memory
 - this reduces the number of arguments and set-up space in dynamics()

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

  ViewVC Help
Powered by ViewVC 1.1.22