|
MOM6
|
Provides integrals of density.
Functions/Subroutines | |
| subroutine, public | 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, dz_neglect, useMassWghtInterp) |
| 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. More... | |
| subroutine, public | 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, dz_neglect, useMassWghtInterp) |
| 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. More... | |
| subroutine, public | 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, dpa, intz_dpa, intx_dpa, inty_dpa, useMassWghtInterp) |
| Compute pressure gradient force integrals by quadrature for the case where T and S are linear profiles. More... | |
| subroutine, public | 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, dpa, intz_dpa, intx_dpa, inty_dpa, useMassWghtInterp) |
| Compute pressure gradient force integrals for layer "k" and the case where T and S are parabolic profiles. More... | |
| subroutine, public | 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, dP_tiny, useMassWghtInterp) |
| 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. More... | |
| subroutine, public | 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, dP_neglect, useMassWghtInterp) |
| 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. More... | |
| subroutine, public | 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, useMassWghtInterp) |
| 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. More... | |
| subroutine, public | 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) |
| Find the depth at which the reconstructed pressure matches P_tgt. More... | |
| real function | frac_dp_at_pos (T_t, T_b, S_t, S_b, z_t, z_b, rho_ref, G_e, pos, EOS) |
| Returns change in anomalous pressure change from top to non-dimensional position pos between z_t and z_b. More... | |
| subroutine, public mom_density_integrals::find_depth_of_pressure_in_cell | ( | real, intent(in) | T_t, |
| real, intent(in) | T_b, | ||
| real, intent(in) | S_t, | ||
| real, intent(in) | S_b, | ||
| real, intent(in) | z_t, | ||
| real, intent(in) | z_b, | ||
| real, intent(in) | P_t, | ||
| real, intent(in) | P_tgt, | ||
| real, intent(in) | rho_ref, | ||
| real, intent(in) | G_e, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, intent(out) | P_b, | ||
| real, intent(out) | z_out, | ||
| real, intent(in), optional | z_tol | ||
| ) |
Find the depth at which the reconstructed pressure matches P_tgt.
| [in] | t_t | Potential temperature at the cell top [degC] |
| [in] | t_b | Potential temperature at the cell bottom [degC] |
| [in] | s_t | Salinity at the cell top [ppt] |
| [in] | s_b | Salinity at the cell bottom [ppt] |
| [in] | z_t | Absolute height of top of cell [Z ~> m] (Boussinesq ????) |
| [in] | z_b | Absolute height of bottom of cell [Z ~> m] |
| [in] | p_t | Anomalous pressure of top of cell, relative to g*rho_ref*z_t [R L2 T-2 ~> Pa] |
| [in] | p_tgt | Target pressure at height z_out, relative to g*rho_ref*z_out [R L2 T-2 ~> Pa] |
| [in] | rho_ref | Reference density with which calculation are anomalous to [R ~> kg m-3] |
| [in] | g_e | Gravitational acceleration [L2 Z-1 T-2 ~> m s-2] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [out] | p_b | Pressure at the bottom of the cell [R L2 T-2 ~> Pa] |
| [out] | z_out | Absolute depth at which anomalous pressure = p_tgt [Z ~> m] |
| [in] | z_tol | The tolerance in finding z_out [Z ~> m] |
Definition at line 1525 of file MOM_density_integrals.F90.
|
private |
Returns change in anomalous pressure change from top to non-dimensional position pos between z_t and z_b.
| [in] | t_t | Potential temperature at the cell top [degC] |
| [in] | t_b | Potential temperature at the cell bottom [degC] |
| [in] | s_t | Salinity at the cell top [ppt] |
| [in] | s_b | Salinity at the cell bottom [ppt] |
| [in] | z_t | The geometric height at the top of the layer [Z ~> m] |
| [in] | z_b | The geometric height at the bottom of the layer [Z ~> m] |
| [in] | rho_ref | A mean density [R ~> kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. |
| [in] | g_e | The Earth's gravitational acceleration [L2 Z-1 T-2 ~> m s-2] |
| [in] | pos | The fractional vertical position, 0 to 1 [nondim] |
| eos | Equation of state structure |
Definition at line 1608 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_density_dz | ( | real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | T, |
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | S, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | z_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | z_b, | ||
| real, intent(in) | rho_ref, | ||
| real, intent(in) | rho_0, | ||
| real, intent(in) | G_e, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout) | dpa, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout), optional | intz_dpa, | ||
| real, dimension(szib_(hi),szj_(hi)), intent(inout), optional | intx_dpa, | ||
| real, dimension(szi_(hi),szjb_(hi)), intent(inout), optional | inty_dpa, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in), optional | bathyT, | ||
| real, intent(in), optional | dz_neglect, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
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.
| [in] | hi | Ocean horizontal index structures for the arrays |
| [in] | t | Potential temperature referenced to the surface [degC] |
| [in] | s | Salinity [ppt] |
| [in] | z_t | Height at the top of the layer in depth units [Z ~> m] |
| [in] | z_b | Height at the bottom of the layer [Z ~> m] |
| [in] | rho_ref | A mean density [R ~> kg m-3] or [kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. |
| [in] | rho_0 | A density [R ~> kg m-3] or [kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state. |
| [in] | g_e | The Earth's gravitational acceleration [L2 Z-1 T-2 ~> m s-2] or [m2 Z-1 s-2 ~> m s-2] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dpa | The change in the pressure anomaly |
| [in,out] | intz_dpa | The integral through the thickness of the |
| [in,out] | intx_dpa | The integral in x of the difference between |
| [in,out] | inty_dpa | The integral in y of the difference between |
| [in] | bathyt | The depth of the bathymetry [Z ~> m] |
| [in] | dz_neglect | A minuscule thickness change [Z ~> m] |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 40 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_density_dz_generic_pcm | ( | real, dimension(szi_(hi),szj_(hi)), intent(in) | T, |
| real, dimension(szi_(hi),szj_(hi)), intent(in) | S, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in) | z_t, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in) | z_b, | ||
| real, intent(in) | rho_ref, | ||
| real, intent(in) | rho_0, | ||
| real, intent(in) | G_e, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout) | dpa, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout), optional | intz_dpa, | ||
| real, dimension( hi %isdb: hi %iedb, hi %jsd: hi %jed), intent(inout), optional | intx_dpa, | ||
| real, dimension( hi %isd: hi %ied, hi %jsdb: hi %jedb), intent(inout), optional | inty_dpa, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in), optional | bathyT, | ||
| real, intent(in), optional | dz_neglect, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
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.
| [in] | hi | Horizontal index type for input variables. |
| [in] | t | Potential temperature of the layer [degC] |
| [in] | s | Salinity of the layer [ppt] |
| [in] | z_t | Height at the top of the layer in depth units [Z ~> m] |
| [in] | z_b | Height at the bottom of the layer [Z ~> m] |
| [in] | rho_ref | A mean density [R ~> kg m-3] or [kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. |
| [in] | rho_0 | A density [R ~> kg m-3] or [kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state. |
| [in] | g_e | The Earth's gravitational acceleration [L2 Z-1 T-2 ~> m s-2] or [m2 Z-1 s-2 ~> m s-2] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dpa | The change in the pressure anomaly |
| [in,out] | intz_dpa | The integral through the thickness of the |
| [in,out] | intx_dpa | The integral in x of the difference between |
| [in,out] | inty_dpa | The integral in y of the difference between |
| [in] | bathyt | The depth of the bathymetry [Z ~> m] |
| [in] | dz_neglect | A minuscule thickness change [Z ~> m] |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 95 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_density_dz_generic_plm | ( | integer, intent(in) | k, |
| type(thermo_var_ptrs), intent(in) | tv, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed, gv %ke), intent(in) | T_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed, gv %ke), intent(in) | T_b, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed, gv %ke), intent(in) | S_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed, gv %ke), intent(in) | S_b, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed, gv %ke+1), intent(in) | e, | ||
| real, intent(in) | rho_ref, | ||
| real, intent(in) | rho_0, | ||
| real, intent(in) | G_e, | ||
| real, intent(in) | dz_subroundoff, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in) | bathyT, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(verticalgrid_type), intent(in) | GV, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout) | dpa, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout), optional | intz_dpa, | ||
| real, dimension(szib_(hi),szj_(hi)), intent(inout), optional | intx_dpa, | ||
| real, dimension(szi_(hi),szjb_(hi)), intent(inout), optional | inty_dpa, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
Compute pressure gradient force integrals by quadrature for the case where T and S are linear profiles.
| [in] | k | Layer index to calculate integrals for |
| [in] | hi | Ocean horizontal index structures for the input arrays |
| [in] | gv | Vertical grid structure |
| [in] | tv | Thermodynamic variables |
| [in] | t_t | Potential temperature at the cell top [degC] |
| [in] | t_b | Potential temperature at the cell bottom [degC] |
| [in] | s_t | Salinity at the cell top [ppt] |
| [in] | s_b | Salinity at the cell bottom [ppt] |
| [in] | e | Height of interfaces [Z ~> m] |
| [in] | rho_ref | A mean density [R ~> kg m-3] or [kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. |
| [in] | rho_0 | A density [R ~> kg m-3] or [kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state. |
| [in] | g_e | The Earth's gravitational acceleration [L2 Z-1 T-2 ~> m s-2] |
| [in] | dz_subroundoff | A minuscule thickness change [Z ~> m] |
| [in] | bathyt | The depth of the bathymetry [Z ~> m] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dpa | The change in the pressure anomaly across the layer [R L2 T-2 ~> Pa] |
| [in,out] | intz_dpa | The integral through the thickness of the layer of |
| [in,out] | intx_dpa | The integral in x of the difference between the |
| [in,out] | inty_dpa | The integral in y of the difference between the |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 298 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_density_dz_generic_ppm | ( | integer, intent(in) | k, |
| type(thermo_var_ptrs), intent(in) | tv, | ||
| real, dimension(szi_(hi),szj_(hi),szk_(gv)), intent(in) | T_t, | ||
| real, dimension(szi_(hi),szj_(hi),szk_(gv)), intent(in) | T_b, | ||
| real, dimension(szi_(hi),szj_(hi),szk_(gv)), intent(in) | S_t, | ||
| real, dimension(szi_(hi),szj_(hi),szk_(gv)), intent(in) | S_b, | ||
| real, dimension(szi_(hi),szj_(hi),szk_(gv)+1), intent(in) | e, | ||
| real, intent(in) | rho_ref, | ||
| real, intent(in) | rho_0, | ||
| real, intent(in) | G_e, | ||
| real, intent(in) | dz_subroundoff, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | bathyT, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(verticalgrid_type), intent(in) | GV, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout) | dpa, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout), optional | intz_dpa, | ||
| real, dimension( hi %isdb: hi %iedb, hi %jsd: hi %jed), intent(inout), optional | intx_dpa, | ||
| real, dimension( hi %isd: hi %ied, hi %jsdb: hi %jedb), intent(inout), optional | inty_dpa, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
Compute pressure gradient force integrals for layer "k" and the case where T and S are parabolic profiles.
| [in] | k | Layer index to calculate integrals for |
| [in] | hi | Ocean horizontal index structures for the input arrays |
| [in] | gv | Vertical grid structure |
| [in] | tv | Thermodynamic variables |
| [in] | t_t | Potential temperature at the cell top [degC] |
| [in] | t_b | Potential temperature at the cell bottom [degC] |
| [in] | s_t | Salinity at the cell top [ppt] |
| [in] | s_b | Salinity at the cell bottom [ppt] |
| [in] | e | Height of interfaces [Z ~> m] |
| [in] | rho_ref | A mean density [R ~> kg m-3] or [kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. |
| [in] | rho_0 | A density [R ~> kg m-3] or [kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state. |
| [in] | g_e | The Earth's gravitational acceleration [m s-2] |
| [in] | dz_subroundoff | A minuscule thickness change [Z ~> m] |
| [in] | bathyt | The depth of the bathymetry [Z ~> m] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dpa | The change in the pressure anomaly across the layer [R L2 T-2 ~> Pa] |
| [in,out] | intz_dpa | The integral through the thickness of the layer of |
| [in,out] | intx_dpa | The integral in x of the difference between the |
| [in,out] | inty_dpa | The integral in y of the difference between the |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 667 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_spec_vol_dp_generic_pcm | ( | real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | T, |
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | S, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | p_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | p_b, | ||
| real, intent(in) | alpha_ref, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout) | dza, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout), optional | intp_dza, | ||
| real, dimension( hi %isdb: hi %iedb, hi %jsd: hi %jed), intent(inout), optional | intx_dza, | ||
| real, dimension( hi %isd: hi %ied, hi %jsdb: hi %jedb), intent(inout), optional | inty_dza, | ||
| integer, intent(in), optional | halo_size, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in), optional | bathyP, | ||
| real, intent(in), optional | dP_neglect, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
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.
| [in] | hi | A horizontal index type structure. |
| [in] | t | Potential temperature of the layer [degC] |
| [in] | s | Salinity of the layer [ppt] |
| [in] | p_t | Pressure atop the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | p_b | Pressure below the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | alpha_ref | 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 | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dza | The change in the geopotential anomaly |
| [in,out] | intp_dza | The integral in pressure through the layer of |
| [in,out] | intx_dza | The integral in x of the difference between |
| [in,out] | inty_dza | The integral in y of the difference between |
| [in] | halo_size | The width of halo points on which to calculate dza. |
| [in] | bathyp | The pressure at the bathymetry [R L2 T-2 ~> Pa] or [Pa] |
| [in] | dp_neglect | A minuscule pressure change with the same units as p_t [R L2 T-2 ~> Pa] or [Pa] |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 1062 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_spec_vol_dp_generic_plm | ( | real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | T_t, |
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | T_b, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | S_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | S_b, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | p_t, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | p_b, | ||
| real, intent(in) | alpha_ref, | ||
| real, intent(in) | dP_neglect, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(in) | bathyP, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout) | dza, | ||
| real, dimension( hi %isd: hi %ied, hi %jsd: hi %jed), intent(inout), optional | intp_dza, | ||
| real, dimension( hi %isdb: hi %iedb, hi %jsd: hi %jed), intent(inout), optional | intx_dza, | ||
| real, dimension( hi %isd: hi %ied, hi %jsdb: hi %jedb), intent(inout), optional | inty_dza, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
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.
| [in] | hi | A horizontal index type structure. |
| [in] | t_t | Potential temperature at the top of the layer [degC] |
| [in] | t_b | Potential temperature at the bottom of the layer [degC] |
| [in] | s_t | Salinity at the top the layer [ppt] |
| [in] | s_b | Salinity at the bottom the layer [ppt] |
| [in] | p_t | Pressure atop the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | p_b | Pressure below the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | alpha_ref | 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. |
| [in] | dp_neglect | !< A miniscule pressure change with the same units as p_t [R L2 T-2 ~> Pa] or [Pa] |
| [in] | bathyp | The pressure at the bathymetry [R L2 T-2 ~> Pa] or [Pa] |
| eos | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dza | The change in the geopotential anomaly |
| [in,out] | intp_dza | The integral in pressure through the layer of |
| [in,out] | intx_dza | The integral in x of the difference between |
| [in,out] | inty_dza | The integral in y of the difference between |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 1277 of file MOM_density_integrals.F90.
| subroutine, public mom_density_integrals::int_specific_vol_dp | ( | real, dimension(szi_(hi),szj_(hi)), intent(in) | T, |
| real, dimension(szi_(hi),szj_(hi)), intent(in) | S, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in) | p_t, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in) | p_b, | ||
| real, intent(in) | alpha_ref, | ||
| type(hor_index_type), intent(in) | HI, | ||
| type(eos_type), pointer | EOS, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout) | dza, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(inout), optional | intp_dza, | ||
| real, dimension(szib_(hi),szj_(hi)), intent(inout), optional | intx_dza, | ||
| real, dimension(szi_(hi),szjb_(hi)), intent(inout), optional | inty_dza, | ||
| integer, intent(in), optional | halo_size, | ||
| real, dimension(szi_(hi),szj_(hi)), intent(in), optional | bathyP, | ||
| real, intent(in), optional | dP_tiny, | ||
| logical, intent(in), optional | useMassWghtInterp | ||
| ) |
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.
| [in] | hi | The horizontal index structure |
| [in] | t | Potential temperature referenced to the surface [degC] |
| [in] | s | Salinity [ppt] |
| [in] | p_t | Pressure at the top of the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | p_b | Pressure at the bottom of the layer [R L2 T-2 ~> Pa] or [Pa] |
| [in] | alpha_ref | 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 | Equation of state structure | |
| [in] | us | A dimensional unit scaling type |
| [in,out] | dza | The change in the geopotential anomaly across |
| [in,out] | intp_dza | The integral in pressure through the layer of the |
| [in,out] | intx_dza | The integral in x of the difference between the |
| [in,out] | inty_dza | The integral in y of the difference between the |
| [in] | halo_size | The width of halo points on which to calculate dza. |
| [in] | bathyp | The pressure at the bathymetry [R L2 T-2 ~> Pa] or [Pa] |
| [in] | dp_tiny | A minuscule pressure change with the same units as p_t [R L2 T-2 ~> Pa] or [Pa] |
| [in] | usemasswghtinterp | If true, uses mass weighting to interpolate T/S for top and bottom integrals. |
Definition at line 1004 of file MOM_density_integrals.F90.