12 #include <gsl/gsl_integration.h>
13 #include <gsl/gsl_multiroots.h>
35 double effective_chemical_potential,
37 const double Ekin = std::sqrt(
p *
p + mass * mass);
39 2 * (1 + statistics * std::exp(-(Ekin - effective_chemical_potential) /
41 const double term2 = (
p *
p) / (temperature * Ekin);
48 gsl_vector *
function) {
52 const double mass = (par->
mass);
54 const double effective_chemical_potential =
58 const double p_radial = gsl_vector_get(roots_array, 0);
60 gsl_vector_set(
function, 0,
62 effective_chemical_potential, statistics));
68 gsl_multiroot_fsolver *solver) {
70 "\n***\nfind_p_at_maximum_of_the_distribution(): iter = %3u \t"
73 iter, gsl_vector_get(solver->x, 0), gsl_vector_get(solver->f, 0));
77 double mass,
double temperature,
double effective_chemical_potential,
78 double statistics,
double p_max_initial_guess,
double solution_precision,
80 const gsl_multiroot_fsolver_type *Solver_name;
81 gsl_multiroot_fsolver *Root_finder;
85 size_t initial_guess_update = 0;
87 const size_t problem_dimension = 1;
90 mass, temperature, effective_chemical_potential, statistics};
92 gsl_multiroot_function MaximumOfDistribution = {
95 double roots_array_initial[1] = {p_max_initial_guess};
97 gsl_vector *roots_array = gsl_vector_alloc(problem_dimension);
98 gsl_vector_set(roots_array, 0, roots_array_initial[0]);
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);
113 status = gsl_multiroot_fsolver_iterate(Root_finder);
121 if (initial_guess_update < 100) {
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,
137 std::cout <<
"\n\nGSL error message:\n"
138 << gsl_strerror(status) << std::endl;
139 logg[LogArea::Distributions::id].warn(
141 "\nfind_p_at_maximum_of_the_distribution\nis stuck!"
142 "\n\nInput parameters:"
144 mass,
"\n temperature [GeV] = ", temperature,
145 "\neffective_chemical_potential = ", effective_chemical_potential,
146 "\n statistics = ", statistics,
147 "\n solution_precision = ", solution_precision,
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 "
162 status = gsl_multiroot_test_residual(Root_finder->f, solution_precision);
164 if (status == GSL_SUCCESS) {
165 p_max[0] = gsl_vector_get(Root_finder->x, 0);
167 }
while (status == GSL_CONTINUE && iter < 100000);
169 gsl_multiroot_fsolver_free(Root_finder);
170 gsl_vector_free(roots_array);
176 double mass,
double temperature,
double effective_chemical_potential,
177 double statistics,
double solution_precision) {
189 double initial_guess_p_max = 0.050;
195 mass, temperature, effective_chemical_potential, statistics,
196 initial_guess_p_max, solution_precision, p_max);
198 double distribution_function_maximum =
199 p_max[0] * p_max[0] *
201 effective_chemical_potential, statistics);
203 return distribution_function_maximum;
211 const std::map<PdgCode, int> &initial_multiplicities,
double volume,
213 : volume_(volume), temperature_(temperature) {
219 constexpr
double solution_precision = 1e-8;
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;
225 const double number_density = number_of_particles / V_in_GeV;
228 const double quantum_statistics = (pdg.
spin() % 2 == 0) ? -1.0 : 1.0;
230 const double particle_mass = ptype.
mass();
231 double chemical_potential = 0.0;
235 spin_degeneracy, particle_mass, number_density,
temperature_,
236 quantum_statistics, solution_precision);
239 particle_mass,
temperature_, chemical_potential, quantum_statistics,
258 const double mass = ptype.
mass();
265 constexpr
double maximum_momentum = 10.0;
266 const double statistics = (pdg.
spin() % 2 == 0) ? -1.0 : 1.0;
267 double sampled_momentum = 0.0, sampled_ratio = 0.0;
271 double distribution_at_sampled_p =
272 sampled_momentum * sampled_momentum *
275 sampled_ratio = distribution_at_sampled_p / distr_max;
278 return sampled_momentum;