mom_density_integrals module reference
Provides integrals of density.
Functions/Subroutines
Calls the appropriate subroutine to calculate analytical and nearly-analytical integrals in z across layers of pressure anomalies, which are required for calculating the finite-volume form pressure accelerations in a Boussinesq model. |
|
Calculates (by numerical quadrature) integrals of pressure anomalies across layers, which are required for calculating the finite-volume form pressure accelerations in a Boussinesq model. |
|
Compute pressure gradient force integrals by quadrature for the case where T and S are linear profiles. |
|
Compute pressure gradient force integrals for layer "k" and the case where T and S are parabolic profiles. |
|
Calls the appropriate subroutine to calculate analytical and nearly-analytical integrals in pressure across layers of geopotential anomalies, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. |
|
This subroutine calculates integrals of specific volume anomalies in pressure across layers, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. |
|
This subroutine calculates integrals of specific volume anomalies in pressure across layers, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. |
|
Diagnose the fractional mass weighting in a layer that might be used with a Boussinesq calculation. |
|
Diagnose the fractional mass weighting in a layer that might be used with a non-Boussinesq calculation. |
|
Find the depth at which the reconstructed pressure matches P_tgt. |
|
Calculate the average in situ specific volume across layers. |
|
Returns change in anomalous pressure change from top to non-dimensional position pos between z_t and z_b [R L2 T-2 ~> Pa]. |
Detailed Description
Provides integrals of density.
Function/Subroutine Documentation
- subroutine mom_density_integrals/int_density_dz(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, EOS, US, dpa, intz_dpa, intx_dpa, inty_dpa, bathyT, SSH, dz_neglect, MassWghtInterp, Z_0p, MassWghtInterpVanOnly, h_nv)
Calls the appropriate subroutine to calculate analytical and nearly-analytical integrals in z across layers of pressure anomalies, which are required for calculating the finite-volume form pressure accelerations in a Boussinesq model.
- Parameters:
hi :: [in] Ocean horizontal index structures for the arrays
t :: [in] Potential temperature referenced to the surface [C ~> degC]
s :: [in] Salinity [S ~> ppt]
z_t ::
z_t[in] Height at the top of the layer in depth units [Z ~> m]z_b ::
z_b[in] Height at the bottom of the layer [Z ~> m]rho_ref ::
rho_ref[in] A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals.rho_0 ::
rho_0[in] A density [R ~> kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state.g_e :: [in] The Earth’s gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
dpa ::
dpa[inout] The change in the pressure anomalyintz_dpa ::
intz_dpa[inout] The integral through the thickness of theintx_dpa ::
intx_dpa[inout] The integral in x of the difference betweeninty_dpa ::
inty_dpa[inout] The integral in y of the difference betweenbathyt :: [in] The depth of the bathymetry [Z ~> m]
ssh :: [in] The sea surface height [Z ~> m]
dz_neglect ::
dz_neglect[in] A minuscule thickness change [Z ~> m]masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
h_nv ::
h_nv[in] Nonvanished height [Z ~> m]z_0p :: [in] The height at which the pressure is 0 [Z ~> m]
- Call to:
- Called from:
- subroutine mom_density_integrals/int_density_dz_generic_pcm(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, EOS, US, dpa, intz_dpa, intx_dpa, inty_dpa, bathyT, SSH, dz_neglect, MassWghtInterp, use_inaccurate_form, Z_0p, MassWghtInterpVanOnly, h_nv)
Calculates (by numerical quadrature) integrals of pressure anomalies across layers, which are required for calculating the finite-volume form pressure accelerations in a Boussinesq model.
- Parameters:
hi :: [in] Horizontal index type for input variables.
t :: [in] Potential temperature of the layer [C ~> degC]
s :: [in] Salinity of the layer [S ~> ppt]
z_t ::
z_t[in] Height at the top of the layer in depth units [Z ~> m]z_b ::
z_b[in] Height at the bottom of the layer [Z ~> m]rho_ref ::
rho_ref[in] A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals.rho_0 ::
rho_0[in] A density [R ~> kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state.g_e :: [in] The Earth’s gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
dpa ::
dpa[inout] The change in the pressure anomalyintz_dpa ::
intz_dpa[inout] The integral through the thickness of theintx_dpa ::
intx_dpa[inout] The integral in x of the difference betweeninty_dpa ::
inty_dpa[inout] The integral in y of the difference betweenbathyt :: [in] The depth of the bathymetry [Z ~> m]
ssh :: [in] The sea surface height [Z ~> m]
dz_neglect ::
dz_neglect[in] A minuscule thickness change [Z ~> m]masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
h_nv ::
h_nv[in] Nonvanished height [Z ~> m]use_inaccurate_form ::
use_inaccurate_form[in] If true, uses an inaccurate form of density anomalies, as was used prior to March 2018.z_0p :: [in] The height at which the pressure is 0 [Z ~> m]
- Call to:
- Called from:
- subroutine mom_density_integrals/int_density_dz_generic_plm(k, tv, T_t, T_b, S_t, S_b, e, rho_ref, rho_0, G_e, dz_subroundoff, bathyT, HI, GV, EOS, US, use_stanley_eos, dpa, intz_dpa, intx_dpa, inty_dpa, MassWghtInterp, use_inaccurate_form, Z_0p, MassWghtInterpVanOnly, h_nv)
Compute pressure gradient force integrals by quadrature for the case where T and S are linear profiles.
- Parameters:
k ::
k[in] Layer index to calculate integrals forhi :: [in] Ocean horizontal index structures for the input arrays
gv :: [in] Vertical grid structure
tv ::
tv[in] Thermodynamic variablest_t :: [in] Potential temperature at the cell top [C ~> degC]
t_b :: [in] Potential temperature at the cell bottom [C ~> degC]
s_t :: [in] Salinity at the cell top [S ~> ppt]
s_b :: [in] Salinity at the cell bottom [S ~> ppt]
e ::
e[in] Height of interfaces [Z ~> m]rho_ref ::
rho_ref[in] A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals.rho_0 ::
rho_0[in] A density [R ~> kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state.g_e :: [in] The Earth’s gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
dz_subroundoff ::
dz_subroundoff[in] A minuscule thickness change [Z ~> m]bathyt :: [in] The depth of the bathymetry [Z ~> m]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
use_stanley_eos ::
use_stanley_eos[in] If true, turn on Stanley SGS T variance parameterizationdpa ::
dpa[inout] The change in the pressure anomaly across the layer [R L2 T-2 ~> Pa]intz_dpa ::
intz_dpa[inout] The integral through the thickness of the layer ofintx_dpa ::
intx_dpa[inout] The integral in x of the difference between theinty_dpa ::
inty_dpa[inout] The integral in y of the difference between themasswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
use_inaccurate_form ::
use_inaccurate_form[in] If true, uses an inaccurate form of density anomalies, as was used prior to March 2018.masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
h_nv ::
h_nv[in] Nonvanished height [Z ~> m]z_0p :: [in] The height at which the pressure is 0 [Z ~> m]
- subroutine mom_density_integrals/int_density_dz_generic_ppm(k, tv, T_t, T_b, S_t, S_b, e, rho_ref, rho_0, G_e, dz_subroundoff, bathyT, HI, GV, EOS, US, use_stanley_eos, dpa, intz_dpa, intx_dpa, inty_dpa, MassWghtInterp, Z_0p, MassWghtInterpVanOnly, h_nv)
Compute pressure gradient force integrals for layer “k” and the case where T and S are parabolic profiles.
- Parameters:
k ::
k[in] Layer index to calculate integrals forhi :: [in] Ocean horizontal index structures for the input arrays
gv :: [in] Vertical grid structure
tv ::
tv[in] Thermodynamic variablest_t :: [in] Potential temperature at the cell top [C ~> degC]
t_b :: [in] Potential temperature at the cell bottom [C ~> degC]
s_t :: [in] Salinity at the cell top [S ~> ppt]
s_b :: [in] Salinity at the cell bottom [S ~> ppt]
e ::
e[in] Height of interfaces [Z ~> m]rho_ref ::
rho_ref[in] A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals.rho_0 ::
rho_0[in] A density [R ~> kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state.g_e :: [in] The Earth’s gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
dz_subroundoff ::
dz_subroundoff[in] A minuscule thickness change [Z ~> m]bathyt :: [in] The depth of the bathymetry [Z ~> m]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
use_stanley_eos ::
use_stanley_eos[in] If true, turn on Stanley SGS T variance parameterizationdpa ::
dpa[inout] The change in the pressure anomaly across the layer [R L2 T-2 ~> Pa]intz_dpa ::
intz_dpa[inout] The integral through the thickness of the layer ofintx_dpa ::
intx_dpa[inout] The integral in x of the difference between theinty_dpa ::
inty_dpa[inout] The integral in y of the difference between themasswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
h_nv ::
h_nv[in] Nonvanished height [Z ~> m]z_0p :: [in] The height at which the pressure is 0 [Z ~> m]
- subroutine mom_density_integrals/int_specific_vol_dp(T, S, p_t, p_b, alpha_ref, HI, EOS, US, dza, intp_dza, intx_dza, inty_dza, halo_size, bathyP, P_surf, dP_tiny, MassWghtInterp, MassWghtInterpVanOnly, p_nv)
Calls the appropriate subroutine to calculate analytical and nearly-analytical integrals in pressure across layers of geopotential anomalies, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. There are essentially no free assumptions, apart from the use of Boole’s rule to do the horizontal integrals, and from a truncation in the series for log(1-eps/1+eps) that assumes that |eps| < 0.34.
- Parameters:
hi :: [in] The horizontal index structure
t :: [in] Potential temperature referenced to the surface [C ~> degC]
s :: [in] Salinity [S ~> ppt]
p_t ::
p_t[in] Pressure at the top of the layer [R L2 T-2 ~> Pa]p_b ::
p_b[in] Pressure at the bottom of the layer [R L2 T-2 ~> Pa]alpha_ref ::
alpha_ref[in] A mean specific volume that is subtracted out to reduce the magnitude of each of the integrals [R-1 ~> m3 kg-1] The calculation is mathematically identical with different values of alpha_ref, but this reduces the effects of roundoff.eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
dza ::
dza[inout] The change in the geopotential anomaly acrossintp_dza ::
intp_dza[inout] The integral in pressure through the layer of theintx_dza ::
intx_dza[inout] The integral in x of the difference between theinty_dza ::
inty_dza[inout] The integral in y of the difference between thehalo_size ::
halo_size[in] The width of halo points on which to calculate dza.bathyp :: [in] The pressure at the bathymetry [R L2 T-2 ~> Pa]
p_surf :: [in] The pressure at the ocean surface [R L2 T-2 ~> Pa]
dp_tiny :: [in] A minuscule pressure change with the same units as p_t [R L2 T-2 ~> Pa]
masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
p_nv ::
p_nv[in] Nonvanished pressure [R L2 T-2 ~> Pa]
- Call to:
- Called from:
mom_interface_heights::dz_to_thickness::dz_to_thickness_eosmom_interface_heights::find_dz_for_etamom_pressureforce_mont::pressureforce_mont_nonbouss
- subroutine mom_density_integrals/int_spec_vol_dp_generic_pcm(T, S, p_t, p_b, alpha_ref, HI, EOS, US, dza, intp_dza, intx_dza, inty_dza, halo_size, bathyP, P_surf, dP_neglect, MassWghtInterp, MassWghtInterpVanOnly, p_nv)
This subroutine calculates integrals of specific volume anomalies in pressure across layers, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. There are essentially no free assumptions, apart from the use of Boole’s rule quadrature to do the integrals.
- Parameters:
hi :: [in] A horizontal index type structure.
t :: [in] Potential temperature of the layer [C ~> degC]
s :: [in] Salinity of the layer [S ~> ppt]
p_t ::
p_t[in] Pressure atop the layer [R L2 T-2 ~> Pa]p_b ::
p_b[in] Pressure below the layer [R L2 T-2 ~> Pa]alpha_ref ::
alpha_ref[in] A mean specific volume that is subtracted out to reduce the magnitude of each of the integrals [R-1 ~> m3 kg-1] The calculation is mathematically identical with different values of alpha_ref, but alpha_ref alters the effects of roundoff, and answers do change.eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
dza ::
dza[inout] The change in the geopotential anomalyintp_dza ::
intp_dza[inout] The integral in pressure through the layer ofintx_dza ::
intx_dza[inout] The integral in x of the difference betweeninty_dza ::
inty_dza[inout] The integral in y of the difference betweenhalo_size ::
halo_size[in] The width of halo points on which to calculate dza.bathyp :: [in] The pressure at the bathymetry [R L2 T-2 ~> Pa]
p_surf :: [in] The pressure at the ocean surface [R L2 T-2 ~> Pa]
dp_neglect :: [in] A minuscule pressure change with the same units as p_t [R L2 T-2 ~> Pa]
masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
p_nv ::
p_nv[in] Nonvanished pressure [R L2 T-2 ~> Pa]
- Call to:
- Called from:
- subroutine mom_density_integrals/int_spec_vol_dp_generic_plm(T_t, T_b, S_t, S_b, p_t, p_b, alpha_ref, dP_neglect, bathyP, HI, EOS, US, dza, intp_dza, intx_dza, inty_dza, P_surf, MassWghtInterp, MassWghtInterpVanOnly, p_nv)
This subroutine calculates integrals of specific volume anomalies in pressure across layers, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. There are essentially no free assumptions, apart from the use of Boole’s rule quadrature to do the integrals.
- Parameters:
hi :: [in] A horizontal index type structure.
t_t :: [in] Potential temperature at the top of the layer [C ~> degC]
t_b :: [in] Potential temperature at the bottom of the layer [C ~> degC]
s_t :: [in] Salinity at the top the layer [S ~> ppt]
s_b :: [in] Salinity at the bottom the layer [S ~> ppt]
p_t ::
p_t[in] Pressure atop the layer [R L2 T-2 ~> Pa]p_b ::
p_b[in] Pressure below the layer [R L2 T-2 ~> Pa]alpha_ref ::
alpha_ref[in] A mean specific volume that is subtracted out to reduce the magnitude of each of the integrals [R-1 ~> m3 kg-1] The calculation is mathematically identical with different values of alpha_ref, but alpha_ref alters the effects of roundoff, and answers do change.dp_neglect :: [in] !< A miniscule pressure change with the same units as p_t [R L2 T-2 ~> Pa]
bathyp :: [in] The pressure at the bathymetry [R L2 T-2 ~> Pa]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
dza ::
dza[inout] The change in the geopotential anomalyintp_dza ::
intp_dza[inout] The integral in pressure through the layer ofintx_dza ::
intx_dza[inout] The integral in x of the difference betweeninty_dza ::
inty_dza[inout] The integral in y of the difference betweenp_surf :: [in] The pressure at the ocean surface [R L2 T-2 ~> Pa]
masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
p_nv ::
p_nv[in] Nonvanished pressure [R L2 T-2 ~> Pa]
- Call to:
- subroutine mom_density_integrals/diagnose_mass_weight_z(z_t, z_b, bathyT, SSH, dz_neglect, MassWghtInterp, HI, MassWt_u, MassWt_v, MassWghtInterpVanOnly, h_nv)
Diagnose the fractional mass weighting in a layer that might be used with a Boussinesq calculation.
- Parameters:
hi :: [in] A horizontal index type structure.
z_t ::
z_t[in] Height at the top of the layer in depth units [Z ~> m]z_b ::
z_b[in] Height at the bottom of the layer [Z ~> m]bathyt :: [in] The depth of the bathymetry [Z ~> m]
ssh :: [in] The sea surface height [Z ~> m]
dz_neglect ::
dz_neglect[in] A minuscule thickness change [Z ~> m]masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswt_u :: [inout] The fractional mass weighting at u-points [nondim]
masswt_v :: [inout] The fractional mass weighting at v-points [nondim]
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
h_nv ::
h_nv[in] Nonvanished height [Z ~> m]
- Called from:
- subroutine mom_density_integrals/diagnose_mass_weight_p(p_t, p_b, bathyP, P_surf, dP_neglect, MassWghtInterp, HI, MassWt_u, MassWt_v, MassWghtInterpVanOnly, p_nv)
Diagnose the fractional mass weighting in a layer that might be used with a non-Boussinesq calculation.
- Parameters:
hi :: [in] A horizontal index type structure.
p_t ::
p_t[in] Pressure atop the layer [R L2 T-2 ~> Pa]p_b ::
p_b[in] Pressure below the layer [R L2 T-2 ~> Pa]dp_neglect :: [in] !< A miniscule pressure change with the same units as p_t [R L2 T-2 ~> Pa]
bathyp :: [in] The pressure at the bathymetry [R L2 T-2 ~> Pa]
p_surf :: [in] The pressure at the ocean surface [R L2 T-2 ~> Pa]
masswghtinterp :: [in] A flag indicating whether and how to use mass weighting to interpolate T/S in integrals
masswt_u :: [inout] The fractional mass weighting at u-points [nondim]
masswt_v :: [inout] The fractional mass weighting at v-points [nondim]
masswghtinterpvanonly :: [in] If true, does not do mass weighting of T/S unless one side smaller than h_nv (i.e. vanished)
p_nv ::
p_nv[in] Nonvanished pressure [R L2 T-2 ~> Pa]
- Called from:
- subroutine mom_density_integrals/find_depth_of_pressure_in_cell(T_t, T_b, S_t, S_b, z_t, z_b, P_t, P_tgt, rho_ref, G_e, EOS, US, P_b, z_out, z_tol, frac_dp_bugfix)
Find the depth at which the reconstructed pressure matches P_tgt.
- Parameters:
t_t :: [in] Potential temperature at the cell top [C ~> degC]
t_b :: [in] Potential temperature at the cell bottom [C ~> degC]
s_t :: [in] Salinity at the cell top [S ~> ppt]
s_b :: [in] Salinity at the cell bottom [S ~> ppt]
z_t ::
z_t[in] Absolute height of top of cell [Z ~> m] (Boussinesq ????)z_b ::
z_b[in] Absolute height of bottom of cell [Z ~> m]p_t :: [in] Anomalous pressure of top of cell, relative to g*rho_ref*z_t [R L2 T-2 ~> Pa]
p_tgt :: [in] Target pressure at height z_out, relative to g*rho_ref*z_out [R L2 T-2 ~> Pa]
rho_ref ::
rho_ref[in] Reference density with which calculation are anomalous to [R ~> kg m-3]g_e :: [in] Gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
eos :: [in] Equation of state structure
us :: [in] A dimensional unit scaling type
p_b :: [out] Pressure at the bottom of the cell [R L2 T-2 ~> Pa]
z_out ::
z_out[out] Absolute depth at which anomalous pressure = p_tgt [Z ~> m]z_tol ::
z_tol[in] The tolerance in finding z_out [Z ~> m]frac_dp_bugfix ::
frac_dp_bugfix[in] If true, use bugfix in frac_dp_at_pos
- Call to:
- Called from:
- subroutine mom_density_integrals/avg_specific_vol(T, S, p_t, dp, HI, EOS, SpV_avg, halo_size)
Calculate the average in situ specific volume across layers.
- Parameters:
hi :: [in] The horizontal index structure
t :: [in] Potential temperature of the layer [C ~> degC]
s :: [in] Salinity of the layer [S ~> ppt]
p_t ::
p_t[in] Pressure at the top of the layer [R L2 T-2 ~> Pa]dp ::
dp[in] Pressure change in the layer [R L2 T-2 ~> Pa]eos :: [in] Equation of state structure
spv_avg :: [inout] The vertical average specific volume
halo_size ::
halo_size[in] The number of halo points in which to work.
- Call to:
- Called from:
- function mom_density_integrals/frac_dp_at_pos(T_t, T_b, S_t, S_b, z_t, z_b, rho_ref, G_e, pos, EOS, frac_dp_bugfix)
Returns change in anomalous pressure change from top to non-dimensional position pos between z_t and z_b [R L2 T-2 ~> Pa].
- Parameters:
t_t :: [in] Potential temperature at the cell top [C ~> degC]
t_b :: [in] Potential temperature at the cell bottom [C ~> degC]
s_t :: [in] Salinity at the cell top [S ~> ppt]
s_b :: [in] Salinity at the cell bottom [S ~> ppt]
z_t ::
z_t[in] The geometric height at the top of the layer [Z ~> m]z_b ::
z_b[in] The geometric height at the bottom of the layer [Z ~> m]rho_ref ::
rho_ref[in] A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals.g_e :: [in] The Earth’s gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
pos ::
pos[in] The fractional vertical position, 0 to 1 [nondim]eos :: [in] Equation of state structure
frac_dp_bugfix ::
frac_dp_bugfix[in] If true, use bugfix in frac_dp_at_pos
- Called from: