12 #include <gsl/gsl_errno.h>
13 #include <gsl/gsl_integration.h>
14 #include <gsl/gsl_multiroots.h>
30 double degeneracy,
double mass,
double temperature,
31 double effective_chemical_potential,
double statistics,
33 const auto integral = (*integrator)(0.0, 1.0, [&](
double x) {
34 return ((1.0 - x) * (1.0 - x) / (x * x * x * x)) *
38 const double integrated_number_density =
39 (degeneracy / (2.0 * M_PI * M_PI)) * integral;
41 return integrated_number_density;
50 double degeneracy,
double mass,
double number_density,
double temperature,
51 double effective_chemical_potential,
double statistics,
53 const double integrated_number_density =
57 return number_density - integrated_number_density;
61 const gsl_vector *roots_array,
void *parameters, gsl_vector *
function) {
66 const double mass = (par->
mass);
74 gsl_vector_set(
function, 0,
76 degeneracy, mass, number_density, temperature,
83 unsigned int iter, gsl_multiroot_fsolver *solver) {
85 "\n***\nfind_effective_chemical_potential(): iter = %3u \t"
88 iter, gsl_vector_get(solver->x, 0), gsl_vector_get(solver->f, 0));
92 double degeneracy,
double mass,
double number_density,
double temperature,
93 double statistics,
double mu_initial_guess,
double solution_precision,
94 Integrator *integrator,
double *effective_chemical_potential) {
95 const gsl_multiroot_fsolver_type *Solver_name;
96 gsl_multiroot_fsolver *Root_finder;
101 const size_t problem_dimension = 1;
104 degeneracy, mass, number_density, temperature, statistics, integrator};
106 gsl_multiroot_function EffectiveChemicalPotential = {
109 problem_dimension, ¶meters};
111 double roots_array_initial[problem_dimension] = {mu_initial_guess};
113 gsl_vector *roots_array = gsl_vector_alloc(problem_dimension);
114 gsl_vector_set(roots_array, 0, roots_array_initial[0]);
116 Solver_name = gsl_multiroot_fsolver_hybrids;
117 Root_finder = gsl_multiroot_fsolver_alloc(Solver_name, problem_dimension);
118 gsl_multiroot_fsolver_set(Root_finder, &EffectiveChemicalPotential,
130 status = gsl_multiroot_fsolver_iterate(Root_finder);
138 std::cout <<
"\nGSL error message:\n"
139 << gsl_strerror(status) <<
"\n\n"
141 const double conversion_factor =
143 logg[LogArea::Main::id].warn(
144 "\n\nThe GSL solver\nfind_effective_chemical_potential\nis stuck!"
145 "\n\nInput parameters:"
147 degeneracy,
"\n mass [GeV] = ", mass,
148 "\nnumber_density [fm^-3] = ", number_density / (conversion_factor),
149 "\n temperature [GeV] = ", temperature,
150 "\n statistics = ", statistics,
151 "\n solution_precision = ", solution_precision,
153 "Initialization cannot continue without calculating the "
154 "chemical potential. \nTry adjusting the initial "
155 "guess or the solution precision.\n\n\n");
156 throw std::runtime_error(
157 "Chemical potential calculation returned no result.\n\n");
161 status = gsl_multiroot_test_residual(Root_finder->f, solution_precision);
163 if (status == GSL_SUCCESS) {
166 }
while (status == GSL_CONTINUE && iter < 100000);
168 gsl_multiroot_fsolver_free(Root_finder);
169 gsl_vector_free(roots_array);
180 double degeneracy,
double mass,
double number_density,
double temperature,
181 double statistics,
double solution_precision) {
182 assert(temperature > 0);
183 assert(degeneracy > 0);
190 double initial_guess = 0.0;
191 if (statistics == -1.0) {
196 initial_guess = 0.999999 * mass;
198 initial_guess = 1.05 * mass;
202 const double integration_precision = 0.1 * solution_precision;
212 double chemical_potential[1];
213 chemical_potential[0] = 0.0;
219 degeneracy, mass, number_density, temperature, statistics, initial_guess,
220 solution_precision, &
integrator_, chemical_potential);
232 const double precision_fm = 1e-3;
234 double nDensityCheck =
238 if (std::abs((number_density - nDensityCheck) / conversion_factor) >
240 logg[LogArea::Main::id].warn(
241 "\n\nThe calculated chemical potential = ", chemical_potential[0],
242 " [GeV] does not reproduce the input number density to within ",
245 "\nnumber_density [GeV^3] = ",
246 number_density,
"\n nDensityCheck [GeV^3] = ", nDensityCheck,
247 "\nnumber_density [fm^-3] = ", number_density / (conversion_factor),
248 "\n nDensityCheck [fm^-3] = ", nDensityCheck / (conversion_factor),
252 degeneracy,
"\n mass [GeV] = ", mass,
253 "\nnumber_density [fm^-3] = ", number_density / (conversion_factor),
254 "\n temperature [GeV] = ", temperature,
255 "\n statistics = ", statistics,
256 "\n solution_precision = ", solution_precision,
258 "Initialization cannot continue without a correct "
259 "calculation of the chemical potential."
260 "\nTry adjusting the initial guess or the solution precision.\n\n\n");
261 throw std::runtime_error(
262 "Chemical potential calculation does not reproduce the "
263 "input number density.\n\n");
272 if ((statistics == -1) && (chemical_potential[0] > mass)) {
273 logg[LogArea::Main::id].warn(
274 "\n\nThe calculated chemical potential indicates the "
275 "input number density is such that Bose-Einstein "
276 "condensate is encountered."
279 degeneracy,
"\n mass [GeV] = ", mass,
280 "\nnumber_density [fm^-3] = ", number_density / (conversion_factor),
281 "\n temperature [GeV] = ", temperature,
282 "\n statistics = ", statistics,
283 "\n solution_precision = ", solution_precision,
285 "Initialization cannot proceed with a pathological "
286 "distribution function leading to the Bose-Einstein "
287 "condensate. \nTry increasing the temperature or "
288 "decreasing the number density.\n\n\n");
289 throw std::runtime_error(
290 "Chemical potential calculation indicates that the Bose-Einstein "
291 "condensate is produced.");
300 const double small_p = 0.005;
302 small_p, mass, temperature, chemical_potential[0], statistics);
303 if (distribution_at_small_p < 0.0) {
304 logg[LogArea::Main::id].warn(
305 "\n\nNegative values of the distribution function at "
306 "small momenta indicate that Bose-Einstein "
307 "condensate is encountered."
310 degeneracy,
"\n mass [GeV] = ", mass,
311 "\nnumber_density [fm^-3] = ", number_density / (conversion_factor),
312 "\n temperature [GeV] = ", temperature,
313 "\n statistics = ", statistics,
314 "\n solution_precision = ", solution_precision,
"\n\nf(p=", small_p,
315 ") = ", distribution_at_small_p,
317 "Initialization cannot proceed with a pathological "
318 "distribution function leading to the Bose-Einstein "
319 "condensate. \nTry increasing the temperature or "
320 "decreasing the number density.\n\n\n");
321 throw std::runtime_error(
322 "Negative values of the distribution function at low momenta "
323 "indicate that the Bose-Einstein condensate is produced.");
326 return chemical_potential[0];
static double root_equation_effective_chemical_potential(double degeneracy, double mass, double number_density, double temperature, double effective_chemical_potential, double statistics, Integrator *integrator)
Root equation for finding the value of the effective chemical potential for one particle species.
Integrator integrator_
A wrapper for gsl numerical integration.
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...
static void print_state_effective_chemical_potential(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 double density_one_species(double degeneracy, double mass, double temperature, double effective_chemical_potential, double statistics, Integrator *integrator)
Vector number density of one particle species, obtained through integrating the Juttner distribution ...
static int root_equation_effective_chemical_potential_for_GSL(const gsl_vector *roots_array, void *parameters, gsl_vector *function)
Root equation for finding the value of the effective chemical potential for one particle species,...
static int find_effective_chemical_potential(double degeneracy, double mass, double number_density, double temperature, double statistics, double mu_initial_guess, double solution_precision, Integrator *integrator, double *effective_chemical_potential)
A GSL root solver for finding the effective chemical potential.
A C++ interface for numerical integration in one dimension with the GSL CQUAD integration functions.
void set_precision(double absolute, double relative)
Set precision for the absolute and relative error of the integration.
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.
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.
Struct, root equations, and procedure for finding the effective chemical potential for a given partic...
double temperature
temperature T of the system [GeV]
double degeneracy
degeneracy g of the particle species
double number_density
number density n of the particle species [GeV^3]
double statistics
statistics quantum statistics of the particles species (+1 for Fermi, -1 for Bose,...
Integrator * integrator
wrapper for gsl numerical integration
double mass
(pole) mass m of the particle species