/[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.60 - (show annotations) (download)
Wed Feb 7 16:28:54 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
Changes since 1.59: +3 -3 lines
Corrected comment about call to OBCS_APPLY_W()

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.59 2001/02/04 14:38:47 cnh 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 phiHyd(i,j,k) = 0. _d 0
190 KappaRU(i,j,k) = 0. _d 0
191 KappaRV(i,j,k) = 0. _d 0
192 sigmaX(i,j,k) = 0. _d 0
193 sigmaY(i,j,k) = 0. _d 0
194 sigmaR(i,j,k) = 0. _d 0
195 ENDDO
196 rhoKM1 (i,j) = 0. _d 0
197 rhok (i,j) = 0. _d 0
198 maskC (i,j) = 0. _d 0
199 ENDDO
200 ENDDO
201
202
203 #ifdef ALLOW_AUTODIFF_TAMC
204 C-- HPF directive to help TAMC
205 CHPF$ INDEPENDENT
206 #endif /* ALLOW_AUTODIFF_TAMC */
207
208 DO bj=myByLo(myThid),myByHi(myThid)
209
210 #ifdef ALLOW_AUTODIFF_TAMC
211 C-- HPF directive to help TAMC
212 CHPF$ INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
213 CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA
214 CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
215 CHPF$& )
216 #endif /* ALLOW_AUTODIFF_TAMC */
217
218 DO bi=myBxLo(myThid),myBxHi(myThid)
219
220 #ifdef ALLOW_AUTODIFF_TAMC
221 act1 = bi - myBxLo(myThid)
222 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
223
224 act2 = bj - myByLo(myThid)
225 max2 = myByHi(myThid) - myByLo(myThid) + 1
226
227 act3 = myThid - 1
228 max3 = nTx*nTy
229
230 act4 = ikey_dynamics - 1
231
232 ikey = (act1 + 1) + act2*max1
233 & + act3*max1*max2
234 & + act4*max1*max2*max3
235 #endif /* ALLOW_AUTODIFF_TAMC */
236
237 C-- Set up work arrays that need valid initial values
238 DO j=1-OLy,sNy+OLy
239 DO i=1-OLx,sNx+OLx
240 rTrans(i,j) = 0. _d 0
241 rVel (i,j,1) = 0. _d 0
242 rVel (i,j,2) = 0. _d 0
243 fVerT (i,j,1) = 0. _d 0
244 fVerT (i,j,2) = 0. _d 0
245 fVerS (i,j,1) = 0. _d 0
246 fVerS (i,j,2) = 0. _d 0
247 fVerU (i,j,1) = 0. _d 0
248 fVerU (i,j,2) = 0. _d 0
249 fVerV (i,j,1) = 0. _d 0
250 fVerV (i,j,2) = 0. _d 0
251 ENDDO
252 ENDDO
253
254 DO k=1,Nr
255 DO j=1-OLy,sNy+OLy
256 DO i=1-OLx,sNx+OLx
257 #ifdef INCLUDE_CONVECT_CALL
258 ConvectCount(i,j,k) = 0.
259 #endif
260 KappaRT(i,j,k) = 0. _d 0
261 KappaRS(i,j,k) = 0. _d 0
262 ENDDO
263 ENDDO
264 ENDDO
265
266 iMin = 1-OLx+1
267 iMax = sNx+OLx
268 jMin = 1-OLy+1
269 jMax = sNy+OLy
270
271
272 C-- Start of diagnostic loop
273 DO k=Nr,1,-1
274
275 #ifdef ALLOW_AUTODIFF_TAMC
276 C? Patrick, is this formula correct now that we change the loop range?
277 C? Do we still need this?
278 kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
279 #endif /* ALLOW_AUTODIFF_TAMC */
280
281 C-- Integrate continuity vertically for vertical velocity
282 CALL INTEGRATE_FOR_W(
283 I bi, bj, k, uVel, vVel,
284 O wVel,
285 I myThid )
286
287 #ifdef ALLOW_OBCS
288 #ifdef ALLOW_NONHYDROSTATIC
289 C-- Apply OBC to W if in N-H mode
290 IF (useOBCS.AND.nonHydrostatic) THEN
291 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
292 ENDIF
293 #endif /* ALLOW_NONHYDROSTATIC */
294 #endif /* ALLOW_OBCS */
295
296 C-- Calculate gradients of potential density for isoneutral
297 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
298 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
299 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
300 CALL FIND_RHO(
301 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
302 I theta, salt,
303 O rhoK,
304 I myThid )
305 IF (k.GT.1) CALL FIND_RHO(
306 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
307 I theta, salt,
308 O rhoKm1,
309 I myThid )
310 CALL GRAD_SIGMA(
311 I bi, bj, iMin, iMax, jMin, jMax, k,
312 I rhoK, rhoKm1, rhoK,
313 O sigmaX, sigmaY, sigmaR,
314 I myThid )
315 ENDIF
316
317 C-- Implicit Vertical Diffusion for Convection
318 c ==> should use sigmaR !!!
319 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
320 CALL CALC_IVDC(
321 I bi, bj, iMin, iMax, jMin, jMax, k,
322 I rhoKm1, rhoK,
323 U ConvectCount, KappaRT, KappaRS,
324 I myTime, myIter, myThid)
325 END IF
326
327 C-- end of diagnostic k loop (Nr:1)
328 ENDDO
329
330 #ifdef ALLOW_OBCS
331 C-- Calculate future values on open boundaries
332 IF (useOBCS) THEN
333 CALL OBCS_CALC( bi, bj, myTime+deltaT,
334 I uVel, vVel, wVel, theta, salt,
335 I myThid )
336 ENDIF
337 #endif /* ALLOW_OBCS */
338
339 C-- Determines forcing terms based on external fields
340 C relaxation terms, etc.
341 CALL EXTERNAL_FORCING_SURF(
342 I bi, bj, iMin, iMax, jMin, jMax,
343 I myThid )
344
345 #ifdef ALLOW_GMREDI
346 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
347 IF (useGMRedi) THEN
348 DO k=1,Nr
349 CALL GMREDI_CALC_TENSOR(
350 I bi, bj, iMin, iMax, jMin, jMax, k,
351 I sigmaX, sigmaY, sigmaR,
352 I myThid )
353 ENDDO
354 #ifdef ALLOW_AUTODIFF_TAMC
355 ELSE
356 DO k=1, Nr
357 CALL GMREDI_CALC_TENSOR_DUMMY(
358 I bi, bj, iMin, iMax, jMin, jMax, k,
359 I sigmaX, sigmaY, sigmaR,
360 I myThid )
361 ENDDO
362 #endif /* ALLOW_AUTODIFF_TAMC */
363 ENDIF
364 #endif /* ALLOW_GMREDI */
365
366 #ifdef ALLOW_KPP
367 C-- Compute KPP mixing coefficients
368 IF (useKPP) THEN
369 CALL KPP_CALC(
370 I bi, bj, myTime, myThid )
371 ENDIF
372 #endif /* ALLOW_KPP */
373
374 #ifdef ALLOW_AUTODIFF_TAMC
375 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
376 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
377 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
378 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
379 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
380 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
381 #endif /* ALLOW_AUTODIFF_TAMC */
382
383 #ifdef ALLOW_AIM
384 C AIM - atmospheric intermediate model, physics package code.
385 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
386 IF ( useAIM ) THEN
387 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
388 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
389 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
390 ENDIF
391 #endif /* ALLOW_AIM */
392
393
394 C-- Start of thermodynamics loop
395 DO k=Nr,1,-1
396
397 C-- km1 Points to level above k (=k-1)
398 C-- kup Cycles through 1,2 to point to layer above
399 C-- kDown Cycles through 2,1 to point to current layer
400
401 km1 = MAX(1,k-1)
402 kup = 1+MOD(k+1,2)
403 kDown= 1+MOD(k,2)
404
405 iMin = 1-OLx+2
406 iMax = sNx+OLx-1
407 jMin = 1-OLy+2
408 jMax = sNy+OLy-1
409
410 #ifdef ALLOW_AUTODIFF_TAMC
411 CPatrick Is this formula correct?
412 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
413 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
414 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
415 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
416 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
417 #endif /* ALLOW_AUTODIFF_TAMC */
418
419 C-- Get temporary terms used by tendency routines
420 CALL CALC_COMMON_FACTORS (
421 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
422 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
423 I myThid)
424
425 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
426 C-- Calculate the total vertical diffusivity
427 CALL CALC_DIFFUSIVITY(
428 I bi,bj,iMin,iMax,jMin,jMax,k,
429 I maskC,maskup,
430 O KappaRT,KappaRS,KappaRU,KappaRV,
431 I myThid)
432 #endif
433
434 C-- Calculate active tracer tendencies (gT,gS,...)
435 C and step forward storing result in gTnm1, gSnm1, etc.
436 IF ( tempStepping ) THEN
437 CALL CALC_GT(
438 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
439 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
440 I KappaRT,
441 U fVerT,
442 I myTime, myThid)
443 CALL TIMESTEP_TRACER(
444 I bi,bj,iMin,iMax,jMin,jMax,k,
445 I theta, gT,
446 U gTnm1,
447 I myIter, myThid)
448 ENDIF
449 IF ( saltStepping ) THEN
450 CALL CALC_GS(
451 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
452 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
453 I KappaRS,
454 U fVerS,
455 I myTime, myThid)
456 CALL TIMESTEP_TRACER(
457 I bi,bj,iMin,iMax,jMin,jMax,k,
458 I salt, gS,
459 U gSnm1,
460 I myIter, myThid)
461 ENDIF
462
463 #ifdef ALLOW_OBCS
464 C-- Apply open boundary conditions
465 IF (useOBCS) THEN
466 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
467 END IF
468 #endif /* ALLOW_OBCS */
469
470 C-- Freeze water
471 IF (allowFreezing) THEN
472 #ifdef ALLOW_AUTODIFF_TAMC
473 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
474 CADJ & , key = kkey, byte = isbyte
475 #endif /* ALLOW_AUTODIFF_TAMC */
476 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
477 END IF
478
479 C-- end of thermodynamic k loop (Nr:1)
480 ENDDO
481
482
483 #ifdef ALLOW_AUTODIFF_TAMC
484 CPatrick? What about this one?
485 maximpl = 6
486 iikey = (ikey-1)*maximpl
487 #endif /* ALLOW_AUTODIFF_TAMC */
488
489 C-- Implicit diffusion
490 IF (implicitDiffusion) THEN
491
492 IF (tempStepping) THEN
493 #ifdef ALLOW_AUTODIFF_TAMC
494 idkey = iikey + 1
495 #endif /* ALLOW_AUTODIFF_TAMC */
496 CALL IMPLDIFF(
497 I bi, bj, iMin, iMax, jMin, jMax,
498 I deltaTtracer, KappaRT, recip_HFacC,
499 U gTNm1,
500 I myThid )
501 ENDIF
502
503 IF (saltStepping) THEN
504 #ifdef ALLOW_AUTODIFF_TAMC
505 idkey = iikey + 2
506 #endif /* ALLOW_AUTODIFF_TAMC */
507 CALL IMPLDIFF(
508 I bi, bj, iMin, iMax, jMin, jMax,
509 I deltaTtracer, KappaRS, recip_HFacC,
510 U gSNm1,
511 I myThid )
512 ENDIF
513
514 #ifdef ALLOW_OBCS
515 C-- Apply open boundary conditions
516 IF (useOBCS) THEN
517 DO K=1,Nr
518 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
519 ENDDO
520 END IF
521 #endif /* ALLOW_OBCS */
522
523 C-- End If implicitDiffusion
524 ENDIF
525
526
527
528 C-- Start of dynamics loop
529 DO k=1,Nr
530
531 C-- km1 Points to level above k (=k-1)
532 C-- kup Cycles through 1,2 to point to layer above
533 C-- kDown Cycles through 2,1 to point to current layer
534
535 km1 = MAX(1,k-1)
536 kup = 1+MOD(k+1,2)
537 kDown= 1+MOD(k,2)
538
539 iMin = 1-OLx+2
540 iMax = sNx+OLx-1
541 jMin = 1-OLy+2
542 jMax = sNy+OLy-1
543
544 C-- Integrate hydrostatic balance for phiHyd with BC of
545 C phiHyd(z=0)=0
546 C distinguishe between Stagger and Non Stagger time stepping
547 IF (staggerTimeStep) THEN
548 CALL CALC_PHI_HYD(
549 I bi,bj,iMin,iMax,jMin,jMax,k,
550 I gTnm1, gSnm1,
551 U phiHyd,
552 I myThid )
553 ELSE
554 CALL CALC_PHI_HYD(
555 I bi,bj,iMin,iMax,jMin,jMax,k,
556 I theta, salt,
557 U phiHyd,
558 I myThid )
559 ENDIF
560
561 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
562 C and step forward storing the result in gUnm1, gVnm1, etc...
563 IF ( momStepping ) THEN
564 CALL CALC_MOM_RHS(
565 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
566 I phiHyd,KappaRU,KappaRV,
567 U fVerU, fVerV,
568 I myTime, myThid)
569 CALL TIMESTEP(
570 I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
571 I myIter, myThid)
572
573 #ifdef ALLOW_OBCS
574 C-- Apply open boundary conditions
575 IF (useOBCS) THEN
576 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
577 END IF
578 #endif /* ALLOW_OBCS */
579
580 #ifdef ALLOW_AUTODIFF_TAMC
581 #ifdef INCLUDE_CD_CODE
582 ELSE
583 DO j=1-OLy,sNy+OLy
584 DO i=1-OLx,sNx+OLx
585 guCD(i,j,k,bi,bj) = 0.0
586 gvCD(i,j,k,bi,bj) = 0.0
587 END DO
588 END DO
589 #endif /* INCLUDE_CD_CODE */
590 #endif /* ALLOW_AUTODIFF_TAMC */
591 ENDIF
592
593
594 C-- end of dynamics k loop (1:Nr)
595 ENDDO
596
597
598
599 C-- Implicit viscosity
600 IF (implicitViscosity.AND.momStepping) THEN
601 #ifdef ALLOW_AUTODIFF_TAMC
602 idkey = iikey + 3
603 #endif /* ALLOW_AUTODIFF_TAMC */
604 CALL IMPLDIFF(
605 I bi, bj, iMin, iMax, jMin, jMax,
606 I deltaTmom, KappaRU,recip_HFacW,
607 U gUNm1,
608 I myThid )
609 #ifdef ALLOW_AUTODIFF_TAMC
610 idkey = iikey + 4
611 #endif /* ALLOW_AUTODIFF_TAMC */
612 CALL IMPLDIFF(
613 I bi, bj, iMin, iMax, jMin, jMax,
614 I deltaTmom, KappaRV,recip_HFacS,
615 U gVNm1,
616 I myThid )
617
618 #ifdef ALLOW_OBCS
619 C-- Apply open boundary conditions
620 IF (useOBCS) THEN
621 DO K=1,Nr
622 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
623 ENDDO
624 END IF
625 #endif /* ALLOW_OBCS */
626
627 #ifdef INCLUDE_CD_CODE
628 #ifdef ALLOW_AUTODIFF_TAMC
629 idkey = iikey + 5
630 #endif /* ALLOW_AUTODIFF_TAMC */
631 CALL IMPLDIFF(
632 I bi, bj, iMin, iMax, jMin, jMax,
633 I deltaTmom, KappaRU,recip_HFacW,
634 U vVelD,
635 I myThid )
636 #ifdef ALLOW_AUTODIFF_TAMC
637 idkey = iikey + 6
638 #endif /* ALLOW_AUTODIFF_TAMC */
639 CALL IMPLDIFF(
640 I bi, bj, iMin, iMax, jMin, jMax,
641 I deltaTmom, KappaRV,recip_HFacS,
642 U uVelD,
643 I myThid )
644 #endif /* INCLUDE_CD_CODE */
645 C-- End If implicitViscosity.AND.momStepping
646 ENDIF
647
648 ENDDO
649 ENDDO
650
651 RETURN
652 END

  ViewVC Help
Powered by ViewVC 1.1.22