/[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.58 - (show annotations) (download)
Fri Feb 2 21:04:48 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
Changes since 1.57: +248 -551 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/model/src/dynamics.F,v 1.57 2001/02/01 19:32:02 heimbach 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 phiHyd(i,j,k) = 0. _d 0
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 #ifdef ALLOW_NONHYDROSTATIC
288 C-- Calculate future values on open boundaries
289 IF (useOBCS.AND.nonHydrostatic) THEN
290 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
291 ENDIF
292 #endif /* ALLOW_NONHYDROSTATIC */
293 #endif /* ALLOW_OBCS */
294
295 C-- Calculate gradients of potential density for isoneutral
296 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
297 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
298 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
299 CALL FIND_RHO(
300 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
301 I theta, salt,
302 O rhoK,
303 I myThid )
304 IF (k.GT.1) CALL FIND_RHO(
305 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
306 I theta, salt,
307 O rhoKm1,
308 I myThid )
309 CALL GRAD_SIGMA(
310 I bi, bj, iMin, iMax, jMin, jMax, k,
311 I rhoK, rhoKm1, rhoK,
312 O sigmaX, sigmaY, sigmaR,
313 I myThid )
314 ENDIF
315
316 C-- Implicit Vertical Diffusion for Convection
317 c ==> should use sigmaR !!!
318 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
319 CALL CALC_IVDC(
320 I bi, bj, iMin, iMax, jMin, jMax, k,
321 I rhoKm1, rhoK,
322 U ConvectCount, KappaRT, KappaRS,
323 I myTime, myIter, myThid)
324 END IF
325
326 C-- end of diagnostic k loop (Nr:1)
327 ENDDO
328
329 #ifdef ALLOW_OBCS
330 C-- Calculate future values on open boundaries
331 IF (useOBCS) THEN
332 CALL OBCS_CALC( bi, bj, myTime+deltaT,
333 I uVel, vVel, wVel, theta, salt,
334 I myThid )
335 ENDIF
336 #endif /* ALLOW_OBCS */
337
338 C-- Determines forcing terms based on external fields
339 C relaxation terms, etc.
340 CALL EXTERNAL_FORCING_SURF(
341 I bi, bj, iMin, iMax, jMin, jMax,
342 I myThid )
343
344 #ifdef ALLOW_GMREDI
345 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
346 IF (useGMRedi) THEN
347 DO k=1,Nr
348 CALL GMREDI_CALC_TENSOR(
349 I bi, bj, iMin, iMax, jMin, jMax, k,
350 I sigmaX, sigmaY, sigmaR,
351 I myThid )
352 ENDDO
353 #ifdef ALLOW_AUTODIFF_TAMC
354 ELSE
355 DO k=1, Nr
356 CALL GMREDI_CALC_TENSOR_DUMMY(
357 I bi, bj, iMin, iMax, jMin, jMax, k,
358 I sigmaX, sigmaY, sigmaR,
359 I myThid )
360 ENDDO
361 #endif /* ALLOW_AUTODIFF_TAMC */
362 ENDIF
363 #endif /* ALLOW_GMREDI */
364
365 #ifdef ALLOW_KPP
366 C-- Compute KPP mixing coefficients
367 IF (useKPP) THEN
368 CALL KPP_CALC(
369 I bi, bj, myTime, myThid )
370 ENDIF
371 #endif /* ALLOW_KPP */
372
373 #ifdef ALLOW_AUTODIFF_TAMC
374 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
375 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
376 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
377 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
378 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
379 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
380 #endif /* ALLOW_AUTODIFF_TAMC */
381
382 #ifdef ALLOW_AIM
383 C AIM - atmospheric intermediate model, physics package code.
384 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
385 IF ( useAIM ) THEN
386 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
387 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
388 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
389 ENDIF
390 #endif /* ALLOW_AIM */
391
392
393 C-- Start of thermodynamics loop
394 DO k=Nr,1,-1
395
396 C-- km1 Points to level above k (=k-1)
397 C-- kup Cycles through 1,2 to point to layer above
398 C-- kDown Cycles through 2,1 to point to current layer
399
400 km1 = MAX(1,k-1)
401 kup = 1+MOD(k+1,2)
402 kDown= 1+MOD(k,2)
403
404 iMin = 1-OLx+2
405 iMax = sNx+OLx-1
406 jMin = 1-OLy+2
407 jMax = sNy+OLy-1
408
409 #ifdef ALLOW_AUTODIFF_TAMC
410 CPatrick Is this formula correct?
411 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
412 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
413 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
414 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
415 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
416 #endif /* ALLOW_AUTODIFF_TAMC */
417
418 C-- Get temporary terms used by tendency routines
419 CALL CALC_COMMON_FACTORS (
420 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
421 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
422 I myThid)
423
424 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
425 C-- Calculate the total vertical diffusivity
426 CALL CALC_DIFFUSIVITY(
427 I bi,bj,iMin,iMax,jMin,jMax,k,
428 I maskC,maskup,
429 O KappaRT,KappaRS,KappaRU,KappaRV,
430 I myThid)
431 #endif
432
433 C-- Calculate active tracer tendencies (gT,gS,...)
434 C and step forward storing result in gTnm1, gSnm1, etc.
435 IF ( tempStepping ) THEN
436 CALL CALC_GT(
437 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
438 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
439 I KappaRT,
440 U fVerT,
441 I myTime, myThid)
442 CALL TIMESTEP_TRACER(
443 I bi,bj,iMin,iMax,jMin,jMax,k,
444 I theta, gT,
445 U gTnm1,
446 I myIter, myThid)
447 ENDIF
448 IF ( saltStepping ) THEN
449 CALL CALC_GS(
450 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
451 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
452 I KappaRS,
453 U fVerS,
454 I myTime, myThid)
455 CALL TIMESTEP_TRACER(
456 I bi,bj,iMin,iMax,jMin,jMax,k,
457 I salt, gS,
458 U gSnm1,
459 I myIter, myThid)
460 ENDIF
461
462 #ifdef ALLOW_OBCS
463 C-- Apply open boundary conditions
464 IF (useOBCS) THEN
465 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
466 END IF
467 #endif /* ALLOW_OBCS */
468
469 C-- Freeze water
470 IF (allowFreezing) THEN
471 #ifdef ALLOW_AUTODIFF_TAMC
472 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
473 CADJ & , key = kkey, byte = isbyte
474 #endif /* ALLOW_AUTODIFF_TAMC */
475 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
476 END IF
477
478 C-- end of thermodynamic k loop (Nr:1)
479 ENDDO
480
481
482 #ifdef ALLOW_AUTODIFF_TAMC
483 CPatrick? What about this one?
484 maximpl = 6
485 iikey = (ikey-1)*maximpl
486 #endif /* ALLOW_AUTODIFF_TAMC */
487
488 C-- Implicit diffusion
489 IF (implicitDiffusion) THEN
490
491 IF (tempStepping) THEN
492 #ifdef ALLOW_AUTODIFF_TAMC
493 idkey = iikey + 1
494 #endif /* ALLOW_AUTODIFF_TAMC */
495 CALL IMPLDIFF(
496 I bi, bj, iMin, iMax, jMin, jMax,
497 I deltaTtracer, KappaRT, recip_HFacC,
498 U gTNm1,
499 I myThid )
500 ENDIF
501
502 IF (saltStepping) THEN
503 #ifdef ALLOW_AUTODIFF_TAMC
504 idkey = iikey + 2
505 #endif /* ALLOW_AUTODIFF_TAMC */
506 CALL IMPLDIFF(
507 I bi, bj, iMin, iMax, jMin, jMax,
508 I deltaTtracer, KappaRS, recip_HFacC,
509 U gSNm1,
510 I myThid )
511 ENDIF
512
513 #ifdef ALLOW_OBCS
514 C-- Apply open boundary conditions
515 IF (useOBCS) THEN
516 DO K=1,Nr
517 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
518 ENDDO
519 END IF
520 #endif /* ALLOW_OBCS */
521
522 C-- End If implicitDiffusion
523 ENDIF
524
525
526
527 C-- Start of dynamics loop
528 DO k=1,Nr
529
530 C-- km1 Points to level above k (=k-1)
531 C-- kup Cycles through 1,2 to point to layer above
532 C-- kDown Cycles through 2,1 to point to current layer
533
534 km1 = MAX(1,k-1)
535 kup = 1+MOD(k+1,2)
536 kDown= 1+MOD(k,2)
537
538 iMin = 1-OLx+2
539 iMax = sNx+OLx-1
540 jMin = 1-OLy+2
541 jMax = sNy+OLy-1
542
543 C-- Integrate hydrostatic balance for phiHyd with BC of
544 C phiHyd(z=0)=0
545 C distinguishe between Stagger and Non Stagger time stepping
546 IF (staggerTimeStep) THEN
547 CALL CALC_PHI_HYD(
548 I bi,bj,iMin,iMax,jMin,jMax,k,
549 I gTnm1, gSnm1,
550 U phiHyd,
551 I myThid )
552 ELSE
553 CALL CALC_PHI_HYD(
554 I bi,bj,iMin,iMax,jMin,jMax,k,
555 I theta, salt,
556 U phiHyd,
557 I myThid )
558 ENDIF
559
560 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
561 C and step forward storing the result in gUnm1, gVnm1, etc...
562 IF ( momStepping ) THEN
563 CALL CALC_MOM_RHS(
564 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
565 I phiHyd,KappaRU,KappaRV,
566 U fVerU, fVerV,
567 I myTime, myThid)
568 CALL TIMESTEP(
569 I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
570 I myIter, myThid)
571
572 #ifdef ALLOW_OBCS
573 C-- Apply open boundary conditions
574 IF (useOBCS) THEN
575 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
576 END IF
577 #endif /* ALLOW_OBCS */
578
579 #ifdef ALLOW_AUTODIFF_TAMC
580 #ifdef INCLUDE_CD_CODE
581 ELSE
582 DO j=1-OLy,sNy+OLy
583 DO i=1-OLx,sNx+OLx
584 guCD(i,j,k,bi,bj) = 0.0
585 gvCD(i,j,k,bi,bj) = 0.0
586 END DO
587 END DO
588 #endif /* INCLUDE_CD_CODE */
589 #endif /* ALLOW_AUTODIFF_TAMC */
590 ENDIF
591
592
593 C-- end of dynamics k loop (1:Nr)
594 ENDDO
595
596
597
598 C-- Implicit viscosity
599 IF (implicitViscosity.AND.momStepping) THEN
600 #ifdef ALLOW_AUTODIFF_TAMC
601 idkey = iikey + 3
602 #endif /* ALLOW_AUTODIFF_TAMC */
603 CALL IMPLDIFF(
604 I bi, bj, iMin, iMax, jMin, jMax,
605 I deltaTmom, KappaRU,recip_HFacW,
606 U gUNm1,
607 I myThid )
608 #ifdef ALLOW_AUTODIFF_TAMC
609 idkey = iikey + 4
610 #endif /* ALLOW_AUTODIFF_TAMC */
611 CALL IMPLDIFF(
612 I bi, bj, iMin, iMax, jMin, jMax,
613 I deltaTmom, KappaRV,recip_HFacS,
614 U gVNm1,
615 I myThid )
616
617 #ifdef ALLOW_OBCS
618 C-- Apply open boundary conditions
619 IF (useOBCS) THEN
620 DO K=1,Nr
621 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
622 ENDDO
623 END IF
624 #endif /* ALLOW_OBCS */
625
626 #ifdef INCLUDE_CD_CODE
627 #ifdef ALLOW_AUTODIFF_TAMC
628 idkey = iikey + 5
629 #endif /* ALLOW_AUTODIFF_TAMC */
630 CALL IMPLDIFF(
631 I bi, bj, iMin, iMax, jMin, jMax,
632 I deltaTmom, KappaRU,recip_HFacW,
633 U vVelD,
634 I myThid )
635 #ifdef ALLOW_AUTODIFF_TAMC
636 idkey = iikey + 6
637 #endif /* ALLOW_AUTODIFF_TAMC */
638 CALL IMPLDIFF(
639 I bi, bj, iMin, iMax, jMin, jMax,
640 I deltaTmom, KappaRV,recip_HFacS,
641 U uVelD,
642 I myThid )
643 #endif /* INCLUDE_CD_CODE */
644 C-- End If implicitViscosity.AND.momStepping
645 ENDIF
646
647 ENDDO
648 ENDDO
649
650 RETURN
651 END

  ViewVC Help
Powered by ViewVC 1.1.22