Skip to content

Commit 5f23058

Browse files
authored
Bugfix for mixing up Rho0 and rho_ref in Boussinesq PGF (#837)
* Remove unused local variable I_Rho in int_density I_Rho = 1.0/rho_0 is never used and therefore removed. * Add runtime flag RHO_PGF_REF_BUG This new parameter addresses a bug in Boussinesq finite volume pressure gradient forces (MOM_PressureForce_FV) that the mean seawater density Rho_0 and reference density Rho_ref are used incorrectly in several instances. It should be noted that by default Rho_ref is Rho_0 and Rho_ref is rarely set to other than Rho_0. * Bugfix: reference height offset for SAL using bpa Recover reference height offset (Z_ref) in calculating bottom pressure for self-attraction and loading in Boussinesq mode.
1 parent 86ed81f commit 5f23058

File tree

2 files changed

+54
-40
lines changed

2 files changed

+54
-40
lines changed

src/core/MOM_PressureForce_FV.F90

Lines changed: 54 additions & 34 deletions
Original file line numberDiff line numberDiff line change
@@ -43,8 +43,9 @@ module MOM_PressureForce_FV
4343
logical :: sal_use_bpa = .false. !< If true, use bottom pressure anomaly instead of SSH
4444
!! to calculate SAL.
4545
logical :: tides = .false. !< If true, apply tidal momentum forcing.
46-
real :: Rho0 !< The density used in the Boussinesq
47-
!! approximation [R ~> kg m-3].
46+
real :: rho_ref !< The reference density that is subtracted off when calculating pressure
47+
!! gradient forces [R ~> kg m-3].
48+
logical :: rho_ref_bug !< If true, recover a bug that mixes GV%Rho0 and CS%rho_ref in Boussinesq mode.
4849
real :: GFS_scale !< A scaling of the surface pressure gradients to
4950
!! allow the use of a reduced gravity model [nondim].
5051
type(time_type), pointer :: Time !< A pointer to the ocean model's clock.
@@ -278,7 +279,7 @@ subroutine PressureForce_FV_nonBouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_
278279

279280
H_to_RL2_T2 = GV%g_Earth*GV%H_to_RZ
280281
dp_neglect = GV%g_Earth*GV%H_to_RZ * GV%H_subroundoff
281-
alpha_ref = 1.0 / CS%Rho0
282+
alpha_ref = 1.0 / CS%rho_ref
282283
I_gEarth = 1.0 / GV%g_Earth
283284
p_nonvanished = GV%g_Earth*GV%H_to_RZ*CS%h_nonvanished
284285

@@ -977,7 +978,10 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
977978
! consistent with what is used in the density integral routines [R L2 T-2 ~> Pa]
978979
real :: p0(SZI_(G)) ! An array of zeros to use for pressure [R L2 T-2 ~> Pa].
979980
real :: dz_geo_sfc ! The change in surface geopotential height between adjacent cells [L2 T-2 ~> m2 s-2]
980-
real :: GxRho ! The gravitational acceleration times density [R L2 Z-1 T-2 ~> Pa m-1]
981+
real :: GxRho0 ! The gravitational acceleration times mean ocean density [R L2 Z-1 T-2 ~> Pa m-1]
982+
real :: GxRho_ref ! The gravitational acceleration times reference density [R L2 Z-1 T-2 ~> Pa m-1]
983+
real :: rho0_int_density ! Rho0 used in int_density_dz_* subroutines [R ~> kg m-3]
984+
real :: rho0_set_pbce ! Rho0 used in set_pbce_Bouss subroutine [R ~> kg m-3]
981985
real :: h_neglect ! A thickness that is so small it is usually lost
982986
! in roundoff and can be neglected [H ~> m].
983987
real :: I_Rho0 ! The inverse of the Boussinesq reference density [R-1 ~> m3 kg-1].
@@ -1029,8 +1033,20 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
10291033
dz_nonvanished = GV%H_to_Z*CS%h_nonvanished
10301034
I_Rho0 = 1.0 / GV%Rho0
10311035
G_Rho0 = GV%g_Earth / GV%Rho0
1032-
GxRho = GV%g_Earth * GV%Rho0
1033-
rho_ref = CS%Rho0
1036+
GxRho0 = GV%g_Earth * GV%Rho0
1037+
rho_ref = CS%rho_ref
1038+
1039+
if (CS%rho_ref_bug) then
1040+
rho0_int_density = rho_ref
1041+
rho0_set_pbce = rho_ref
1042+
GxRho_ref = GxRho0
1043+
I_g_rho = 1.0 / (rho_ref * GV%g_Earth)
1044+
else
1045+
rho0_int_density = GV%Rho0
1046+
rho0_set_pbce = GV%Rho0
1047+
GxRho_ref = GV%g_Earth * rho_ref
1048+
I_g_rho = 1.0 / (GV%rho0 * GV%g_Earth)
1049+
endif
10341050

10351051
if ((CS%id_MassWt_u > 0) .or. (CS%id_MassWt_v > 0)) then
10361052
MassWt_u(:,:,:) = 0.0 ; MassWt_v(:,:,:) = 0.0
@@ -1141,17 +1157,16 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
11411157
if (use_p_atm) then
11421158
!$OMP parallel do default(shared)
11431159
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
1144-
pa(i,j,1) = GxRho*(e(i,j,1) - G%Z_ref) + p_atm(i,j)
1160+
pa(i,j,1) = GxRho_ref * (e(i,j,1) - G%Z_ref) + p_atm(i,j)
11451161
enddo ; enddo
11461162
else
11471163
!$OMP parallel do default(shared)
11481164
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
1149-
pa(i,j,1) = GxRho*(e(i,j,1) - G%Z_ref)
1165+
pa(i,j,1) = GxRho_ref * (e(i,j,1) - G%Z_ref)
11501166
enddo ; enddo
11511167
endif
11521168

11531169
if (CS%use_SSH_in_Z0p .and. use_p_atm) then
1154-
I_g_rho = 1.0 / (CS%rho0*GV%g_Earth)
11551170
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
11561171
Z_0p(i,j) = e(i,j,1) + p_atm(i,j) * I_g_rho
11571172
enddo ; enddo
@@ -1177,23 +1192,23 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
11771192
if ( use_ALE .and. CS%Recon_Scheme > 0 ) then
11781193
if ( CS%Recon_Scheme == 1 ) then
11791194
call int_density_dz_generic_plm(k, tv, T_t, T_b, S_t, S_b, e, &
1180-
rho_ref, CS%Rho0, GV%g_Earth, dz_neglect, G%bathyT, &
1195+
rho_ref, rho0_int_density, GV%g_Earth, dz_neglect, G%bathyT, &
11811196
G%HI, GV, tv%eqn_of_state, US, CS%use_stanley_pgf, dpa(:,:,k), intz_dpa(:,:,k), &
11821197
intx_dpa(:,:,k), inty_dpa(:,:,k), &
11831198
MassWghtInterp=CS%MassWghtInterp, &
11841199
use_inaccurate_form=CS%use_inaccurate_pgf_rho_anom, Z_0p=Z_0p, &
11851200
MassWghtInterpVanOnly=CS%MassWghtInterpVanOnly, h_nv=dz_nonvanished)
11861201
elseif ( CS%Recon_Scheme == 2 ) then
11871202
call int_density_dz_generic_ppm(k, tv, T_t, T_b, S_t, S_b, e, &
1188-
rho_ref, CS%Rho0, GV%g_Earth, dz_neglect, G%bathyT, &
1203+
rho_ref, rho0_int_density, GV%g_Earth, dz_neglect, G%bathyT, &
11891204
G%HI, GV, tv%eqn_of_state, US, CS%use_stanley_pgf, dpa(:,:,k), intz_dpa(:,:,k), &
11901205
intx_dpa(:,:,k), inty_dpa(:,:,k), &
11911206
MassWghtInterp=CS%MassWghtInterp, Z_0p=Z_0p, &
11921207
MassWghtInterpVanOnly=CS%MassWghtInterpVanOnly, h_nv=dz_nonvanished)
11931208
endif
11941209
else
11951210
call int_density_dz(tv_tmp%T(:,:,k), tv_tmp%S(:,:,k), e(:,:,K), e(:,:,K+1), &
1196-
rho_ref, CS%Rho0, GV%g_Earth, G%HI, tv%eqn_of_state, US, dpa(:,:,k), &
1211+
rho_ref, rho0_int_density, GV%g_Earth, G%HI, tv%eqn_of_state, US, dpa(:,:,k), &
11971212
intz_dpa(:,:,k), intx_dpa(:,:,k), inty_dpa(:,:,k), G%bathyT, e(:,:,1), dz_neglect, &
11981213
CS%MassWghtInterp, Z_0p=Z_0p, &
11991214
MassWghtInterpVanOnly=CS%MassWghtInterpVanOnly, h_nv=dz_nonvanished)
@@ -1239,7 +1254,7 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
12391254
if (CS%sal_use_bpa) then
12401255
!$OMP parallel do default(shared)
12411256
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
1242-
pbot(i,j) = pa(i,j,nz+1) - GxRho * e(i,j,nz+1)
1257+
pbot(i,j) = pa(i,j,nz+1) - GxRho_ref * (e(i,j,nz+1) - G%Z_ref)
12431258
enddo ; enddo
12441259
call calc_SAL(pbot, e_sal, G, CS%SAL_CSp, tmp_scale=US%Z_to_m)
12451260
else
@@ -1253,7 +1268,7 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
12531268
!$OMP parallel do default(shared)
12541269
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
12551270
e(i,j,K) = e(i,j,K) - e_sal(i,j)
1256-
pa(i,j,K) = pa(i,j,K) - GxRho * e_sal(i,j)
1271+
pa(i,j,K) = pa(i,j,K) - GxRho_ref * e_sal(i,j)
12571272
enddo ; enddo
12581273
enddo ; endif
12591274
endif
@@ -1265,7 +1280,7 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
12651280
!$OMP parallel do default(shared)
12661281
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
12671282
e(i,j,K) = e(i,j,K) - (e_tidal_eq(i,j) + e_tidal_sal(i,j))
1268-
pa(i,j,K) = pa(i,j,K) - GxRho * (e_tidal_eq(i,j) + e_tidal_sal(i,j))
1283+
pa(i,j,K) = pa(i,j,K) - GxRho_ref * (e_tidal_eq(i,j) + e_tidal_sal(i,j))
12691284
enddo ; enddo
12701285
enddo ; endif
12711286
endif
@@ -1288,7 +1303,7 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
12881303
!$OMP parallel do default(shared) private(p_surf_EOS)
12891304
do j=Jsq,Jeq+1
12901305
! P_surf_EOS here is consistent with the pressure that is used in the int_density_dz routines.
1291-
do i=Isq,Ieq+1 ; p_surf_EOS(i) = -GxRho*(e(i,j,1) - Z_0p(i,j)) ; enddo
1306+
do i=Isq,Ieq+1 ; p_surf_EOS(i) = -GxRho0*(e(i,j,1) - Z_0p(i,j)) ; enddo
12921307
call calculate_density(T_top(:,j), S_top(:,j), p_surf_EOS, rho_top(:,j), &
12931308
tv%eqn_of_state, EOSdom, rho_ref=rho_ref)
12941309
enddo
@@ -1316,8 +1331,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
13161331
S5(1) = S_top(i,j) ; S5(5) = S_top(i+1,j)
13171332
pa5(1) = pa(i,j,1) ; pa5(5) = pa(i+1,j,1)
13181333
! Pressure input to density EOS is consistent with the pressure used in the int_density_dz routines.
1319-
p5(1) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1320-
p5(5) = -GxRho*(e(i+1,j,1) - Z_0p(i,j))
1334+
p5(1) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1335+
p5(5) = -GxRho0*(e(i+1,j,1) - Z_0p(i,j))
13211336
do m=2,4
13221337
wt_R = 0.25*real(m-1)
13231338
T5(m) = T5(1) + (T5(5)-T5(1))*wt_R
@@ -1351,8 +1366,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
13511366
S5(1) = S_top(i,j) ; S5(5) = S_top(i,j+1)
13521367
pa5(1) = pa(i,j,1) ; pa5(5) = pa(i,j+1,1)
13531368
! Pressure input to density EOS is consistent with the pressure used in the int_density_dz routines.
1354-
p5(1) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1355-
p5(5) = -GxRho*(e(i,j+1,1) - Z_0p(i,j))
1369+
p5(1) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1370+
p5(5) = -GxRho0*(e(i,j+1,1) - Z_0p(i,j))
13561371

13571372
do m=2,4
13581373
wt_R = 0.25*real(m-1)
@@ -1464,8 +1479,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
14641479
! This is a typical case in the open ocean, so use the topmost interface.
14651480
T_int_W(I,j) = T_top(i,j) ; T_int_E(I,j) = T_top(i+1,j)
14661481
S_int_W(I,j) = S_top(i,j) ; S_int_E(I,j) = S_top(i+1,j)
1467-
p_int_W(I,j) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1468-
p_int_E(I,j) = -GxRho*(e(i+1,j,1) - Z_0p(i,j))
1482+
p_int_W(I,j) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1483+
p_int_E(I,j) = -GxRho0*(e(i+1,j,1) - Z_0p(i,j))
14691484
intx_pa_nonlin(I,j) = intx_pa(I,j,1) - 0.5*(pa(i,j,1) + pa(i+1,j,1))
14701485
dgeo_x(I,j) = GV%g_Earth * (e(i+1,j,1)-e(i,j,1))
14711486
seek_x_cor(I,j) = .false.
@@ -1485,8 +1500,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
14851500
S_int_W(I,j) = S_b(i,j,k) ; S_int_E(I,j) = S_b(i+1,j,k)
14861501
! These pressures are only used for the equation of state, and are only a function of
14871502
! height, consistent with the expressions in the int_density_dz routines.
1488-
p_int_W(I,j) = -GxRho*(e(i,j,K+1) - Z_0p(i,j))
1489-
p_int_E(I,j) = -GxRho*(e(i+1,j,K+1) - Z_0p(i,j))
1503+
p_int_W(I,j) = -GxRho0*(e(i,j,K+1) - Z_0p(i,j))
1504+
p_int_E(I,j) = -GxRho0*(e(i+1,j,K+1) - Z_0p(i,j))
14901505

14911506
intx_pa_nonlin(I,j) = intx_pa(I,j,K+1) - 0.5*(pa(i,j,K+1) + pa(i+1,j,K+1))
14921507
dgeo_x(I,j) = GV%g_Earth * (e(i+1,j,K+1)-e(i,j,K+1))
@@ -1503,8 +1518,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
15031518
do j=js,je ; do I=Isq,Ieq ; if (seek_x_cor(I,j)) then
15041519
T_int_W(I,j) = T_top(i,j) ; T_int_E(I,j) = T_top(i+1,j)
15051520
S_int_W(I,j) = S_top(i,j) ; S_int_E(I,j) = S_top(i+1,j)
1506-
p_int_W(I,j) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1507-
p_int_E(I,j) = -GxRho*(e(i+1,j,1) - Z_0p(i,j))
1521+
p_int_W(I,j) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1522+
p_int_E(I,j) = -GxRho0*(e(i+1,j,1) - Z_0p(i,j))
15081523
intx_pa_nonlin(I,j) = intx_pa(I,j,1) - 0.5*(pa(i,j,1) + pa(i+1,j,1))
15091524
dgeo_x(I,j) = GV%g_Earth * (e(i+1,j,1)-e(i,j,1))
15101525
seek_x_cor(I,j) = .false.
@@ -1546,8 +1561,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
15461561
! This is a typical case in the open ocean, so use the topmost interface.
15471562
T_int_S(i,J) = T_top(i,j) ; T_int_N(i,J) = T_top(i,j+1)
15481563
S_int_S(i,J) = S_top(i,j) ; S_int_N(i,J) = S_top(i,j+1)
1549-
p_int_S(i,J) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1550-
p_int_N(i,J) = -GxRho*(e(i,j+1,1) - Z_0p(i,j))
1564+
p_int_S(i,J) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1565+
p_int_N(i,J) = -GxRho0*(e(i,j+1,1) - Z_0p(i,j))
15511566
inty_pa_nonlin(i,J) = inty_pa(i,J,1) - 0.5*(pa(i,j,1) + pa(i,j+1,1))
15521567
dgeo_y(i,J) = GV%g_Earth * (e(i,j+1,1)-e(i,j,1))
15531568
seek_y_cor(i,J) = .false.
@@ -1567,8 +1582,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
15671582
S_int_S(i,J) = S_b(i,j,k) ; S_int_N(i,J) = S_b(i,j+1,k)
15681583
! These pressures are only used for the equation of state, and are only a function of
15691584
! height, consistent with the expressions in the int_density_dz routines.
1570-
p_int_S(i,J) = -GxRho*(e(i,j,K+1) - Z_0p(i,j))
1571-
p_int_N(i,J) = -GxRho*(e(i,j+1,K+1) - Z_0p(i,j))
1585+
p_int_S(i,J) = -GxRho0*(e(i,j,K+1) - Z_0p(i,j))
1586+
p_int_N(i,J) = -GxRho0*(e(i,j+1,K+1) - Z_0p(i,j))
15721587
inty_pa_nonlin(i,J) = inty_pa(i,J,K+1) - 0.5*(pa(i,j,K+1) + pa(i,j+1,K+1))
15731588
dgeo_y(i,J) = GV%g_Earth * (e(i,j+1,K+1)-e(i,j,K+1))
15741589
seek_y_cor(i,J) = .false.
@@ -1584,8 +1599,8 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
15841599
do J=Jsq,Jeq ; do i=is,ie ; if (seek_y_cor(i,J)) then
15851600
T_int_S(i,J) = T_top(i,j) ; T_int_N(i,J) = T_top(i,j+1)
15861601
S_int_S(i,J) = S_top(i,j) ; S_int_N(i,J) = S_top(i,j+1)
1587-
p_int_S(i,J) = -GxRho*(e(i,j,1) - Z_0p(i,j))
1588-
p_int_N(i,J) = -GxRho*(e(i,j+1,1) - Z_0p(i,j))
1602+
p_int_S(i,J) = -GxRho0*(e(i,j,1) - Z_0p(i,j))
1603+
p_int_N(i,J) = -GxRho0*(e(i,j+1,1) - Z_0p(i,j))
15891604
inty_pa_nonlin(i,J) = inty_pa(i,J,1) - 0.5*(pa(i,j,1) + pa(i,j+1,1))
15901605
dgeo_y(i,J) = GV%g_Earth * (e(i,j+1,1)-e(i,j,1))
15911606
seek_y_cor(i,J) = .false.
@@ -1709,7 +1724,7 @@ subroutine PressureForce_FV_Bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm
17091724
endif
17101725

17111726
if (present(pbce)) then
1712-
call set_pbce_Bouss(e, tv_tmp, G, GV, US, CS%Rho0, CS%GFS_scale, pbce)
1727+
call set_pbce_Bouss(e, tv_tmp, G, GV, US, rho0_set_pbce, CS%GFS_scale, pbce)
17131728
endif
17141729

17151730
if (present(eta)) then
@@ -1836,11 +1851,16 @@ subroutine PressureForce_FV_init(Time, G, GV, US, param_file, diag, CS, SAL_CSp,
18361851
call get_param(param_file, mdl, "DEBUG", CS%debug, &
18371852
"If true, write out verbose debugging data.", &
18381853
default=.false., debuggingParam=.true., do_not_log=.true.)
1839-
call get_param(param_file, mdl, "RHO_PGF_REF", CS%Rho0, &
1854+
call get_param(param_file, mdl, "RHO_PGF_REF", CS%rho_ref, &
18401855
"The reference density that is subtracted off when calculating pressure "//&
18411856
"gradient forces. Its inverse is subtracted off of specific volumes when "//&
18421857
"in non-Boussinesq mode. The default is RHO_0.", &
18431858
units="kg m-3", default=GV%Rho0*US%R_to_kg_m3, scale=US%kg_m3_to_R)
1859+
call get_param(param_file, mdl, "RHO_PGF_REF_BUG", CS%rho_ref_bug, &
1860+
"If true, recover a bug that RHO_0 (the mean seawater density in Boussinesq mode) "//&
1861+
"and RHO_PGF_REF (the subtracted reference density in finite volume pressure "//&
1862+
"gradient forces) are incorrectly interchanged in several instances in Boussinesq mode.", &
1863+
default=.true.)
18441864
call get_param(param_file, mdl, "TIDES", CS%tides, &
18451865
"If true, apply tidal momentum forcing.", default=.false.)
18461866
call get_param(param_file, '', "DEFAULT_ANSWER_DATE", default_answer_date, default=99991231)

src/core/MOM_density_integrals.F90

Lines changed: 0 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -170,7 +170,6 @@ subroutine int_density_dz_generic_pcm(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, &
170170
real :: rho_anom ! The depth averaged density anomaly [R ~> kg m-3]
171171
real, parameter :: C1_90 = 1.0/90.0 ! A rational constant [nondim]
172172
real :: GxRho ! The product of the gravitational acceleration and reference density [R L2 Z-1 T-2 ~> Pa m-1]
173-
real :: I_Rho ! The inverse of the Boussinesq density [R-1 ~> m3 kg-1]
174173
real :: dz ! The layer thickness [Z ~> m]
175174
real :: dz_x(5,HI%iscB:HI%iecB) ! Layer thicknesses along an x-line of subgrid locations [Z ~> m]
176175
real :: dz_y(5,HI%isc:HI%iec) ! Layer thicknesses along a y-line of subgrid locations [Z ~> m]
@@ -204,7 +203,6 @@ subroutine int_density_dz_generic_pcm(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, &
204203
js = HI%jsc ; je = HI%jec
205204

206205
GxRho = G_e * rho_0
207-
I_Rho = 1.0 / rho_0
208206
if (present(Z_0p)) then
209207
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
210208
z0pres(i,j) = Z_0p(i,j)
@@ -511,7 +509,6 @@ subroutine int_density_dz_generic_plm(k, tv, T_t, T_b, S_t, S_b, e, rho_ref, &
511509
! with height at the 5 sub-column locations [R L2 T-2 ~> Pa]
512510
real, parameter :: C1_90 = 1.0/90.0 ! A rational constant [nondim]
513511
real :: GxRho ! The product of the gravitational acceleration and reference density [R L2 Z-1 T-2 ~> Pa m-1]
514-
real :: I_Rho ! The inverse of the Boussinesq density [R-1 ~> m3 kg-1]
515512
real :: dz(HI%iscB:HI%iecB+1) ! Layer thicknesses at tracer points [Z ~> m]
516513
real :: dz_x(5,HI%iscB:HI%iecB) ! Layer thicknesses along an x-line of subgrid locations [Z ~> m]
517514
real :: dz_y(5,HI%isc:HI%iec) ! Layer thicknesses along a y-line of subgrid locations [Z ~> m]
@@ -538,7 +535,6 @@ subroutine int_density_dz_generic_plm(k, tv, T_t, T_b, S_t, S_b, e, rho_ref, &
538535
Isq = HI%IscB ; Ieq = HI%IecB ; Jsq = HI%JscB ; Jeq = HI%JecB
539536

540537
GxRho = G_e * rho_0
541-
I_Rho = 1.0 / rho_0
542538
if (present(Z_0p)) then
543539
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
544540
z0pres(i,j) = Z_0p(i,j)
@@ -966,7 +962,6 @@ subroutine int_density_dz_generic_ppm(k, tv, T_t, T_b, S_t, S_b, e, &
966962
! with height at the 5 sub-column locations [R L2 T-2 ~> Pa]
967963
real, parameter :: C1_90 = 1.0/90.0 ! A rational constant [nondim]
968964
real :: GxRho ! The gravitational acceleration times density [R L2 Z-1 T-2 ~> kg m-2 s-2]
969-
real :: I_Rho ! The inverse of the Boussinesq density [R-1 ~> m3 kg-1]
970965
real :: dz ! Layer thicknesses at tracer points [Z ~> m]
971966
real :: dz_x(5,HI%iscB:HI%iecB) ! Layer thicknesses along an x-line of subgrid locations [Z ~> m]
972967
real :: dz_y(5,HI%isc:HI%iec) ! Layer thicknesses along a y-line of subgrid locations [Z ~> m]
@@ -996,7 +991,6 @@ subroutine int_density_dz_generic_ppm(k, tv, T_t, T_b, S_t, S_b, e, &
996991
Isq = HI%IscB ; Ieq = HI%IecB ; Jsq = HI%JscB ; Jeq = HI%JecB
997992

998993
GxRho = G_e * rho_0
999-
I_Rho = 1.0 / rho_0
1000994
if (present(Z_0p)) then
1001995
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
1002996
z0pres(i,j) = Z_0p(i,j)

0 commit comments

Comments
 (0)