/[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.64 - (show annotations) (download)
Tue Mar 6 16:59:44 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.63: +11 -14 lines
separate the state variable "eta" from the 2D solver solution cg2d_x
    change Time-Average routines names (new package)

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.63 2001/02/20 15:06:21 jmc 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 "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 #ifdef ALLOW_TIMEAVE
43 #include "TIMEAVE_STATV.h"
44 #endif
45
46 C == Routine arguments ==
47 C myTime - Current time in simulation
48 C myIter - Current iteration number in simulation
49 C myThid - Thread number for this instance of the routine.
50 _RL myTime
51 INTEGER myIter
52 INTEGER myThid
53
54 C == Local variables
55 C xA, yA - Per block temporaries holding face areas
56 C uTrans, vTrans, rTrans - Per block temporaries holding flow
57 C transport
58 C o uTrans: Zonal transport
59 C o vTrans: Meridional transport
60 C o rTrans: Vertical transport
61 C maskC,maskUp o maskC: land/water mask for tracer cells
62 C o maskUp: land/water mask for W points
63 C fVer[STUV] o fVer: Vertical flux term - note fVer
64 C is "pipelined" in the vertical
65 C so we need an fVer for each
66 C variable.
67 C rhoK, rhoKM1 - Density at current level, and level above
68 C phiHyd - Hydrostatic part of the potential phiHydi.
69 C In z coords phiHydiHyd is the hydrostatic
70 C pressure anomaly
71 C In p coords phiHydiHyd is the geopotential
72 C surface height
73 C anomaly.
74 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
75 C phiSurfY or geopotentiel (atmos) in X and Y direction
76 C KappaRT, - Total diffusion in vertical for T and S.
77 C KappaRS (background + spatially varying, isopycnal term).
78 C iMin, iMax - Ranges and sub-block indices on which calculations
79 C jMin, jMax are applied.
80 C bi, bj
81 C k, kup, - Index for layer above and below. kup and kDown
82 C kDown, km1 are switched with layer to be the appropriate
83 C index into fVerTerm.
84 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
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 phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
101 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
102 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
103 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
104 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
105 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
106 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
107
108 C This is currently used by IVDC and Diagnostics
109 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
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 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
118 c CHARACTER*(MAX_LEN_MBUF) suff
119 c LOGICAL DIFFERENT_MULTIPLE
120 c EXTERNAL DIFFERENT_MULTIPLE
121 Cjmc(end)
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 rho = rho ( theta[n], salt[n] )
148 C b = b(rho, theta)
149 C K31 = K31 ( rho )
150 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
151 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
152 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
153 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
154 C
155 C "Time-stepping" or "Prediction"
156 C ================================
157 C The models variables are stepped forward with the appropriate
158 C time-stepping scheme (currently we use Adams-Bashforth II)
159 C - For momentum, the result is always *only* a "prediction"
160 C in that the flow may be divergent and will be "corrected"
161 C later with a surface pressure gradient.
162 C - Normally for tracers the result is the new field at time
163 C level [n+1} *BUT* in the case of implicit diffusion the result
164 C is also *only* a prediction.
165 C - We denote "predictors" with an asterisk (*).
166 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
167 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
168 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
169 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
170 C With implicit diffusion:
171 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
172 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
173 C (1 + dt * K * d_zz) theta[n] = theta*
174 C (1 + dt * K * d_zz) salt[n] = salt*
175 C---
176
177 #ifdef ALLOW_AUTODIFF_TAMC
178 C-- dummy statement to end declaration part
179 ikey = 1
180 #endif /* ALLOW_AUTODIFF_TAMC */
181
182 C-- Set up work arrays with valid (i.e. not NaN) values
183 C These inital values do not alter the numerical results. They
184 C just ensure that all memory references are to valid floating
185 C point numbers. This prevents spurious hardware signals due to
186 C uninitialised but inert locations.
187 DO j=1-OLy,sNy+OLy
188 DO i=1-OLx,sNx+OLx
189 xA(i,j) = 0. _d 0
190 yA(i,j) = 0. _d 0
191 uTrans(i,j) = 0. _d 0
192 vTrans(i,j) = 0. _d 0
193 DO k=1,Nr
194 phiHyd(i,j,k) = 0. _d 0
195 KappaRU(i,j,k) = 0. _d 0
196 KappaRV(i,j,k) = 0. _d 0
197 sigmaX(i,j,k) = 0. _d 0
198 sigmaY(i,j,k) = 0. _d 0
199 sigmaR(i,j,k) = 0. _d 0
200 ENDDO
201 rhoKM1 (i,j) = 0. _d 0
202 rhok (i,j) = 0. _d 0
203 maskC (i,j) = 0. _d 0
204 phiSurfX(i,j) = 0. _d 0
205 phiSurfY(i,j) = 0. _d 0
206 ENDDO
207 ENDDO
208
209
210 #ifdef ALLOW_AUTODIFF_TAMC
211 C-- HPF directive to help TAMC
212 CHPF$ INDEPENDENT
213 #endif /* ALLOW_AUTODIFF_TAMC */
214
215 DO bj=myByLo(myThid),myByHi(myThid)
216
217 #ifdef ALLOW_AUTODIFF_TAMC
218 C-- HPF directive to help TAMC
219 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV
220 CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA
221 CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
222 CHPF$& )
223 #endif /* ALLOW_AUTODIFF_TAMC */
224
225 DO bi=myBxLo(myThid),myBxHi(myThid)
226
227 #ifdef ALLOW_AUTODIFF_TAMC
228 act1 = bi - myBxLo(myThid)
229 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
230
231 act2 = bj - myByLo(myThid)
232 max2 = myByHi(myThid) - myByLo(myThid) + 1
233
234 act3 = myThid - 1
235 max3 = nTx*nTy
236
237 act4 = ikey_dynamics - 1
238
239 ikey = (act1 + 1) + act2*max1
240 & + act3*max1*max2
241 & + act4*max1*max2*max3
242 #endif /* ALLOW_AUTODIFF_TAMC */
243
244 C-- Set up work arrays that need valid initial values
245 DO j=1-OLy,sNy+OLy
246 DO i=1-OLx,sNx+OLx
247 rTrans(i,j) = 0. _d 0
248 fVerT (i,j,1) = 0. _d 0
249 fVerT (i,j,2) = 0. _d 0
250 fVerS (i,j,1) = 0. _d 0
251 fVerS (i,j,2) = 0. _d 0
252 fVerU (i,j,1) = 0. _d 0
253 fVerU (i,j,2) = 0. _d 0
254 fVerV (i,j,1) = 0. _d 0
255 fVerV (i,j,2) = 0. _d 0
256 ENDDO
257 ENDDO
258
259 DO k=1,Nr
260 DO j=1-OLy,sNy+OLy
261 DO i=1-OLx,sNx+OLx
262 C This is currently also used by IVDC and Diagnostics
263 ConvectCount(i,j,k) = 0.
264 KappaRT(i,j,k) = 0. _d 0
265 KappaRS(i,j,k) = 0. _d 0
266 ENDDO
267 ENDDO
268 ENDDO
269
270 iMin = 1-OLx+1
271 iMax = sNx+OLx
272 jMin = 1-OLy+1
273 jMax = sNy+OLy
274
275
276 C-- Start of diagnostic loop
277 DO k=Nr,1,-1
278
279 #ifdef ALLOW_AUTODIFF_TAMC
280 C? Patrick, is this formula correct now that we change the loop range?
281 C? Do we still need this?
282 kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
283 #endif /* ALLOW_AUTODIFF_TAMC */
284
285 C-- Integrate continuity vertically for vertical velocity
286 CALL INTEGRATE_FOR_W(
287 I bi, bj, k, uVel, vVel,
288 O wVel,
289 I myThid )
290
291 #ifdef ALLOW_OBCS
292 #ifdef ALLOW_NONHYDROSTATIC
293 C-- Apply OBC to W if in N-H mode
294 IF (useOBCS.AND.nonHydrostatic) THEN
295 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
296 ENDIF
297 #endif /* ALLOW_NONHYDROSTATIC */
298 #endif /* ALLOW_OBCS */
299
300 C-- Calculate gradients of potential density for isoneutral
301 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
302 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
303 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
304 CALL FIND_RHO(
305 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
306 I theta, salt,
307 O rhoK,
308 I myThid )
309 IF (k.GT.1) CALL FIND_RHO(
310 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
311 I theta, salt,
312 O rhoKm1,
313 I myThid )
314 CALL GRAD_SIGMA(
315 I bi, bj, iMin, iMax, jMin, jMax, k,
316 I rhoK, rhoKm1, rhoK,
317 O sigmaX, sigmaY, sigmaR,
318 I myThid )
319 ENDIF
320
321 C-- Implicit Vertical Diffusion for Convection
322 c ==> should use sigmaR !!!
323 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
324 CALL CALC_IVDC(
325 I bi, bj, iMin, iMax, jMin, jMax, k,
326 I rhoKm1, rhoK,
327 U ConvectCount, KappaRT, KappaRS,
328 I myTime, myIter, myThid)
329 ENDIF
330
331 C-- end of diagnostic k loop (Nr:1)
332 ENDDO
333
334 #ifdef ALLOW_OBCS
335 C-- Calculate future values on open boundaries
336 IF (useOBCS) THEN
337 CALL OBCS_CALC( bi, bj, myTime+deltaT,
338 I uVel, vVel, wVel, theta, salt,
339 I myThid )
340 ENDIF
341 #endif /* ALLOW_OBCS */
342
343 C-- Determines forcing terms based on external fields
344 C relaxation terms, etc.
345 CALL EXTERNAL_FORCING_SURF(
346 I bi, bj, iMin, iMax, jMin, jMax,
347 I myThid )
348
349 #ifdef ALLOW_GMREDI
350 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
351 IF (useGMRedi) THEN
352 DO k=1,Nr
353 CALL GMREDI_CALC_TENSOR(
354 I bi, bj, iMin, iMax, jMin, jMax, k,
355 I sigmaX, sigmaY, sigmaR,
356 I myThid )
357 ENDDO
358 #ifdef ALLOW_AUTODIFF_TAMC
359 ELSE
360 DO k=1, Nr
361 CALL GMREDI_CALC_TENSOR_DUMMY(
362 I bi, bj, iMin, iMax, jMin, jMax, k,
363 I sigmaX, sigmaY, sigmaR,
364 I myThid )
365 ENDDO
366 #endif /* ALLOW_AUTODIFF_TAMC */
367 ENDIF
368 #endif /* ALLOW_GMREDI */
369
370 #ifdef ALLOW_KPP
371 C-- Compute KPP mixing coefficients
372 IF (useKPP) THEN
373 CALL KPP_CALC(
374 I bi, bj, myTime, myThid )
375 ENDIF
376 #endif /* ALLOW_KPP */
377
378 #ifdef ALLOW_AUTODIFF_TAMC
379 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
380 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
381 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
382 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
383 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
384 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
385 #endif /* ALLOW_AUTODIFF_TAMC */
386
387 #ifdef ALLOW_AIM
388 C AIM - atmospheric intermediate model, physics package code.
389 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
390 IF ( useAIM ) THEN
391 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
392 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
393 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
394 ENDIF
395 #endif /* ALLOW_AIM */
396
397
398 C-- Start of thermodynamics loop
399 DO k=Nr,1,-1
400
401 C-- km1 Points to level above k (=k-1)
402 C-- kup Cycles through 1,2 to point to layer above
403 C-- kDown Cycles through 2,1 to point to current layer
404
405 km1 = MAX(1,k-1)
406 kup = 1+MOD(k+1,2)
407 kDown= 1+MOD(k,2)
408
409 iMin = 1-OLx+2
410 iMax = sNx+OLx-1
411 jMin = 1-OLy+2
412 jMax = sNy+OLy-1
413
414 #ifdef ALLOW_AUTODIFF_TAMC
415 CPatrick Is this formula correct?
416 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
417 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
418 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
419 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
420 #endif /* ALLOW_AUTODIFF_TAMC */
421
422 C-- Get temporary terms used by tendency routines
423 CALL CALC_COMMON_FACTORS (
424 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
425 O xA,yA,uTrans,vTrans,rTrans,maskC,maskUp,
426 I myThid)
427
428 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
429 C-- Calculate the total vertical diffusivity
430 CALL CALC_DIFFUSIVITY(
431 I bi,bj,iMin,iMax,jMin,jMax,k,
432 I maskC,maskup,
433 O KappaRT,KappaRS,KappaRU,KappaRV,
434 I myThid)
435 #endif
436
437 C-- Calculate active tracer tendencies (gT,gS,...)
438 C and step forward storing result in gTnm1, gSnm1, etc.
439 IF ( tempStepping ) THEN
440 CALL CALC_GT(
441 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
442 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
443 I KappaRT,
444 U fVerT,
445 I myTime, myThid)
446 CALL TIMESTEP_TRACER(
447 I bi,bj,iMin,iMax,jMin,jMax,k,
448 I theta, gT,
449 U gTnm1,
450 I myIter, myThid)
451 ENDIF
452 IF ( saltStepping ) THEN
453 CALL CALC_GS(
454 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
455 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
456 I KappaRS,
457 U fVerS,
458 I myTime, myThid)
459 CALL TIMESTEP_TRACER(
460 I bi,bj,iMin,iMax,jMin,jMax,k,
461 I salt, gS,
462 U gSnm1,
463 I myIter, myThid)
464 ENDIF
465
466 #ifdef ALLOW_OBCS
467 C-- Apply open boundary conditions
468 IF (useOBCS) THEN
469 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
470 END IF
471 #endif /* ALLOW_OBCS */
472
473 C-- Freeze water
474 IF (allowFreezing) THEN
475 #ifdef ALLOW_AUTODIFF_TAMC
476 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
477 CADJ & , key = kkey, byte = isbyte
478 #endif /* ALLOW_AUTODIFF_TAMC */
479 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
480 END IF
481
482 C-- end of thermodynamic k loop (Nr:1)
483 ENDDO
484
485
486 #ifdef ALLOW_AUTODIFF_TAMC
487 CPatrick? What about this one?
488 maximpl = 6
489 iikey = (ikey-1)*maximpl
490 #endif /* ALLOW_AUTODIFF_TAMC */
491
492 C-- Implicit diffusion
493 IF (implicitDiffusion) THEN
494
495 IF (tempStepping) THEN
496 #ifdef ALLOW_AUTODIFF_TAMC
497 idkey = iikey + 1
498 #endif /* ALLOW_AUTODIFF_TAMC */
499 CALL IMPLDIFF(
500 I bi, bj, iMin, iMax, jMin, jMax,
501 I deltaTtracer, KappaRT, recip_HFacC,
502 U gTNm1,
503 I myThid )
504 ENDIF
505
506 IF (saltStepping) THEN
507 #ifdef ALLOW_AUTODIFF_TAMC
508 idkey = iikey + 2
509 #endif /* ALLOW_AUTODIFF_TAMC */
510 CALL IMPLDIFF(
511 I bi, bj, iMin, iMax, jMin, jMax,
512 I deltaTtracer, KappaRS, recip_HFacC,
513 U gSNm1,
514 I myThid )
515 ENDIF
516
517 #ifdef ALLOW_OBCS
518 C-- Apply open boundary conditions
519 IF (useOBCS) THEN
520 DO K=1,Nr
521 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
522 ENDDO
523 END IF
524 #endif /* ALLOW_OBCS */
525
526 C-- End If implicitDiffusion
527 ENDIF
528
529 C-- Start computation of dynamics
530 iMin = 1-OLx+2
531 iMax = sNx+OLx-1
532 jMin = 1-OLy+2
533 jMax = sNy+OLy-1
534
535 C-- Explicit part of the Surface Pressure Gradient (add in TIMESTEP)
536 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
537 IF (implicSurfPress.NE.1.) THEN
538 DO j=jMin,jMax
539 DO i=iMin,iMax
540 phiSurfX(i,j) = _recip_dxC(i,j,bi,bj)*gBaro
541 & *(etaN(i,j,bi,bj)-etaN(i-1,j,bi,bj))
542 phiSurfY(i,j) = _recip_dyC(i,j,bi,bj)*gBaro
543 & *(etaN(i,j,bi,bj)-etaN(i,j-1,bi,bj))
544 ENDDO
545 ENDDO
546 ENDIF
547
548 C-- Start of dynamics loop
549 DO k=1,Nr
550
551 C-- km1 Points to level above k (=k-1)
552 C-- kup Cycles through 1,2 to point to layer above
553 C-- kDown Cycles through 2,1 to point to current layer
554
555 km1 = MAX(1,k-1)
556 kup = 1+MOD(k+1,2)
557 kDown= 1+MOD(k,2)
558
559 C-- Integrate hydrostatic balance for phiHyd with BC of
560 C phiHyd(z=0)=0
561 C distinguishe between Stagger and Non Stagger time stepping
562 IF (staggerTimeStep) THEN
563 CALL CALC_PHI_HYD(
564 I bi,bj,iMin,iMax,jMin,jMax,k,
565 I gTnm1, gSnm1,
566 U phiHyd,
567 I myThid )
568 ELSE
569 CALL CALC_PHI_HYD(
570 I bi,bj,iMin,iMax,jMin,jMax,k,
571 I theta, salt,
572 U phiHyd,
573 I myThid )
574 ENDIF
575
576 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
577 C and step forward storing the result in gUnm1, gVnm1, etc...
578 IF ( momStepping ) THEN
579 CALL CALC_MOM_RHS(
580 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
581 I phiHyd,KappaRU,KappaRV,
582 U fVerU, fVerV,
583 I myTime, myThid)
584 CALL TIMESTEP(
585 I bi,bj,iMin,iMax,jMin,jMax,k,
586 I phiHyd, phiSurfX, phiSurfY,
587 I myIter, myThid)
588
589 #ifdef ALLOW_OBCS
590 C-- Apply open boundary conditions
591 IF (useOBCS) THEN
592 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
593 END IF
594 #endif /* ALLOW_OBCS */
595
596 #ifdef ALLOW_AUTODIFF_TAMC
597 #ifdef INCLUDE_CD_CODE
598 ELSE
599 DO j=1-OLy,sNy+OLy
600 DO i=1-OLx,sNx+OLx
601 guCD(i,j,k,bi,bj) = 0.0
602 gvCD(i,j,k,bi,bj) = 0.0
603 END DO
604 END DO
605 #endif /* INCLUDE_CD_CODE */
606 #endif /* ALLOW_AUTODIFF_TAMC */
607 ENDIF
608
609
610 C-- end of dynamics k loop (1:Nr)
611 ENDDO
612
613
614
615 C-- Implicit viscosity
616 IF (implicitViscosity.AND.momStepping) THEN
617 #ifdef ALLOW_AUTODIFF_TAMC
618 idkey = iikey + 3
619 #endif /* ALLOW_AUTODIFF_TAMC */
620 CALL IMPLDIFF(
621 I bi, bj, iMin, iMax, jMin, jMax,
622 I deltaTmom, KappaRU,recip_HFacW,
623 U gUNm1,
624 I myThid )
625 #ifdef ALLOW_AUTODIFF_TAMC
626 idkey = iikey + 4
627 #endif /* ALLOW_AUTODIFF_TAMC */
628 CALL IMPLDIFF(
629 I bi, bj, iMin, iMax, jMin, jMax,
630 I deltaTmom, KappaRV,recip_HFacS,
631 U gVNm1,
632 I myThid )
633
634 #ifdef ALLOW_OBCS
635 C-- Apply open boundary conditions
636 IF (useOBCS) THEN
637 DO K=1,Nr
638 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
639 ENDDO
640 END IF
641 #endif /* ALLOW_OBCS */
642
643 #ifdef INCLUDE_CD_CODE
644 #ifdef ALLOW_AUTODIFF_TAMC
645 idkey = iikey + 5
646 #endif /* ALLOW_AUTODIFF_TAMC */
647 CALL IMPLDIFF(
648 I bi, bj, iMin, iMax, jMin, jMax,
649 I deltaTmom, KappaRU,recip_HFacW,
650 U vVelD,
651 I myThid )
652 #ifdef ALLOW_AUTODIFF_TAMC
653 idkey = iikey + 6
654 #endif /* ALLOW_AUTODIFF_TAMC */
655 CALL IMPLDIFF(
656 I bi, bj, iMin, iMax, jMin, jMax,
657 I deltaTmom, KappaRV,recip_HFacS,
658 U uVelD,
659 I myThid )
660 #endif /* INCLUDE_CD_CODE */
661 C-- End If implicitViscosity.AND.momStepping
662 ENDIF
663
664 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
665 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
666 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
667 c WRITE(suff,'(I10.10)') myIter+1
668 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
669 c ENDIF
670 Cjmc(end)
671
672 #ifdef ALLOW_TIMEAVE
673 IF (taveFreq.GT.0.) THEN
674 CALL TIMEAVE_CUMULATE(phiHydtave, phiHyd, Nr,
675 I deltaTclock, bi, bj, myThid)
676 IF (ivdc_kappa.NE.0.) THEN
677 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
678 I deltaTclock, bi, bj, myThid)
679 ENDIF
680 ENDIF
681 #endif /* ALLOW_TIMEAVE */
682
683 ENDDO
684 ENDDO
685
686 RETURN
687 END

  ViewVC Help
Powered by ViewVC 1.1.22