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