/[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.32 - (show annotations) (download)
Mon Jan 3 03:04:37 2005 UTC (19 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint57g_post, checkpoint57y_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint58, checkpoint57n_post, checkpoint57z_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57f_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57y_pre, checkpoint57c_post, checkpoint57c_pre, checkpoint57e_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint57h_done, checkpoint57j_post, checkpoint57f_pre, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint57x_post
Changes since 1.31: +2 -1 lines
really add diagnostics of RHO (include PACKAGES_CONFIG.h was missing)

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

  ViewVC Help
Powered by ViewVC 1.1.22