jesterTOV.eos module

The equation of state (EOS) module provides classes and functions for modeling neutron star matter equations of state.

class jesterTOV.eos.Interpolate_EOS_model[source]

Bases: object

Base class to interpolate EOS data.

interpolate_eos(n: Float[Array, 'n_points'], p: Float[Array, 'n_points'], e: Float[Array, 'n_points'])[source]

Given n, p and e, interpolate to obtain necessary auxiliary quantities.

Parameters:
  • n (Float[Array, n_points]) – Number densities. Expected units are n[fm^-3]

  • p (Float[Array, n_points]) – Pressure values. Expected units are p[MeV / fm^3]

  • e (Float[Array, n_points]) – Energy densities. Expected units are e[MeV / fm^3]

Returns:

Interpolated values of n, p, hs (enthalpy) e, and dloge_dlogps.

Return type:

tuple

class jesterTOV.eos.MetaModel_EOS_model(kappas: tuple[Float, Float, Float, Float, Float, Float] = (0.0, 0.0, 0.0, 0.0, 0.0, 0.0), v_nq: list[float] = [0.0, 0.0, 0.0, 0.0, 0.0], b_sat: Float = 17.0, b_sym: Float = 25.0, nsat: Float = 0.16, nmin_MM_nsat: Float = 0.75, nmax_nsat: Float = 12, ndat: Int = 200, crust_name: str = 'DH', max_n_crust_nsat: Float = 0.5, ndat_spline: Int = 10, proton_fraction: bool | float | None = None)[source]

Bases: Interpolate_EOS_model

MetaModel_EOS_model is a class to interpolate EOS data with a meta-model.

Parameters:

Interpolate_EOS_model (object) – Base class of interpolation EOS data.

compute_b(delta: Array | float)[source]
compute_cs2(n: Array, p: Array, e: Array, x: Array, delta: Array | float, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array)[source]
compute_energy(x: Array, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array) Array[source]
compute_f_1(delta: Array | float)[source]
compute_f_star(delta: Array | float)[source]
compute_f_star2(delta: Array | float)[source]
compute_f_star3(delta: Array | float)[source]
compute_pressure(x: Array, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array) Array[source]
compute_proton_fraction(coefficient_sym: list, n: Array) Float[Array, 'n_points'][source]

Computes the proton fraction for a given number density.

Parameters:

n (Float[Array, "n_points"]) – Number density in fm^-3.

Returns:

Proton fraction as a function of the number density.

Return type:

Float[Array, “n_points”]

compute_v(v_sat: Array, v_sym2: Array, delta: Array | float) Array[source]
compute_x(n: Array)[source]
construct_eos(NEP_dict: dict) tuple[source]

Construct the EOS.

Parameters:

NEP_dict (dict) – Dictionary with the NEP keys to be passed to the metamodel EOS class.

Returns:

EOS quantities (see Interpolate_EOS_model), as well as the chemical potential and speed of sound.

Return type:

tuple

esym(coefficient_sym: list, x: Array)[source]
u(x: Array, b: Array, alpha: Int)[source]
class jesterTOV.eos.MetaModel_with_CSE_EOS_model(nsat: Float = 0.16, nmin_MM_nsat: Float = 0.75, nmax_nsat: Float = 12, ndat_metamodel: Int = 100, ndat_CSE: Int = 100, **metamodel_kwargs)[source]

Bases: Interpolate_EOS_model

MetaModel_with_CSE_EOS_model is a class to interpolate EOS data with a meta-model and using the CSE.

construct_eos(NEP_dict: dict, ngrids: Float[Array, 'n_grid_point'], cs2grids: Float[Array, 'n_grid_point']) tuple[source]

Construct the EOS

Parameters:
  • NEP_dict (dict) – Dictionary with the NEP keys to be passed to the metamodel EOS class.

  • ngrids (Float[Array, n_grid_point]) – Density grid points of densities for the CSE part of the EOS.

  • cs2grids (Float[Array, n_grid_point]) – Speed-of-sound squared grid points of densities for the CSE part of the EOS.

Returns:

EOS quantities (see Interpolate_EOS_model), as well as the chemical potential and speed of sound.

Return type:

tuple

class jesterTOV.eos.MetaModel_with_peakCSE_EOS_model(nsat: Float = 0.16, nmin_MM_nsat: Float = 0.75, nmax_nsat: Float = 12, ndat_metamodel: Int = 100, ndat_CSE: Int = 100, **metamodel_kwargs)[source]

Bases: Interpolate_EOS_model

MetaModel_with_peakCSE_EOS_model is a class to interpolate EOS data with a meta-model and using the CSE. The parametrization of the CSE is based on the peakCSE model, which is a Gaussian peak with a logit growth rate, in order to guarantee consistency with pQCD at the highest densities.

Parameters:

Interpolate_EOS_model (object) – Base class of interpolation EOS data.

construct_eos(NEP_dict: dict, peakCSE_dict: dict)[source]
jesterTOV.eos.construct_family(eos: tuple, ndat: Int = 50, min_nsat: Float = 2) tuple[Float[Array, 'ndat'], Float[Array, 'ndat'], Float[Array, 'ndat'], Float[Array, 'ndat']][source]

Solve the TOV equations and generate the M, R and Lambda curves for the given EOS.

Parameters:
  • eos (tuple) – Tuple of the EOS data (ns, ps, hs, es).

  • ndat (int, optional) – Number of datapoints used when constructing the central pressure grid. Defaults to 50.

  • min_nsat (int, optional) – Starting density for central pressure in numbers of nsat (assumed to be 0.16 fm^-3). Defaults to 2.

Returns:

log(pcs), masses in solar masses, radii in km, and dimensionless tidal deformabilities

Return type:

tuple[Float[Array, “ndat”], Float[Array, “ndat”], Float[Array, “ndat”], Float[Array, “ndat”]]

jesterTOV.eos.construct_family_nonGR(eos: tuple, ndat: Int = 50, min_nsat: Float = 2) tuple[Float[Array, 'ndat'], Float[Array, 'ndat'], Float[Array, 'ndat'], Float[Array, 'ndat']][source]

Solve the post-TOV equations and generate the M, R and Lambda curves.

Parameters:
  • eos (tuple) – Tuple of the EOS data (ns, ps, hs, es).

  • ndat (int, optional) – Number of datapoints used when constructing the central pressure grid. Defaults to 50.

  • min_nsat (int, optional) – Starting density for central pressure in numbers of nsat (assumed to be 0.16 fm^-3). Defaults to 2.

Returns:

log(pcs), masses in solar masses, radii in km, and dimensionless tidal deformabilities

Return type:

tuple[Float[Array, “ndat”], Float[Array, “ndat”], Float[Array, “ndat”], Float[Array, “ndat”]]

jesterTOV.eos.load_crust(name: str) tuple[Array, Array, Array][source]

Load a crust file from the default directory.

Parameters:

name (str) – Name of the crust to load, or a filename if a file outside of jose is supplied.

Returns:

Number densities [fm^-3], pressures [MeV / fm^-3], and energy densities [MeV / fm^-3] of the crust.

Return type:

tuple[Array, Array, Array]

jesterTOV.eos.locate_lowest_non_causal_point(cs2)[source]

Key Classes

class jesterTOV.eos.MetaModel_EOS_model(kappas: tuple[Float, Float, Float, Float, Float, Float] = (0.0, 0.0, 0.0, 0.0, 0.0, 0.0), v_nq: list[float] = [0.0, 0.0, 0.0, 0.0, 0.0], b_sat: Float = 17.0, b_sym: Float = 25.0, nsat: Float = 0.16, nmin_MM_nsat: Float = 0.75, nmax_nsat: Float = 12, ndat: Int = 200, crust_name: str = 'DH', max_n_crust_nsat: Float = 0.5, ndat_spline: Int = 10, proton_fraction: bool | float | None = None)[source]

Bases: Interpolate_EOS_model

MetaModel_EOS_model is a class to interpolate EOS data with a meta-model.

Parameters:

Interpolate_EOS_model (object) – Base class of interpolation EOS data.

compute_b(delta: Array | float)[source]
compute_cs2(n: Array, p: Array, e: Array, x: Array, delta: Array | float, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array)[source]
compute_energy(x: Array, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array) Array[source]
compute_f_1(delta: Array | float)[source]
compute_f_star(delta: Array | float)[source]
compute_f_star2(delta: Array | float)[source]
compute_f_star3(delta: Array | float)[source]
compute_pressure(x: Array, f_1: Array, f_star: Array, f_star2: Array, f_star3: Array, b: Array, v: Array) Array[source]
compute_proton_fraction(coefficient_sym: list, n: Array) Float[Array, 'n_points'][source]

Computes the proton fraction for a given number density.

Parameters:

n (Float[Array, "n_points"]) – Number density in fm^-3.

Returns:

Proton fraction as a function of the number density.

Return type:

Float[Array, “n_points”]

compute_v(v_sat: Array, v_sym2: Array, delta: Array | float) Array[source]
compute_x(n: Array)[source]
construct_eos(NEP_dict: dict) tuple[source]

Construct the EOS.

Parameters:

NEP_dict (dict) – Dictionary with the NEP keys to be passed to the metamodel EOS class.

Returns:

EOS quantities (see Interpolate_EOS_model), as well as the chemical potential and speed of sound.

Return type:

tuple

esym(coefficient_sym: list, x: Array)[source]
u(x: Array, b: Array, alpha: Int)[source]

Mathematical Background

The metamodel equation of state is based on the formulation from Margueron et al. (2021), which provides a flexible framework for modeling nuclear matter at various densities.

The total energy density is given by:

\[\varepsilon(n, \delta) = \varepsilon_{\text{kinetic}}(n, \delta) + \varepsilon_{\text{potential}}(n, \delta)\]

where \(n\) is the baryon number density and \(\delta = 1 - 2Y_p\) is the isospin asymmetry parameter.