Version: SMASH-2.2
quantumsampling.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2020-
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #include "smash/quantumsampling.h"
11 
12 #include <gsl/gsl_integration.h>
13 #include <gsl/gsl_multiroots.h>
14 #include <cstdlib>
15 #include <iostream>
16 #include <random>
17 
19 #include "smash/constants.h"
20 #include "smash/distributions.h"
21 #include "smash/logging.h"
22 #include "smash/particletype.h"
23 
24 namespace smash {
25 
26 /*
27  * Root equations and GSL procedure for finding the momentum for which the
28  * maximum of a given Juttner distribution occurs. This is needed for a method
29  * of sampling the distribution function in which one samples uniformly below
30  * the maximum of the distribution.
31  */
32 
33 double QuantumSampling::p_max_root_equation(double p, double mass,
34  double temperature,
35  double effective_chemical_potential,
36  double statistics) {
37  const double Ekin = std::sqrt(p * p + mass * mass);
38  const double term1 =
39  2 * (1 + statistics * std::exp(-(Ekin - effective_chemical_potential) /
40  temperature));
41  const double term2 = (p * p) / (temperature * Ekin);
42 
43  return term1 - term2;
44 }
45 
46 int QuantumSampling::p_max_root_equation_for_GSL(const gsl_vector *roots_array,
47  void *parameters,
48  gsl_vector *function) {
50  static_cast<struct ParametersForMaximumMomentumRootFinder *>(parameters);
51 
52  const double mass = (par->mass);
53  const double temperature = (par->temperature);
54  const double effective_chemical_potential =
56  const double statistics = (par->statistics);
57 
58  const double p_radial = gsl_vector_get(roots_array, 0);
59 
60  gsl_vector_set(function, 0,
61  p_max_root_equation(p_radial, mass, temperature,
62  effective_chemical_potential, statistics));
63 
64  return GSL_SUCCESS;
65 }
66 
67 void QuantumSampling::print_state_p_max(unsigned int iter,
68  gsl_multiroot_fsolver *solver) {
69  printf(
70  "\n***\nfind_p_at_maximum_of_the_distribution(): iter = %3u \t"
71  "x = % .3f \t"
72  "f(x) = % .3e \n",
73  iter, gsl_vector_get(solver->x, 0), gsl_vector_get(solver->f, 0));
74 }
75 
77  double mass, double temperature, double effective_chemical_potential,
78  double statistics, double p_max_initial_guess, double solution_precision,
79  double *p_max) {
80  const gsl_multiroot_fsolver_type *Solver_name;
81  gsl_multiroot_fsolver *Root_finder;
82 
83  int status;
84  size_t iter = 0;
85  size_t initial_guess_update = 0;
86 
87  const size_t problem_dimension = 1;
88 
89  struct ParametersForMaximumMomentumRootFinder parameters = {
90  mass, temperature, effective_chemical_potential, statistics};
91 
92  gsl_multiroot_function MaximumOfDistribution = {
93  &p_max_root_equation_for_GSL, problem_dimension, &parameters};
94 
95  double roots_array_initial[1] = {p_max_initial_guess};
96 
97  gsl_vector *roots_array = gsl_vector_alloc(problem_dimension);
98  gsl_vector_set(roots_array, 0, roots_array_initial[0]);
99 
100  Solver_name = gsl_multiroot_fsolver_hybrids;
101  Root_finder = gsl_multiroot_fsolver_alloc(Solver_name, problem_dimension);
102  gsl_multiroot_fsolver_set(Root_finder, &MaximumOfDistribution, roots_array);
103 
104  // print_state_p_max (iter, Root_finder);
105 
106  do {
107  iter++;
108 
109  /*
110  * gsl_multiroot_fsolver_iterate returns either 0 for a correct behavior,
111  * or an error code (a positive integer) when the solver is stuck
112  */
113  status = gsl_multiroot_fsolver_iterate(Root_finder);
114 
115  // print_state_p_max (iter, Root_finder);
116 
117  /*
118  * Check whether the solver is stuck
119  */
120  if (status) {
121  if (initial_guess_update < 100) {
122  /*
123  * In case the solution is not found for the (somewhat small) default
124  * initial guess of p_max_initial_guess = 0.05 [GeV], the value of the
125  * p_max_initial_guess is increased and the root solving procedure is
126  * restarted. This can take place up to a 100 times, with the largest
127  * possible initial guess of p_max_initial_guess = 5.05 [GeV].
128  */
129  p_max_initial_guess += 0.05;
130  initial_guess_update++;
131  roots_array_initial[0] = p_max_initial_guess;
132  gsl_vector_set(roots_array, 0, roots_array_initial[0]);
133  gsl_multiroot_fsolver_set(Root_finder, &MaximumOfDistribution,
134  roots_array);
135  iter = 0;
136  } else {
137  std::cout << "\n\nGSL error message:\n"
138  << gsl_strerror(status) << std::endl;
139  logg[LogArea::Distributions::id].warn(
140  "\n\nThe GSL solver"
141  "\nfind_p_at_maximum_of_the_distribution\nis stuck!"
142  "\n\nInput parameters:"
143  "\n mass [GeV] = ",
144  mass, "\n temperature [GeV] = ", temperature,
145  "\neffective_chemical_potential = ", effective_chemical_potential,
146  "\n statistics = ", statistics,
147  "\n solution_precision = ", solution_precision,
148  "\n\n"
149  "Initialization cannot sample the momenta without "
150  "calculating the distribution maximum."
151  "\nTry adjusting the initial guess (which is "
152  "looped over in the GSL procedure) or the "
153  "solution precision."
154  "\nUncomment print_state_p_max to check solver progress.\n\n\n");
155  throw std::runtime_error(
156  "QuantumSampling::find_p_at_maximum_of_the_distribution returned "
157  "no result.\n\n");
158  continue;
159  }
160  }
161 
162  status = gsl_multiroot_test_residual(Root_finder->f, solution_precision);
163 
164  if (status == GSL_SUCCESS) {
165  p_max[0] = gsl_vector_get(Root_finder->x, 0);
166  }
167  } while (status == GSL_CONTINUE && iter < 100000);
168 
169  gsl_multiroot_fsolver_free(Root_finder);
170  gsl_vector_free(roots_array);
171 
172  return 0;
173 }
174 
176  double mass, double temperature, double effective_chemical_potential,
177  double statistics, double solution_precision) {
178  /*
179  * Momentum at which the distribution function has its maximum.
180  */
181  double p_max[1];
182  p_max[0] = 0.0;
183 
184  /*
185  * Initial guess for the value of p_max, in GeV. This value is
186  * looped over within the GSL solver, so that many initial guesses
187  * are probed if the solution is not found.
188  */
189  double initial_guess_p_max = 0.050;
190 
191  /*
192  * Calling the GSL distribution maximum finder
193  */
195  mass, temperature, effective_chemical_potential, statistics,
196  initial_guess_p_max, solution_precision, p_max);
197 
198  double distribution_function_maximum =
199  p_max[0] * p_max[0] *
200  juttner_distribution_func(p_max[0], mass, temperature,
201  effective_chemical_potential, statistics);
202 
203  return distribution_function_maximum;
204 }
205 
206 /*
207  * Initializing the QuantumSampling object triggers calculation of the
208  * chemical potential and distribution maximum for all species present.
209  */
211  const std::map<PdgCode, int> &initial_multiplicities, double volume,
212  double temperature)
213  : volume_(volume), temperature_(temperature) {
214  /*
215  * This is the precision which we expect from the solution; note that
216  * solution precision also goes into the precision of calculating the
217  * integrals involved etc. Recommended precision is at least 1e-7.
218  */
219  constexpr double solution_precision = 1e-8;
220 
221  for (const auto &pdg_and_mult : initial_multiplicities) {
222  const PdgCode pdg = pdg_and_mult.first;
223  const int number_of_particles = pdg_and_mult.second;
224  const double V_in_GeV = volume_ / (hbarc * hbarc * hbarc);
225  const double number_density = number_of_particles / V_in_GeV;
226  const double spin_degeneracy = pdg.spin_degeneracy();
227  // '+' for fermions, '-' for bosons
228  const double quantum_statistics = (pdg.spin() % 2 == 0) ? -1.0 : 1.0;
229  const ParticleType &ptype = ParticleType::find(pdg);
230  const double particle_mass = ptype.mass();
231  double chemical_potential = 0.0;
232  ChemicalPotentialSolver mu_solver;
233  // Calling the wrapper for the GSL chemical potential finder
234  chemical_potential = mu_solver.effective_chemical_potential(
235  spin_degeneracy, particle_mass, number_density, temperature_,
236  quantum_statistics, solution_precision);
237  effective_chemical_potentials_[pdg] = chemical_potential;
238  const double distribution_function_maximum = maximum_of_the_distribution(
239  particle_mass, temperature_, chemical_potential, quantum_statistics,
240  solution_precision);
241  distribution_function_maximums_[pdg] = distribution_function_maximum;
242  }
243 }
244 
245 /*
246  * Sampling radial momenta of given particle species from Bose, Boltzmann, or
247  * Fermi distribution. The choice between the distributions is made based on
248  * the <statistics> variable:
249  * -1 for Bose
250  * 0 fot Boltzmann
251  * +1 for Fermi
252  *
253  * This sampler is the simplest implementation of sampling based on sampling
254  * from a uniform distribution.
255  */
256 double QuantumSampling::sample(const PdgCode pdg) {
257  const ParticleType &ptype = ParticleType::find(pdg);
258  const double mass = ptype.mass();
259  const double mu = effective_chemical_potentials_.find(pdg)->second;
260  const double distr_max = distribution_function_maximums_.find(pdg)->second;
261  /*
262  * The variable maximum_momentum denotes the "far right" boundary of the
263  * sampled region; we assume that no particle has momentum larger than 10 GeV
264  */
265  constexpr double maximum_momentum = 10.0; // in [GeV]
266  const double statistics = (pdg.spin() % 2 == 0) ? -1.0 : 1.0;
267  double sampled_momentum = 0.0, sampled_ratio = 0.0;
268 
269  do {
270  sampled_momentum = random::uniform(0.0, maximum_momentum);
271  double distribution_at_sampled_p =
272  sampled_momentum * sampled_momentum *
273  juttner_distribution_func(sampled_momentum, mass, temperature_, mu,
274  statistics);
275  sampled_ratio = distribution_at_sampled_p / distr_max;
276  } while (random::canonical() > sampled_ratio);
277 
278  return sampled_momentum;
279 }
280 
281 } // namespace smash
A class which encapsulates a GSL algorithm for finding the effective chemical potential and supportin...
double effective_chemical_potential(double degeneracy, double mass, double number_density, double temperature, double statistics, double solution_precision)
Convenience wrapper for finding the effective chemical potential for a given particle species and per...
Particle type contains the static properties of a particle species.
Definition: particletype.h:97
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
Definition: particletype.cc:99
double mass() const
Definition: particletype.h:144
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:108
unsigned int spin() const
Definition: pdgcode.h:529
unsigned int spin_degeneracy() const
Definition: pdgcode.h:565
QuantumSampling(const std::map< PdgCode, int > &initial_multiplicities, double volume, double temperature)
Constructor of a QuantumSampling object.
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...
const double temperature_
Temperature [GeV].
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...
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.
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 procedur...
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...
double sample(const PdgCode pdg)
Sampling radial momenta of given particle species from Boltzmann, Bose, or Fermi distribution.
Collection of useful constants that are known at compile time.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
constexpr int p
Proton.
T uniform(T min, T max)
Definition: random.h:88
T canonical()
Definition: random.h:113
Definition: action.h:24
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,...
constexpr double hbarc
GeV <-> fm conversion factor.
Definition: constants.h:25
Struct object that holds the parameters relevant to finding the momentum for which the maximum of the...
double statistics
quantum statistics of the particles species (+1 for Fermi, -1 for Bose, 0 for Boltzmann)
double effective_chemical_potential
effective chemical potential mu^* of the particle species