/[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.2.1 - (show annotations) (download)
Fri Mar 30 22:56:00 2001 UTC (23 years, 1 month ago) by jmc
Branch: pre38
Changes since 1.65: +8 -11 lines
use the 3D global center-cell maskC instead of a local one.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.65 2001/03/08 20:25:01 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 maskUp 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 Potential (=pressure/rho0) anomaly
70 C In p coords phiHydiHyd is the geopotential
71 C surface height anomaly.
72 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
73 C phiSurfY or geopotentiel (atmos) in X and Y direction
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 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
89 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
90 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
91 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
92 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
93 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL phiSurfY(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 used by IVDC and Diagnostics
106 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
107
108 INTEGER iMin, iMax
109 INTEGER jMin, jMax
110 INTEGER bi, bj
111 INTEGER i, j
112 INTEGER k, km1, kup, kDown
113
114 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
115 c CHARACTER*(MAX_LEN_MBUF) suff
116 c LOGICAL DIFFERENT_MULTIPLE
117 c EXTERNAL DIFFERENT_MULTIPLE
118 Cjmc(end)
119
120 #ifdef ALLOW_AUTODIFF_TAMC
121 INTEGER isbyte
122 PARAMETER( isbyte = 4 )
123
124 INTEGER act1, act2, act3, act4
125 INTEGER max1, max2, max3
126 INTEGER iikey, kkey
127 INTEGER maximpl
128 #endif /* ALLOW_AUTODIFF_TAMC */
129
130 C--- The algorithm...
131 C
132 C "Correction Step"
133 C =================
134 C Here we update the horizontal velocities with the surface
135 C pressure such that the resulting flow is either consistent
136 C with the free-surface evolution or the rigid-lid:
137 C U[n] = U* + dt x d/dx P
138 C V[n] = V* + dt x d/dy P
139 C
140 C "Calculation of Gs"
141 C ===================
142 C This is where all the accelerations and tendencies (ie.
143 C physics, parameterizations etc...) are calculated
144 C rho = rho ( theta[n], salt[n] )
145 C b = b(rho, theta)
146 C K31 = K31 ( rho )
147 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
148 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
149 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
150 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
151 C
152 C "Time-stepping" or "Prediction"
153 C ================================
154 C The models variables are stepped forward with the appropriate
155 C time-stepping scheme (currently we use Adams-Bashforth II)
156 C - For momentum, the result is always *only* a "prediction"
157 C in that the flow may be divergent and will be "corrected"
158 C later with a surface pressure gradient.
159 C - Normally for tracers the result is the new field at time
160 C level [n+1} *BUT* in the case of implicit diffusion the result
161 C is also *only* a prediction.
162 C - We denote "predictors" with an asterisk (*).
163 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
164 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
165 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
166 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
167 C With implicit diffusion:
168 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
169 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
170 C (1 + dt * K * d_zz) theta[n] = theta*
171 C (1 + dt * K * d_zz) salt[n] = salt*
172 C---
173
174 #ifdef ALLOW_AUTODIFF_TAMC
175 C-- dummy statement to end declaration part
176 ikey = 1
177 #endif /* ALLOW_AUTODIFF_TAMC */
178
179 C-- Set up work arrays with valid (i.e. not NaN) values
180 C These inital values do not alter the numerical results. They
181 C just ensure that all memory references are to valid floating
182 C point numbers. This prevents spurious hardware signals due to
183 C uninitialised but inert locations.
184 DO j=1-OLy,sNy+OLy
185 DO i=1-OLx,sNx+OLx
186 xA(i,j) = 0. _d 0
187 yA(i,j) = 0. _d 0
188 uTrans(i,j) = 0. _d 0
189 vTrans(i,j) = 0. _d 0
190 DO k=1,Nr
191 phiHyd(i,j,k) = 0. _d 0
192 KappaRU(i,j,k) = 0. _d 0
193 KappaRV(i,j,k) = 0. _d 0
194 sigmaX(i,j,k) = 0. _d 0
195 sigmaY(i,j,k) = 0. _d 0
196 sigmaR(i,j,k) = 0. _d 0
197 ENDDO
198 rhoKM1 (i,j) = 0. _d 0
199 rhok (i,j) = 0. _d 0
200 phiSurfX(i,j) = 0. _d 0
201 phiSurfY(i,j) = 0. _d 0
202 ENDDO
203 ENDDO
204
205
206 #ifdef ALLOW_AUTODIFF_TAMC
207 C-- HPF directive to help TAMC
208 CHPF$ INDEPENDENT
209 #endif /* ALLOW_AUTODIFF_TAMC */
210
211 DO bj=myByLo(myThid),myByHi(myThid)
212
213 #ifdef ALLOW_AUTODIFF_TAMC
214 C-- HPF directive to help TAMC
215 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV
216 CHPF$& ,phiHyd,utrans,vtrans,xA,yA
217 CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
218 CHPF$& )
219 #endif /* ALLOW_AUTODIFF_TAMC */
220
221 DO bi=myBxLo(myThid),myBxHi(myThid)
222
223 #ifdef ALLOW_AUTODIFF_TAMC
224 act1 = bi - myBxLo(myThid)
225 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
226
227 act2 = bj - myByLo(myThid)
228 max2 = myByHi(myThid) - myByLo(myThid) + 1
229
230 act3 = myThid - 1
231 max3 = nTx*nTy
232
233 act4 = ikey_dynamics - 1
234
235 ikey = (act1 + 1) + act2*max1
236 & + act3*max1*max2
237 & + act4*max1*max2*max3
238 #endif /* ALLOW_AUTODIFF_TAMC */
239
240 C-- Set up work arrays that need valid initial values
241 DO j=1-OLy,sNy+OLy
242 DO i=1-OLx,sNx+OLx
243 rTrans(i,j) = 0. _d 0
244 fVerT (i,j,1) = 0. _d 0
245 fVerT (i,j,2) = 0. _d 0
246 fVerS (i,j,1) = 0. _d 0
247 fVerS (i,j,2) = 0. _d 0
248 fVerU (i,j,1) = 0. _d 0
249 fVerU (i,j,2) = 0. _d 0
250 fVerV (i,j,1) = 0. _d 0
251 fVerV (i,j,2) = 0. _d 0
252 ENDDO
253 ENDDO
254
255 DO k=1,Nr
256 DO j=1-OLy,sNy+OLy
257 DO i=1-OLx,sNx+OLx
258 C This is currently also used by IVDC and Diagnostics
259 ConvectCount(i,j,k) = 0.
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 ENDIF
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 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,
421 O xA,yA,uTrans,vTrans,rTrans,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 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,
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,
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 C-- Start computation of dynamics
526 iMin = 1-OLx+2
527 iMax = sNx+OLx-1
528 jMin = 1-OLy+2
529 jMax = sNy+OLy-1
530
531 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
532 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
533 IF (implicSurfPress.NE.1.) THEN
534 CALL CALC_GRAD_PHI_SURF(
535 I bi,bj,iMin,iMax,jMin,jMax,
536 I etaN,
537 O phiSurfX,phiSurfY,
538 I myThid )
539 ENDIF
540
541 C-- Start of dynamics loop
542 DO k=1,Nr
543
544 C-- km1 Points to level above k (=k-1)
545 C-- kup Cycles through 1,2 to point to layer above
546 C-- kDown Cycles through 2,1 to point to current layer
547
548 km1 = MAX(1,k-1)
549 kup = 1+MOD(k+1,2)
550 kDown= 1+MOD(k,2)
551
552 C-- Integrate hydrostatic balance for phiHyd with BC of
553 C phiHyd(z=0)=0
554 C distinguishe between Stagger and Non Stagger time stepping
555 IF (staggerTimeStep) THEN
556 CALL CALC_PHI_HYD(
557 I bi,bj,iMin,iMax,jMin,jMax,k,
558 I gTnm1, gSnm1,
559 U phiHyd,
560 I myThid )
561 ELSE
562 CALL CALC_PHI_HYD(
563 I bi,bj,iMin,iMax,jMin,jMax,k,
564 I theta, salt,
565 U phiHyd,
566 I myThid )
567 ENDIF
568
569 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
570 C and step forward storing the result in gUnm1, gVnm1, etc...
571 IF ( momStepping ) THEN
572 CALL CALC_MOM_RHS(
573 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
574 I phiHyd,KappaRU,KappaRV,
575 U fVerU, fVerV,
576 I myTime, myThid)
577 CALL TIMESTEP(
578 I bi,bj,iMin,iMax,jMin,jMax,k,
579 I phiHyd, phiSurfX, phiSurfY,
580 I myIter, myThid)
581
582 #ifdef ALLOW_OBCS
583 C-- Apply open boundary conditions
584 IF (useOBCS) THEN
585 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
586 END IF
587 #endif /* ALLOW_OBCS */
588
589 #ifdef ALLOW_AUTODIFF_TAMC
590 #ifdef INCLUDE_CD_CODE
591 ELSE
592 DO j=1-OLy,sNy+OLy
593 DO i=1-OLx,sNx+OLx
594 guCD(i,j,k,bi,bj) = 0.0
595 gvCD(i,j,k,bi,bj) = 0.0
596 END DO
597 END DO
598 #endif /* INCLUDE_CD_CODE */
599 #endif /* ALLOW_AUTODIFF_TAMC */
600 ENDIF
601
602
603 C-- end of dynamics k loop (1:Nr)
604 ENDDO
605
606
607
608 C-- Implicit viscosity
609 IF (implicitViscosity.AND.momStepping) THEN
610 #ifdef ALLOW_AUTODIFF_TAMC
611 idkey = iikey + 3
612 #endif /* ALLOW_AUTODIFF_TAMC */
613 CALL IMPLDIFF(
614 I bi, bj, iMin, iMax, jMin, jMax,
615 I deltaTmom, KappaRU,recip_HFacW,
616 U gUNm1,
617 I myThid )
618 #ifdef ALLOW_AUTODIFF_TAMC
619 idkey = iikey + 4
620 #endif /* ALLOW_AUTODIFF_TAMC */
621 CALL IMPLDIFF(
622 I bi, bj, iMin, iMax, jMin, jMax,
623 I deltaTmom, KappaRV,recip_HFacS,
624 U gVNm1,
625 I myThid )
626
627 #ifdef ALLOW_OBCS
628 C-- Apply open boundary conditions
629 IF (useOBCS) THEN
630 DO K=1,Nr
631 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
632 ENDDO
633 END IF
634 #endif /* ALLOW_OBCS */
635
636 #ifdef INCLUDE_CD_CODE
637 #ifdef ALLOW_AUTODIFF_TAMC
638 idkey = iikey + 5
639 #endif /* ALLOW_AUTODIFF_TAMC */
640 CALL IMPLDIFF(
641 I bi, bj, iMin, iMax, jMin, jMax,
642 I deltaTmom, KappaRU,recip_HFacW,
643 U vVelD,
644 I myThid )
645 #ifdef ALLOW_AUTODIFF_TAMC
646 idkey = iikey + 6
647 #endif /* ALLOW_AUTODIFF_TAMC */
648 CALL IMPLDIFF(
649 I bi, bj, iMin, iMax, jMin, jMax,
650 I deltaTmom, KappaRV,recip_HFacS,
651 U uVelD,
652 I myThid )
653 #endif /* INCLUDE_CD_CODE */
654 C-- End If implicitViscosity.AND.momStepping
655 ENDIF
656
657 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
658 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
659 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
660 c WRITE(suff,'(I10.10)') myIter+1
661 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
662 c ENDIF
663 Cjmc(end)
664
665 #ifdef ALLOW_TIMEAVE
666 IF (taveFreq.GT.0.) THEN
667 CALL TIMEAVE_CUMULATE(phiHydtave, phiHyd, Nr,
668 I deltaTclock, bi, bj, myThid)
669 IF (ivdc_kappa.NE.0.) THEN
670 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
671 I deltaTclock, bi, bj, myThid)
672 ENDIF
673 ENDIF
674 #endif /* ALLOW_TIMEAVE */
675
676 ENDDO
677 ENDDO
678
679 RETURN
680 END

  ViewVC Help
Powered by ViewVC 1.1.22