/[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.37 - (show annotations) (download)
Mon Sep 8 21:45:13 2008 UTC (15 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.36: +15 -6 lines
fix for useShelfIce & useDOWN_SLOPE & useDynP_inEos_Zc

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

  ViewVC Help
Powered by ViewVC 1.1.22