/[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.27 - (show annotations) (download)
Sun Feb 9 02:58:39 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
Changes since 1.26: +32 -48 lines
o New S/R for diagnostic of bottom pressure
  (phi0surf contribution was missing in checkpoint48d_post)

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_phi_hyd.F,v 1.26 2003/02/09 02:00:50 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: CALC_PHI_HYD
8 C !INTERFACE:
9 SUBROUTINE CALC_PHI_HYD(
10 I bi, bj, iMin, iMax, jMin, jMax, K,
11 I tFld, sFld,
12 U phiHyd,
13 O dPhiHydX, dPhiHydY,
14 I myTime, myIter, myThid)
15 C !DESCRIPTION: \bv
16 C *==========================================================*
17 C | SUBROUTINE CALC_PHI_HYD |
18 C | o Integrate the hydrostatic relation to find the Hydros. |
19 C *==========================================================*
20 C | Potential (ocean: Pressure/rho ; atmos = geopotential)|
21 C | On entry: |
22 C | tFld,sFld are the current thermodynamics quantities|
23 C | (unchanged on exit) |
24 C | phiHyd(i,j,1:k-1) is the hydrostatic Potential |
25 C | at cell centers (tracer points) |
26 C | - 1:k-1 layers are valid |
27 C | - k:Nr layers are invalid |
28 C | phiHyd(i,j,k) is the hydrostatic Potential |
29 C | (ocean only_^) at cell the interface k (w point above) |
30 C | On exit: |
31 C | phiHyd(i,j,1:k) is the hydrostatic Potential |
32 C | at cell centers (tracer points) |
33 C | - 1:k layers are valid |
34 C | - k+1:Nr layers are invalid |
35 C | phiHyd(i,j,k+1) is the hydrostatic Potential (P/rho) |
36 C | (ocean only-^) at cell the interface k+1 (w point below)|
37 C | Atmosphere: |
38 C | integr_GeoPot allows to select one integration method |
39 C | (see the list below) |
40 C *==========================================================*
41 C \ev
42 C !USES:
43 IMPLICIT NONE
44 C == Global variables ==
45 #include "SIZE.h"
46 #include "GRID.h"
47 #include "EEPARAMS.h"
48 #include "PARAMS.h"
49 c #include "FFIELDS.h"
50 #ifdef ALLOW_AUTODIFF_TAMC
51 #include "tamc.h"
52 #include "tamc_keys.h"
53 #endif /* ALLOW_AUTODIFF_TAMC */
54 #include "SURFACE.h"
55 #include "DYNVARS.h"
56
57 C !INPUT/OUTPUT PARAMETERS:
58 C == Routine arguments ==
59 INTEGER bi,bj,iMin,iMax,jMin,jMax,K
60 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
61 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
62 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
63 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
64 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
65 _RL myTime
66 INTEGER myIter, myThid
67
68 #ifdef INCLUDE_PHIHYD_CALCULATION_CODE
69
70 C !LOCAL VARIABLES:
71 C == Local variables ==
72 INTEGER i,j, Kp1
73 _RL zero, one, half
74 _RL alphaRho(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL dRloc,dRlocKp1,locAlpha
76 _RL ddPI, ddPIm, ddPIp, ratioRp, ratioRm
77 INTEGER iMnLoc,jMnLoc
78 PARAMETER ( zero= 0. _d 0 , one= 1. _d 0 , half= .5 _d 0 )
79 LOGICAL useDiagPhiRlow
80 CEOP
81 useDiagPhiRlow = .TRUE.
82
83 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
84 C Atmosphere:
85 C integr_GeoPot => select one option for the integration of the Geopotential:
86 C = 0 : Energy Conserving Form, No hFac ;
87 C = 1 : Finite Volume Form, with hFac, linear in P by Half level;
88 C =2,3: Finite Difference Form, with hFac, linear in P between 2 Tracer levels
89 C 2 : case Tracer level at the middle of InterFace_W;
90 C 3 : case InterFace_W at the middle of Tracer levels;
91 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
92
93 #ifdef ALLOW_AUTODIFF_TAMC
94 act1 = bi - myBxLo(myThid)
95 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
96
97 act2 = bj - myByLo(myThid)
98 max2 = myByHi(myThid) - myByLo(myThid) + 1
99
100 act3 = myThid - 1
101 max3 = nTx*nTy
102
103 act4 = ikey_dynamics - 1
104
105 ikey = (act1 + 1) + act2*max1
106 & + act3*max1*max2
107 & + act4*max1*max2*max3
108 #endif /* ALLOW_AUTODIFF_TAMC */
109
110
111 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
112 IF ( buoyancyRelation .eq. 'OCEANIC' ) THEN
113 C This is the hydrostatic pressure calculation for the Ocean
114 C which uses the FIND_RHO() routine to calculate density
115 C before integrating g*rho over the current layer/interface
116 #ifdef ALLOW_AUTODIFF_TAMC
117 CADJ GENERAL
118 #endif /* ALLOW_AUTODIFF_TAMC */
119
120 dRloc=drC(k)
121 IF (k.EQ.1) dRloc=drF(1)
122 IF (k.EQ.Nr) THEN
123 dRlocKp1=0.
124 ELSE
125 dRlocKp1=drC(k+1)
126 ENDIF
127
128 C-- If this is the top layer we impose the boundary condition
129 C P(z=eta) = P(atmospheric_loading)
130 IF (k.EQ.1) THEN
131 DO j=jMin,jMax
132 DO i=iMin,iMax
133 c phiHyd(i,j,k) = phi0surf(i,j,bi,bj)
134 phiHyd(i,j,k) = 0.
135 ENDDO
136 ENDDO
137 ENDIF
138
139 C Calculate density
140 #ifdef ALLOW_AUTODIFF_TAMC
141 kkey = (ikey-1)*Nr + k
142 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
143 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
144 #endif /* ALLOW_AUTODIFF_TAMC */
145 CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, k,
146 & tFld, sFld,
147 & alphaRho, myThid)
148
149 C Quasi-hydrostatic terms are added in as if they modify the buoyancy
150 IF (quasiHydrostatic) THEN
151 CALL QUASIHYDROSTATICTERMS(bi,bj,k,alphaRho,myThid)
152 ENDIF
153
154 C--- Diagnose Hydrostatic pressure at the bottom:
155 IF (useDiagPhiRlow) THEN
156 CALL DIAGS_PHI_RLOW(
157 I k, bi, bj, iMin,iMax, jMin,jMax,
158 I phiHyd, alphaRho, tFld, sFld,
159 I myTime, myIter, myThid)
160 ENDIF
161
162 C--- Hydrostatic pressure at cell centers
163
164 IF (integr_GeoPot.EQ.1) THEN
165 C -- Finite Volume Form
166
167 DO j=jMin,jMax
168 DO i=iMin,iMax
169
170 C---------- This discretization is the "finite volume" form
171 C which has not been used to date since it does not
172 C conserve KE+PE exactly even though it is more natural
173 C
174 IF (k.LT.Nr) phiHyd(i,j,k+1)=phiHyd(i,j,k)
175 & + drF(K)*gravity*alphaRho(i,j)*recip_rhoConst
176 phiHyd(i,j,k)=phiHyd(i,j,k)+
177 & + half*drF(K)*gravity*alphaRho(i,j)*recip_rhoConst
178
179 ENDDO
180 ENDDO
181
182 ELSE
183 C -- Finite Difference Form
184
185 DO j=jMin,jMax
186 DO i=iMin,iMax
187
188 C---------- This discretization is the "energy conserving" form
189 C which has been used since at least Adcroft et al., MWR 1997
190 C
191 phiHyd(i,j,k)=phiHyd(i,j,k)
192 & +half*dRloc*gravity*alphaRho(i,j)*recip_rhoConst
193 IF (k.LT.Nr) phiHyd(i,j,k+1)=phiHyd(i,j,k)
194 & +half*dRlocKp1*gravity*alphaRho(i,j)*recip_rhoConst
195
196 ENDDO
197 ENDDO
198
199 C -- end if integr_GeoPot = ...
200 ENDIF
201
202 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
203 ELSEIF ( buoyancyRelation .eq. 'OCEANICP' ) THEN
204 C This is the hydrostatic pressure calculation for the Ocean
205 C which uses the FIND_RHO() routine to calculate density
206 C before integrating (1/rho)'*dp over the current layer/interface
207 #ifdef ALLOW_AUTODIFF_TAMC
208 CADJ GENERAL
209 #endif /* ALLOW_AUTODIFF_TAMC */
210
211 dRloc=drC(k)
212 IF (k.EQ.1) dRloc=drF(1)
213 IF (k.EQ.Nr) THEN
214 dRlocKp1=0.
215 ELSE
216 dRlocKp1=drC(k+1)
217 ENDIF
218
219 IF (k.EQ.1) THEN
220 DO j=jMin,jMax
221 DO i=iMin,iMax
222 c phiHyd(i,j,k) = phi0surf(i,j,bi,bj)
223 phiHyd(i,j,k) = 0.
224 ENDDO
225 ENDDO
226 ENDIF
227
228 C-- Calculate density
229 #ifdef ALLOW_AUTODIFF_TAMC
230 kkey = (ikey-1)*Nr + k
231 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
232 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
233 #endif /* ALLOW_AUTODIFF_TAMC */
234 CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, k,
235 & tFld, sFld,
236 & alphaRho, myThid)
237 #ifdef ALLOW_AUTODIFF_TAMC
238 CADJ STORE alphaRho (:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
239 #endif /* ALLOW_AUTODIFF_TAMC */
240
241 C-- Calculate specific volume anomaly : alpha' = 1/rho - alpha_Cst
242 DO j=jMin,jMax
243 DO i=iMin,iMax
244 locAlpha=alphaRho(i,j)+rhoConst
245 alphaRho(i,j)=maskC(i,j,k,bi,bj)*
246 & (one/locAlpha - recip_rhoConst)
247 ENDDO
248 ENDDO
249
250 C--- Diagnose Sea-surface height (Hydrostatic geopotential at r=Rlow):
251 IF (useDiagPhiRlow) THEN
252 CALL DIAGS_PHI_RLOW(
253 I k, bi, bj, iMin,iMax, jMin,jMax,
254 I phiHyd, alphaRho, tFld, sFld,
255 I myTime, myIter, myThid)
256 ENDIF
257
258 C---- Hydrostatic pressure at cell centers
259
260 IF (integr_GeoPot.EQ.1) THEN
261 C -- Finite Volume Form
262
263 DO j=jMin,jMax
264 DO i=iMin,iMax
265
266 C---------- This discretization is the "finite volume" form
267 C which has not been used to date since it does not
268 C conserve KE+PE exactly even though it is more natural
269 C
270 IF (k.LT.Nr) phiHyd(i,j,k+1)=phiHyd(i,j,k)
271 & + hFacC(i,j,k,bi,bj)*drF(K)*alphaRho(i,j)
272 phiHyd(i,j,k)=phiHyd(i,j,k)
273 & +(hFacC(i,j,k,bi,bj)-half)*drF(K)*alphaRho(i,j)
274
275 ENDDO
276 ENDDO
277
278 ELSE
279 C -- Finite Difference Form
280
281 DO j=jMin,jMax
282 DO i=iMin,iMax
283
284 C---------- This discretization is the "energy conserving" form
285
286 phiHyd(i,j,k)=phiHyd(i,j,k)
287 & + half*dRloc*alphaRho(i,j)
288 IF (k.LT.Nr) phiHyd(i,j,k+1)=phiHyd(i,j,k)
289 & + half*dRlocKp1*alphaRho(i,j)
290
291 ENDDO
292 ENDDO
293
294 C -- end if integr_GeoPot = ...
295 ENDIF
296
297 ELSEIF ( buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
298 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
299 C This is the hydrostatic geopotential calculation for the Atmosphere
300 C The ideal gas law is used implicitly here rather than calculating
301 C the specific volume, analogous to the oceanic case.
302
303 C Integrate d Phi / d pi
304
305 IF (integr_GeoPot.EQ.0) THEN
306 C -- Energy Conserving Form, No hFac --
307 C------------ The integration for the first level phi(k=1) is the same
308 C for both the "finite volume" and energy conserving methods.
309 Ci *NOTE* o Working with geopotential Anomaly, the geopotential boundary
310 C condition is simply Phi-prime(Ro_surf)=0.
311 C o convention ddPI > 0 (same as drF & drC)
312 C-----------------------------------------------------------------------
313 IF (K.EQ.1) THEN
314 ddPIp=atm_Cp*( ((rF(K)/atm_Po)**atm_kappa)
315 & -((rC(K)/atm_Po)**atm_kappa) )
316 DO j=jMin,jMax
317 DO i=iMin,iMax
318 c phiHyd(i,j,K)= phi0surf(i,j,bi,bj)+
319 phiHyd(i,j,K)=
320 & ddPIp*maskC(i,j,K,bi,bj)
321 & *(tFld(I,J,K,bi,bj)-tRef(K))
322 ENDDO
323 ENDDO
324 ELSE
325 C-------- This discretization is the energy conserving form
326 ddPI=atm_Cp*( ((rC(K-1)/atm_Po)**atm_kappa)
327 & -((rC( K )/atm_Po)**atm_kappa) )*0.5
328 DO j=jMin,jMax
329 DO i=iMin,iMax
330 phiHyd(i,j,K)=phiHyd(i,j,K-1)
331 & +ddPI*maskC(i,j,K-1,bi,bj)
332 & *(tFld(I,J,K-1,bi,bj)-tRef(K-1))
333 & +ddPI*maskC(i,j, K ,bi,bj)
334 & *(tFld(I,J, K ,bi,bj)-tRef( K ))
335 C Old code (atmos-exact) looked like this
336 Cold phiHyd(i,j,K)=phiHyd(i,j,K-1) - ddPI*
337 Cold & (tFld(I,J,K-1,bi,bj)+tFld(I,J,K,bi,bj)-2.*tRef(K))
338 ENDDO
339 ENDDO
340 ENDIF
341 C end: Energy Conserving Form, No hFac --
342 C-----------------------------------------------------------------------
343
344 ELSEIF (integr_GeoPot.EQ.1) THEN
345 C -- Finite Volume Form, with hFac, linear in P by Half level --
346 C---------
347 C Finite Volume formulation consistent with Partial Cell, linear in p by piece
348 C Note: a true Finite Volume form should be linear between 2 Interf_W :
349 C phi_C = (phi_W_k+ phi_W_k+1)/2 ; but not accurate in Stratosphere (low p)
350 C also: if Interface_W at the middle between tracer levels, this form
351 C is close to the Energy Cons. form in the Interior, except for the
352 C non-linearity in PI(p)
353 C---------
354 IF (K.EQ.1) THEN
355 ddPIp=atm_Cp*( ((rF(K)/atm_Po)**atm_kappa)
356 & -((rC(K)/atm_Po)**atm_kappa) )
357 DO j=jMin,jMax
358 DO i=iMin,iMax
359 c phiHyd(i,j,K)= phi0surf(i,j,bi,bj)+
360 phiHyd(i,j,K)=
361 & ddPIp*_hFacC(I,J, K ,bi,bj)
362 & *(tFld(I,J, K ,bi,bj)-tRef( K ))
363 ENDDO
364 ENDDO
365 ELSE
366 ddPIm=atm_Cp*( ((rC(K-1)/atm_Po)**atm_kappa)
367 & -((rF( K )/atm_Po)**atm_kappa) )
368 ddPIp=atm_Cp*( ((rF( K )/atm_Po)**atm_kappa)
369 & -((rC( K )/atm_Po)**atm_kappa) )
370 DO j=jMin,jMax
371 DO i=iMin,iMax
372 phiHyd(i,j,K) = phiHyd(i,j,K-1)
373 & +ddPIm*_hFacC(I,J,K-1,bi,bj)
374 & *(tFld(I,J,K-1,bi,bj)-tRef(K-1))
375 & +ddPIp*_hFacC(I,J, K ,bi,bj)
376 & *(tFld(I,J, K ,bi,bj)-tRef( K ))
377 ENDDO
378 ENDDO
379 ENDIF
380 C end: Finite Volume Form, with hFac, linear in P by Half level --
381 C-----------------------------------------------------------------------
382
383 ELSEIF (integr_GeoPot.EQ.2) THEN
384 C -- Finite Difference Form, with hFac, Tracer Lev. = middle --
385 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
386 C Finite Difference formulation consistent with Partial Cell,
387 C case Tracer level at the middle of InterFace_W
388 C linear between 2 Tracer levels ; conserve energy in the Interior
389 C---------
390 Kp1 = min(Nr,K+1)
391 IF (K.EQ.1) THEN
392 ddPIm=atm_Cp*( ((rF( K )/atm_Po)**atm_kappa)
393 & -((rC( K )/atm_Po)**atm_kappa) ) * 2. _d 0
394 ddPIp=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa)
395 & -((rC(Kp1)/atm_Po)**atm_kappa) )
396 DO j=jMin,jMax
397 DO i=iMin,iMax
398 c phiHyd(i,j,K)= phi0surf(i,j,bi,bj)+
399 phiHyd(i,j,K)=
400 & ( ddPIm*max(zero, _hFacC(i,j,K,bi,bj)-half)
401 & +ddPIp*min(zero, _hFacC(i,j,K,bi,bj)-half) )
402 & *(tFld(i,j, K ,bi,bj)-tRef( K ))
403 & * maskC(i,j, K ,bi,bj)
404 ENDDO
405 ENDDO
406 ELSE
407 ddPIm=atm_Cp*( ((rC(K-1)/atm_Po)**atm_kappa)
408 & -((rC( K )/atm_Po)**atm_kappa) )
409 ddPIp=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa)
410 & -((rC(Kp1)/atm_Po)**atm_kappa) )
411 DO j=jMin,jMax
412 DO i=iMin,iMax
413 phiHyd(i,j,K) = phiHyd(i,j,K-1)
414 & + ddPIm*0.5
415 & *(tFld(i,j,K-1,bi,bj)-tRef(K-1))
416 & * maskC(i,j,K-1,bi,bj)
417 & +(ddPIm*max(zero, _hFacC(i,j,K,bi,bj)-half)
418 & +ddPIp*min(zero, _hFacC(i,j,K,bi,bj)-half) )
419 & *(tFld(i,j, K ,bi,bj)-tRef( K ))
420 & * maskC(i,j, K ,bi,bj)
421 ENDDO
422 ENDDO
423 ENDIF
424 C end: Finite Difference Form, with hFac, Tracer Lev. = middle --
425 C-----------------------------------------------------------------------
426
427 ELSEIF (integr_GeoPot.EQ.3) THEN
428 C -- Finite Difference Form, with hFac, Interface_W = middle --
429 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
430 C Finite Difference formulation consistent with Partial Cell,
431 C Valid & accurate if Interface_W at middle between tracer levels
432 C linear in p between 2 Tracer levels ; conserve energy in the Interior
433 C---------
434 Kp1 = min(Nr,K+1)
435 IF (K.EQ.1) THEN
436 ratioRm=0.5*drF(K)/(rF(k)-rC(K))
437 ratioRp=drF(K)*recip_drC(Kp1)
438 ddPIm=atm_Cp*( ((rF( K )/atm_Po)**atm_kappa)
439 & -((rC( K )/atm_Po)**atm_kappa) ) * 2. _d 0
440 ddPIp=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa)
441 & -((rC(Kp1)/atm_Po)**atm_kappa) )
442 DO j=jMin,jMax
443 DO i=iMin,iMax
444 c phiHyd(i,j,K)= phi0surf(i,j,bi,bj)+
445 phiHyd(i,j,K)=
446 & ( ddPIm*max(zero,(_hFacC(i,j,K,bi,bj)-one)*ratioRm+half)
447 & +ddPIp*min(zero, _hFacC(i,j,K,bi,bj)*ratioRp -half) )
448 & *(tFld(i,j, K ,bi,bj)-tRef( K ))
449 & * maskC(i,j, K ,bi,bj)
450 ENDDO
451 ENDDO
452 ELSE
453 ratioRm=drF(K)*recip_drC(K)
454 ratioRp=drF(K)*recip_drC(Kp1)
455 ddPIm=atm_Cp*( ((rC(K-1)/atm_Po)**atm_kappa)
456 & -((rC( K )/atm_Po)**atm_kappa) )
457 ddPIp=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa)
458 & -((rC(Kp1)/atm_Po)**atm_kappa) )
459 DO j=jMin,jMax
460 DO i=iMin,iMax
461 phiHyd(i,j,K) = phiHyd(i,j,K-1)
462 & + ddPIm*0.5
463 & *(tFld(i,j,K-1,bi,bj)-tRef(K-1))
464 & * maskC(i,j,K-1,bi,bj)
465 & +(ddPIm*max(zero,(_hFacC(i,j,K,bi,bj)-one)*ratioRm+half)
466 & +ddPIp*min(zero, _hFacC(i,j,K,bi,bj)*ratioRp -half) )
467 & *(tFld(i,j, K ,bi,bj)-tRef( K ))
468 & * maskC(i,j, K ,bi,bj)
469 ENDDO
470 ENDDO
471 ENDIF
472 C end: Finite Difference Form, with hFac, Interface_W = middle --
473 C-----------------------------------------------------------------------
474
475 ELSE
476 STOP 'CALC_PHI_HYD: Bad integr_GeoPot option !'
477 ENDIF
478
479 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
480 ELSE
481 STOP 'CALC_PHI_HYD: Bad value of buoyancyRelation !'
482 ENDIF
483
484 IF (momPressureForcing) THEN
485 iMnLoc = MAX(1-Olx+1,iMin)
486 jMnLoc = MAX(1-Oly+1,jMin)
487 CALL CALC_GRAD_PHI_HYD(
488 I k, bi, bj, iMnLoc,iMax, jMnLoc,jMax,
489 I phiHyd, alphaRho, tFld, sFld,
490 O dPhiHydX, dPhiHydY,
491 I myTime, myIter, myThid)
492 ENDIF
493
494 #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
495
496 RETURN
497 END

  ViewVC Help
Powered by ViewVC 1.1.22