/[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.36 - (show annotations) (download)
Mon Aug 11 22:25:52 2008 UTC (15 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61c
Changes since 1.35: +34 -29 lines
replace calls to "FIND_RHO" with recent version "FIND_RHO_2D"

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

  ViewVC Help
Powered by ViewVC 1.1.22