/[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.38 - (show annotations) (download)
Mon Sep 22 17:55:16 2008 UTC (15 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61f, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61h, checkpoint61i
Changes since 1.37: +34 -26 lines
store 3-D (in-situ) density in commom block (this save 1 rho computation);
 re-arange rho diagnostics.

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_phi_hyd.F,v 1.37 2008/09/08 21:45:13 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 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 STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=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 (k.EQ.1 .AND. addSurfPhiAnom) THEN
193 DO j=jMin,jMax
194 DO i=iMin,iMax
195 phiHydF(i,j) = surfPhiFac*etaH(i,j,bi,bj)
196 & *gravity*alphaRho(i,j)*recip_rhoConst
197 ENDDO
198 ENDDO
199 ENDIF
200 #endif /* NONLIN_FRSURF */
201
202 C---- Hydrostatic pressure at cell centers
203
204 IF (integr_GeoPot.EQ.1) THEN
205 C -- Finite Volume Form
206
207 DO j=jMin,jMax
208 DO i=iMin,iMax
209
210 C---------- This discretization is the "finite volume" form
211 C which has not been used to date since it does not
212 C conserve KE+PE exactly even though it is more natural
213 C
214 phiHydC(i,j)=phiHydF(i,j)
215 & + half*drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
216 phiHydF(i,j)=phiHydF(i,j)
217 & + drF(k)*gravity*alphaRho(i,j)*recip_rhoConst
218 ENDDO
219 ENDDO
220
221 ELSE
222 C -- Finite Difference Form
223
224 dRlocM=half*drC(k)
225 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
226 IF (k.EQ.Nr) THEN
227 dRlocP=rC(k)-rF(k+1)
228 ELSE
229 dRlocP=half*drC(k+1)
230 ENDIF
231
232 DO j=jMin,jMax
233 DO i=iMin,iMax
234
235 C---------- This discretization is the "energy conserving" form
236 C which has been used since at least Adcroft et al., MWR 1997
237 C
238 phiHydC(i,j)=phiHydF(i,j)
239 & +dRlocM*gravity*alphaRho(i,j)*recip_rhoConst
240 phiHydF(i,j)=phiHydC(i,j)
241 & +dRlocP*gravity*alphaRho(i,j)*recip_rhoConst
242 ENDDO
243 ENDDO
244
245 C -- end if integr_GeoPot = ...
246 ENDIF
247
248 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
249 ELSEIF ( buoyancyRelation .EQ. 'OCEANICP' ) THEN
250 C This is the hydrostatic pressure calculation for the Ocean
251 C which uses the FIND_RHO() routine to calculate density
252 C before integrating (1/rho)'*dp over the current layer/interface
253 #ifdef ALLOW_AUTODIFF_TAMC
254 CADJ GENERAL
255 #endif /* ALLOW_AUTODIFF_TAMC */
256
257 IF ( implicitIntGravWave .OR. myIter.LT.0 ) THEN
258 C-- Calculate density
259 #ifdef ALLOW_AUTODIFF_TAMC
260 kkey = (ikey-1)*Nr + k
261 CADJ STORE tFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
262 CADJ STORE sFld (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
263 #endif /* ALLOW_AUTODIFF_TAMC */
264 CALL FIND_RHO_2D(
265 I iMin, iMax, jMin, jMax, k,
266 I tFld(1-OLx,1-OLy,k,bi,bj),
267 I sFld(1-OLx,1-OLy,k,bi,bj),
268 O alphaRho,
269 I k, bi, bj, myThid )
270 #ifdef ALLOW_AUTODIFF_TAMC
271 CADJ STORE alphaRho (:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
272 #endif /* ALLOW_AUTODIFF_TAMC */
273 ELSE
274 DO j=jMin,jMax
275 DO i=iMin,iMax
276 alphaRho(i,j) = rhoInSitu(i,j,k,bi,bj)
277 ENDDO
278 ENDDO
279 ENDIF
280
281 C-- Calculate specific volume anomaly : alpha' = 1/rho - alpha_Cst
282 DO j=jMin,jMax
283 DO i=iMin,iMax
284 locAlpha=alphaRho(i,j)+rhoConst
285 alphaRho(i,j)=maskC(i,j,k,bi,bj)*
286 & (one/locAlpha - recip_rhoConst)
287 ENDDO
288 ENDDO
289
290 C---- Hydrostatic pressure at cell centers
291
292 IF (integr_GeoPot.EQ.1) THEN
293 C -- Finite Volume Form
294
295 DO j=jMin,jMax
296 DO i=iMin,iMax
297
298 C---------- This discretization is the "finite volume" form
299 C which has not been used to date since it does not
300 C conserve KE+PE exactly even though it is more natural
301 C
302 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
303 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
304 #ifdef NONLIN_FRSURF
305 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
306 #endif
307 phiHydC(i,j) = ddRloc*alphaRho(i,j)
308 c--to reproduce results of c48d_post: uncomment those 4+1 lines
309 c phiHydC(i,j)=phiHydF(i,j)
310 c & +(hFacC(i,j,k,bi,bj)-half)*drF(k)*alphaRho(i,j)
311 c phiHydF(i,j)=phiHydF(i,j)
312 c & + hFacC(i,j,k,bi,bj)*drF(k)*alphaRho(i,j)
313 ELSE
314 phiHydC(i,j) = phiHydF(i,j) + half*drF(k)*alphaRho(i,j)
315 c phiHydF(i,j) = phiHydF(i,j) + drF(k)*alphaRho(i,j)
316 ENDIF
317 c-- and comment this last one:
318 phiHydF(i,j) = phiHydC(i,j) + half*drF(k)*alphaRho(i,j)
319 c-----
320 ENDDO
321 ENDDO
322
323 ELSE
324 C -- Finite Difference Form, with Part-Cell Bathy
325
326 dRlocM=half*drC(k)
327 IF (k.EQ.1) dRlocM=rF(k)-rC(k)
328 IF (k.EQ.Nr) THEN
329 dRlocP=rC(k)-rF(k+1)
330 ELSE
331 dRlocP=half*drC(k+1)
332 ENDIF
333 rec_dRm = one/(rF(k)-rC(k))
334 rec_dRp = one/(rC(k)-rF(k+1))
335
336 DO j=jMin,jMax
337 DO i=iMin,iMax
338
339 C---------- This discretization is the "energy conserving" form
340
341 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
342 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
343 #ifdef NONLIN_FRSURF
344 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
345 #endif
346 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*dRlocM
347 & +MIN(zero,ddRloc)*rec_dRp*dRlocP
348 & )*alphaRho(i,j)
349 ELSE
350 phiHydC(i,j) = phiHydF(i,j) + dRlocM*alphaRho(i,j)
351 ENDIF
352 phiHydF(i,j) = phiHydC(i,j) + dRlocP*alphaRho(i,j)
353 ENDDO
354 ENDDO
355
356 C -- end if integr_GeoPot = ...
357 ENDIF
358
359 ELSEIF ( buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
360 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
361 C This is the hydrostatic geopotential calculation for the Atmosphere
362 C The ideal gas law is used implicitly here rather than calculating
363 C the specific volume, analogous to the oceanic case.
364
365 C-- virtual potential temperature anomaly (including water vapour effect)
366 DO j=jMin,jMax
367 DO i=iMin,iMax
368 alphaRho(i,j)=maskC(i,j,k,bi,bj)
369 & *( tFld(i,j,k,bi,bj)*(sFld(i,j,k,bi,bj)*atm_Rq+one)
370 & -tRef(k) )
371 ENDDO
372 ENDDO
373
374 C--- Integrate d Phi / d pi
375
376 IF (integr_GeoPot.EQ.0) THEN
377 C -- Energy Conserving Form, accurate with Full cell topo --
378 C------------ The integration for the first level phi(k=1) is the same
379 C for both the "finite volume" and energy conserving methods.
380 C *NOTE* o Working with geopotential Anomaly, the geopotential boundary
381 C condition is simply Phi-prime(Ro_surf)=0.
382 C o convention ddPI > 0 (same as drF & drC)
383 C-----------------------------------------------------------------------
384 IF (k.EQ.1) THEN
385 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
386 & -((rC( k )/atm_Po)**atm_kappa) )
387 ELSE
388 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
389 & -((rC( k )/atm_Po)**atm_kappa) )*half
390 ENDIF
391 IF (k.EQ.Nr) THEN
392 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
393 & -((rF(k+1)/atm_Po)**atm_kappa) )
394 ELSE
395 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
396 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
397 ENDIF
398 C-------- This discretization is the energy conserving form
399 DO j=jMin,jMax
400 DO i=iMin,iMax
401 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
402 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
403 ENDDO
404 ENDDO
405 C end: Energy Conserving Form, No hFac --
406 C-----------------------------------------------------------------------
407
408 ELSEIF (integr_GeoPot.EQ.1) THEN
409 C -- Finite Volume Form, with Part-Cell Topo, linear in P by Half level
410 C---------
411 C Finite Volume formulation consistent with Partial Cell, linear in p by piece
412 C Note: a true Finite Volume form should be linear between 2 Interf_W :
413 C phi_C = (phi_W_k+ phi_W_k+1)/2 ; but not accurate in Stratosphere (low p)
414 C also: if Interface_W at the middle between tracer levels, this form
415 C is close to the Energy Cons. form in the Interior, except for the
416 C non-linearity in PI(p)
417 C---------
418 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
419 & -((rC( k )/atm_Po)**atm_kappa) )
420 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
421 & -((rF(k+1)/atm_Po)**atm_kappa) )
422 DO j=jMin,jMax
423 DO i=iMin,iMax
424 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
425 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
426 #ifdef NONLIN_FRSURF
427 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
428 #endif
429 phiHydC(i,j) = ddRloc*recip_drF(k)*2. _d 0
430 & *ddPIm*alphaRho(i,j)
431 ELSE
432 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
433 ENDIF
434 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
435 ENDDO
436 ENDDO
437 C end: Finite Volume Form, with Part-Cell Topo, linear in P by Half level
438 C-----------------------------------------------------------------------
439
440 ELSEIF ( integr_GeoPot.EQ.2
441 & .OR. integr_GeoPot.EQ.3 ) THEN
442 C -- Finite Difference Form, with Part-Cell Topo,
443 C works with Interface_W at the middle between 2.Tracer_Level
444 C and with Tracer_Level at the middle between 2.Interface_W.
445 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
446 C Finite Difference formulation consistent with Partial Cell,
447 C Valid & accurate if Interface_W at middle between tracer levels
448 C linear in p between 2 Tracer levels ; conserve energy in the Interior
449 C---------
450 IF (k.EQ.1) THEN
451 ddPIm=atm_Cp*( ((rF( k )/atm_Po)**atm_kappa)
452 & -((rC( k )/atm_Po)**atm_kappa) )
453 ELSE
454 ddPIm=atm_Cp*( ((rC(k-1)/atm_Po)**atm_kappa)
455 & -((rC( k )/atm_Po)**atm_kappa) )*half
456 ENDIF
457 IF (k.EQ.Nr) THEN
458 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
459 & -((rF(k+1)/atm_Po)**atm_kappa) )
460 ELSE
461 ddPIp=atm_Cp*( ((rC( k )/atm_Po)**atm_kappa)
462 & -((rC(k+1)/atm_Po)**atm_kappa) )*half
463 ENDIF
464 rec_dRm = one/(rF(k)-rC(k))
465 rec_dRp = one/(rC(k)-rF(k+1))
466 DO j=jMin,jMax
467 DO i=iMin,iMax
468 IF (k.EQ.kSurfC(i,j,bi,bj)) THEN
469 ddRloc = Ro_surf(i,j,bi,bj)-rC(k)
470 #ifdef NONLIN_FRSURF
471 ddRloc = ddRloc + surfPhiFac*etaH(i,j,bi,bj)
472 #endif
473 phiHydC(i,j) =( MAX(zero,ddRloc)*rec_dRm*ddPIm
474 & +MIN(zero,ddRloc)*rec_dRp*ddPIp )
475 & *alphaRho(i,j)
476 ELSE
477 phiHydC(i,j) = phiHydF(i,j) +ddPIm*alphaRho(i,j)
478 ENDIF
479 phiHydF(i,j) = phiHydC(i,j) +ddPIp*alphaRho(i,j)
480 ENDDO
481 ENDDO
482 C end: Finite Difference Form, with Part-Cell Topo
483 C-----------------------------------------------------------------------
484
485 ELSE
486 STOP 'CALC_PHI_HYD: Bad integr_GeoPot option !'
487 ENDIF
488
489 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
490 ELSE
491 STOP 'CALC_PHI_HYD: Bad value of buoyancyRelation !'
492 ENDIF
493
494 C--- Diagnose Phi at boundary r=R_low :
495 C = Ocean bottom pressure (Ocean, Z-coord.)
496 C = Sea-surface height (Ocean, P-coord.)
497 C = Top atmosphere height (Atmos, P-coord.)
498 IF (useDiagPhiRlow) THEN
499 CALL DIAGS_PHI_RLOW(
500 I k, bi, bj, iMin,iMax, jMin,jMax,
501 I phiHydF, phiHydC, alphaRho, tFld, sFld,
502 I myTime, myIter, myThid)
503 ENDIF
504
505 C--- Diagnose Full Hydrostatic Potential at cell center level
506 CALL DIAGS_PHI_HYD(
507 I k, bi, bj, iMin,iMax, jMin,jMax,
508 I phiHydC,
509 I myTime, myIter, myThid)
510
511 IF (momPressureForcing) THEN
512 CALL CALC_GRAD_PHI_HYD(
513 I k, bi, bj, iMin,iMax, jMin,jMax,
514 I phiHydC, alphaRho, tFld, sFld,
515 O dPhiHydX, dPhiHydY,
516 I myTime, myIter, myThid)
517 ENDIF
518
519 #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
520
521 RETURN
522 END

  ViewVC Help
Powered by ViewVC 1.1.22