/[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.34 - (show annotations) (download)
Mon Mar 20 14:22:26 2006 UTC (18 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58f_post, checkpoint58d_post, checkpoint58t_post, checkpoint58m_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, mitgcm_mapl_00, checkpoint58r_post, checkpoint58n_post, checkpoint58k_post, checkpoint58l_post, checkpoint58g_post, checkpoint58h_post, checkpoint58j_post, checkpoint58i_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.33: +4 -2 lines
move quasiquasihydrostaticterms.F to pkg/mom_common/mom_quasihydrostatic.F
 - fix bug (rhoConst was missing).
 - deal with curvilinear (spherical) grid.

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

  ViewVC Help
Powered by ViewVC 1.1.22