Version: SMASH-3.1
smash::QuantumSampling Class Reference

#include <quantumsampling.h>

This class:

  • Calculates chemical potentials given density of particle species
  • Calculates maxima of a Juttner distribution for these chemical potentials
  • Samples Juttner distribution. This is the main intent of this class, while previous points are auxiliary calculations for it.

Definition at line 31 of file quantumsampling.h.

Classes

struct  ParametersForMaximumMomentumRootFinder
 Struct object that holds the parameters relevant to finding the momentum for which the maximum of the distribution occurs. More...
 

Public Member Functions

 QuantumSampling (const std::map< PdgCode, int > &initial_multiplicities, double volume, double temperature)
 Constructor of a QuantumSampling object. More...
 
double sample (const PdgCode pdg)
 Sampling radial momenta of given particle species from Boltzmann, Bose, or Fermi distribution. More...
 

Static Public Member Functions

static double p_max_root_equation (double p, double mass, double temperature, double effective_chemical_potential, double statistics)
 Root equation for finding the radial momentum at which the Juttner distribution function has its maximum. More...
 
static int p_max_root_equation_for_GSL (const gsl_vector *roots_array, void *parameters, gsl_vector *function)
 Root equation for finding the radial momentum at which the Juttner distribution function has its maximum, suited for the GSL root finding procedure. More...
 
static void print_state_p_max (unsigned int iter, gsl_multiroot_fsolver *solver)
 A GSL utility which allows for printing out the status of the solver during the root finding procedure. More...
 
static int find_p_at_maximum_of_the_distribution (double mass, double temperature, double effective_chemical_potential, double statistics, double p_max_initial_guess, double solution_precision, double *p_max)
 A GSL root solver for finding the radial momentum value at which the maximum of the given Juttner distribution function occurs. More...
 
static double maximum_of_the_distribution (double mass, double temperature, double effective_chemical_potential, double statistics, double solution_precision)
 A convenience wrapper for finding the maximum value of the Juttner distribution, returning the value of the distribution for the momentum at which the maximum occurs, identified by find_p_at_maximum_of_the_distribution(). More...
 

Private Attributes

std::map< PdgCode, double > effective_chemical_potentials_
 Tabulated effective chemical potentials for every particle species. More...
 
std::map< PdgCode, double > distribution_function_maximums_
 Tabulated distribution function maxima for every particle species. More...
 
const double volume_
 Volume [fm^3] in which particles sre sampled. More...
 
const double temperature_
 Temperature [GeV]. More...
 

Constructor & Destructor Documentation

◆ QuantumSampling()

smash::QuantumSampling::QuantumSampling ( const std::map< PdgCode, int > &  initial_multiplicities,
double  volume,
double  temperature 
)

Constructor of a QuantumSampling object.

Parameters
[in]initial_multiplicitiesa map of pdg codes of samples particle species and corresponding multiplicities
[in]volumevolume V in which the particles are sampled [fm^3], needed to calculate the density of the species
[in]temperaturetemperature T of the system [GeV]

Definition at line 211 of file quantumsampling.cc.

214  : volume_(volume), temperature_(temperature) {
215  /*
216  * This is the precision which we expect from the solution; note that
217  * solution precision also goes into the precision of calculating the
218  * integrals involved etc. Recommended precision is at least 1e-7.
219  */
220  constexpr double solution_precision = 1e-8;
221 
222  for (const auto &pdg_and_mult : initial_multiplicities) {
223  const PdgCode pdg = pdg_and_mult.first;
224  const int number_of_particles = pdg_and_mult.second;
225  const double V_in_GeV = volume_ / (hbarc * hbarc * hbarc);
226  const double number_density = number_of_particles / V_in_GeV;
227  const double spin_degeneracy = pdg.spin_degeneracy();
228  // '+' for fermions, '-' for bosons
229  const double quantum_statistics = (pdg.spin() % 2 == 0) ? -1.0 : 1.0;
230  const ParticleType &ptype = ParticleType::find(pdg);
231  const double particle_mass = ptype.mass();
232  double chemical_potential = 0.0;
233  ChemicalPotentialSolver mu_solver;
234  // Calling the wrapper for the GSL chemical potential finder
235  chemical_potential = mu_solver.effective_chemical_potential(
236  spin_degeneracy, particle_mass, number_density, temperature_,
237  quantum_statistics, solution_precision);
238  effective_chemical_potentials_[pdg] = chemical_potential;
239  const double distribution_function_maximum = maximum_of_the_distribution(
240  particle_mass, temperature_, chemical_potential, quantum_statistics,
241  solution_precision);
242  distribution_function_maximums_[pdg] = distribution_function_maximum;
243  }
244 }
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
Definition: particletype.cc:99
const double temperature_
Temperature [GeV].
std::map< PdgCode, double > distribution_function_maximums_
Tabulated distribution function maxima for every particle species.
static double maximum_of_the_distribution(double mass, double temperature, double effective_chemical_potential, double statistics, double solution_precision)
A convenience wrapper for finding the maximum value of the Juttner distribution, returning the value ...
std::map< PdgCode, double > effective_chemical_potentials_
Tabulated effective chemical potentials for every particle species.
const double volume_
Volume [fm^3] in which particles sre sampled.
constexpr double hbarc
GeV <-> fm conversion factor.
Definition: constants.h:25

Member Function Documentation

◆ p_max_root_equation()

double smash::QuantumSampling::p_max_root_equation ( double  p,
double  mass,
double  temperature,
double  effective_chemical_potential,
double  statistics 
)
static

Root equation for finding the radial momentum at which the Juttner distribution function has its maximum.

Parameters
[in]pradial momentum, i.e., length of the momentum vector [GeV]
[in]mass(pole) mass m of the particle species [GeV]
[in]temperaturetemperature T of the system [GeV]
[in]effective_chemical_potentialeffective chemical potential mu of the system [GeV]
[in]statisticsquantum statistics of the particles species (+1 for Fermi, -1 for Bose, 0 for Boltzmann)
Returns
the extremum equation for the maximum of the Juttner distribution

Definition at line 34 of file quantumsampling.cc.

37  {
38  const double Ekin = std::sqrt(p * p + mass * mass);
39  const double term1 =
40  2 * (1 + statistics * std::exp(-(Ekin - effective_chemical_potential) /
41  temperature));
42  const double term2 = (p * p) / (temperature * Ekin);
43 
44  return term1 - term2;
45 }
constexpr int p
Proton.

◆ p_max_root_equation_for_GSL()

int smash::QuantumSampling::p_max_root_equation_for_GSL ( const gsl_vector *  roots_array,
void *  parameters,
gsl_vector *  function 
)
static

Root equation for finding the radial momentum at which the Juttner distribution function has its maximum, suited for the GSL root finding procedure.

Parameters
[in]roots_arrayan array holding the current best estimate of the roots of the solved equation
[in]parametersrefers to the parameters as provided in the GSL root solving procedure
[in]functionrefers to the root equation(s) as provided in the GSL root solving procedure (where it's called "function")
Returns
the root equation suited for GSL root finding procedure

Definition at line 47 of file quantumsampling.cc.

49  {
50  struct ParametersForMaximumMomentumRootFinder *par =
51  static_cast<struct ParametersForMaximumMomentumRootFinder *>(parameters);
52 
53  const double mass = (par->mass);
54  const double temperature = (par->temperature);
55  const double effective_chemical_potential =
56  (par->effective_chemical_potential);
57  const double statistics = (par->statistics);
58 
59  const double p_radial = gsl_vector_get(roots_array, 0);
60 
61  gsl_vector_set(function, 0,
62  p_max_root_equation(p_radial, mass, temperature,
63  effective_chemical_potential, statistics));
64 
65  return GSL_SUCCESS;
66 }
static double p_max_root_equation(double p, double mass, double temperature, double effective_chemical_potential, double statistics)
Root equation for finding the radial momentum at which the Juttner distribution function has its maxi...

◆ print_state_p_max()

void smash::QuantumSampling::print_state_p_max ( unsigned int  iter,
gsl_multiroot_fsolver *  solver 
)
static

A GSL utility which allows for printing out the status of the solver during the root finding procedure.

Parameters
[in]itervariable keeping track of how many steps in the root solving procedure have been taken
[in]solverGSL solver object, which has acces to the current best estimate of the roots and the corresponding function values
Returns
message about the current state of the solver

Definition at line 68 of file quantumsampling.cc.

69  {
70  std::printf(
71  "\n***\nfind_p_at_maximum_of_the_distribution(): iter = %3u \t"
72  "x = % .3f \t"
73  "f(x) = % .3e \n",
74  iter, gsl_vector_get(solver->x, 0), gsl_vector_get(solver->f, 0));
75 }

◆ find_p_at_maximum_of_the_distribution()

int smash::QuantumSampling::find_p_at_maximum_of_the_distribution ( double  mass,
double  temperature,
double  effective_chemical_potential,
double  statistics,
double  p_max_initial_guess,
double  solution_precision,
double *  p_max 
)
static

A GSL root solver for finding the radial momentum value at which the maximum of the given Juttner distribution function occurs.

For the value of the distribution at the maximum, one shoud use the function maximum_of_the_distribution().

Parameters
[in]mass(pole) mass m of the particle species [GeV]
[in]temperaturetemperature T of the system [GeV]
[in]effective_chemical_potentialeffective chemical potential mu of the system [GeV]
[in]statisticsquantum statistics of the particles species (+1 for Fermi, -1 for Bose, 0 for Boltzmann)
[in]p_max_initial_guessthe initial guess for the value of the solution [GeV]
[in]solution_precisionprecision with which the solution is found
[out]p_maxthe solution (momentum for which the distribution takes on the maximum value) stored in an array object [GeV]

Definition at line 77 of file quantumsampling.cc.

80  {
81  const gsl_multiroot_fsolver_type *Solver_name;
82  gsl_multiroot_fsolver *Root_finder;
83 
84  int status;
85  size_t iter = 0;
86  size_t initial_guess_update = 0;
87 
88  const size_t problem_dimension = 1;
89 
90  struct ParametersForMaximumMomentumRootFinder parameters = {
91  mass, temperature, effective_chemical_potential, statistics};
92 
93  gsl_multiroot_function MaximumOfDistribution = {
94  &p_max_root_equation_for_GSL, problem_dimension, &parameters};
95 
96  double roots_array_initial[1] = {p_max_initial_guess};
97 
98  gsl_vector *roots_array = gsl_vector_alloc(problem_dimension);
99  gsl_vector_set(roots_array, 0, roots_array_initial[0]);
100 
101  Solver_name = gsl_multiroot_fsolver_hybrids;
102  Root_finder = gsl_multiroot_fsolver_alloc(Solver_name, problem_dimension);
103  gsl_multiroot_fsolver_set(Root_finder, &MaximumOfDistribution, roots_array);
104 
105  // print_state_p_max (iter, Root_finder);
106 
107  do {
108  iter++;
109 
110  /*
111  * gsl_multiroot_fsolver_iterate returns either 0 for a correct behavior,
112  * or an error code (a positive integer) when the solver is stuck
113  */
114  status = gsl_multiroot_fsolver_iterate(Root_finder);
115 
116  // print_state_p_max (iter, Root_finder);
117 
118  /*
119  * Check whether the solver is stuck
120  */
121  if (status) {
122  if (initial_guess_update < 100) {
123  /*
124  * In case the solution is not found for the (somewhat small) default
125  * initial guess of p_max_initial_guess = 0.05 [GeV], the value of the
126  * p_max_initial_guess is increased and the root solving procedure is
127  * restarted. This can take place up to a 100 times, with the largest
128  * possible initial guess of p_max_initial_guess = 5.05 [GeV].
129  */
130  p_max_initial_guess += 0.05;
131  initial_guess_update++;
132  roots_array_initial[0] = p_max_initial_guess;
133  gsl_vector_set(roots_array, 0, roots_array_initial[0]);
134  gsl_multiroot_fsolver_set(Root_finder, &MaximumOfDistribution,
135  roots_array);
136  iter = 0;
137  } else {
138  std::cout << "\n\nGSL error message:\n"
139  << gsl_strerror(status) << std::endl;
140  logg[LogArea::Distributions::id].warn(
141  "\n\nThe GSL solver"
142  "\nfind_p_at_maximum_of_the_distribution\nis stuck!"
143  "\n\nInput parameters:"
144  "\n mass [GeV] = ",
145  mass, "\n temperature [GeV] = ", temperature,
146  "\neffective_chemical_potential = ", effective_chemical_potential,
147  "\n statistics = ", statistics,
148  "\n solution_precision = ", solution_precision,
149  "\n\n"
150  "Initialization cannot sample the momenta without "
151  "calculating the distribution maximum."
152  "\nTry adjusting the initial guess (which is "
153  "looped over in the GSL procedure) or the "
154  "solution precision."
155  "\nUncomment print_state_p_max to check solver progress.\n\n\n");
156  throw std::runtime_error(
157  "QuantumSampling::find_p_at_maximum_of_the_distribution returned "
158  "no result.\n\n");
159  continue;
160  }
161  }
162 
163  status = gsl_multiroot_test_residual(Root_finder->f, solution_precision);
164 
165  if (status == GSL_SUCCESS) {
166  p_max[0] = gsl_vector_get(Root_finder->x, 0);
167  }
168  } while (status == GSL_CONTINUE && iter < 100000);
169 
170  gsl_multiroot_fsolver_free(Root_finder);
171  gsl_vector_free(roots_array);
172 
173  return 0;
174 }
static int p_max_root_equation_for_GSL(const gsl_vector *roots_array, void *parameters, gsl_vector *function)
Root equation for finding the radial momentum at which the Juttner distribution function has its maxi...
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39

◆ maximum_of_the_distribution()

double smash::QuantumSampling::maximum_of_the_distribution ( double  mass,
double  temperature,
double  effective_chemical_potential,
double  statistics,
double  solution_precision 
)
static

A convenience wrapper for finding the maximum value of the Juttner distribution, returning the value of the distribution for the momentum at which the maximum occurs, identified by find_p_at_maximum_of_the_distribution().

Parameters
[in]mass(pole) mass m of the particle species [GeV]
[in]temperaturetemperature T of the system [GeV]
[in]effective_chemical_potentialeffective chemical potential mu of the system [GeV]
[in]statisticsquantum statistics of the particles species (+1 for Fermi, -1 for Bose, 0 for Boltzmann)
[in]solution_precisionprecision with which the solution is found

Definition at line 176 of file quantumsampling.cc.

178  {
179  /*
180  * Momentum at which the distribution function has its maximum.
181  */
182  double p_max[1];
183  p_max[0] = 0.0;
184 
185  /*
186  * Initial guess for the value of p_max, in GeV. This value is
187  * looped over within the GSL solver, so that many initial guesses
188  * are probed if the solution is not found.
189  */
190  double initial_guess_p_max = 0.050;
191 
192  /*
193  * Calling the GSL distribution maximum finder
194  */
196  mass, temperature, effective_chemical_potential, statistics,
197  initial_guess_p_max, solution_precision, p_max);
198 
199  double distribution_function_maximum =
200  p_max[0] * p_max[0] *
201  juttner_distribution_func(p_max[0], mass, temperature,
202  effective_chemical_potential, statistics);
203 
204  return distribution_function_maximum;
205 }
static int find_p_at_maximum_of_the_distribution(double mass, double temperature, double effective_chemical_potential, double statistics, double p_max_initial_guess, double solution_precision, double *p_max)
A GSL root solver for finding the radial momentum value at which the maximum of the given Juttner dis...
double juttner_distribution_func(double momentum_radial, double mass, double temperature, double effective_chemical_potential, double statistics)
Relativistic Juttner distribution function is just a convenience wrapper for displaying Fermi,...

◆ sample()

double smash::QuantumSampling::sample ( const PdgCode  pdg)

Sampling radial momenta of given particle species from Boltzmann, Bose, or Fermi distribution.

This sampler uses the simplest rejection sampling.

Parameters
[in]pdgthe pdg code of the sampled particle species return the sampled momentum [GeV]

Definition at line 257 of file quantumsampling.cc.

257  {
258  const ParticleType &ptype = ParticleType::find(pdg);
259  const double mass = ptype.mass();
260  const double mu = effective_chemical_potentials_.find(pdg)->second;
261  const double distr_max = distribution_function_maximums_.find(pdg)->second;
262  /*
263  * The variable maximum_momentum denotes the "far right" boundary of the
264  * sampled region; we assume that no particle has momentum larger than 10 GeV
265  */
266  constexpr double maximum_momentum = 10.0; // in [GeV]
267  const double statistics = (pdg.spin() % 2 == 0) ? -1.0 : 1.0;
268  double sampled_momentum = 0.0, sampled_ratio = 0.0;
269 
270  do {
271  sampled_momentum = random::uniform(0.0, maximum_momentum);
272  double distribution_at_sampled_p =
273  sampled_momentum * sampled_momentum *
274  juttner_distribution_func(sampled_momentum, mass, temperature_, mu,
275  statistics);
276  sampled_ratio = distribution_at_sampled_p / distr_max;
277  } while (random::canonical() > sampled_ratio);
278 
279  return sampled_momentum;
280 }
T uniform(T min, T max)
Definition: random.h:88
T canonical()
Definition: random.h:113

Member Data Documentation

◆ effective_chemical_potentials_

std::map<PdgCode, double> smash::QuantumSampling::effective_chemical_potentials_
private

Tabulated effective chemical potentials for every particle species.

Definition at line 156 of file quantumsampling.h.

◆ distribution_function_maximums_

std::map<PdgCode, double> smash::QuantumSampling::distribution_function_maximums_
private

Tabulated distribution function maxima for every particle species.

Definition at line 158 of file quantumsampling.h.

◆ volume_

const double smash::QuantumSampling::volume_
private

Volume [fm^3] in which particles sre sampled.

Definition at line 160 of file quantumsampling.h.

◆ temperature_

const double smash::QuantumSampling::temperature_
private

Temperature [GeV].

Definition at line 162 of file quantumsampling.h.


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