/[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.63 - (show annotations) (download)
Tue Feb 20 15:06:21 2001 UTC (23 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint36
Changes since 1.62: +28 -11 lines
implement a Crank-Nickelson barotropic time-stepping

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

  ViewVC Help
Powered by ViewVC 1.1.22