Version: SMASH-3.1
quantumsampling.h
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2020,2022
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #ifndef SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
11 #define SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
12 
13 #include <map>
14 
15 #include "gsl/gsl_multiroots.h"
16 #include "gsl/gsl_vector.h"
17 
18 #include "smash/pdgcode.h"
19 #include "smash/random.h"
20 
21 namespace smash {
22 
32  public:
41  QuantumSampling(const std::map<PdgCode, int>& initial_multiplicities,
42  double volume, double temperature);
43 
50  double mass;
52  double temperature;
59  double statistics;
60  };
61 
74  static double p_max_root_equation(double p, double mass, double temperature,
75  double effective_chemical_potential,
76  double statistics);
77 
90  static int p_max_root_equation_for_GSL(const gsl_vector* roots_array,
91  void* parameters,
92  gsl_vector* function);
93 
103  static void print_state_p_max(unsigned int iter,
104  gsl_multiroot_fsolver* solver);
105 
124  double mass, double temperature, double effective_chemical_potential,
125  double statistics, double p_max_initial_guess, double solution_precision,
126  double* p_max);
127 
141  static double maximum_of_the_distribution(double mass, double temperature,
142  double effective_chemical_potential,
143  double statistics,
144  double solution_precision);
145 
152  double sample(const PdgCode pdg);
153 
154  private:
156  std::map<PdgCode, double> effective_chemical_potentials_;
158  std::map<PdgCode, double> distribution_function_maximums_;
160  const double volume_;
162  const double temperature_;
163 };
164 
165 } // namespace smash
166 
167 #endif // SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:111
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.
constexpr int p
Proton.
Definition: action.h:24
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