/[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.65 - (show annotations) (download)
Thu Mar 8 20:25:01 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint37
Branch point for: pre38
Changes since 1.64: +9 -13 lines
all potentials (cg2d_x, cg3d_x, phiHyd) have units of P/rho in ocean AND atmos

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

  ViewVC Help
Powered by ViewVC 1.1.22