Version: SMASH-3.1
smash::HadronGasEos Class Reference

#include <hadgas_eos.h>

Class to handle the equation of state (EoS) of the hadron gas, consisting of all hadrons included in SMASH.

This implementation deals with an ideal Boltzmann gas and allows to compute:

  • energy density \(\epsilon\), pressure \(p\), density \(n\), net baryon density \(n_B\), net strangeness \(n_S\) and net charge density \(n_Q\) as a function of temperature \(T\), baryon chemical potential \(\mu_B\), strange chemical potential \(\mu_S\) and charge chemical potential \(\mu_Q\).
  • Temperature and chemical potentials given energy-, net baryon- and net strangeness density. This requires solving a system of nonlinear equations.

Definition at line 125 of file hadgas_eos.h.

Classes

struct  eparams
 Another structure for passing energy density to the gnu library. More...
 
struct  rparams
 A structure for passing equation parameters to the gnu library. More...
 

Public Member Functions

 HadronGasEos (bool tabulate, bool account_for_widths)
 Constructor of HadronGasEos. More...
 
 ~HadronGasEos ()
 
std::array< double, 4 > solve_eos (double e, double nb, double ns, double nq, std::array< double, 4 > initial_approximation)
 Compute temperature and chemical potentials given energy-, net baryon-, net strangeness- and net charge density and an inital approximation. More...
 
std::array< double, 4 > solve_eos (double e, double nb, double ns, double nq)
 Compute temperature and chemical potentials given energy-, net baryon-, net strangeness- and net charge density without an inital approximation. More...
 
std::array< double, 4 > solve_eos_initial_approximation (double e, double nb, double nq)
 Compute a reasonable initial approximation for solve_eos. More...
 
void from_table (EosTable::table_element &res, double e, double nb, double nq) const
 Get the element of eos table. More...
 
bool is_tabulated () const
 Create an EoS table or not? More...
 
bool account_for_resonance_widths () const
 If resonance spectral functions are taken into account. More...
 

Static Public Member Functions

static double energy_density (double T, double mub, double mus, double muq)
 Compute energy density. More...
 
static double density (double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute particle number density. More...
 
static double pressure (double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute pressure \( p = n T \). More...
 
static double net_baryon_density (double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute net baryon density. More...
 
static double net_strange_density (double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute net strangeness density. More...
 
static double net_charge_density (double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute net charge density. More...
 
static double partial_density (const ParticleType &ptype, double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
 Compute partial density of one hadron sort. More...
 
static double sample_mass_thermal (const ParticleType &ptype, double beta)
 Sample resonance mass in a thermal medium. More...
 
static double mus_net_strangeness0 (double T, double mub, double muq)
 Compute strangeness chemical potential, requiring that net strangeness = 0. More...
 
static bool is_eos_particle (const ParticleType &ptype)
 Check if a particle belongs to the EoS. More...
 

Private Member Functions

std::string print_solver_state (size_t iter) const
 Helpful printout, useful for debugging if gnu equation solving goes crazy. More...
 

Static Private Member Functions

static double scaled_partial_density_auxiliary (double m_over_T, double mu_over_T)
 Function used to avoid duplications in density calculations. More...
 
static double scaled_partial_density (const ParticleType &ptype, double beta, double mub, double mus, double muq, bool account_for_width=false)
 Compute (unnormalized) density of one hadron sort - helper functions used to reduce code duplication. More...
 
static int set_eos_solver_equations (const gsl_vector *x, void *params, gsl_vector *f)
 Interface EoS equations to be solved to gnu library. More...
 
static double e_equation (double T, void *params)
 

Private Attributes

EosTable eos_table_ = EosTable(1.e-1, 1.e-1, 1.e-1, 90, 90, 90)
 EOS Table to be used. More...
 
gsl_vector * x_
 Variables used by gnu equation solver. More...
 
gsl_multiroot_fsolver * solver_
 
const bool tabulate_
 Create an EoS table or not? More...
 
const bool account_for_resonance_widths_
 Use pole masses of resonances or integrate over spectral functions. More...
 

Static Private Attributes

static constexpr double prefactor_
 Constant factor, that appears in front of many thermodyn. expressions. More...
 
static constexpr double tolerance_ = 1.e-8
 Precision of equation solving. More...
 
static constexpr size_t n_equations_ = 4
 Number of equations in the system of equations to be solved. More...
 

Constructor & Destructor Documentation

◆ HadronGasEos()

smash::HadronGasEos::HadronGasEos ( bool  tabulate,
bool  account_for_widths 
)

Constructor of HadronGasEos.

Parameters
[in]tabulateWhether the equation of state should be tabulated Tabulation takes time once (typically around 5 minutes), but makes the further usage of the class much faster. Tabulated values are saved in a file and loaded at the next run.
[in]account_for_widthsWhether equation of state should account for resonance spectral functions. Normally one wants to do it, if HadronGasEos is used for density calculations, for example in the box initialization. However, it is not recommended to account for spectral functions if EoS is tabulated (tabulate = true), because this makes tabulation incredibly slow. Therefore, for HadronGasEos to be used in thermalizer, this option has to be false. Also note that presently width account is not implemented for energy density calculation.

Definition at line 197 of file hadgas_eos.cc.

198  : x_(gsl_vector_alloc(n_equations_)),
199  tabulate_(tabulate),
200  account_for_resonance_widths_(account_for_width) {
201  const gsl_multiroot_fsolver_type *solver_type;
202  solver_type = gsl_multiroot_fsolver_hybrid;
203  solver_ = gsl_multiroot_fsolver_alloc(solver_type, n_equations_);
205  logg[LResonances].error(
206  "Compilation of hadron gas EoS table requested with"
207  " account of resonance spectral functions. This is not "
208  "advised, as it will likely take a few days to finish."
209  " Besides, the effect of resonance widths is currently not "
210  "implemented for energy density computation, so the computed"
211  " table will be inconsistent anyways.");
212  }
213  if (tabulate_) {
214  eos_table_.compile_table(*this);
215  }
216 }
void compile_table(HadronGasEos &eos, const std::string &eos_savefile_name="hadgas_eos.dat")
Computes the actual content of the table (for EosTable description see documentation of the construct...
Definition: hadgas_eos.cc:34
gsl_multiroot_fsolver * solver_
Definition: hadgas_eos.h:452
EosTable eos_table_
EOS Table to be used.
Definition: hadgas_eos.h:441
gsl_vector * x_
Variables used by gnu equation solver.
Definition: hadgas_eos.h:449
const bool tabulate_
Create an EoS table or not?
Definition: hadgas_eos.h:455
const bool account_for_resonance_widths_
Use pole masses of resonances or integrate over spectral functions.
Definition: hadgas_eos.h:458
static constexpr size_t n_equations_
Number of equations in the system of equations to be solved.
Definition: hadgas_eos.h:438
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
static constexpr int LResonances

◆ ~HadronGasEos()

smash::HadronGasEos::~HadronGasEos ( )

Definition at line 218 of file hadgas_eos.cc.

218  {
219  gsl_multiroot_fsolver_free(solver_);
220  gsl_vector_free(x_);
221 }

Member Function Documentation

◆ energy_density()

double smash::HadronGasEos::energy_density ( double  T,
double  mub,
double  mus,
double  muq 
)
static

Compute energy density.

Grand-canonical Boltzmann ideal gas, consisting of all hadrons in SMASH:

\[ \epsilon = \sum \frac{g_i m_i^2 T^2}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B_i + \mu_S S_i + \mu_Q Q_i}{T} \right) \times \left[ 3 K_2\left( \frac{m_i}{T}\right) + \frac{m_i}{T} K_1\left( \frac{m_i}{T}\right)\right] \]

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
Returns
energy density e [GeV/fm \(^3\)]

Definition at line 281 of file hadgas_eos.cc.

282  {
283  if (T < really_small) {
284  return 0.0;
285  }
286  const double beta = 1.0 / T;
287  double e = 0.0;
288  for (const ParticleType &ptype : ParticleType::list_all()) {
289  if (!is_eos_particle(ptype)) {
290  continue;
291  }
292  const double z = ptype.mass() * beta;
293  double x = beta * (mub * ptype.baryon_number() + mus * ptype.strangeness() +
294  muq * ptype.charge() - ptype.mass());
295  if (x < -500.0) {
296  return 0.0;
297  }
298  x = std::exp(x);
299  const size_t g = ptype.spin() + 1;
300  // Small mass case, z*z*K_2(z) -> 2, z*z*z*K_1(z) -> 0 at z->0
301  e += (z < really_small) ? 3.0 * g * x
302  : z * z * g * x *
303  (3.0 * gsl_sf_bessel_Kn_scaled(2, z) +
304  z * gsl_sf_bessel_K1_scaled(z));
305  }
306  e *= prefactor_ * T * T * T * T;
307  return e;
308 }
static constexpr double prefactor_
Constant factor, that appears in front of many thermodyn. expressions.
Definition: hadgas_eos.h:431
static bool is_eos_particle(const ParticleType &ptype)
Check if a particle belongs to the EoS.
Definition: hadgas_eos.h:355
static const ParticleTypeList & list_all()
Definition: particletype.cc:51
T beta(T a, T b)
Draws a random number from a beta-distribution, where probability density of is .
Definition: random.h:329
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37

◆ density()

double smash::HadronGasEos::density ( double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
static

Compute particle number density.

Grand-canonical Boltzmann ideal gas, consisting of all hadrons in SMASH:

\[ n = \sum \frac{g_i m_i^2 T}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B_i + \mu_S S_i + \mu_Q Q_i}{T} \right) K_2\left( \frac{m_i}{T}\right) \]

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
particle number density n [fm \(^{-3}\)]

Definition at line 310 of file hadgas_eos.cc.

311  {
312  if (T < really_small) {
313  return 0.0;
314  }
315  const double beta = 1.0 / T;
316  double rho = 0.0;
317  for (const ParticleType &ptype : ParticleType::list_all()) {
318  if (!is_eos_particle(ptype)) {
319  continue;
320  }
321  rho +=
322  scaled_partial_density(ptype, beta, mub, mus, muq, account_for_width);
323  }
324  rho *= prefactor_ * T * T * T;
325  return rho;
326 }
static double scaled_partial_density(const ParticleType &ptype, double beta, double mub, double mus, double muq, bool account_for_width=false)
Compute (unnormalized) density of one hadron sort - helper functions used to reduce code duplication.
Definition: hadgas_eos.cc:237

◆ pressure()

static double smash::HadronGasEos::pressure ( double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
inlinestatic

Compute pressure \( p = n T \).

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
pressure p [GeV/fm \(^{-3}\)]

Definition at line 195 of file hadgas_eos.h.

196  {
197  return T * density(T, mub, mus, muq, account_for_resonance_widths);
198  }
bool account_for_resonance_widths() const
If resonance spectral functions are taken into account.
Definition: hadgas_eos.h:363
static double density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute particle number density.
Definition: hadgas_eos.cc:310

◆ net_baryon_density()

double smash::HadronGasEos::net_baryon_density ( double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
static

Compute net baryon density.

Grand-canonical Boltzmann ideal gas, consisting of all hadrons in SMASH:

\[ n_B = \sum B_i \frac{g_i m_i^2 T}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B_i + \mu_S S_i + \mu_Q Q_i}{T} \right) K_2\left( \frac{m_i}{T}\right) \]

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
net baryon density \(n_B\) [fm \(^{-3}\)]

Definition at line 328 of file hadgas_eos.cc.

329  {
330  if (T < really_small) {
331  return 0.0;
332  }
333  const double beta = 1.0 / T;
334  double rho = 0.0;
335  for (const ParticleType &ptype : ParticleType::list_all()) {
336  if (!ptype.is_baryon() || !is_eos_particle(ptype)) {
337  continue;
338  }
339  rho +=
340  scaled_partial_density(ptype, beta, mub, mus, muq, account_for_width) *
341  ptype.baryon_number();
342  }
343  rho *= prefactor_ * T * T * T;
344  return rho;
345 }

◆ net_strange_density()

double smash::HadronGasEos::net_strange_density ( double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
static

Compute net strangeness density.

Grand-canonical Boltzmann ideal gas, consisting of all hadrons in SMASH:

\[ n_S = \sum S_i \frac{g_i m_i^2 T}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B_i + \mu_S S_i + \mu_Q Q_i}{T} \right) K_2\left( \frac{m_i}{T}\right) \]

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
net strangeness density density \(n_S\) [fm \(^{-3}\)]

Definition at line 347 of file hadgas_eos.cc.

348  {
349  if (T < really_small) {
350  return 0.0;
351  }
352  const double beta = 1.0 / T;
353  double rho = 0.0;
354  for (const ParticleType &ptype : ParticleType::list_all()) {
355  if (ptype.strangeness() == 0 || !is_eos_particle(ptype)) {
356  continue;
357  }
358  rho +=
359  scaled_partial_density(ptype, beta, mub, mus, muq, account_for_width) *
360  ptype.strangeness();
361  }
362  rho *= prefactor_ * T * T * T;
363  return rho;
364 }

◆ net_charge_density()

double smash::HadronGasEos::net_charge_density ( double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
static

Compute net charge density.

Grand-canonical Boltzmann ideal gas, consisting of all hadrons in SMASH:

\[ n_Q = \sum Q_i \frac{g_i m_i^2 T}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B_i + \mu_S S_i + \mu_Q Q_i}{T} \right) K_2\left( \frac{m_i}{T}\right) \]

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
net charge density density \(n_S\) [fm \(^{-3}\)]

Definition at line 366 of file hadgas_eos.cc.

367  {
368  if (T < really_small) {
369  return 0.0;
370  }
371  const double beta = 1.0 / T;
372  double rho = 0.0;
373  for (const ParticleType &ptype : ParticleType::list_all()) {
374  if (ptype.charge() == 0 || !is_eos_particle(ptype)) {
375  continue;
376  }
377  rho +=
378  scaled_partial_density(ptype, beta, mub, mus, muq, account_for_width) *
379  ptype.charge();
380  }
381  rho *= prefactor_ * T * T * T;
382  return rho;
383 }

◆ partial_density()

double smash::HadronGasEos::partial_density ( const ParticleType ptype,
double  T,
double  mub,
double  mus,
double  muq,
bool  account_for_resonance_widths = false 
)
static

Compute partial density of one hadron sort.

Grand-canonical Boltzmann ideal gas:

\[ n = \frac{g m^2 T}{2\pi^2(\hbar c)^3} exp \left(\frac{\mu_B B + \mu_S S + \mu_Q Q_i}{T} \right) K_2\left( \frac{m}{T}\right) \]

Parameters
[in]ptypethe hadron sort, for which partial density is computed
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_resonance_widthsif false, pole masses are used; if true, then integration over spectral function is included
Returns
partial density of the given hadron sort \(n\) [fm \(^{-3}\)]

Definition at line 270 of file hadgas_eos.cc.

272  {
273  if (T < really_small) {
274  return 0.0;
275  }
276  return prefactor_ * T * T * T *
277  scaled_partial_density(ptype, 1.0 / T, mub, mus, muq,
278  account_for_width);
279 }

◆ sample_mass_thermal()

double smash::HadronGasEos::sample_mass_thermal ( const ParticleType ptype,
double  beta 
)
static

Sample resonance mass in a thermal medium.

Samples mass from the distribution

\[ dN/dm \sim A(m) m^2 K_2\left( \frac{m}{T}\right) \]

For stable particles always returns pole mass.

Parameters
[in]ptypethe hadron sort, for which mass is sampled
[in]betainverse temperature 1/T [1/GeV]
Returns
sampled mass

Definition at line 385 of file hadgas_eos.cc.

386  {
387  if (ptype.is_stable()) {
388  return ptype.mass();
389  }
390  // Sampling mass m from A(m) x^2 BesselK_2(x), where x = beta m.
391  // Strategy employs the idea of importance sampling:
392  // -- Sample mass from the simple Breit-Wigner first, then
393  // reject by the ratio.
394  // -- Note that f(x) = x^2 BesselK_2(x) monotonously decreases
395  // and has maximum at x = xmin. That is why instead of f(x)
396  // the ratio f(x)/f(xmin) is used.
397 
398  const double max_mass = 5.0; // GeV
399  double m, q;
400  {
401  // Allow underflows in exponentials
402  DisableFloatTraps guard(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
403  const double w0 = ptype.width_at_pole();
404  const double mth = ptype.min_mass_spectral();
405  const double m0 = ptype.mass();
406  double max_ratio =
407  m0 * m0 * std::exp(-beta * m0) * gsl_sf_bessel_Kn_scaled(2, m0 * beta) *
408  ptype.spectral_function(m0) / ptype.spectral_function_simple(m0);
409  // Heuristic adaptive maximum search to find max_ratio
410  constexpr int npoints = 31;
411  double m_lower = mth, m_upper = max_mass, m_where_max = m0;
412 
413  for (size_t n_iterations = 0; n_iterations < 2; n_iterations++) {
414  const double dm = (m_upper - m_lower) / npoints;
415  for (size_t i = 1; i < npoints; i++) {
416  m = m_lower + dm * i;
417  const double thermal_factor =
418  m * m * std::exp(-beta * m) * gsl_sf_bessel_Kn_scaled(2, m * beta);
419  q = ptype.spectral_function(m) * thermal_factor /
420  ptype.spectral_function_simple(m);
421  if (q > max_ratio) {
422  max_ratio = q;
423  m_where_max = m;
424  }
425  }
426  m_lower = m_where_max - (m_where_max - m_lower) * 0.1;
427  m_upper = m_where_max + (m_upper - m_where_max) * 0.1;
428  }
429  // Safety factor
430  max_ratio *= 1.5;
431 
432  do {
433  // sample mass from A(m)
434  do {
435  m = random::cauchy(m0, 0.5 * w0, mth, max_mass);
436  const double thermal_factor =
437  m * m * std::exp(-beta * m) * gsl_sf_bessel_Kn_scaled(2, m * beta);
438  q = ptype.spectral_function(m) * thermal_factor /
439  ptype.spectral_function_simple(m);
440  } while (q < random::uniform(0., max_ratio));
441  if (q > max_ratio) {
442  logg[LResonances].warn(
443  ptype.name(), " - maximum increased in",
444  " sample_mass_thermal from ", max_ratio, " to ", q, ", mass = ", m,
445  " previously assumed maximum at m = ", m_where_max);
446  max_ratio = q;
447  m_where_max = m;
448  } else {
449  break;
450  }
451  } while (true);
452  }
453  return m;
454 }
T uniform(T min, T max)
Definition: random.h:88
T cauchy(T pole, T width, T min, T max)
Draws a random number from a Cauchy distribution (sometimes also called Lorentz or non-relativistic B...
Definition: random.h:307

◆ solve_eos() [1/2]

std::array< double, 4 > smash::HadronGasEos::solve_eos ( double  e,
double  nb,
double  ns,
double  nq,
std::array< double, 4 >  initial_approximation 
)

Compute temperature and chemical potentials given energy-, net baryon-, net strangeness- and net charge density and an inital approximation.

Parameters
[in]eenergy density [GeV/fm \(^3\)]
[in]nbnet baryon density [fm \(^{-3}\)]
[in]nsnet strangeness density [fm \(^{-3}\)]
[in]nqnet charge density [fm \(^{-3}\)]
[in]initial_approximation(T [GeV], mub [GeV], mus [GeV]) to use as starting point
Returns
array of 4 values: temperature, baryon chemical potential, strange chemical potential and charge chemical potential

Definition at line 586 of file hadgas_eos.cc.

588  {
589  int residual_status = GSL_SUCCESS;
590  size_t iter = 0;
591 
592  struct rparams p = {e, nb, ns, nq, account_for_resonance_widths_};
593  gsl_multiroot_function f = {&HadronGasEos::set_eos_solver_equations,
594  n_equations_, &p};
595 
596  gsl_vector_set(x_, 0, initial_approximation[0]);
597  gsl_vector_set(x_, 1, initial_approximation[1]);
598  gsl_vector_set(x_, 2, initial_approximation[2]);
599  gsl_vector_set(x_, 3, initial_approximation[3]);
600 
601  gsl_multiroot_fsolver_set(solver_, &f, x_);
602  do {
603  iter++;
604  const auto iterate_status = gsl_multiroot_fsolver_iterate(solver_);
605  // std::cout << print_solver_state(iter);
606 
607  // Avoiding too low temperature
608  if (gsl_vector_get(solver_->x, 0) < 0.015) {
609  return {0.0, 0.0, 0.0, 0.0};
610  }
611 
612  // check if solver is stuck
613  if (iterate_status) {
614  break;
615  }
616  residual_status = gsl_multiroot_test_residual(solver_->f, tolerance_);
617  } while (residual_status == GSL_CONTINUE && iter < 1000);
618 
619  if (residual_status != GSL_SUCCESS) {
620  std::stringstream solver_parameters;
621  solver_parameters << "\nSolver run with "
622  << "e = " << e << ", nb = " << nb << ", ns = " << ns
623  << ", nq = " << nq
624  << ", init. approx.: " << initial_approximation[0] << " "
625  << initial_approximation[1] << " "
626  << initial_approximation[2] << " "
627  << initial_approximation[3] << std::endl;
628  logg[LResonances].warn(gsl_strerror(residual_status) +
629  solver_parameters.str() + print_solver_state(iter));
630  }
631 
632  return {gsl_vector_get(solver_->x, 0), gsl_vector_get(solver_->x, 1),
633  gsl_vector_get(solver_->x, 2), gsl_vector_get(solver_->x, 3)};
634 }
static int set_eos_solver_equations(const gsl_vector *x, void *params, gsl_vector *f)
Interface EoS equations to be solved to gnu library.
Definition: hadgas_eos.cc:480
std::string print_solver_state(size_t iter) const
Helpful printout, useful for debugging if gnu equation solving goes crazy.
Definition: hadgas_eos.cc:636
static constexpr double tolerance_
Precision of equation solving.
Definition: hadgas_eos.h:435
constexpr int p
Proton.

◆ solve_eos() [2/2]

std::array<double, 4> smash::HadronGasEos::solve_eos ( double  e,
double  nb,
double  ns,
double  nq 
)
inline

Compute temperature and chemical potentials given energy-, net baryon-, net strangeness- and net charge density without an inital approximation.

Parameters
[in]eenergy density [GeV/fm \(^3\)]
[in]nbnet baryon density [fm \(^{-3}\)]
[in]nsnet strangeness density [fm \(^{-3}\)]
[in]nqnet charge density [fm \(^{-3}\)]
Returns
array of 4 values: temperature, baryon chemical potential and strange chemical potential and charge

Definition at line 322 of file hadgas_eos.h.

322  {
323  return solve_eos(e, nb, ns, nq, solve_eos_initial_approximation(e, nb, nq));
324  }
std::array< double, 4 > solve_eos(double e, double nb, double ns, double nq, std::array< double, 4 > initial_approximation)
Compute temperature and chemical potentials given energy-, net baryon-, net strangeness- and net char...
Definition: hadgas_eos.cc:586
std::array< double, 4 > solve_eos_initial_approximation(double e, double nb, double nq)
Compute a reasonable initial approximation for solve_eos.
Definition: hadgas_eos.cc:506

◆ solve_eos_initial_approximation()

std::array< double, 4 > smash::HadronGasEos::solve_eos_initial_approximation ( double  e,
double  nb,
double  nq 
)

Compute a reasonable initial approximation for solve_eos.

Parameters
[in]eenergy density [GeV/fm \(^3\)]
[in]nbnet baryon density [fm \(^{-3}\)]
[in]nqnet charge density [fm \(^{-3}\)]
Returns
array of 3 values: temperature, baryon chemical potential and strange chemical potential

Definition at line 506 of file hadgas_eos.cc.

508  {
509  assert(e >= 0.0);
510  // 1. Get temperature from energy density assuming zero chemical potentials
511  int degeneracies_sum = 0.0;
512  for (const ParticleType &ptype : ParticleType::list_all()) {
513  if (is_eos_particle(ptype)) {
514  degeneracies_sum += ptype.spin() + 1;
515  }
516  }
517  // Temperature in case of massless gas. For massive it should be larger.
518  const double T_min = std::pow(e / prefactor_ / 6 / degeneracies_sum, 1. / 4.);
519  // Simply assume that the temperature is not higher than 2 GeV.
520  const double T_max = 2.0;
521 
522  struct eparams parameters = {e};
523  gsl_function F = {&e_equation, &parameters};
524  const gsl_root_fsolver_type *T = gsl_root_fsolver_brent;
525  gsl_root_fsolver *e_solver;
526  e_solver = gsl_root_fsolver_alloc(T);
527  gsl_root_fsolver_set(e_solver, &F, T_min, T_max);
528 
529  int iter = 0, status, max_iter = 100;
530  double T_init = 0.0;
531 
532  do {
533  iter++;
534  status = gsl_root_fsolver_iterate(e_solver);
535  if (status != GSL_SUCCESS) {
536  break;
537  }
538  T_init = gsl_root_fsolver_root(e_solver);
539  double x_lo = gsl_root_fsolver_x_lower(e_solver);
540  double x_hi = gsl_root_fsolver_x_upper(e_solver);
541  status = gsl_root_test_interval(x_lo, x_hi, 0.0, 0.001);
542  } while (status == GSL_CONTINUE && iter < max_iter);
543 
544  if (status != GSL_SUCCESS) {
545  std::stringstream err_msg;
546  err_msg << "Solver of equation for temperature with e = " << e
547  << " failed to converge. Maybe Tmax = " << T_max << " is too small?"
548  << std::endl;
549  throw std::runtime_error(gsl_strerror(status) + err_msg.str());
550  }
551 
552  gsl_root_fsolver_free(e_solver);
553 
554  // 2. Get the baryon chemical potential for muS = muQ = 0 with previously
555  // obtained T
556  double n_only_baryons = 0.0;
557  for (const ParticleType &ptype : ParticleType::list_all()) {
558  if (is_eos_particle(ptype) && ptype.baryon_number() == 1) {
559  n_only_baryons +=
560  scaled_partial_density(ptype, 1.0 / T_init, 0.0, 0.0, 0.0);
561  }
562  }
563  const double nb_scaled = nb / prefactor_ / (T_init * T_init * T_init);
564  double mub_init = T_init * std::asinh(nb_scaled / n_only_baryons / 2.0);
565 
566  // 3. Get the charge chemical potential assuming muB = muS = 0 with previously
567  // obtained T
568  double n_only_charge_1_particles = 0.0;
569  for (const ParticleType &ptype : ParticleType::list_all()) {
570  if (is_eos_particle(ptype) && ptype.charge() == 1) {
571  n_only_charge_1_particles +=
572  scaled_partial_density(ptype, 1.0 / T_init, 0.0, 0.0, 0.0);
573  }
574  }
575  const double q_scaled = nq / prefactor_ / (T_init * T_init * T_init);
576  double muq_init =
577  T_init * std::asinh(q_scaled / n_only_charge_1_particles / 2.0);
578 
579  // 4. Get the strange chemical potential, where mus = 0 is typically a good
580  // initial approximation
581  std::array<double, 4> initial_approximation = {T_init, mub_init, 0.0,
582  muq_init};
583  return initial_approximation;
584 }
static double e_equation(double T, void *params)
Definition: hadgas_eos.cc:501

◆ mus_net_strangeness0()

double smash::HadronGasEos::mus_net_strangeness0 ( double  T,
double  mub,
double  muq 
)
static

Compute strangeness chemical potential, requiring that net strangeness = 0.

Parameters
[in]Ttemperature [GeV]
[in]mubbaryon chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
Returns
strangeness chemical potential [GeV]

Definition at line 456 of file hadgas_eos.cc.

456  {
457  // Binary search
458  double mus_u = mub + T;
459  double mus_l = 0.0;
460  double mus, rhos;
461  size_t iteration = 0;
462  // 50 iterations should give precision 2^-50 ~ 10^-15
463  const size_t max_iteration = 50;
464  do {
465  mus = 0.5 * (mus_u + mus_l);
466  rhos = net_strange_density(T, mub, mus, muq);
467  if (rhos > 0.0) {
468  mus_u = mus;
469  } else {
470  mus_l = mus;
471  }
472  iteration++;
473  } while (std::abs(rhos) > tolerance_ && iteration < max_iteration);
474  if (iteration == max_iteration) {
475  throw std::runtime_error("Solving rho_s = 0: too many iterations.");
476  }
477  return mus;
478 }
static double net_strange_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net strangeness density.
Definition: hadgas_eos.cc:347

◆ from_table()

void smash::HadronGasEos::from_table ( EosTable::table_element res,
double  e,
double  nb,
double  nq 
) const
inline

Get the element of eos table.

Definition at line 349 of file hadgas_eos.h.

350  {
351  eos_table_.get(res, e, nb, nq);
352  }
void get(table_element &res, double e, double nb, double nq) const
Obtain interpolated p/T/muB/muS/muQ from the tabulated equation of state given energy density,...
Definition: hadgas_eos.cc:162

◆ is_eos_particle()

static bool smash::HadronGasEos::is_eos_particle ( const ParticleType ptype)
inlinestatic

Check if a particle belongs to the EoS.

Definition at line 355 of file hadgas_eos.h.

355  {
356  return ptype.is_hadron() && ptype.pdgcode().charmness() == 0;
357  }

◆ is_tabulated()

bool smash::HadronGasEos::is_tabulated ( ) const
inline

Create an EoS table or not?

Definition at line 360 of file hadgas_eos.h.

360 { return tabulate_; }

◆ account_for_resonance_widths()

bool smash::HadronGasEos::account_for_resonance_widths ( ) const
inline

If resonance spectral functions are taken into account.

Definition at line 363 of file hadgas_eos.h.

363  {
365  }

◆ scaled_partial_density_auxiliary()

double smash::HadronGasEos::scaled_partial_density_auxiliary ( double  m_over_T,
double  mu_over_T 
)
staticprivate

Function used to avoid duplications in density calculations.

Parameters
[in]m_over_Tmass to temperature ratio \( m/T \)
[in]mu_over_Tchemical potential to temperature ratio \( \mu/T \)
Returns
calculated \( (m/T)^2 exp(\mu/T) K_2(m/T) \)

Definition at line 223 of file hadgas_eos.cc.

224  {
225  double x = mu_over_T - m_over_T;
226  if (x < -500.0) {
227  return 0.0;
228  }
229  x = std::exp(x);
230  // In the case of small masses: K_n(z) -> (n-1)!/2 *(2/z)^n, z -> 0,
231  // z*z*K_2(z) -> 2
232  return (m_over_T < really_small)
233  ? 2.0 * x
234  : m_over_T * m_over_T * x * gsl_sf_bessel_Kn_scaled(2, m_over_T);
235 }

◆ scaled_partial_density()

double smash::HadronGasEos::scaled_partial_density ( const ParticleType ptype,
double  beta,
double  mub,
double  mus,
double  muq,
bool  account_for_width = false 
)
staticprivate

Compute (unnormalized) density of one hadron sort - helper functions used to reduce code duplication.

Parameters
[in]ptypethe hadron sort, for which partial density is computed
[in]betainverse temperature [1/GeV]
[in]mubbaryon chemical potential [GeV]
[in]musstrangeness chemical potential [GeV]
[in]muqcharge chemical potential [GeV]
[in]account_for_widthTake hadron spectral functions into account or not. When taken into account, they result in a considerable slow down.
Returns
partial (unnormalized) density of the given hadron sort \(n\) [fm \(^{-3}\)]

Definition at line 237 of file hadgas_eos.cc.

240  {
241  const double m_over_T = ptype.mass() * beta;
242  double mu_over_T = beta * (ptype.baryon_number() * mub +
243  ptype.strangeness() * mus + ptype.charge() * muq);
244  const double g = ptype.spin() + 1;
245  if (ptype.is_stable() || !account_for_width) {
246  return g * scaled_partial_density_auxiliary(m_over_T, mu_over_T);
247  } else {
248  // Integral \int_{threshold}^{\infty} A(m) N_{thermal}(m) dm,
249  // where A(m) is the spectral function of the resonance.
250  const double m0 = ptype.mass();
251  const double w0 = ptype.width_at_pole();
252  const double mth = ptype.min_mass_spectral();
253  const double u_min = std::atan(2.0 * (mth - m0) / w0);
254  const double u_max = 0.5 * M_PI;
255  Integrator integrate;
256  const double result =
257  g * integrate(u_min, u_max, [&](double u) {
258  // One of many possible variable substitutions. Not clear if it has
259  // any advantages, except transforming (m_th, inf) to finite interval.
260  const double tanu = std::tan(u);
261  const double m = m0 + 0.5 * w0 * tanu;
262  const double jacobian = 0.5 * w0 * (1.0 + tanu * tanu);
263  return ptype.spectral_function(m) * jacobian *
264  scaled_partial_density_auxiliary(m * beta, mu_over_T);
265  });
266  return result;
267  }
268 }
static double scaled_partial_density_auxiliary(double m_over_T, double mu_over_T)
Function used to avoid duplications in density calculations.
Definition: hadgas_eos.cc:223
static Integrator integrate
Definition: decaytype.cc:143

◆ set_eos_solver_equations()

int smash::HadronGasEos::set_eos_solver_equations ( const gsl_vector *  x,
void *  params,
gsl_vector *  f 
)
staticprivate

Interface EoS equations to be solved to gnu library.

Definition at line 480 of file hadgas_eos.cc.

481  {
482  double e = reinterpret_cast<struct rparams *>(params)->e;
483  double nb = reinterpret_cast<struct rparams *>(params)->nb;
484  double ns = reinterpret_cast<struct rparams *>(params)->ns;
485  double nq = reinterpret_cast<struct rparams *>(params)->nq;
486  bool w = reinterpret_cast<struct rparams *>(params)->account_for_width;
487 
488  const double T = gsl_vector_get(x, 0);
489  const double mub = gsl_vector_get(x, 1);
490  const double mus = gsl_vector_get(x, 2);
491  const double muq = gsl_vector_get(x, 3);
492 
493  gsl_vector_set(f, 0, energy_density(T, mub, mus, muq) - e);
494  gsl_vector_set(f, 1, net_baryon_density(T, mub, mus, muq, w) - nb);
495  gsl_vector_set(f, 2, net_strange_density(T, mub, mus, muq, w) - ns);
496  gsl_vector_set(f, 3, net_charge_density(T, mub, mus, muq, w) - nq);
497 
498  return GSL_SUCCESS;
499 }
static double net_charge_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net charge density.
Definition: hadgas_eos.cc:366
static double net_baryon_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net baryon density.
Definition: hadgas_eos.cc:328
static double energy_density(double T, double mub, double mus, double muq)
Compute energy density.
Definition: hadgas_eos.cc:281

◆ e_equation()

double smash::HadronGasEos::e_equation ( double  T,
void *  params 
)
staticprivate
See also
set_eos_solver_equations()

Definition at line 501 of file hadgas_eos.cc.

501  {
502  const double edens = reinterpret_cast<struct eparams *>(params)->edens;
503  return edens - energy_density(T, 0.0, 0.0, 0.0);
504 }

◆ print_solver_state()

std::string smash::HadronGasEos::print_solver_state ( size_t  iter) const
private

Helpful printout, useful for debugging if gnu equation solving goes crazy.

Parameters
[in]itercurrent value of iterator
Returns
debug output string with iter, x and f(x) from solver

Definition at line 636 of file hadgas_eos.cc.

636  {
637  std::stringstream s;
638  // clang-format off
639  s << "iter = " << iter << ","
640  << " x = " << gsl_vector_get(solver_->x, 0) << " "
641  << gsl_vector_get(solver_->x, 1) << " "
642  << gsl_vector_get(solver_->x, 2) << ", "
643  << gsl_vector_get(solver_->x, 3) << ", "
644  << "f(x) = " << gsl_vector_get(solver_->f, 0) << " "
645  << gsl_vector_get(solver_->f, 1) << " "
646  << gsl_vector_get(solver_->f, 2) << " "
647  << gsl_vector_get(solver_->f, 3) << std::endl;
648  // clang-format on
649  return s.str();
650 }

Member Data Documentation

◆ prefactor_

constexpr double smash::HadronGasEos::prefactor_
staticconstexprprivate
Initial value:
=
0.5 * M_1_PI * M_1_PI / (hbarc * hbarc * hbarc)
constexpr double hbarc
GeV <-> fm conversion factor.
Definition: constants.h:25

Constant factor, that appears in front of many thermodyn. expressions.

Definition at line 431 of file hadgas_eos.h.

◆ tolerance_

constexpr double smash::HadronGasEos::tolerance_ = 1.e-8
staticconstexprprivate

Precision of equation solving.

Definition at line 435 of file hadgas_eos.h.

◆ n_equations_

constexpr size_t smash::HadronGasEos::n_equations_ = 4
staticconstexprprivate

Number of equations in the system of equations to be solved.

Definition at line 438 of file hadgas_eos.h.

◆ eos_table_

EosTable smash::HadronGasEos::eos_table_ = EosTable(1.e-1, 1.e-1, 1.e-1, 90, 90, 90)
private

EOS Table to be used.

Definition at line 441 of file hadgas_eos.h.

◆ x_

gsl_vector* smash::HadronGasEos::x_
private

Variables used by gnu equation solver.

They are stored here to allocate and deallocate memory for them only once. It is expected that this class will be used for solving the EoS many times, so multiple allocations and frees are unwanted.

Definition at line 449 of file hadgas_eos.h.

◆ solver_

gsl_multiroot_fsolver* smash::HadronGasEos::solver_
private
See also
x_

Definition at line 452 of file hadgas_eos.h.

◆ tabulate_

const bool smash::HadronGasEos::tabulate_
private

Create an EoS table or not?

Definition at line 455 of file hadgas_eos.h.

◆ account_for_resonance_widths_

const bool smash::HadronGasEos::account_for_resonance_widths_
private

Use pole masses of resonances or integrate over spectral functions.

Definition at line 458 of file hadgas_eos.h.


The documentation for this class was generated from the following files: