/[MITgcm]/MITgcm/model/src/calc_phi_hyd.F
ViewVC logotype

Contents of /MITgcm/model/src/calc_phi_hyd.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.42 - (show annotations) (download)
Tue Dec 18 01:16:53 2012 UTC (11 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64b
Changes since 1.41: +49 -37 lines
quasi-hydrostatic extended to p & p* coordinate

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_phi_hyd.F,v 1.41 2012/04/11 04:02:05 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: CALC_PHI_HYD
9 C !INTERFACE:
10 SUBROUTINE CALC_PHI_HYD(
11 I bi, bj, iMin, iMax, jMin, jMax, k,
12 I tFld, sFld,
13 U phiHydF,
14 O phiHydC, dPhiHydX, dPhiHydY,
15 I myTime, myIter, myThid )
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE CALC_PHI_HYD |
19 C | o Integrate the hydrostatic relation to find the Hydros. |
20 C *==========================================================*
21 C | Potential (ocean: Pressure/rho ; atmos = geopotential)
22 C | On entry:
23 C | tFld,sFld are the current thermodynamics quantities
24 C | (unchanged on exit)
25 C | phiHydF(i,j) is the hydrostatic Potential anomaly
26 C | at middle between tracer points k-1,k
27 C | On exit:
28 C | phiHydC(i,j) is the hydrostatic Potential anomaly
29 C | at cell centers (tracer points), level k
30 C | phiHydF(i,j) is the hydrostatic Potential anomaly
31 C | at middle between tracer points k,k+1
32 C | dPhiHydX,Y hydrostatic Potential gradient (X&Y dir)
33 C | at cell centers (tracer points), level k
34 C | integr_GeoPot allows to select one integration method
35 C | 1= Finite volume form ; else= Finite difference form
36 C *==========================================================*
37 C \ev
38 C !USES:
39 IMPLICIT NONE
40 C == Global variables ==
41 #include "SIZE.h"
42 #include "GRID.h"
43 #include "EEPARAMS.h"
44 #include "PARAMS.h"
45 #ifdef ALLOW_AUTODIFF_TAMC
46 #include "tamc.h"
47 #include "tamc_keys.h"
48 #endif /* ALLOW_AUTODIFF_TAMC */
49 #include "SURFACE.h"
50 #include "DYNVARS.h"
51
52 C !INPUT/OUTPUT PARAMETERS:
53 C == Routine arguments ==
54 C bi, bj, k :: tile and level indices
55 C iMin,iMax,jMin,jMax :: computational domain
56 C tFld :: potential temperature
57 C sFld :: salinity
58 C phiHydF :: hydrostatic potential anomaly at middle between
59 C 2 centers (entry: Interf_k ; output: Interf_k+1)
60 C phiHydC :: hydrostatic potential anomaly at cell center
61 C dPhiHydX,Y :: gradient (X & Y dir.) of hydrostatic potential anom.
62 C myTime :: current time
63 C myIter :: current iteration number
64 C myThid :: thread number for this instance of the routine.
65 INTEGER bi,bj,iMin,iMax,jMin,jMax,k
66 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
67 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
68 c _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
69 _RL phiHydF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70 _RL phiHydC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71 _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 _RL myTime
74 INTEGER myIter, myThid
75
76 #ifdef INCLUDE_PHIHYD_CALCULATION_CODE
77
78 C !LOCAL VARIABLES:
79 C == Local variables ==
80 INTEGER i,j
81 _RL alphaRho(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL dRlocM,dRlocP, ddRloc, locAlpha
83 _RL ddPIm, ddPIp, rec_dRm, rec_dRp
84 _RL surfPhiFac
85 LOGICAL useDiagPhiRlow, addSurfPhiAnom
86 CEOP
87 useDiagPhiRlow = .TRUE.
88 addSurfPhiAnom = select_rStar.EQ.0 .AND. nonlinFreeSurf.GE.4
89 surfPhiFac = 0.
90 IF (addSurfPhiAnom) surfPhiFac = 1.
91
92 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
93 C Atmosphere:
94 C integr_GeoPot => select one option for the integration of the Geopotential:
95 C = 0 : Energy Conserving Form, accurate with Topo full cell;
96 C = 1 : Finite Volume Form, with Part-Cell, linear in P by Half level;
97 C =2,3: Finite Difference Form, with Part-Cell,
98 C linear in P between 2 Tracer levels.
99 C can handle both cases: Tracer lev at the middle of InterFace_W
100 C and InterFace_W at the middle of Tracer lev;
101 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
102
103 #ifdef ALLOW_AUTODIFF_TAMC
104 act1 = bi - myBxLo(myThid)
105 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
106
107 act2 = bj - myByLo(myThid)
108 max2 = myByHi(myThid) - myByLo(myThid) + 1
109
110 act3 = myThid - 1
111 max3 = nTx*nTy
112
113 act4 = ikey_dynamics - 1
114
115 ikey = (act1 + 1) + act2*max1
116 & + act3*max1*max2
117 & + act4*max1*max2*max3
118 #endif /* ALLOW_AUTODIFF_TAMC */
119
120 C-- Initialize phiHydF to zero :
121 C note: atmospheric_loading or Phi_topo anomaly are incorporated
122 C later in S/R calc_grad_phi_hyd
123 IF (k.EQ.1) THEN
124 DO j=1-OLy,sNy+OLy
125 DO i=1-OLx,sNx+OLx
126 phiHydF(i,j) = 0.
127 ENDDO
128 ENDDO
129 ENDIF
130
131 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
132 IF ( buoyancyRelation .EQ. 'OCEANIC' ) THEN
133 C This is the hydrostatic pressure calculation for the Ocean
134 C which uses the FIND_RHO() routine to calculate density
135 C before integrating g*rho over the current layer/interface
136 #ifdef ALLOW_AUTODIFF_TAMC
137 CADJ GENERAL
138 #endif /* ALLOW_AUTODIFF_TAMC */
139
140 IF ( implicitIntGravWave .OR. myIter.LT.0 ) THEN
141 C--- Calculate density
142 #ifdef ALLOW_AUTODIFF_TAMC
143 kkey = (ikey-1)*Nr + k
144 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
145 CADJ & kind = isbyte
146 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
147 CADJ & kind = isbyte
148 #endif /* ALLOW_AUTODIFF_TAMC */
149 CALL FIND_RHO_2D(
150 I iMin, iMax, jMin, jMax, k,
151 I tFld(1-OLx,1-OLy,k,bi,bj),
152 I sFld(1-OLx,1-OLy,k,bi,bj),
153 O alphaRho,
154 I k, bi, bj, myThid )
155 ELSE
156 DO j=jMin,jMax
157 DO i=iMin,iMax
158 alphaRho(i,j) = rhoInSitu(i,j,k,bi,bj)
159 ENDDO
160 ENDDO
161 ENDIF
162
163 #ifdef ALLOW_SHELFICE
164 C mask rho, so that there is no contribution of phiHyd from
165 C overlying shelfice (whose density we do not know)
166 IF ( useShelfIce .AND. useDOWN_SLOPE ) THEN
167 C- note: does not work for down_slope pkg which needs rho below the bottom.
168 C setting rho=0 above the ice-shelf base is enough (and works in both cases)
169 C but might be slower (--> keep original masking if not using down_slope pkg)
170 DO j=jMin,jMax
171 DO i=iMin,iMax
172 IF ( k.LT.kSurfC(i,j,bi,bj) ) alphaRho(i,j) = 0. _d 0
173 ENDDO
174 ENDDO
175 ELSEIF ( useShelfIce ) THEN
176 DO j=jMin,jMax
177 DO i=iMin,iMax
178 alphaRho(i,j) = alphaRho(i,j)*maskC(i,j,k,bi,bj)
179 ENDDO
180 ENDDO
181 ENDIF
182 #endif /* ALLOW_SHELFICE */
183
184 #ifdef ALLOW_MOM_COMMON
185 C-- Quasi-hydrostatic terms are added in as if they modify the buoyancy
186 IF (quasiHydrostatic) THEN
187 CALL MOM_QUASIHYDROSTATIC(bi,bj,k,uVel,vVel,alphaRho,myThid)
188 ENDIF
189 #endif /* ALLOW_MOM_COMMON */
190
191 #ifdef NONLIN_FRSURF
192 IF ( addSurfPhiAnom .AND.
193 & uniformFreeSurfLev .AND. k.EQ.1 ) THEN
194 DO j=jMin,jMax
195 DO i=iMin,iMax
196 phiHydF(i,j) = surfPhiFac*etaH(i,j,bi,bj)
197 & *gravity*alphaRho(i,j)*recip_rhoConst
198 ENDDO
199 ENDDO
200 ENDIF
201 #endif /* NONLIN_FRSURF */
202
203 C---- Hydrostatic pressure at cell centers
204
205 IF (integr_GeoPot.EQ.1) THEN
206 C -- Finite Volume Form
207
208 C---------- This discretization is the "finite volume" form
209 C which has not been used to date since it does not
210 C conserve KE+PE exactly even though it is more natural
211
212 IF ( uniformFreeSurfLev ) THEN
213 DO j=jMin,jMax
214 DO i=iMin,iMax
215 phiHydC(i,j) = phiHydF(i,j)
216 & + halfRL*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
217 phiHydF(i,j) = phiHydF(i,j)
218 & + drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
219 ENDDO
220 ENDDO
221 ELSE
222 DO j=jMin,jMax
223 DO i=iMin,iMax
224 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
225 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
226 #ifdef NONLIN_FRSURF
227 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
228 #endif
229 phiHydC(i,j) = ddRloc*gravity*alphaRho(i,j)*recip_rhoConst
230 ELSE
231 phiHydC(i,j) = phiHydF(i,j)
232 & + halfRL*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
233 ENDIF
234 phiHydF(i,j) = phiHydC(i,j)
235 & + halfRL*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
236 ENDDO
237 ENDDO
238 ENDIF
239
240 ELSE
241 C -- Finite Difference Form
242
243 C---------- This discretization is the "energy conserving" form
244 C which has been used since at least Adcroft et al., MWR 1997
245
246 dRlocM = halfRL*drC(k)
247 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
248 IF (k.EQ.Nr) THEN
249 dRlocP=rC(k)-rF(k+1)
250 ELSE
251 dRlocP=halfRL*drC(k+1)
252 ENDIF
253 IF ( uniformFreeSurfLev ) THEN
254 DO j=jMin,jMax
255 DO i=iMin,iMax
256 phiHydC(i,j) = phiHydF(i,j)
257 & +dRlocM*gravity*alphaRho(i,j)*recip_rhoConst
258 phiHydF(i,j) = phiHydC(i,j)
259 & +dRlocP*gravity*alphaRho(i,j)*recip_rhoConst
260 ENDDO
261 ENDDO
262 ELSE
263 rec_dRm = oneRL/(rF(k)-rC(k))
264 rec_dRp = oneRL/(rC(k)-rF(k+1))
265 DO j=jMin,jMax
266 DO i=iMin,iMax
267 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
268 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
269 #ifdef NONLIN_FRSURF
270 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
271 #endif
272 phiHydC(i,j) =( MAX(zeroRL,ddRloc)*rec_dRm*dRlocM
273 & +MIN(zeroRL,ddRloc)*rec_dRp*dRlocP
274 & )*gravity*alphaRho(i,j)*recip_rhoConst
275 ELSE
276 phiHydC(i,j) = phiHydF(i,j)
277 & +dRlocM*gravity*alphaRho(i,j)*recip_rhoConst
278 ENDIF
279 phiHydF(i,j) = phiHydC(i,j)
280 & +dRlocP*gravity*alphaRho(i,j)*recip_rhoConst
281 ENDDO
282 ENDDO
283 ENDIF
284
285 C -- end if integr_GeoPot = ...
286 ENDIF
287
288 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
289 ELSEIF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
290 C This is the hydrostatic pressure calculation for the Ocean
291 C which uses the FIND_RHO() routine to calculate density before
292 C integrating (1/rho)_prime*dp over the current layer/interface
293 #ifdef ALLOW_AUTODIFF_TAMC
294 CADJ GENERAL
295 #endif /* ALLOW_AUTODIFF_TAMC */
296
297 IF ( implicitIntGravWave .OR. myIter.LT.0 ) THEN
298 C-- Calculate density
299 #ifdef ALLOW_AUTODIFF_TAMC
300 kkey = (ikey-1)*Nr + k
301 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
302 CADJ & kind = isbyte
303 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
304 CADJ & kind = isbyte
305 #endif /* ALLOW_AUTODIFF_TAMC */
306 CALL FIND_RHO_2D(
307 I iMin, iMax, jMin, jMax, k,
308 I tFld(1-OLx,1-OLy,k,bi,bj),
309 I sFld(1-OLx,1-OLy,k,bi,bj),
310 O alphaRho,
311 I k, bi, bj, myThid )
312 #ifdef ALLOW_AUTODIFF_TAMC
313 CADJ STORE alphaRho (:,:) = comlev1_bibj_k, key=kkey, byte=isbyte,
314 CADJ & kind = isbyte
315 #endif /* ALLOW_AUTODIFF_TAMC */
316 ELSE
317 DO j=jMin,jMax
318 DO i=iMin,iMax
319 alphaRho(i,j) = rhoInSitu(i,j,k,bi,bj)
320 ENDDO
321 ENDDO
322 ENDIF
323
324 C-- Calculate specific volume anomaly : alpha_prime = 1/rho - alpha_Cst
325 DO j=jMin,jMax
326 DO i=iMin,iMax
327 locAlpha=alphaRho(i,j)+rhoConst
328 alphaRho(i,j)=maskC(i,j,k,bi,bj)*
329 & (oneRL/locAlpha - recip_rhoConst)
330 ENDDO
331 ENDDO
332
333 #ifdef ALLOW_MOM_COMMON
334 C-- Quasi-hydrostatic terms are added as if they modify the specific-volume
335 IF (quasiHydrostatic) THEN
336 CALL MOM_QUASIHYDROSTATIC(bi,bj,k,uVel,vVel,alphaRho,myThid)
337 ENDIF
338 #endif /* ALLOW_MOM_COMMON */
339
340 C---- Hydrostatic pressure at cell centers
341
342 IF (integr_GeoPot.EQ.1) THEN
343 C -- Finite Volume Form
344
345 DO j=jMin,jMax
346 DO i=iMin,iMax
347
348 C---------- This discretization is the "finite volume" form
349 C which has not been used to date since it does not
350 C conserve KE+PE exactly even though it is more natural
351
352 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
353 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
354 #ifdef NONLIN_FRSURF
355 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
356 #endif
357 phiHydC(i,j) = ddRloc*alphaRho(i,j)
358 c--to reproduce results of c48d_post: uncomment those 4+1 lines
359 c phiHydC(i,j)=phiHydF(i,j)
360 c & +(hFacC(i,j,k,bi,bj)-halfRL)*drF(k)*alphaRho(i,j)
361 c phiHydF(i,j)=phiHydF(i,j)
362 c & + hFacC(i,j,k,bi,bj)*drF(k)*alphaRho(i,j)
363 ELSE
364 phiHydC(i,j) = phiHydF(i,j) + halfRL*drF(k)*alphaRho(i,j)
365 c phiHydF(i,j) = phiHydF(i,j) + drF(k)*alphaRho(i,j)
366 ENDIF
367 c-- and comment this last one:
368 phiHydF(i,j) = phiHydC(i,j) + halfRL*drF(k)*alphaRho(i,j)
369 c-----
370 ENDDO
371 ENDDO
372
373 ELSE
374 C -- Finite Difference Form, with Part-Cell Bathy
375
376 dRlocM = halfRL*drC(k)
377 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
378 IF (k.EQ.Nr) THEN
379 dRlocP=rC(k)-rF(k+1)
380 ELSE
381 dRlocP=halfRL*drC(k+1)
382 ENDIF
383 rec_dRm = oneRL/(rF(k)-rC(k))
384 rec_dRp = oneRL/(rC(k)-rF(k+1))
385
386 DO j=jMin,jMax
387 DO i=iMin,iMax
388
389 C---------- This discretization is the "energy conserving" form
390
391 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
392 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
393 #ifdef NONLIN_FRSURF
394 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
395 #endif
396 phiHydC(i,j) =( MAX(zeroRL,ddRloc)*rec_dRm*dRlocM
397 & +MIN(zeroRL,ddRloc)*rec_dRp*dRlocP
398 & )*alphaRho(i,j)
399 ELSE
400 phiHydC(i,j) = phiHydF(i,j) + dRlocM*alphaRho(i,j)
401 ENDIF
402 phiHydF(i,j) = phiHydC(i,j) + dRlocP*alphaRho(i,j)
403 ENDDO
404 ENDDO
405
406 C -- end if integr_GeoPot = ...
407 ENDIF
408
409 ELSEIF ( buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
410 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
411 C This is the hydrostatic geopotential calculation for the Atmosphere
412 C The ideal gas law is used implicitly here rather than calculating
413 C the specific volume, analogous to the oceanic case.
414
415 C-- virtual potential temperature anomaly (including water vapour effect)
416 DO j=jMin,jMax
417 DO i=iMin,iMax
418 alphaRho(i,j) = ( tFld(i,j,k,bi,bj)
419 & *( sFld(i,j,k,bi,bj)*atm_Rq + oneRL )
420 & - tRef(k) )*maskC(i,j,k,bi,bj)
421 ENDDO
422 ENDDO
423
424 #ifdef ALLOW_MOM_COMMON
425 C-- Quasi-hydrostatic terms are added in as if they modify the Pot.Temp
426 IF (quasiHydrostatic) THEN
427 CALL MOM_QUASIHYDROSTATIC(bi,bj,k,uVel,vVel,alphaRho,myThid)
428 ENDIF
429 #endif /* ALLOW_MOM_COMMON */
430
431 C--- Integrate d Phi / d pi
432
433 IF (integr_GeoPot.EQ.0) THEN
434 C -- Energy Conserving Form, accurate with Full cell topo --
435 C------------ The integration for the first level phi(k=1) is the same
436 C for both the "finite volume" and energy conserving methods.
437 C *NOTE* o Working with geopotential Anomaly, the geopotential boundary
438 C condition is simply Phi-prime(Ro_surf)=0.
439 C o convention ddPI > 0 (same as drF & drC)
440 C-----------------------------------------------------------------------
441 IF (k.EQ.1) THEN
442 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
443 & -((rC( k )/atm_Po)**atm_kappa) )
444 ELSE
445 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
446 & -((rC( k )/atm_Po)**atm_kappa) )*halfRL
447 ENDIF
448 IF (k.EQ.Nr) THEN
449 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
450 & -((rF(k+1)/atm_Po)**atm_kappa) )
451 ELSE
452 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
453 & -((rC(k+1)/atm_Po)**atm_kappa) )*halfRL
454 ENDIF
455 C-------- This discretization is the energy conserving form
456 DO j=jMin,jMax
457 DO i=iMin,iMax
458 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
459 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
460 ENDDO
461 ENDDO
462 C end: Energy Conserving Form, No hFac --
463 C-----------------------------------------------------------------------
464
465 ELSEIF (integr_GeoPot.EQ.1) THEN
466 C -- Finite Volume Form, with Part-Cell Topo, linear in P by Half level
467 C---------
468 C Finite Volume formulation consistent with Partial Cell, linear in p by piece
469 C Note: a true Finite Volume form should be linear between 2 Interf_W :
470 C phi_C = (phi_W_k+ phi_W_k+1)/2 ; but not accurate in Stratosphere (low p)
471 C also: if Interface_W at the middle between tracer levels, this form
472 C is close to the Energy Cons. form in the Interior, except for the
473 C non-linearity in PI(p)
474 C---------
475 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
476 & -((rC( k )/atm_Po)**atm_kappa) )
477 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
478 & -((rF(k+1)/atm_Po)**atm_kappa) )
479 DO j=jMin,jMax
480 DO i=iMin,iMax
481 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
482 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
483 #ifdef NONLIN_FRSURF
484 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
485 #endif
486 phiHydC(i,j) = ddRloc*recip_drF(k)*2. _d 0
487 & *ddPIm*alphaRho(i,j)
488 ELSE
489 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
490 ENDIF
491 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
492 ENDDO
493 ENDDO
494 C end: Finite Volume Form, with Part-Cell Topo, linear in P by Half level
495 C-----------------------------------------------------------------------
496
497 ELSEIF ( integr_GeoPot.EQ.2
498 & .OR. integr_GeoPot.EQ.3 ) THEN
499 C -- Finite Difference Form, with Part-Cell Topo,
500 C works with Interface_W at the middle between 2.Tracer_Level
501 C and with Tracer_Level at the middle between 2.Interface_W.
502 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
503 C Finite Difference formulation consistent with Partial Cell,
504 C Valid & accurate if Interface_W at middle between tracer levels
505 C linear in p between 2 Tracer levels ; conserve energy in the Interior
506 C---------
507 IF (k.EQ.1) THEN
508 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
509 & -((rC( k )/atm_Po)**atm_kappa) )
510 ELSE
511 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
512 & -((rC( k )/atm_Po)**atm_kappa) )*halfRL
513 ENDIF
514 IF (k.EQ.Nr) THEN
515 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
516 & -((rF(k+1)/atm_Po)**atm_kappa) )
517 ELSE
518 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
519 & -((rC(k+1)/atm_Po)**atm_kappa) )*halfRL
520 ENDIF
521 rec_dRm = oneRL/(rF(k)-rC(k))
522 rec_dRp = oneRL/(rC(k)-rF(k+1))
523 DO j=jMin,jMax
524 DO i=iMin,iMax
525 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
526 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
527 #ifdef NONLIN_FRSURF
528 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
529 #endif
530 phiHydC(i,j) =( MAX(zeroRL,ddRloc)*rec_dRm*ddPIm
531 & +MIN(zeroRL,ddRloc)*rec_dRp*ddPIp
532 & )*alphaRho(i,j)
533 ELSE
534 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
535 ENDIF
536 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
537 ENDDO
538 ENDDO
539 C end: Finite Difference Form, with Part-Cell Topo
540 C-----------------------------------------------------------------------
541
542 ELSE
543 STOP 'CALC_PHI_HYD: Bad integr_GeoPot option !'
544 ENDIF
545
546 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
547 ELSE
548 STOP 'CALC_PHI_HYD: Bad value of buoyancyRelation !'
549 ENDIF
550
551 C--- Diagnose Phi at boundary r=R_low :
552 C = Ocean bottom pressure (Ocean, Z-coord.)
553 C = Sea-surface height (Ocean, P-coord.)
554 C = Top atmosphere height (Atmos, P-coord.)
555 IF (useDiagPhiRlow) THEN
556 CALL DIAGS_PHI_RLOW(
557 I k, bi, bj, iMin,iMax, jMin,jMax,
558 I phiHydF, phiHydC, alphaRho, tFld, sFld,
559 I myTime, myIter, myThid)
560 ENDIF
561
562 C--- Diagnose Full Hydrostatic Potential at cell center level
563 CALL DIAGS_PHI_HYD(
564 I k, bi, bj, iMin,iMax, jMin,jMax,
565 I phiHydC,
566 I myTime, myIter, myThid)
567
568 IF (momPressureForcing) THEN
569 CALL CALC_GRAD_PHI_HYD(
570 I k, bi, bj, iMin,iMax, jMin,jMax,
571 I phiHydC, alphaRho, tFld, sFld,
572 O dPhiHydX, dPhiHydY,
573 I myTime, myIter, myThid)
574 ENDIF
575
576 #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
577
578 RETURN
579 END

  ViewVC Help
Powered by ViewVC 1.1.22