10 #include <gsl/gsl_sf_bessel.h>
16 #include <boost/filesystem.hpp>
31 : de_(de), dnb_(dnb), n_e_(n_e), n_nb_(n_nb) {
36 const std::string &eos_savefile_name) {
37 bool table_read_success =
false, table_consistency =
true;
38 if (boost::filesystem::exists(eos_savefile_name)) {
39 std::cout <<
"Reading table from file " << eos_savefile_name << std::endl;
41 file.open(eos_savefile_name, std::ios::in);
45 for (
size_t ie = 0; ie <
n_e_; ie++) {
46 for (
size_t inb = 0; inb <
n_nb_; inb++) {
47 double p, T, mub, mus;
48 file >>
p >> T >> mub >> mus;
52 table_read_success =
true;
53 std::cout <<
"Table consumed successfully." << std::endl;
56 if (table_read_success) {
58 std::cout <<
"Checking consistency of the table... " << std::endl;
59 constexpr
size_t number_of_steps = 50;
60 const size_t ie_step = 1 +
n_e_ / number_of_steps;
61 const size_t inb_step = 1 +
n_nb_ / number_of_steps;
62 for (
size_t ie = 0; ie <
n_e_; ie += ie_step) {
63 for (
size_t inb = 0; inb <
n_nb_; inb += inb_step) {
71 const double eps = 1.e-3;
73 if ((std::abs(
de_ * ie - e_comp) > eps ||
74 std::abs(
dnb_ * inb - nb_comp) > eps || std::abs(ns_comp) > eps ||
75 std::abs(x.
p - p_comp) > eps) &&
77 std::cout <<
"discrepancy: " <<
de_ * ie <<
" = " << e_comp <<
", "
78 <<
dnb_ * inb <<
" = " << nb_comp <<
", " << x.
p <<
" = "
79 << p_comp <<
", 0 = " << ns_comp << std::endl;
80 table_consistency =
false;
81 goto finish_consistency_check;
86 finish_consistency_check:
88 if (!table_read_success || !table_consistency) {
89 std::cout <<
"Compiling an EoS table..." << std::endl;
90 const double ns = 0.0;
91 for (
size_t ie = 0; ie <
n_e_; ie++) {
92 std::cout << ie <<
"/" <<
n_e_ <<
"\r" << std::flush;
93 const double e =
de_ * ie;
94 for (
size_t inb = 0; inb <
n_nb_; inb++) {
95 const double nb =
dnb_ * inb;
103 std::array<double, 3> init_approx;
107 init_approx = {2.0 * x.
T - y.
T, 2.0 * x.
mub - y.
mub,
112 const std::array<double, 3> res = eos.
solve_eos(e, nb, ns, init_approx);
113 const double T = res[0];
114 const double mub = res[1];
115 const double mus = res[2];
121 std::cout <<
"Saving table to file " << eos_savefile_name << std::endl;
123 file.open(eos_savefile_name, std::ios::out);
124 file <<
de_ <<
" " <<
dnb_ << std::endl;
125 file <<
n_e_ <<
" " <<
n_nb_ << std::endl;
126 file << std::setprecision(7);
128 for (
size_t ie = 0; ie <
n_e_; ie++) {
129 for (
size_t inb = 0; inb <
n_nb_; inb++) {
131 file << x.
p <<
" " << x.
T <<
" " << x.
mub <<
" " << x.
mus << std::endl;
138 const size_t ie = static_cast<size_t>(std::floor(e /
de_));
139 const size_t inb = static_cast<size_t>(std::floor(nb /
dnb_));
141 if (ie >=
n_e_ - 1 || inb >=
n_nb_ - 1) {
142 res = {-1.0, -1.0, -1.0, -1.0};
145 const double ae = e /
de_ - ie;
146 const double an = nb /
dnb_ - inb;
151 res.
p = ae * (an * s4.
p + (1.0 - an) * s2.
p) +
152 (1.0 - ae) * (an * s3.
p + (1.0 - an) * s1.
p);
153 res.
T = ae * (an * s4.
T + (1.0 - an) * s2.
T) +
154 (1.0 - ae) * (an * s3.
T + (1.0 - an) * s1.
T);
155 res.
mub = ae * (an * s4.
mub + (1.0 - an) * s2.
mub) +
156 (1.0 - ae) * (an * s3.
mub + (1.0 - an) * s1.
mub);
157 res.
mus = ae * (an * s4.
mus + (1.0 - an) * s2.
mus) +
158 (1.0 - ae) * (an * s3.
mus + (1.0 - an) * s1.
mus);
163 : x_(gsl_vector_alloc(n_equations_)),
165 account_for_resonance_widths_(account_for_width) {
166 const gsl_multiroot_fsolver_type *solver_type;
167 solver_type = gsl_multiroot_fsolver_hybrid;
171 "Compilation of hadron gas EoS table requested with"
172 " account of resonance spectral functions. This is not "
173 "advised, as it will likely take a few days to finish."
174 " Besides, the effect of resonance widths is currently not "
175 "implemented for energy density computation, so the computed"
176 " table will be inconsistent anyways.");
184 gsl_multiroot_fsolver_free(
solver_);
190 double x = mu_over_T - m_over_T;
199 : m_over_T * m_over_T * x * gsl_sf_bessel_Kn_scaled(2, m_over_T);
203 double beta,
double mub,
double mus,
204 bool account_for_width) {
205 const double m_over_T = ptype.
mass() *
beta;
208 const double g = ptype.
spin() + 1;
209 if (ptype.
is_stable() || !account_for_width) {
214 const double m0 = ptype.
mass();
217 const double u_min = std::atan(2.0 * (mth - m0) / w0);
218 const double u_max = 0.5 * M_PI;
220 const double result =
221 g *
integrate(u_min, u_max, [&](
double u) {
224 const double tanu = std::tan(u);
225 const double m = m0 + 0.5 * w0 * tanu;
226 const double jacobian = 0.5 * w0 * (1.0 + tanu * tanu);
235 double mub,
double mus,
236 bool account_for_width) {
248 const double beta = 1.0 / T;
254 const double z = ptype.mass() *
beta;
255 double x =
beta * (mub * ptype.baryon_number() + mus * ptype.strangeness() -
261 const size_t g = ptype.spin() + 1;
265 (3.0 * gsl_sf_bessel_Kn_scaled(2, z) +
266 z * gsl_sf_bessel_K1_scaled(z));
273 bool account_for_width) {
277 const double beta = 1.0 / T;
290 bool account_for_width) {
294 const double beta = 1.0 / T;
301 ptype.baryon_number();
308 bool account_for_width) {
312 const double beta = 1.0 / T;
338 const double max_mass = 5.0;
345 const double m0 = ptype.
mass();
347 m0 * m0 * std::exp(-
beta * m0) * gsl_sf_bessel_Kn_scaled(2, m0 *
beta) *
350 constexpr
int npoints = 31;
351 double m_lower = mth, m_upper = max_mass, m_where_max = m0;
353 for (
size_t n_iterations = 0; n_iterations < 2; n_iterations++) {
354 const double dm = (m_upper - m_lower) / npoints;
355 for (
size_t i = 1; i < npoints; i++) {
356 m = m_lower + dm * i;
357 const double thermal_factor =
358 m * m * std::exp(-
beta * m) * gsl_sf_bessel_Kn_scaled(2, m *
beta);
366 m_lower = m_where_max - (m_where_max - m_lower) * 0.1;
367 m_upper = m_where_max + (m_upper - m_where_max) * 0.1;
376 const double thermal_factor =
377 m * m * std::exp(-
beta * m) * gsl_sf_bessel_Kn_scaled(2, m *
beta);
383 ptype.
name(),
" - maximum increased in",
384 " sample_mass_thermal from ", max_ratio,
" to ", q,
", mass = ", m,
385 " previously assumed maximum at m = ", m_where_max);
398 double mus_u = mub + T;
401 size_t iteration = 0;
403 const size_t max_iteration = 50;
405 mus = 0.5 * (mus_u + mus_l);
413 }
while (std::abs(rhos) >
tolerance_ && iteration < max_iteration);
414 if (iteration == max_iteration) {
415 throw std::runtime_error(
"Solving rho_s = 0: too many iterations.");
422 double e = reinterpret_cast<struct rparams *>(params)->e;
423 double nb = reinterpret_cast<struct rparams *>(params)->nb;
424 double ns = reinterpret_cast<struct rparams *>(params)->ns;
425 bool w = reinterpret_cast<struct rparams *>(params)->account_for_width;
427 const double T = gsl_vector_get(x, 0);
428 const double mub = gsl_vector_get(x, 1);
429 const double mus = gsl_vector_get(x, 2);
439 const double edens = reinterpret_cast<struct eparams *>(params)->edens;
447 int degeneracies_sum = 0.0;
450 degeneracies_sum += ptype.spin() + 1;
454 const double T_min = std::pow(e /
prefactor_ / 6 / degeneracies_sum, 1. / 4.);
456 const double T_max = 2.0;
458 struct eparams parameters = {e};
460 const gsl_root_fsolver_type *T = gsl_root_fsolver_brent;
461 gsl_root_fsolver *e_solver;
462 e_solver = gsl_root_fsolver_alloc(T);
463 gsl_root_fsolver_set(e_solver, &F, T_min, T_max);
465 int iter = 0, status, max_iter = 100;
470 status = gsl_root_fsolver_iterate(e_solver);
471 if (status != GSL_SUCCESS) {
474 T_init = gsl_root_fsolver_root(e_solver);
475 double x_lo = gsl_root_fsolver_x_lower(e_solver);
476 double x_hi = gsl_root_fsolver_x_upper(e_solver);
477 status = gsl_root_test_interval(x_lo, x_hi, 0.0, 0.001);
478 }
while (status == GSL_CONTINUE && iter < max_iter);
480 if (status != GSL_SUCCESS) {
481 std::stringstream err_msg;
482 err_msg <<
"Solver of equation for temperature with e = " << e
483 <<
" failed to converge. Maybe Tmax = " << T_max <<
" is too small?"
485 throw std::runtime_error(gsl_strerror(status) + err_msg.str());
488 gsl_root_fsolver_free(e_solver);
491 double n_only_baryons = 0.0;
497 const double nb_scaled = nb /
prefactor_ / (T_init * T_init * T_init);
498 double mub_init = T_init * std::asinh(nb_scaled / n_only_baryons / 2.0);
502 std::array<double, 3> initial_approximation = {T_init, mub_init, 0.0};
503 return initial_approximation;
507 double e,
double nb,
double ns,
508 std::array<double, 3> initial_approximation) {
509 int residual_status = GSL_SUCCESS;
516 gsl_vector_set(
x_, 0, initial_approximation[0]);
517 gsl_vector_set(
x_, 1, initial_approximation[1]);
518 gsl_vector_set(
x_, 2, initial_approximation[2]);
520 gsl_multiroot_fsolver_set(
solver_, &f,
x_);
523 const auto iterate_status = gsl_multiroot_fsolver_iterate(
solver_);
527 if (gsl_vector_get(
solver_->x, 0) < 0.015) {
528 return {0.0, 0.0, 0.0};
532 if (iterate_status) {
536 }
while (residual_status == GSL_CONTINUE && iter < 1000);
538 if (residual_status != GSL_SUCCESS) {
539 std::stringstream solver_parameters;
540 solver_parameters <<
"\nSolver run with "
541 <<
"e = " << e <<
", nb = " << nb <<
", ns = " << ns
542 <<
", init. approx.: " << initial_approximation[0] <<
" "
543 << initial_approximation[1] <<
" "
544 << initial_approximation[2] << std::endl;
549 return {gsl_vector_get(
solver_->x, 0), gsl_vector_get(
solver_->x, 1),
550 gsl_vector_get(
solver_->x, 2)};
556 s <<
"iter = " << iter <<
","
557 <<
" x = " << gsl_vector_get(
solver_->x, 0) <<
" "
558 << gsl_vector_get(
solver_->x, 1) <<
" "
559 << gsl_vector_get(
solver_->x, 2) <<
", "
560 <<
"f(x) = " << gsl_vector_get(
solver_->f, 0) <<
" "
561 << gsl_vector_get(
solver_->f, 1) <<
" "
562 << gsl_vector_get(
solver_->f, 2) << std::endl;