Physical Constants

O2scl

Constant contents

Constant introduction

In order to avoid confusing numerical differences when using multiprecision arithmetic, physical constants are template functions which return a value given a user-specified floating-point type. Physical constants are promoted to higher precision by adding zeros in the base-10 representation. The constants are in the namespace o2scl_const. The numerical values are periodically updated with CODATA, the Particle Data Book, and other databases. A handful of constants have a short name for their values in double precision, see below.

The find_constants class contains a simple constant database which can be searched at compiled time and also provides the constant database to acol -constants (see The acol Command Line Utility).

These physical constants can also be used to create unit conversion factors, described in Unit Conversions.

Namespace o2scl_const

Top

namespace o2scl_const

Constants.

CODATA 2014 values were from [Mohr16] and previous versions contained constants from [Luzum11] and [Mohr12] .

CODATA 2018 values are from physics.nist.gov/constants. IAU 2015 values are the nominal values from arXiv:1510.07674 and arXiv:1605.09788 .

Unit prefixes

const double quetta = 1e30
const double ronna = 1e27
const double yotta = 1e24
const double zetta = 1e21
const double exa = 1e18
const double peta = 1e15
const double tera = 1e12
const double giga = 1e9
const double mega = 1e6
const double kilo = 1e3
const double milli = 1e-3
const double micro = 1e-6
const double nano = 1e-9
const double pico = 1e-12
const double femto = 1e-15
const double atto = 1e-18
const double zepto = 1e-21
const double yocto = 1e-24
const double ronto = 1e-27
const double quecto = 1e-30

Unit systems

static const size_t o2scl_mks = 1

MKS units.

static const size_t o2scl_cgs = 2

CGS units.

A few double-precision constants for convenience

const double pi = boost::math::constants::pi<double>()

\( \pi \)

const double pi2 = boost::math::constants::pi_sqr<double>()

\( \pi^2 \)

const double root_pi = boost::math::constants::root_pi<double>()

\( \sqrt{\pi} \)

const double hc_mev_fm = hc_mev_fm_f<double>()

\( \hbar c \) in MeV fm (exact)

const double hc_mev_cm = hc_mev_fm * 1.0e-13

\( \hbar c \) in MeV cm (exact)

Mathematical constants

template<class fp_t>
fp_t pi_f()

\( \pi \)

template<class fp_t>
fp_t pi2_f()

\( \pi^2 \)

template<class fp_t>
fp_t root_pi_f()

\( \sqrt{\pi} \)

template<class fp_t>
fp_t zeta2_f()

\( \zeta(2) \)

template<class fp_t>
fp_t zeta3_f()

\( \zeta(3) \)

template<class fp_t>
fp_t euler_f()

The Euler-Mascheroni constant.

Physical constants

template<class fp_t>
fp_t fine_structure_f()

Fine structure constant (CODATA 2018 value)

template<class fp_t>
fp_t avogadro_f()

Avogadro’s number (exact)

template<class fp_t>
fp_t speed_of_light_f(size_t system = o2scl_mks)

Speed of light (exact)

template<class fp_t>
fp_t planck_f(size_t system = o2scl_mks)

Planck constant (exact)

template<class fp_t>
fp_t hbar_f(size_t system = o2scl_mks)

Reduced Planck constant (derived; exact)

template<class fp_t>
fp_t hbarc_f(size_t system = o2scl_mks)

Reduced Planck’s constant times speed of light \( \hbar c \) (derived; exact)

template<class fp_t>
fp_t boltzmann_f(size_t system = o2scl_mks)

Boltzmann’s constant (exact)

template<class fp_t>
fp_t gravitational_constant_f(size_t system = o2scl_mks)

Gravitational constant (CODATA 2018 value)

template<class fp_t>
fp_t elem_charge_f()

Elementary charge in C (exact)

template<class fp_t>
fp_t electron_volt_f(size_t system = o2scl_mks)

Electron volt (exact)

template<class fp_t>
fp_t hc_mev_fm_f()

Reduced Planck constant times speed of light in \( \mathrm{MeV}~\mathrm{fm} \) (derived; exact)

template<class fp_t>
fp_t stefan_boltz_cons_f(size_t system = o2scl_mks)

Stefan-Boltzmann constant ( \( \mathrm{kg} / \mathrm{K}^4 \mathrm{s}^3 \) in MKS and \( \mathrm{g} / \mathrm{K}^4 \mathrm{s}^3 \) in CGS; derived; exact)

template<class fp_t>
fp_t sin2_theta_weak_f()

\( \sin^2 \theta_W \) (unitless; PDG 2020 value)

template<class fp_t>
fp_t gfermi_gev2_f()

Fermi coupling constant in \( \mathrm{GeV}^{-2} \) (CODATA 2018 value)

template<class fp_t>
fp_t gfermi_f(size_t system = o2scl_mks)

Fermi coupling constant in ( \( \mathrm{s}^4 / \mathrm{m}^4 \mathrm{kg}^2 \) in MKS and \( \mathrm{s}^4 / \mathrm{cm}^4 \mathrm{g}^2 \) in CGS; derived from gfermi_gev2_f())

Astrophysical constants

template<class fp_t>
fp_t astronomical_unit_f(size_t system = o2scl_mks)

Astronomical unit (IAU 2009 value; now exact)

template<class fp_t>
fp_t parsec_f(size_t system = o2scl_mks)

Parsec (derived; exact)

template<class fp_t>
fp_t grav_accel_f(size_t system = o2scl_mks)

Acccleration due to gravity ( \( \mathrm{m} / \mathrm{s}^2 \) in MKS and \( \mathrm{cm} / \mathrm{s}^2 \) in CGS ) (CODATA 2018; now exact)

template<class fp_t>
fp_t sidereal_year_f(size_t system = o2scl_mks)

Sidereal year in s (from https://pdg.lbl.gov/2021/reviews/contents_sports.html)

template<class fp_t>
fp_t tropical_year_f(size_t system = o2scl_mks)

Tropical year in s (from https://pdg.lbl.gov/2021/reviews/contents_sports.html)

template<class fp_t>
fp_t julian_year_f(size_t system = o2scl_mks)

Julian year in s (exact)

template<class fp_t>
fp_t light_year_f(size_t system = o2scl_mks)

Light year in \( \mathrm{cm} \) (derived; exact)

Particle masses

template<class fp_t>
fp_t mass_neutron_f(size_t system = o2scl_mks)

Neutron mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_proton_f(size_t system = o2scl_mks)

Proton mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_proton_amu_f()

Proton mass in AMU (CODATA 2018 value)

template<class fp_t>
fp_t mass_electron_f(size_t system = o2scl_mks)

Electron mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_muon_f(size_t system = o2scl_mks)

Muon mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_tau_f(size_t system = o2scl_mks)

Tau mass (CODATA 2018 value)

Nuclear masses

template<class fp_t>
fp_t mass_deuteron_f(size_t system = o2scl_mks)

Deuteron mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_triton_f(size_t system = o2scl_mks)

Triton mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_helion_f(size_t system = o2scl_mks)

Helion mass (CODATA 2018 value)

template<class fp_t>
fp_t mass_alpha_f(size_t system = o2scl_mks)

Alpha mass (CODATA 2018 value)

template<class fp_t>
fp_t unified_atomic_mass_f(size_t system = o2scl_mks)

Deuteron mass (CODATA 2018 value)

template<class fp_t>
fp_t bohr_radius_f(size_t system = o2scl_mks)

Bohr radius (CODATA 2018 value)

template<class fp_t>
fp_t thomson_csec_f(size_t system = o2scl_mks)

Thomson cross section (MKS in \( m^2 \) and CGS in \( cm^2 \) ; CODATA 2018 value)

Chemical constants

template<class fp_t>
fp_t rydberg_f(size_t system = o2scl_mks)

Rydberg constant ( \( \mathrm{kg} \mathrm{m}^2 / \mathrm{s}^2 \) in MKS and \( \mathrm{g} \mathrm{cm}^2 / \mathrm{s}^2 \) in CGS; CODATA 2018 value)

This is referred to as the “Rydberg constant” by GSL and is equal to \( h c R_{\infty} \). The value \( R_{\infty} \) is the inverse Rydberg wavelength (also sometimes referred to as the Rydberg constant).

template<class fp_t>
fp_t molar_gas_f(size_t system = o2scl_mks)

Molar gas constant, \( R\) , ( \( \mathrm{kg} \mathrm{m}^2 / \mathrm{K} \mathrm{mol} \mathrm{s}^2 \) in MKS \( \mathrm{g} \mathrm{cm}^2 / \mathrm{K} \mathrm{mol} \mathrm{s}^2 \) in CGS; derived; exact)

template<class fp_t>
fp_t std_gas_volume_f(size_t system = o2scl_mks)

Molar volume of ideal gas at standard T and P ( \( \mathrm{m}/\mathrm{mol} \) in MKS and ( \( \mathrm{cm}/\mathrm{mol} \) in CGS; CODATA 2018 value)

Electromagnetic constants

template<class fp_t>
fp_t electron_mag_mom_f(size_t system = o2scl_mks)

Electron magnetic moment ( \( \mathrm{J}/\mathrm{T} \) in MKS; CODATA 2018 value)

Note

This quantity is defined without a minus sign.

template<class fp_t>
fp_t proton_mag_mom_f(size_t system = o2scl_mks)

Proton magnetic moment ( \( \mathrm{J}/\mathrm{T} \) in MKS; CODATA 2018 value)

template<class fp_t>
fp_t roentgen_f(size_t system = o2scl_mks)

Roentgen units?

template<class fp_t>
fp_t bohr_magneton_f(size_t system = o2scl_mks)

Bohr magneton ( \( \mathrm{J}/\mathrm{T} \) in MKS; CODATA 2018 value)

template<class fp_t>
fp_t nuclear_magneton_f(size_t system = o2scl_mks)

Nuclear magneton ( \( \mathrm{J}/\mathrm{T} \) in MKS; CODATA 2018 value)

template<class fp_t>
fp_t faraday_f(size_t system = o2scl_mks)

Faraday constant in A s / mol (derived; exact)

template<class fp_t>
fp_t vacuum_permittivity_f(size_t system = o2scl_mks)

Vacuum permittivity in MKS units (CODATA 2018)

template<class fp_t>
fp_t vacuum_permeability_f(size_t system = o2scl_mks)

Vacuum permeability in kg m / A^2 s^2.

template<class fp_t>
fp_t e2_gaussian_f()

Electron charge squared in Gaussian units (derived)

In Gaussian Units:

\[\begin{split}\begin{eqnarray*} &\vec{\nabla} \cdot \vec{E} = 4 \pi \rho \, , \quad \vec{E}=-\vec{\nabla} \Phi \, , \quad \nabla^2 \Phi = - 4 \pi \rho \, , &\\& F=\frac{q_1 q_2}{r^2} \, , \quad W=\frac{1}{2} \int \rho V d^3 x =\frac{1}{8 \pi} \int | \vec{E} |^2 d^3 x \, , \quad \alpha=\frac{e^2}{\hbar c}=\frac{1}{137}& \end{eqnarray*}\end{split}\]

template<class fp_t>
fp_t e2_hlorentz_f()

Electron charge sqaured in Heaviside-Lorentz units where \(\hbar=c=1\) (derived)

In Heaviside-Lorentz units:

\[\begin{split}\begin{eqnarray*} &\vec{\nabla} \cdot \vec{E} = \rho \, , \quad \vec{E}=-\vec{\nabla} \Phi \, , \quad \nabla^2 \Phi = - \rho \, , &\\& F=\frac{q_1 q_2}{4 \pi r^2} \, , \quad W=\frac{1}{2} \int \rho V d^3 x =\frac{1}{2} \int | \vec{E} |^2 d^3 x \, , \quad \alpha=\frac{e^2}{4 \pi}=\frac{1}{137}& \end{eqnarray*}\end{split}\]

template<class fp_t>
fp_t elem_charge_squared_f()

Electron charge squared in SI(MKS) units (derived)

In MKS units:

\[\begin{split}\begin{eqnarray*} &\vec{\nabla} \cdot \vec{E} = \rho \, , \quad \vec{E}=-\vec{\nabla} \Phi \, , \quad \nabla^2 \Phi = - \rho \, , &\\& F=\frac{1}{4 \pi \varepsilon_0}\frac{q_1 q_2}{r^2} \, , \quad W=\frac{1}{2} \int \rho V d^3 x =\frac{\varepsilon_0}{2} \int | \vec{E} |^2 d^3 x \, , \quad \alpha=\frac{e^2}{4 \pi \varepsilon_0 \hbar c}=\frac{1}{137}& \end{eqnarray*}\end{split}\]

Note the conversion formulas

\[ q_HL=\sqrt{4 \pi} q_G = \frac{1}{\sqrt{\varepsilon_0}} q_{SI} \]
as mentioned, e.g. in pg. 13 of D. Griffiths Intro to Elem. Particles.

template<class fp_t>
fp_t ec_gauss_fm2_f()

1 \(\mathrm{Gauss}\) times the electron charge in Gaussian units in \(\mathrm{fm}^{-2}\)

template<class fp_t>
fp_t gauss2_fm4_f()

Conversion factor from \( \mathrm{Gauss}^2 \) to \(\mathrm{fm}^{-4}\) in Gaussian units.

This is useful, e.g. in converting magnetic field squared to an energy density.

Time measurements (s)

template<class fp_t>
fp_t minute_f(size_t system = o2scl_mks)

Minute.

template<class fp_t>
fp_t hour_f(size_t system = o2scl_mks)

Hour.

template<class fp_t>
fp_t day_f(size_t system = o2scl_mks)

Day.

template<class fp_t>
fp_t week_f(size_t system = o2scl_mks)

Week.

Length measurements (MKS in m and CGS in cm)

template<class fp_t>
fp_t inch_f(size_t system = o2scl_mks)

Inch.

template<class fp_t>
fp_t foot_f(size_t system = o2scl_mks)

Foot.

template<class fp_t>
fp_t yard_f(size_t system = o2scl_mks)

Yard.

template<class fp_t>
fp_t mile_f(size_t system = o2scl_mks)

Mile.

template<class fp_t>
fp_t nautical_mile_f(size_t system = o2scl_mks)

Nautical mile.

template<class fp_t>
fp_t fathom_f(size_t system = o2scl_mks)

Fathom.

template<class fp_t>
fp_t mil_f(size_t system = o2scl_mks)

Mil.

template<class fp_t>
fp_t point_f(size_t system = o2scl_mks)

Point.

template<class fp_t>
fp_t texpoint_f(size_t system = o2scl_mks)

Texpoint.

template<class fp_t>
fp_t micron_f(size_t system = o2scl_mks)

Micron.

template<class fp_t>
fp_t angstrom_f(size_t system = o2scl_mks)

Angstrom.

Particle masses from PDG 2020

template<class fp_t>
fp_t mass_lambda_MeV_f()

\( \Lambda \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_sigma_minus_MeV_f()

\( \Sigma^{-} \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_sigma_zero_MeV_f()

\( \Sigma^{0} \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_sigma_plus_MeV_f()

\( \Sigma^{+} \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_cascade_zero_MeV_f()

\( \Xi^{0} \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_cascade_minus_MeV_f()

\( \Xi^{-} \) hyperon mass in \( \mathrm{MeV} \) (used value labeled “OUR FIT”)

template<class fp_t>
fp_t mass_up_MeV_f()

Up quark mass in \( \mathrm{MeV} \) (used value labeled “OUR EVALUATION”)

template<class fp_t>
fp_t mass_down_MeV_f()

Down quark mass in \( \mathrm{MeV} \) (used value labeled “OUR EVALUATION”)

template<class fp_t>
fp_t mass_strange_MeV_f()

Strange quark mass in \( \mathrm{MeV} \) (used value labeled “OUR EVALUATION”)

Solar system properties

template<class fp_t>
fp_t solar_mass_param_f(size_t system = o2scl_mks)

Solar mass times gravitational constant ( \( \mathrm{m}^3 / \mathrm{s}^2 \) in MKS and \( \mathrm{cm}^3 / \mathrm{s}^2 \) in CGS; IAU 2015 value, see https://arxiv.org/abs/1510.07674)

Note that this value differs slightly in Barycentric Coordinate Time and Barycentric Dynamical Time. This is the IAU’s nominal value.

template<class fp_t>
fp_t solar_mass_f(size_t system = o2scl_mks)

Mass of the sun in g (derived)

template<class fp_t>
fp_t schwarzchild_radius_f(size_t system = o2scl_mks)

Schwarzchild radius (MKS in m and CGS in cm) (derived)

template<class fp_t>
fp_t solar_radius_f(size_t system = o2scl_mks)

Radius of the sun (MKS in m and CGS in cm) (IAU 2015 nominal value)

template<class fp_t>
fp_t solar_temperature_f(size_t system = o2scl_mks)

Temperature of the sun’s photosphere in K (IAU 2015 nominal value)

template<class fp_t>
fp_t solar_luminosity_f(size_t system = o2scl_mks)

Luminosity of sun in erg/s (IAU 2015 nominal value)

template<class fp_t>
fp_t solar_irradiance_f()

Solar total irradiance in W/m^2 (IAU 2015 nominal value)

template<class fp_t>
fp_t earth_mass_param_f(size_t system = o2scl_mks)

Earth mass times gravitational constant ( \( \mathrm{m}^3 / \mathrm{s}^2 \) in MKS and \( \mathrm{cm}^3 / \mathrm{s}^2 \) in CGS; IAU 2015 nominal values)

template<class fp_t>
fp_t earth_mass_f(size_t system = o2scl_mks)

Mass of the earth in g (derived)

template<class fp_t>
fp_t earth_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of earth (MKS in m and CGS in cm) (IAU 2015 value)

template<class fp_t>
fp_t earth_radius_pol_f(size_t system = o2scl_mks)

Polar radius of earth (MKS in m and CGS in cm) (IAU 2015 value)

template<class fp_t>
fp_t jupiter_mass_param_f(size_t system = o2scl_mks)

Jupter mass times gravitational constant ( \( \mathrm{m}^3 / \mathrm{s}^2 \) in MKS and \( \mathrm{cm}^3 / \mathrm{s}^2 \) in CGS; IAU 2015 nominal values)

template<class fp_t>
fp_t jupiter_mass_f(size_t system = o2scl_mks)

Mass of jupiter in g (derived)

template<class fp_t>
fp_t jupiter_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of jupiter (MKS in m and CGS in cm) (IAU 2015 value)

template<class fp_t>
fp_t jupiter_radius_pol_f(size_t system = o2scl_mks)

Polar radius of jupiter (MKS in m and CGS in cm) (IAU 2015 value)

template<class fp_t>
fp_t mercury_mass_f(size_t system = o2scl_mks)

Mass of mercury in g.

template<class fp_t>
fp_t mercury_radius_f(size_t system = o2scl_mks)

Radius of mercury (MKS in m and CGS in cm)

template<class fp_t>
fp_t venus_mass_f(size_t system = o2scl_mks)

Mass of venus in g.

template<class fp_t>
fp_t venus_radius_f(size_t system = o2scl_mks)

Radius of venus (MKS in m and CGS in cm)

template<class fp_t>
fp_t mars_mass_f(size_t system = o2scl_mks)

Mass of mars in g.

template<class fp_t>
fp_t mars_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of mars (MKS in m and CGS in cm)

template<class fp_t>
fp_t mars_radius_pol_f(size_t system = o2scl_mks)

Polar radius of mars (MKS in m and CGS in cm)

template<class fp_t>
fp_t saturn_mass_f(size_t system = o2scl_mks)

Mass of saturn in g.

template<class fp_t>
fp_t saturn_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of saturn (MKS in m and CGS in cm)

template<class fp_t>
fp_t saturn_radius_pol_f(size_t system = o2scl_mks)

Polar radius of saturn (MKS in m and CGS in cm)

template<class fp_t>
fp_t uranus_mass_f(size_t system = o2scl_mks)

Mass of uranus in g.

template<class fp_t>
fp_t uranus_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of uranus (MKS in m and CGS in cm)

template<class fp_t>
fp_t uranus_radius_pol_f(size_t system = o2scl_mks)

Polar radius of uranus (MKS in m and CGS in cm)

template<class fp_t>
fp_t neptune_mass_f(size_t system = o2scl_mks)

Mass of neptune in g.

template<class fp_t>
fp_t neptune_radius_eq_f(size_t system = o2scl_mks)

Equatorial radius of neptune (MKS in m and CGS in cm)

template<class fp_t>
fp_t neptune_radius_pol_f(size_t system = o2scl_mks)

Polar radius of neptune (MKS in m and CGS in cm)

template<class fp_t>
fp_t pluto_mass_f(size_t system = o2scl_mks)

Mass of pluto in g.

template<class fp_t>
fp_t pluto_radius_f(size_t system = o2scl_mks)

Radius of pluto (MKS in m and CGS in cm)