Other Classes

Class eos

class eos

Phenomenological EOS for homogeneous nucleonic matter.

Subclassed by eos_nuclei

Grid specification

size_t n_nB2
size_t n_Ye2
size_t n_T2
size_t n_S2
std::vector<double> nB_grid2
std::vector<double> Ye_grid2
std::vector<double> T_grid2
std::vector<double> S_grid2
std::string nB_grid_spec

The function for default baryon density grid.

This parameter is used by the new_table() function, and the check-virial and eos-deriv commands.

std::string Ye_grid_spec

The function for default electron fraction grid.

This parameter is used by the new_table() function, and the check-virial and eos-deriv commands.

std::string T_grid_spec

The function for default temperature grid.

This parameter is used by the new_table() function, and the check-virial and eos-deriv commands.

std::string S_grid_spec

The function for default strangeness grid.

This parameter is used by the new_table() function, and the check-virial and eos-deriv commands.

Particle objects [protected]

o2scl::fermion neutron

Neutron.

o2scl::fermion proton

Proton.

o2scl::eos_leptons elep

New lepton object.

o2scl::fermion electron

Electron/positron.

o2scl::fermion muon

Muon/anti-muon.

o2scl::boson photon

Photon.

o2scl::fermion n_chiral

Neutron for chiral part.

o2scl::fermion p_chiral

Proton for chiral part.

o2scl::fermion neutrino

Neutrino.

Settings [public]

bool old_version

If true, use the EOS from the Du et al. (2019) paper instead of the Du et al. (2022) update (default false)

bool use_alt_eos

Use an alternate EOS rather than the Du et al. combined EOS (default false)

bool strange_axis

If true, strangeness is included as an additional tensor rank (default false)

bool test_ns_cs2

If true, test the neutron star speed of sound (default true)

bool ns_record

If true, save the results of the neutron star fit to a file, and immediately exit (default false)

bool old_ns_fit

If true, use the old neutron star fit (default true)

This defaults to true because the old fit performs a bit better than the new one. The new fit was never used in a publication.

int verbose

Generic verbosity parameter (default 0)

See also cs2_verbose and function_verbose .

bool output_files

If true, create output files for individual EOSs.

double a_virial

Coefficient for modulation of virial EOS (default 10.0)

double b_virial

Coefficient for modulation of virial EOS (default 10.0)

bool include_muons

If true, include muons (default false)

bool select_cs2_test

If true, test cs2 in the select_internal() function (default true)

Main EOS parameters [protected]

double qmc_alpha

The first exponent for density in the QMC EOS (unitless)

double qmc_a

The first coefficient for the QMC EOS (in MeV)

double phi

The speed of sound in neutron star matter at \( n_B=2~\mathrm{fm}^{-3} \).

double eos_S

The symmetry energy.

double eos_L

The slope of the symmetry energy.

int i_ns

The index of the neutron star model.

int i_skyrme

The index of the Skyrme model.

Internal variables [protected]

o2scl::table_units nstar_tab

The table which stores the neutron star EOS results.

o2scl::table_units UNEDF_tab

The table which stores the Skyrme fits.

bool model_selected

If true, a model has been selected (default false)

o2scl::rng rng

Random number generator.

EOS outputs

double f_deg

The free energy of degenerate matter.

double f_virial

The virial free energy.

double s_virial

The virial entropy.

double Lambda_bar_14

The value of \( \bar{\Lambda} \) for a 1.4 solar mass neutron star.

The fit to the neutron star EOS [protected]

std::vector<double> ns_fit_parms

Parameters for the function which fits the neutron star EOS.

double chi2_ns

The chi-squared for the neutron star fit.

double ns_nb_max

The maximum baryon density at which the neutron star EOS is causal.

This quantity is determined by ns_fit()

double energy_density_ns(double nn)

Compute the energy density (in \( \mathrm{fm}^{-4} \)) of neutron matter at high density from the neutron star data using the most recent fit (without the rest mass contribution)

Note

Currently this just returns the value of ed_fit() .

double fit_fun(size_t np, const std::vector<double> &parms, double nb)

The fit function for the energy per particle in units of MeV as a function of the baryon density (in \( \mathrm{fm}^{-3} \) )

Note this function does not include the rest mass energy density for the nucleons.

double ed_fit(double nb)

The energy density (in \( \mathrm{fm}^{-4} \) ) as a function of baryon density (in \( \mathrm{fm}^{-3} \) )

Note this function does not include the rest mass energy density for the nucleons.

double dmudn_fit(double nb)

The inverse susceptibility (in \( \mathrm{fm}^{2} \) ) as a function of baryon density (in \( \mathrm{fm}^{-3} \) )

double cs2_fit(double nb)

The speed of sound as a function of baryon density (in \( \mathrm{fm}^{-3} \) )

void min_max_cs2(double &cs2_min, double &cs2_max)

Compute the minimum and maximum speed of sound between 0.08 and ns_nb_max.

void ns_fit(int row)

Fit neutron star data from Bamr to an analytical expression.

double mu_fit(double nb)

The baryon number chemical potential (in \( \mathrm{fm}^{-1} \) ) as a function of number density (in \( \mathrm{fm}^{-3} \) )

Note this function does not include the rest mass for the nucleons.

Parameter objects

o2scl::cli::parameter_int p_verbose
o2scl::cli::parameter_bool p_old_ns_fit
o2scl::cli::parameter_bool p_ns_record
o2scl::cli::parameter_bool p_include_muons
o2scl::cli::parameter_bool p_select_cs2_test
o2scl::cli::parameter_bool p_test_ns_cs2
o2scl::cli::parameter_bool p_use_alt_eos
o2scl::cli::parameter_double p_a_virial
o2scl::cli::parameter_double p_b_virial
o2scl::cli::parameter_int p_cs2_verbose
o2scl::cli::parameter_string p_nB_grid_spec
o2scl::cli::parameter_string p_Ye_grid_spec
o2scl::cli::parameter_string p_T_grid_spec
o2scl::cli::parameter_string p_S_grid_spec
bool rmf_fields

If true, then RMF fields are included.

Base physics objects [protected]

o2scl::eos_had_virial vsd

The virial equation solver (now part of O2scl)

o2scl::fermion_rel relf

Object for computing electron/positron thermodynamic integrals.

o2scl::thermo th2

Thermodynamic quantities.

o2scl::thermo th_chiral

Thermodynamic quantities for chiral part.

o2scl::eos_had_skyrme sk

Base EOS model.

o2scl::eos_had_rmf rmf

Base RMF model.

o2scl::eos_had_rmf_hyp rmf_hyp

Base RMF model.

o2scl::eos_had_skyrme sk_Tcorr

Skyrme model for finite-temperature correction.

o2scl::eos_had_temp_eden_base *eos_Tcorr

Pointer to EOS for finite-temperature corrections.

eos_crust_virial_v2 ecv

The virial EOS.

o2scl::eos_had_skyrme sk_alt

Alternative skryme model.

o2scl::eos_had_temp_base *eosp_alt

Pointer to alternative model.

std::string alt_name

Name of the alternate EOS.

The parameters for the QMC energy density [protected]

double qmc_beta

The second exponent for density in the QMC EOS (unitless)

double qmc_b

The second coefficient for the QMC EOS (in MeV)

double qmc_n0

The saturation density of the QMC EOS, equal to \( 0.16~\mathrm{fm}^{-3} \).

Output saturation properties [protected]

double eos_EoA

The binding energy per particle.

double eos_K

The incompressibility.

double eos_n0

The saturation density.

Basic EOS functions [protected]

double free_energy_density(o2scl::fermion &n, o2scl::fermion &p, double T, o2scl::thermo &th)

Return the total free energy density of matter (without the rest mass contribution for the nucleons)

double free_energy_density_detail(o2scl::fermion &n, o2scl::fermion &p, double T, o2scl::thermo &th, std::map<std::string, double> &vdet)

Compute the free energy returning several details as output parameters.

f1 is g*f_virial f2 is (1-g)*f_skyrme f3 is (1-g)*delta^2*esym f4 is (1-g)*delta f_hot

so that the total homogeneous free energy is f1+f2+f3+f4

virtual double free_energy_density_virial(o2scl::fermion &n, o2scl::fermion &p, double T, o2scl::thermo &th, double &dmundnn, double &dmundnp, double &dmupdnn, double &dmupdnp, double &dmundT, double &dmupdT)

Compute the free energy density using the virial expansion including derivative information.

inline virtual double free_energy_density_virial(o2scl::fermion &n, o2scl::fermion &p, double T, o2scl::thermo &th)

Compute the free energy density using the virial expansion.

double free_energy_density_alt(o2scl::fermion &n, o2scl::fermion &p, double nn, double np, double T, o2scl::thermo &th)

Alternate form of free_energy_density() for computing derivatives.

This function does not include electrons or photons.

double free_energy_density_ep(double nn, double np, double T)

Alternate form of free_energy_density() which includes electrons, positrons, and photons.

double entropy(o2scl::fermion &n, o2scl::fermion &p, double nn, double np, double T, o2scl::thermo &th)

Compute the entropy density including photons, electrons, and positrons.

This function is used in cs2_func() .

double ed(o2scl::fermion &n, o2scl::fermion &p, double nn, double np, double T, o2scl::thermo &th)

Compute energy density including photons and electons (without the rest mass energy density for the nucleons)

double cs2_func(o2scl::fermion &n, o2scl::fermion &p, double T, o2scl::thermo &th)

Compute the squared speed of sound.

The temperature should be in \( 1/\mathrm{fm} \).

Other EOS functions [protected]

double energy_density_qmc(double nn, double pn)

Compute the energy density (in \( \mathrm{fm}^{-4} \)) of neutron matter from quantum Monte Carlo (without the rest mass contribution)

int new_ns_eos(double nb, o2scl::fermion &n, double &e_ns, double &densdnn)

Construct a new neutron star EOS which ensures causality at high densities.

double dfdnn_total(o2scl::fermion &n, o2scl::fermion &p, double nn, double pn, double T, o2scl::thermo &th)

Compute dfdnn including photons and electons.

This function is used in cs2_func() .

double dfdnp_total(o2scl::fermion &n, o2scl::fermion &p, double nn, double pn, double T, o2scl::thermo &th)

Compute dfdnp including photons and electons.

This function is used in cs2_func() .

int solve_Ye(size_t nv, const ubvector &x, ubvector &y, double nb, double T, double muL)

Solve for Ye to ensure a specified value of muL at fixed T.

int solve_coeff_big(size_t nv, const ubvector &x, ubvector &y, double nb_last, double cs_ns_2, double cs_ns_last)

solve for a1 and a2 when cs_ns(2.0)>cs_ns(1.28)

int solve_coeff_small(size_t nv, const ubvector &x, ubvector &y, double nb_last, double cs_ns_2, double cs_ns_last)

solve for a1 and a2 when cs_ns(2.0)<cs_ns(1.28)

int select_internal(int i_ns_loc, int i_skyrme_loc, double qmc_alpha_loc, double qmc_a_loc, double eos_L_loc, double eos_S_loc, double phi_loc)

Internal select function.

Constructor and destructor

eos()
inline virtual ~eos()

Command-line interface functions [public]

int table_Ye(std::vector<std::string> &sv, bool itive_com)

Construct a table at fixed electron fraction.

<filename> <Ye>

Constructs the EOS as a table3d object and outputs to <filename>. Currently stores Fint, Pint, Sint, g, msn, and msp.

int hrg_load(std::vector<std::string> &sv, bool itive_com)

Load the HRG table.

<filename>

Loads a list of resonances from a text file.

int alt_model(std::vector<std::string> &sv, bool itive_com)

Use alternate, rather than the Du et al. EOS.

<”Skyrme”> <name> or <”RMF”> <name>, etc.

For a Skyrme model, the first argument should be the word Skyrme, and the second should be the name of the desired Skyrme model. Similarly for an RMF model. (Hyperons are not yet fully supported.)

int table_nB(std::vector<std::string> &sv, bool itive_com)

Construct a table at fixed baryon density.

<filename> <nB>

Constructs the EOS as a table3d object and outputs to <filename>. Currently stores Fint, Pint, Sint, g, msn, and msp.

int pns_eos(std::vector<std::string> &sv, bool itive_com)

Construct the PNS EOS and the M-R curve.

<entropy per baryon> <lepton fraction> <output filename>

Use YL=0 for beta equilibrium. Currently always uses a cold crust.

int table_full(std::vector<std::string> &sv, bool itive_com)

Construct a full 3D EOS table without nuclei.

<filename>

This constructs a full 3D EOS table without nuclei using the specified model. The resulting file has several tensor_grid objects including Fint, Eint, Pint, Sint, mun, mup, cs2, mue, F, E, P, and S. This function does not yet support muons or strangeness.

int test_deriv(std::vector<std::string> &sv, bool itive_com)

Test the first derivatives of the free energy (no nuclei)

(no parameters)

This function tests the first derivatives of the homogeneous matter EOS without nuclei. The model must be selected before running this function. It tests derivatives over a range of densities for three temperatures and two electron fractions.

int select_model(std::vector<std::string> &sv, bool itive_com)

Select an EOS model.

<i_ns> <i_skyrme> <alpha> <L>

int vir_comp(std::vector<std::string> &sv, bool itive_com)

Compare the virial and full EOS.

Params.

Help.

int point(std::vector<std::string> &sv, bool itive_com)

Evaluate the EOS at one (nB,Ye,T) point.

Params.

Help.

int random(std::vector<std::string> &sv, bool itive_com)

Select a random EOS model.

(No arguments.)

Select a random EOS, checking several physical constraints and re-selecting a new random EOS until all the constraints are met.

int comp_figs(std::vector<std::string> &sv, bool itive_com)

Compute the data for the comparison figures.

Params.

Help.

int mcarlo_data(std::vector<std::string> &sv, bool itive_com)

Compute the data for the Monte Carlo figures.

Params.

Help.

int vir_fit(std::vector<std::string> &sv, bool itive_com)

Fit the virial EOS.

<no parameters>

Fit the virial EOS table to the functional form specified in Du et al. (2019) using eos_crust_virial_v2::fit() .

int test_eg(std::vector<std::string> &sv, bool itive_com)

Test the electron and photon EOS.

[filename]

This function tests the electron and photon EOS to ensure that it does not call the error handler (it does not test accuracy). It uses a larger grid than the default EOS grid and stores the results in tensor_grid objects. If a file is specified, these tensor_grid objects are output to the specified file.

This function does not require an EOS table or any hadronic model specification.

int eos_sn(std::vector<std::string> &sv, bool itive_com)

Compare to other EOSs?

Params.

Help.

double free_energy_density_s(o2scl::fermion &n, o2scl::fermion &p, double Y_s, double T, o2scl::thermo &th)

Free energy density as a function of the strangeness fraction.

double free_energy_density_detail_s(o2scl::fermion &n, o2scl::fermion &p, double Y_s, double T, o2scl::thermo &th, std::map<std::string, double> &vdet)

Free energy density as a function of the strangeness fraction.

Miscellaneous functions [public]

int solve_fixed_sonb_YL(size_t nv, const ubvector &x, ubvector &y, double nB, double sonb, double YL)

Solve for fixed entropy per baryon and fixed lepton fraction.

int solve_T(size_t nv, const ubvector &x, ubvector &y, double nb, double Ye, double sonb)

Solve for T to ensure a specified value of sonb at fixed Ye.

virtual void setup_cli(o2scl::cli &cl, bool read_docs = true)

Setup the command-line interface with commands and parameters.

int comm_set(std::vector<std::string> &sv, bool itive_com)

Public Members

int cs2_verbose

If true, increase the verbosity for cs2 (default false)

std::vector<o2scl::part_pdg_db::pdg_entry> part_db

Particle database.

std::vector<o2scl::fermion> res_f

Fermionic resonances.

std::vector<o2scl::boson> res_b

Bosonic resonances.

Protected Functions

int process_grid_spec()

Process the grid specification strings and store the values in the arrays.

Class partition_func

class partition_func

Compute partition functions using Fowler prescription.

Public Functions

double delta_small_iand(double E)

Integrand for partition function when \( \delta \) is smaller than \( E_d \).

From eq. (25) & (27) in Shen 10.

double delta_small_iand_prime(double E)

Integrand for derivative of partition function with respect to temperature when \( \delta \) is smaller than \( E_d \).

From eq. (26) & (27) in Shen 10.

double delta_large_iand(double E)

Integrand when \( \delta \) is greater than \( E_d \).

From eq. (25) and (30) in Shen 10.

double delta_large_iand_prime(double E)

Integrand for temperature derivative when \( \delta \) is greater than \( E_d \).

From eq. (26) and (30) in Shen 10.

Public Members

double T_MeV

Temperature in MeV.

double a

Nuclear level spacing in 1/MeV.

double delta

Backshift parameter in MeV.

double C

Coefficient for large \( \delta \) in 1/MeV.

double Tc

Critical temperature connecting low and high energies in MeV.

Class eos_crust_virial_v2

class eos_crust_virial_v2 : public o2scl::eos_crust_virial

An updated version of o2scl::eos_crust_virial with a better fit for the virial coefficients.

Constructor and destructor

eos_crust_virial_v2()
inline virtual ~eos_crust_virial_v2()

Public Functions

double bn0_free()

Free nucleon virial coefficient.

double bpn0_free()

Free isospin nucleon virial coefficient.

double bn1_free()

Free spin nucleon virial coefficient.

double bpn1_free()

Free spin-isospin nucleon virial coefficient.

inline double bn0(double T)

Nucleon virial coefficient.

inline double bn1(double T)

Spin nucleon virial coefficient.

inline double bpn0(double T)

Isopin nucleon virial coefficient.

inline double bpn1(double T)

Spin-isospin nucleon virial coefficient.

double bna(double T)

Nucleon-alpha virial coefficient.

The temperature must be specified in MeV

double bpna(double T)

Nucleon-alpha isospin virial coefficient.

The temperature must be specified in MeV

double f0(double lambda, double T)

Fermi-liquid parameter \( F_0 \).

The value of lambda should be in 1/MeV and the temperature should be specified in MeV. The result is returned in units of 1/MeV^2.

double f0p(double lambda, double T)

Fermi-liquid parameter \( F_0^{\prime} \).

The value of lambda should be in 1/MeV and the temperature should be specified in MeV. The result is returned in units of 1/MeV^2.

double g0(double lambda, double T)

Fermi-liquid parameter \( G_0 \).

The value of lambda should be in 1/MeV and the temperature should be specified in MeV. The result is returned in units of 1/MeV^2.

double g0p(double lambda, double T)

Fermi-liquid parameter \( G_0^{\prime} \).

The value of lambda should be in 1/MeV and the temperature should be specified in MeV. The result is returned in units of 1/MeV^2.

double bn_func(size_t np, const std::vector<double> &par, double T)

The neutron-neutron virial coefficient given the function parameters specified in par.

double bpn_func(size_t np, const std::vector<double> &par, double T)

The neutron-proton virial coefficient given the function parameters specified in par.

double bn_f(double T)

The neutron-neutron virial coefficient.

The temperature should be specified in MeV.

double bpn_f(double T)

The neutron-proton virial coefficient.

The temperature should be specified in MeV.

double dbndT_f(double T)

The temperature derivative of the neutron-neutron virial coefficient.

The temperature should be specified in MeV.

double dbpndT_f(double T)

The temperature derivative of the neutron-proton virial coefficient.

The temperature should be specified in MeV.

virtual void fit(bool show_fit = false)

Perform the fit to the scattering data.

Public Members

bool include_deuteron

If true, include the deuteron contribution in the virial coefficients.

std::vector<double> ba_T

Temperature grid for alpha-nucleon virial coefficients.

std::vector<double> vba

Alpha-nucleon virial coefficient.

std::vector<double> vbpna

Isospin alpha-nucleon virial coefficient.

o2scl::interp_vec<std::vector<double>> iv_ba

Interpolator for alpha-nucleon virial coefficient.

o2scl::interp_vec<std::vector<double>> iv_bpna

Interpolator for isospin alpha-nucleon virial coefficient.

std::vector<double> bn_params

The current neutron-neutron virial coefficient parameters.

std::vector<double> bpn_params

The current neutron-proton virial coefficient parameters.

Public Static Attributes

static const size_t bn_np = 10

The number of neutron-neutron virial coefficient parameters.

static const size_t bpn_np = 6

The number of neutron-proton virial coefficient parameters.

Class eos_had_skyrme_ext

class eos_had_skyrme_ext : public o2scl::eos_had_skyrme

Extended Skyrme hadronic equation of state.

This is the modified Skyrme model proposed by Holt and Lim.

The new parameters

double t4
double x4
double t5
double x5
double alpha2
double alpha3

Basic usage

eos_had_skyrme_ext()

Create a blank Skyrme EOS.

inline virtual ~eos_had_skyrme_ext()

Destructor.

virtual int calc_temp_e(o2scl::fermion &ne, o2scl::fermion &pr, double temper, o2scl::thermo &th)

Equation of state as a function of densities at finite temperature.

virtual int calc_e(o2scl::fermion &ne, o2scl::fermion &pr, o2scl::thermo &lt)

Equation of state as a function of densities at zero temperature.

inline virtual const char *type()

Return string denoting type (“eos_had_skyrme_ext”)

Protected Functions

template<class fermion_t>
inline void base_thermo(fermion_t &ne, fermion_t &pr, double ltemper, o2scl::thermo &locth, double term, double term2, double ham1, double ham2, double ham3, double ham4, double ham5, double ham6)

Compute the base thermodynamic quantities.

This function computes the energy density, pressure, entropy, and chemical potentials.