/[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.31 - (show annotations) (download)
Mon Jan 3 02:34:01 2005 UTC (19 years, 4 months ago) by jmc
Branch: MAIN
Changes since 1.30: +11 -1 lines
add diagnostics of RHO (-rhoConst)

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_phi_hyd.F,v 1.30 2003/08/01 04:03:54 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 phiHydF,
13 O phiHydC, 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 | phiHydF(i,j) is the hydrostatic Potential anomaly
25 C | at middle between tracer points k-1,k
26 C | On exit:
27 C | phiHydC(i,j) is the hydrostatic Potential anomaly
28 C | at cell centers (tracer points), level k
29 C | phiHydF(i,j) is the hydrostatic Potential anomaly
30 C | at middle between tracer points k,k+1
31 C | dPhiHydX,Y hydrostatic Potential gradient (X&Y dir)
32 C | at cell centers (tracer points), level k
33 C | integr_GeoPot allows to select one integration method
34 C | 1= Finite volume form ; else= Finite difference form
35 C *==========================================================*
36 C \ev
37 C !USES:
38 IMPLICIT NONE
39 C == Global variables ==
40 #include "SIZE.h"
41 #include "GRID.h"
42 #include "EEPARAMS.h"
43 #include "PARAMS.h"
44 #ifdef ALLOW_AUTODIFF_TAMC
45 #include "tamc.h"
46 #include "tamc_keys.h"
47 #endif /* ALLOW_AUTODIFF_TAMC */
48 #include "SURFACE.h"
49 #include "DYNVARS.h"
50
51 C !INPUT/OUTPUT PARAMETERS:
52 C == Routine arguments ==
53 C bi, bj, k :: tile and level indices
54 C iMin,iMax,jMin,jMax :: computational domain
55 C tFld :: potential temperature
56 C sFld :: salinity
57 C phiHydF :: hydrostatic potential anomaly at middle between
58 C 2 centers (entry: Interf_k ; output: Interf_k+1)
59 C phiHydC :: hydrostatic potential anomaly at cell center
60 C dPhiHydX,Y :: gradient (X & Y dir.) of hydrostatic potential anom.
61 C myTime :: current time
62 C myIter :: current iteration number
63 C myThid :: thread number for this instance of the routine.
64 INTEGER bi,bj,iMin,iMax,jMin,jMax,k
65 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
66 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
67 c _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
68 _RL phiHydF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69 _RL phiHydC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
71 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
72 _RL myTime
73 INTEGER myIter, myThid
74
75 #ifdef INCLUDE_PHIHYD_CALCULATION_CODE
76
77 C !LOCAL VARIABLES:
78 C == Local variables ==
79 INTEGER i,j
80 _RL zero, one, half
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 INTEGER iMnLoc,jMnLoc
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 C--- Calculate density
143 #ifdef ALLOW_AUTODIFF_TAMC
144 kkey = (ikey-1)*Nr + k
145 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
146 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
147 #endif /* ALLOW_AUTODIFF_TAMC */
148 CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, k,
149 & tFld, sFld,
150 & alphaRho, myThid)
151
152 #ifdef ALLOW_DIAGNOSTICS
153 IF ( useDiagnostics )
154 & CALL DIAGNOSTICS_FILL(alphaRho,'RHOAnoma',k,1,2,bi,bj,myThid)
155 #endif
156
157 C Quasi-hydrostatic terms are added in as if they modify the buoyancy
158 IF (quasiHydrostatic) THEN
159 CALL QUASIHYDROSTATICTERMS(bi,bj,k,alphaRho,myThid)
160 ENDIF
161
162 #ifdef NONLIN_FRSURF
163 IF (k.EQ.1 .AND. addSurfPhiAnom) THEN
164 DO j=jMin,jMax
165 DO i=iMin,iMax
166 phiHydF(i,j) = surfPhiFac*etaH(i,j,bi,bj)
167 & *gravity*alphaRho(i,j)*recip_rhoConst
168 ENDDO
169 ENDDO
170 ENDIF
171 #endif /* NONLIN_FRSURF */
172
173 C---- Hydrostatic pressure at cell centers
174
175 IF (integr_GeoPot.EQ.1) THEN
176 C -- Finite Volume Form
177
178 DO j=jMin,jMax
179 DO i=iMin,iMax
180
181 C---------- This discretization is the "finite volume" form
182 C which has not been used to date since it does not
183 C conserve KE+PE exactly even though it is more natural
184 C
185 phiHydC(i,j)=phiHydF(i,j)
186 & + half*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
187 phiHydF(i,j)=phiHydF(i,j)
188 & + drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
189 ENDDO
190 ENDDO
191
192 ELSE
193 C -- Finite Difference Form
194
195 dRlocM=half*drC(k)
196 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
197 IF (k.EQ.Nr) THEN
198 dRlocP=rC(k)-rF(k+1)
199 ELSE
200 dRlocP=half*drC(k+1)
201 ENDIF
202
203 DO j=jMin,jMax
204 DO i=iMin,iMax
205
206 C---------- This discretization is the "energy conserving" form
207 C which has been used since at least Adcroft et al., MWR 1997
208 C
209 phiHydC(i,j)=phiHydF(i,j)
210 & +dRlocM*gravity*alphaRho(i,j)*recip_rhoConst
211 phiHydF(i,j)=phiHydC(i,j)
212 & +dRlocP*gravity*alphaRho(i,j)*recip_rhoConst
213 ENDDO
214 ENDDO
215
216 C -- end if integr_GeoPot = ...
217 ENDIF
218
219 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
220 ELSEIF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
221 C This is the hydrostatic pressure calculation for the Ocean
222 C which uses the FIND_RHO() routine to calculate density
223 C before integrating (1/rho)'*dp over the current layer/interface
224 #ifdef ALLOW_AUTODIFF_TAMC
225 CADJ GENERAL
226 #endif /* ALLOW_AUTODIFF_TAMC */
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 #ifdef ALLOW_DIAGNOSTICS
242 IF ( useDiagnostics )
243 & CALL DIAGNOSTICS_FILL(alphaRho,'RHOAnoma',k,1,2,bi,bj,myThid)
244 #endif
245
246 C-- Calculate specific volume anomaly : alpha' = 1/rho - alpha_Cst
247 DO j=jMin,jMax
248 DO i=iMin,iMax
249 locAlpha=alphaRho(i,j)+rhoConst
250 alphaRho(i,j)=maskC(i,j,k,bi,bj)*
251 & (one/locAlpha - recip_rhoConst)
252 ENDDO
253 ENDDO
254
255 C---- Hydrostatic pressure at cell centers
256
257 IF (integr_GeoPot.EQ.1) THEN
258 C -- Finite Volume Form
259
260 DO j=jMin,jMax
261 DO i=iMin,iMax
262
263 C---------- This discretization is the "finite volume" form
264 C which has not been used to date since it does not
265 C conserve KE+PE exactly even though it is more natural
266 C
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) = ddRloc*alphaRho(i,j)
273 c--to reproduce results of c48d_post: uncomment those 4+1 lines
274 c phiHydC(i,j)=phiHydF(i,j)
275 c & +(hFacC(i,j,k,bi,bj)-half)*drF(k)*alphaRho(i,j)
276 c phiHydF(i,j)=phiHydF(i,j)
277 c & + hFacC(i,j,k,bi,bj)*drF(k)*alphaRho(i,j)
278 ELSE
279 phiHydC(i,j) = phiHydF(i,j) + half*drF(k)*alphaRho(i,j)
280 c phiHydF(i,j) = phiHydF(i,j) + drF(k)*alphaRho(i,j)
281 ENDIF
282 c-- and comment this last one:
283 phiHydF(i,j) = phiHydC(i,j) + half*drF(k)*alphaRho(i,j)
284 c-----
285 ENDDO
286 ENDDO
287
288 ELSE
289 C -- Finite Difference Form, with Part-Cell Bathy
290
291 dRlocM=half*drC(k)
292 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
293 IF (k.EQ.Nr) THEN
294 dRlocP=rC(k)-rF(k+1)
295 ELSE
296 dRlocP=half*drC(k+1)
297 ENDIF
298 rec_dRm = one/(rF(k)-rC(k))
299 rec_dRp = one/(rC(k)-rF(k+1))
300
301 DO j=jMin,jMax
302 DO i=iMin,iMax
303
304 C---------- This discretization is the "energy conserving" form
305
306 IF (k.EQ.ksurfC(i,j,bi,bj)) THEN
307 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
308 #ifdef NONLIN_FRSURF
309 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
310 #endif
311 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*dRlocM
312 & +MIN(zero,ddRloc)*rec_dRp*dRlocP
313 & )*alphaRho(i,j)
314 ELSE
315 phiHydC(i,j) = phiHydF(i,j) + dRlocM*alphaRho(i,j)
316 ENDIF
317 phiHydF(i,j) = phiHydC(i,j) + dRlocP*alphaRho(i,j)
318 ENDDO
319 ENDDO
320
321 C -- end if integr_GeoPot = ...
322 ENDIF
323
324 ELSEIF ( buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
325 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
326 C This is the hydrostatic geopotential calculation for the Atmosphere
327 C The ideal gas law is used implicitly here rather than calculating
328 C the specific volume, analogous to the oceanic case.
329
330 C-- virtual potential temperature anomaly (including water vapour effect)
331 DO j=jMin,jMax
332 DO i=iMin,iMax
333 alphaRho(i,j)=maskC(i,j,k,bi,bj)
334 & *( tFld(i,j,k,bi,bj)*(sFld(i,j,k,bi,bj)*atm_Rq+one)
335 & -tRef(k) )
336 ENDDO
337 ENDDO
338
339 C--- Integrate d Phi / d pi
340
341 IF (integr_GeoPot.EQ.0) THEN
342 C -- Energy Conserving Form, accurate with Full cell topo --
343 C------------ The integration for the first level phi(k=1) is the same
344 C for both the "finite volume" and energy conserving methods.
345 C *NOTE* o Working with geopotential Anomaly, the geopotential boundary
346 C condition is simply Phi-prime(Ro_surf)=0.
347 C o convention ddPI > 0 (same as drF & drC)
348 C-----------------------------------------------------------------------
349 IF (k.EQ.1) THEN
350 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
351 & -((rC( k )/atm_Po)**atm_kappa) )
352 ELSE
353 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
354 & -((rC( k )/atm_Po)**atm_kappa) )*half
355 ENDIF
356 IF (k.EQ.Nr) THEN
357 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
358 & -((rF(k+1)/atm_Po)**atm_kappa) )
359 ELSE
360 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
361 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
362 ENDIF
363 C-------- This discretization is the energy conserving form
364 DO j=jMin,jMax
365 DO i=iMin,iMax
366 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
367 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
368 ENDDO
369 ENDDO
370 C end: Energy Conserving Form, No hFac --
371 C-----------------------------------------------------------------------
372
373 ELSEIF (integr_GeoPot.EQ.1) THEN
374 C -- Finite Volume Form, with Part-Cell Topo, linear in P by Half level
375 C---------
376 C Finite Volume formulation consistent with Partial Cell, linear in p by piece
377 C Note: a true Finite Volume form should be linear between 2 Interf_W :
378 C phi_C = (phi_W_k+ phi_W_k+1)/2 ; but not accurate in Stratosphere (low p)
379 C also: if Interface_W at the middle between tracer levels, this form
380 C is close to the Energy Cons. form in the Interior, except for the
381 C non-linearity in PI(p)
382 C---------
383 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
384 & -((rC( k )/atm_Po)**atm_kappa) )
385 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
386 & -((rF(k+1)/atm_Po)**atm_kappa) )
387 DO j=jMin,jMax
388 DO i=iMin,iMax
389 IF (k.EQ.ksurfC(i,j,bi,bj)) THEN
390 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
391 #ifdef NONLIN_FRSURF
392 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
393 #endif
394 phiHydC(i,j) = ddRloc*recip_drF(k)*2. _d 0
395 & *ddPIm*alphaRho(i,j)
396 ELSE
397 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
398 ENDIF
399 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
400 ENDDO
401 ENDDO
402 C end: Finite Volume Form, with Part-Cell Topo, linear in P by Half level
403 C-----------------------------------------------------------------------
404
405 ELSEIF ( integr_GeoPot.EQ.2
406 & .OR. integr_GeoPot.EQ.3 ) THEN
407 C -- Finite Difference Form, with Part-Cell Topo,
408 C works with Interface_W at the middle between 2.Tracer_Level
409 C and with Tracer_Level at the middle between 2.Interface_W.
410 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
411 C Finite Difference formulation consistent with Partial Cell,
412 C Valid & accurate if Interface_W at middle between tracer levels
413 C linear in p between 2 Tracer levels ; conserve energy in the Interior
414 C---------
415 IF (k.EQ.1) THEN
416 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
417 & -((rC( k )/atm_Po)**atm_kappa) )
418 ELSE
419 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
420 & -((rC( k )/atm_Po)**atm_kappa) )*half
421 ENDIF
422 IF (k.EQ.Nr) THEN
423 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
424 & -((rF(k+1)/atm_Po)**atm_kappa) )
425 ELSE
426 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
427 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
428 ENDIF
429 rec_dRm = one/(rF(k)-rC(k))
430 rec_dRp = one/(rC(k)-rF(k+1))
431 DO j=jMin,jMax
432 DO i=iMin,iMax
433 IF (k.EQ.ksurfC(i,j,bi,bj)) THEN
434 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
435 #ifdef NONLIN_FRSURF
436 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
437 #endif
438 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*ddPIm
439 & +MIN(zero,ddRloc)*rec_dRp*ddPIp )
440 & *alphaRho(i,j)
441 ELSE
442 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
443 ENDIF
444 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
445 ENDDO
446 ENDDO
447 C end: Finite Difference Form, with Part-Cell Topo
448 C-----------------------------------------------------------------------
449
450 ELSE
451 STOP 'CALC_PHI_HYD: Bad integr_GeoPot option !'
452 ENDIF
453
454 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
455 ELSE
456 STOP 'CALC_PHI_HYD: Bad value of buoyancyRelation !'
457 ENDIF
458
459 C--- Diagnose Phi at boundary r=R_low :
460 C = Ocean bottom pressure (Ocean, Z-coord.)
461 C = Sea-surface height (Ocean, P-coord.)
462 C = Top atmosphere height (Atmos, P-coord.)
463 IF (useDiagPhiRlow) THEN
464 CALL DIAGS_PHI_RLOW(
465 I k, bi, bj, iMin,iMax, jMin,jMax,
466 I phiHydF, phiHydC, alphaRho, tFld, sFld,
467 I myTime, myIter, myThid)
468 ENDIF
469
470 C--- Diagnose Full Hydrostatic Potential at cell center level
471 CALL DIAGS_PHI_HYD(
472 I k, bi, bj, iMin,iMax, jMin,jMax,
473 I phiHydC,
474 I myTime, myIter, myThid)
475
476 IF (momPressureForcing) THEN
477 iMnLoc = MAX(1-Olx+1,iMin)
478 jMnLoc = MAX(1-Oly+1,jMin)
479 CALL CALC_GRAD_PHI_HYD(
480 I k, bi, bj, iMnLoc,iMax, jMnLoc,jMax,
481 I phiHydC, alphaRho, tFld, sFld,
482 O dPhiHydX, dPhiHydY,
483 I myTime, myIter, myThid)
484 ENDIF
485
486 #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
487
488 RETURN
489 END

  ViewVC Help
Powered by ViewVC 1.1.22