Version: SMASH-2.2
quantumsampling.h
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 #ifndef SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
11 #define SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
12 
13 #include <gsl/gsl_multiroots.h>
14 #include <gsl/gsl_vector.h>
15 #include <map>
16 
17 #include "smash/pdgcode.h"
18 #include "smash/random.h"
19 
20 namespace smash {
21 
31  public:
40  QuantumSampling(const std::map<PdgCode, int>& initial_multiplicities,
41  double volume, double temperature);
42 
49  double mass;
51  double temperature;
58  double statistics;
59  };
60 
73  static double p_max_root_equation(double p, double mass, double temperature,
74  double effective_chemical_potential,
75  double statistics);
76 
89  static int p_max_root_equation_for_GSL(const gsl_vector* roots_array,
90  void* parameters,
91  gsl_vector* function);
92 
102  static void print_state_p_max(unsigned int iter,
103  gsl_multiroot_fsolver* solver);
104 
123  double mass, double temperature, double effective_chemical_potential,
124  double statistics, double p_max_initial_guess, double solution_precision,
125  double* p_max);
126 
140  static double maximum_of_the_distribution(double mass, double temperature,
141  double effective_chemical_potential,
142  double statistics,
143  double solution_precision);
144 
151  double sample(const PdgCode pdg);
152 
153  private:
155  std::map<PdgCode, double> effective_chemical_potentials_;
157  std::map<PdgCode, double> distribution_function_maximums_;
159  const double volume_;
161  const double temperature_;
162 };
163 
164 } // namespace smash
165 
166 #endif // SRC_INCLUDE_SMASH_QUANTUMSAMPLING_H_
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:108
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