/[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.40 - (show annotations) (download)
Tue Mar 16 00:08:27 2010 UTC (14 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63g, checkpoint63, checkpoint63l, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.39: +4 -4 lines
avoid unbalanced quote (single or double) in commented line

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_phi_hyd.F,v 1.39 2009/02/13 21:56:48 heimbach 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 zero, one, half
82 _RL alphaRho(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL dRlocM,dRlocP, ddRloc, locAlpha
84 _RL ddPIm, ddPIp, rec_dRm, rec_dRp
85 _RL surfPhiFac
86 PARAMETER ( zero= 0. _d 0 , one= 1. _d 0 , half= .5 _d 0 )
87 LOGICAL useDiagPhiRlow, addSurfPhiAnom
88 CEOP
89 useDiagPhiRlow = .TRUE.
90 addSurfPhiAnom = select_rStar.EQ.0 .AND. nonlinFreeSurf.GT.3
91 surfPhiFac = 0.
92 IF (addSurfPhiAnom) surfPhiFac = 1.
93
94 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
95 C Atmosphere:
96 C integr_GeoPot => select one option for the integration of the Geopotential:
97 C = 0 : Energy Conserving Form, accurate with Topo full cell;
98 C = 1 : Finite Volume Form, with Part-Cell, linear in P by Half level;
99 C =2,3: Finite Difference Form, with Part-Cell,
100 C linear in P between 2 Tracer levels.
101 C can handle both cases: Tracer lev at the middle of InterFace_W
102 C and InterFace_W at the middle of Tracer lev;
103 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
104
105 #ifdef ALLOW_AUTODIFF_TAMC
106 act1 = bi - myBxLo(myThid)
107 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
108
109 act2 = bj - myByLo(myThid)
110 max2 = myByHi(myThid) - myByLo(myThid) + 1
111
112 act3 = myThid - 1
113 max3 = nTx*nTy
114
115 act4 = ikey_dynamics - 1
116
117 ikey = (act1 + 1) + act2*max1
118 & + act3*max1*max2
119 & + act4*max1*max2*max3
120 #endif /* ALLOW_AUTODIFF_TAMC */
121
122 C-- Initialize phiHydF to zero :
123 C note: atmospheric_loading or Phi_topo anomaly are incorporated
124 C later in S/R calc_grad_phi_hyd
125 IF (k.EQ.1) THEN
126 DO j=1-Oly,sNy+Oly
127 DO i=1-Olx,sNx+Olx
128 phiHydF(i,j) = 0.
129 ENDDO
130 ENDDO
131 ENDIF
132
133 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
134 IF ( buoyancyRelation .EQ. 'OCEANIC' ) THEN
135 C This is the hydrostatic pressure calculation for the Ocean
136 C which uses the FIND_RHO() routine to calculate density
137 C before integrating g*rho over the current layer/interface
138 #ifdef ALLOW_AUTODIFF_TAMC
139 CADJ GENERAL
140 #endif /* ALLOW_AUTODIFF_TAMC */
141
142 IF ( implicitIntGravWave .OR. myIter.LT.0 ) THEN
143 C--- Calculate density
144 #ifdef ALLOW_AUTODIFF_TAMC
145 kkey = (ikey-1)*Nr + k
146 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
147 CADJ & kind = isbyte
148 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
149 CADJ & kind = isbyte
150 #endif /* ALLOW_AUTODIFF_TAMC */
151 CALL FIND_RHO_2D(
152 I iMin, iMax, jMin, jMax, k,
153 I tFld(1-OLx,1-OLy,k,bi,bj),
154 I sFld(1-OLx,1-OLy,k,bi,bj),
155 O alphaRho,
156 I k, bi, bj, myThid )
157 ELSE
158 DO j=jMin,jMax
159 DO i=iMin,iMax
160 alphaRho(i,j) = rhoInSitu(i,j,k,bi,bj)
161 ENDDO
162 ENDDO
163 ENDIF
164
165 #ifdef ALLOW_SHELFICE
166 C mask rho, so that there is no contribution of phiHyd from
167 C overlying shelfice (whose density we do not know)
168 IF ( useShelfIce .AND. useDOWN_SLOPE ) THEN
169 C- note: does not work for down_slope pkg which needs rho below the bottom.
170 C setting rho=0 above the ice-shelf base is enough (and works in both cases)
171 C but might be slower (--> keep original masking if not using down_slope pkg)
172 DO j=jMin,jMax
173 DO i=iMin,iMax
174 IF ( k.LT.kSurfC(i,j,bi,bj) ) alphaRho(i,j) = 0. _d 0
175 ENDDO
176 ENDDO
177 ELSEIF ( useShelfIce ) THEN
178 DO j=jMin,jMax
179 DO i=iMin,iMax
180 alphaRho(i,j) = alphaRho(i,j)*maskC(i,j,k,bi,bj)
181 ENDDO
182 ENDDO
183 ENDIF
184 #endif /* ALLOW_SHELFICE */
185
186 #ifdef ALLOW_MOM_COMMON
187 C Quasi-hydrostatic terms are added in as if they modify the buoyancy
188 IF (quasiHydrostatic) THEN
189 CALL MOM_QUASIHYDROSTATIC(bi,bj,k,uVel,vVel,alphaRho,myThid)
190 ENDIF
191 #endif /* ALLOW_MOM_COMMON */
192
193 #ifdef NONLIN_FRSURF
194 IF (k.EQ.1 .AND. addSurfPhiAnom) THEN
195 DO j=jMin,jMax
196 DO i=iMin,iMax
197 phiHydF(i,j) = surfPhiFac*etaH(i,j,bi,bj)
198 & *gravity*alphaRho(i,j)*recip_rhoConst
199 ENDDO
200 ENDDO
201 ENDIF
202 #endif /* NONLIN_FRSURF */
203
204 C---- Hydrostatic pressure at cell centers
205
206 IF (integr_GeoPot.EQ.1) THEN
207 C -- Finite Volume Form
208
209 DO j=jMin,jMax
210 DO i=iMin,iMax
211
212 C---------- This discretization is the "finite volume" form
213 C which has not been used to date since it does not
214 C conserve KE+PE exactly even though it is more natural
215 C
216 phiHydC(i,j)=phiHydF(i,j)
217 & + half*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
218 phiHydF(i,j)=phiHydF(i,j)
219 & + drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
220 ENDDO
221 ENDDO
222
223 ELSE
224 C -- Finite Difference Form
225
226 dRlocM=half*drC(k)
227 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
228 IF (k.EQ.Nr) THEN
229 dRlocP=rC(k)-rF(k+1)
230 ELSE
231 dRlocP=half*drC(k+1)
232 ENDIF
233
234 DO j=jMin,jMax
235 DO i=iMin,iMax
236
237 C---------- This discretization is the "energy conserving" form
238 C which has been used since at least Adcroft et al., MWR 1997
239 C
240 phiHydC(i,j)=phiHydF(i,j)
241 & +dRlocM*gravity*alphaRho(i,j)*recip_rhoConst
242 phiHydF(i,j)=phiHydC(i,j)
243 & +dRlocP*gravity*alphaRho(i,j)*recip_rhoConst
244 ENDDO
245 ENDDO
246
247 C -- end if integr_GeoPot = ...
248 ENDIF
249
250 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
251 ELSEIF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
252 C This is the hydrostatic pressure calculation for the Ocean
253 C which uses the FIND_RHO() routine to calculate density before
254 C integrating (1/rho)_prime*dp over the current layer/interface
255 #ifdef ALLOW_AUTODIFF_TAMC
256 CADJ GENERAL
257 #endif /* ALLOW_AUTODIFF_TAMC */
258
259 IF ( implicitIntGravWave .OR. myIter.LT.0 ) THEN
260 C-- Calculate density
261 #ifdef ALLOW_AUTODIFF_TAMC
262 kkey = (ikey-1)*Nr + k
263 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
264 CADJ & kind = isbyte
265 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte,
266 CADJ & kind = isbyte
267 #endif /* ALLOW_AUTODIFF_TAMC */
268 CALL FIND_RHO_2D(
269 I iMin, iMax, jMin, jMax, k,
270 I tFld(1-OLx,1-OLy,k,bi,bj),
271 I sFld(1-OLx,1-OLy,k,bi,bj),
272 O alphaRho,
273 I k, bi, bj, myThid )
274 #ifdef ALLOW_AUTODIFF_TAMC
275 CADJ STORE alphaRho (:,:) = comlev1_bibj_k, key=kkey, byte=isbyte,
276 CADJ & kind = isbyte
277 #endif /* ALLOW_AUTODIFF_TAMC */
278 ELSE
279 DO j=jMin,jMax
280 DO i=iMin,iMax
281 alphaRho(i,j) = rhoInSitu(i,j,k,bi,bj)
282 ENDDO
283 ENDDO
284 ENDIF
285
286 C-- Calculate specific volume anomaly : alpha_prime = 1/rho - alpha_Cst
287 DO j=jMin,jMax
288 DO i=iMin,iMax
289 locAlpha=alphaRho(i,j)+rhoConst
290 alphaRho(i,j)=maskC(i,j,k,bi,bj)*
291 & (one/locAlpha - recip_rhoConst)
292 ENDDO
293 ENDDO
294
295 C---- Hydrostatic pressure at cell centers
296
297 IF (integr_GeoPot.EQ.1) THEN
298 C -- Finite Volume Form
299
300 DO j=jMin,jMax
301 DO i=iMin,iMax
302
303 C---------- This discretization is the "finite volume" form
304 C which has not been used to date since it does not
305 C conserve KE+PE exactly even though it is more natural
306 C
307 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
308 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
309 #ifdef NONLIN_FRSURF
310 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
311 #endif
312 phiHydC(i,j) = ddRloc*alphaRho(i,j)
313 c--to reproduce results of c48d_post: uncomment those 4+1 lines
314 c phiHydC(i,j)=phiHydF(i,j)
315 c & +(hFacC(i,j,k,bi,bj)-half)*drF(k)*alphaRho(i,j)
316 c phiHydF(i,j)=phiHydF(i,j)
317 c & + hFacC(i,j,k,bi,bj)*drF(k)*alphaRho(i,j)
318 ELSE
319 phiHydC(i,j) = phiHydF(i,j) + half*drF(k)*alphaRho(i,j)
320 c phiHydF(i,j) = phiHydF(i,j) + drF(k)*alphaRho(i,j)
321 ENDIF
322 c-- and comment this last one:
323 phiHydF(i,j) = phiHydC(i,j) + half*drF(k)*alphaRho(i,j)
324 c-----
325 ENDDO
326 ENDDO
327
328 ELSE
329 C -- Finite Difference Form, with Part-Cell Bathy
330
331 dRlocM=half*drC(k)
332 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
333 IF (k.EQ.Nr) THEN
334 dRlocP=rC(k)-rF(k+1)
335 ELSE
336 dRlocP=half*drC(k+1)
337 ENDIF
338 rec_dRm = one/(rF(k)-rC(k))
339 rec_dRp = one/(rC(k)-rF(k+1))
340
341 DO j=jMin,jMax
342 DO i=iMin,iMax
343
344 C---------- This discretization is the "energy conserving" form
345
346 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
347 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
348 #ifdef NONLIN_FRSURF
349 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
350 #endif
351 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*dRlocM
352 & +MIN(zero,ddRloc)*rec_dRp*dRlocP
353 & )*alphaRho(i,j)
354 ELSE
355 phiHydC(i,j) = phiHydF(i,j) + dRlocM*alphaRho(i,j)
356 ENDIF
357 phiHydF(i,j) = phiHydC(i,j) + dRlocP*alphaRho(i,j)
358 ENDDO
359 ENDDO
360
361 C -- end if integr_GeoPot = ...
362 ENDIF
363
364 ELSEIF ( buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
365 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
366 C This is the hydrostatic geopotential calculation for the Atmosphere
367 C The ideal gas law is used implicitly here rather than calculating
368 C the specific volume, analogous to the oceanic case.
369
370 C-- virtual potential temperature anomaly (including water vapour effect)
371 DO j=jMin,jMax
372 DO i=iMin,iMax
373 alphaRho(i,j)=maskC(i,j,k,bi,bj)
374 & *( tFld(i,j,k,bi,bj)*(sFld(i,j,k,bi,bj)*atm_Rq+one)
375 & -tRef(k) )
376 ENDDO
377 ENDDO
378
379 C--- Integrate d Phi / d pi
380
381 IF (integr_GeoPot.EQ.0) THEN
382 C -- Energy Conserving Form, accurate with Full cell topo --
383 C------------ The integration for the first level phi(k=1) is the same
384 C for both the "finite volume" and energy conserving methods.
385 C *NOTE* o Working with geopotential Anomaly, the geopotential boundary
386 C condition is simply Phi-prime(Ro_surf)=0.
387 C o convention ddPI > 0 (same as drF & drC)
388 C-----------------------------------------------------------------------
389 IF (k.EQ.1) THEN
390 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
391 & -((rC( k )/atm_Po)**atm_kappa) )
392 ELSE
393 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
394 & -((rC( k )/atm_Po)**atm_kappa) )*half
395 ENDIF
396 IF (k.EQ.Nr) THEN
397 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
398 & -((rF(k+1)/atm_Po)**atm_kappa) )
399 ELSE
400 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
401 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
402 ENDIF
403 C-------- This discretization is the energy conserving form
404 DO j=jMin,jMax
405 DO i=iMin,iMax
406 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
407 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
408 ENDDO
409 ENDDO
410 C end: Energy Conserving Form, No hFac --
411 C-----------------------------------------------------------------------
412
413 ELSEIF (integr_GeoPot.EQ.1) THEN
414 C -- Finite Volume Form, with Part-Cell Topo, linear in P by Half level
415 C---------
416 C Finite Volume formulation consistent with Partial Cell, linear in p by piece
417 C Note: a true Finite Volume form should be linear between 2 Interf_W :
418 C phi_C = (phi_W_k+ phi_W_k+1)/2 ; but not accurate in Stratosphere (low p)
419 C also: if Interface_W at the middle between tracer levels, this form
420 C is close to the Energy Cons. form in the Interior, except for the
421 C non-linearity in PI(p)
422 C---------
423 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
424 & -((rC( k )/atm_Po)**atm_kappa) )
425 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
426 & -((rF(k+1)/atm_Po)**atm_kappa) )
427 DO j=jMin,jMax
428 DO i=iMin,iMax
429 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
430 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
431 #ifdef NONLIN_FRSURF
432 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
433 #endif
434 phiHydC(i,j) = ddRloc*recip_drF(k)*2. _d 0
435 & *ddPIm*alphaRho(i,j)
436 ELSE
437 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
438 ENDIF
439 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
440 ENDDO
441 ENDDO
442 C end: Finite Volume Form, with Part-Cell Topo, linear in P by Half level
443 C-----------------------------------------------------------------------
444
445 ELSEIF ( integr_GeoPot.EQ.2
446 & .OR. integr_GeoPot.EQ.3 ) THEN
447 C -- Finite Difference Form, with Part-Cell Topo,
448 C works with Interface_W at the middle between 2.Tracer_Level
449 C and with Tracer_Level at the middle between 2.Interface_W.
450 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
451 C Finite Difference formulation consistent with Partial Cell,
452 C Valid & accurate if Interface_W at middle between tracer levels
453 C linear in p between 2 Tracer levels ; conserve energy in the Interior
454 C---------
455 IF (k.EQ.1) THEN
456 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
457 & -((rC( k )/atm_Po)**atm_kappa) )
458 ELSE
459 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
460 & -((rC( k )/atm_Po)**atm_kappa) )*half
461 ENDIF
462 IF (k.EQ.Nr) THEN
463 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
464 & -((rF(k+1)/atm_Po)**atm_kappa) )
465 ELSE
466 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
467 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
468 ENDIF
469 rec_dRm = one/(rF(k)-rC(k))
470 rec_dRp = one/(rC(k)-rF(k+1))
471 DO j=jMin,jMax
472 DO i=iMin,iMax
473 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
474 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
475 #ifdef NONLIN_FRSURF
476 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
477 #endif
478 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*ddPIm
479 & +MIN(zero,ddRloc)*rec_dRp*ddPIp )
480 & *alphaRho(i,j)
481 ELSE
482 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
483 ENDIF
484 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
485 ENDDO
486 ENDDO
487 C end: Finite Difference Form, with Part-Cell Topo
488 C-----------------------------------------------------------------------
489
490 ELSE
491 STOP 'CALC_PHI_HYD: Bad integr_GeoPot option !'
492 ENDIF
493
494 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
495 ELSE
496 STOP 'CALC_PHI_HYD: Bad value of buoyancyRelation !'
497 ENDIF
498
499 C--- Diagnose Phi at boundary r=R_low :
500 C = Ocean bottom pressure (Ocean, Z-coord.)
501 C = Sea-surface height (Ocean, P-coord.)
502 C = Top atmosphere height (Atmos, P-coord.)
503 IF (useDiagPhiRlow) THEN
504 CALL DIAGS_PHI_RLOW(
505 I k, bi, bj, iMin,iMax, jMin,jMax,
506 I phiHydF, phiHydC, alphaRho, tFld, sFld,
507 I myTime, myIter, myThid)
508 ENDIF
509
510 C--- Diagnose Full Hydrostatic Potential at cell center level
511 CALL DIAGS_PHI_HYD(
512 I k, bi, bj, iMin,iMax, jMin,jMax,
513 I phiHydC,
514 I myTime, myIter, myThid)
515
516 IF (momPressureForcing) THEN
517 CALL CALC_GRAD_PHI_HYD(
518 I k, bi, bj, iMin,iMax, jMin,jMax,
519 I phiHydC, alphaRho, tFld, sFld,
520 O dPhiHydX, dPhiHydY,
521 I myTime, myIter, myThid)
522 ENDIF
523
524 #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
525
526 RETURN
527 END

  ViewVC Help
Powered by ViewVC 1.1.22