7 #ifndef SRC_INCLUDE_GRANDCAN_THERMALIZER_H_ 8 #define SRC_INCLUDE_GRANDCAN_THERMALIZER_H_ 94 double nb()
const {
return nb_; }
96 double ns()
const {
return ns_; }
98 double e()
const {
return e_; }
100 double p()
const {
return p_; }
104 double T()
const {
return T_; }
236 const std::array<int, 3> n_cells,
237 const std::array<double, 3> origin,
bool periodicity,
238 double e_critical,
double t_start,
double delta_t,
242 const std::array<double, 3> lat_sizes,
243 const std::array<double, 3> origin,
bool periodicity)
245 lat_sizes, conf.take({
"Cell_Number"}), origin, periodicity,
246 conf.take({
"Critical_Edens"}), conf.take({
"Start_Time"}),
247 conf.take({
"Timestep"}),
249 conf.take({
"Microcanonical"},
false)) {}
256 const double t = clock->current_time();
257 const int n =
static_cast<int>(std::floor((t - t_start_) / period_));
258 return (t > t_start_ &&
259 t < t_start_ + n * period_ + clock->timestep_duration());
268 void update_thermalizer_lattice(
const Particles& particles,
270 bool ignore_cells_under_threshold =
true);
280 void renormalize_momenta(ParticleList& plist,
293 void sample_multinomial(
HadronClass particle_class,
int N);
305 void sample_in_random_cell_BF_algo(ParticleList& plist,
const double time,
319 void thermalize_BF_algo(
QuantumNumbers& conserved_initial,
double time,
328 template <
typename F>
331 N_total_in_cells_ = 0.0;
332 for (
auto cell_index : cells_to_sample_) {
334 const double gamma = 1.0 / std::sqrt(1.0 - cell.
v().
sqr());
337 if (condition(i->strangeness(), i->baryon_number(), i->charge())) {
339 N_tot += cell_volume_ * gamma *
344 N_in_cells_.push_back(N_tot);
345 N_total_in_cells_ += N_tot;
359 template <
typename F>
364 double partial_sum = 0.0;
365 int index_only_thermalized = -1;
366 while (partial_sum < r) {
367 index_only_thermalized++;
368 partial_sum += N_in_cells_[index_only_thermalized];
370 const int cell_index = cells_to_sample_[index_only_thermalized];
372 const ThreeVector cell_center = lat_->cell_center(cell_index);
373 const double gamma = 1.0 / std::sqrt(1.0 - cell.
v().
sqr());
374 const double N_in_cell = N_in_cells_[index_only_thermalized];
380 if (!condition(i->strangeness(), i->baryon_number(), i->charge())) {
384 cell_volume_ * gamma *
393 const double m = type_to_sample->
mass();
400 particle.
set_4momentum(m, phitheta.threevec() * momentum_radial);
414 void thermalize_mode_algo(
QuantumNumbers& conserved_initial,
double time);
422 void thermalize(
const Particles& particles,
double time,
int ntest);
430 void print_statistics(
const Clock& clock)
const;
434 double e_crit()
const {
return e_crit_; }
446 ParticleTypePtrList res;
449 res.push_back(&ptype);
459 const int B = eos_typelist_[typelist_index]->baryon_number();
460 const int S = eos_typelist_[typelist_index]->strangeness();
461 const int ch = eos_typelist_[typelist_index]->charge();
474 return mult_classes_[
static_cast<size_t>(cl)];
483 std::unique_ptr<RectangularLattice<ThermLatticeNode>>
lat_;
529 #endif // SRC_INCLUDE_GRANDCAN_THERMALIZER_H_ const double e_crit_
Critical energy density above which cells are thermalized.
void add_particle(const ParticleData &p, double factor)
Add particle contribution to Tmu0, nb and ns May look like unused at first glance, but it is actually used by update_lattice, where the node type of the lattice is templated.
A class to pre-calculate and store parameters relevant for density calculation.
ThermalizationAlgorithm
Defines the algorithm used for the forced thermalization.
HadronClass get_class(size_t typelist_index) const
Defines the class of hadrons by quantum numbers.
const size_t N_sorts_
Number of different species to be sampled.
The ThreeVector class represents a physical three-vector with the components .
double p() const
Get pressure in the rest frame.
ThermLatticeNode()
Default constructor of thermal quantities on the lattice returning thermodynamic quantities in comput...
Class to handle the equation of state (EoS) of the hadron gas, consisting of all hadrons included int...
FourVector Tmu0_
Four-momentum flow of the cell.
void add_particle_for_derivatives(const ParticleData &, double, ThreeVector)
dummy function for update_lattice
double mub() const
Get the net baryon chemical potential.
double ns() const
Get net strangeness density of the cell in the computational frame.
ParticleList particles_to_insert() const
List of newly created particles to be inserted in the simulation.
double e() const
Get energy density in the rest frame.
const double period_
Defines periodicity of the lattice in fm.
Interface to the SMASH configuration files.
static bool is_eos_particle(const ParticleType &ptype)
Check if a particle belongs to the EoS.
double mus() const
Get the net strangeness chemical potential.
double T() const
Get the temperature.
const ThermalizationAlgorithm algorithm_
Algorithm to choose for sampling of particles obeying conservation laws.
static double partial_density(const ParticleType &ptype, double T, double mub, double mus, bool account_for_resonance_widths=false)
Compute partial density of one hadron sort.
RectangularLattice< ThermLatticeNode > & lattice() const
Getter function for the lattice.
void set_formation_time(const double &form_time)
Set the absolute formation time.
void set_rest_frame_quantities(double T0, double mub0, double mus0, const ThreeVector v0)
Set all the rest frame quantities to some values, this is useful for testing.
A container class to hold all the arrays on the lattice and access them.
ThreeVector v_
Velocity of the rest frame.
double p_
Pressure in the rest frame.
static const ParticleTypeList & list_all()
double ns_
Net strangeness density of the cell in the computational frame.
const bool BF_enforce_microcanonical_
Enforce energy conservation as part of BF sampling algorithm or not.
double mub_
Net baryon chemical potential.
GrandCanThermalizer(Configuration &conf, const std::array< double, 3 > lat_sizes, const std::array< double, 3 > origin, bool periodicity)
ParticleList sampled_list_
Newly generated particles by thermalizer.
double e_crit() const
Get the critical energy density.
ThreeVector v() const
Get 3-velocity of the rest frame.
Particle type contains the static properties of a particle species.
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
Non-strange mesons (S = 0) with electric cherge Q < 0.
A container for storing conserved values.
void compute_rest_frame_quantities(HadronGasEos &eos)
Temperature, chemical potentials and rest frame velocity are calculated given the hadron gas equation...
double sample_momenta_from_thermal(const double temperature, const double mass)
Samples a momentum from the Maxwell-Boltzmann (thermal) distribution in a faster way, given by Scott Pratt (see Pratt:2014vja) APPENDIX: ALGORITHM FOR GENERATING PARTICLES math trick: for distribution, sample x by: where are uniform random numbers between [0,1) for : , where is used as rejection weight.
void compute_N_in_cells_mode_algo(F &&condition)
Computes average number of particles in each cell for the mode algorithm.
std::vector< size_t > cells_to_sample_
Cells above critical energy density.
Mesons with strangeness S < 0.
double N_total_in_cells_
Total number of particles over all cells in thermalization region.
The GrandCanThermalizer class implements the following functionality:
Clock tracks the time in the simulation.
Neutral non-strange mesons.
Non-strange mesons (S = 0) with electric cherge Q > 0.
ParticleData sample_in_random_cell_mode_algo(const double time, F &&condition)
Samples one particle and the species, cell, momentum and coordinate are chosen from the corresponding...
std::vector< double > mult_sort_
Real number multiplicity for each particle type.
ParticleList to_remove_
Particles to be removed after this thermalization step.
Mesons with strangeness S > 0.
std::vector< double > N_in_cells_
Number of particles to be sampled in one cell.
double nb_
Net baryon density of the cell in the computational frame.
ParticleList particles_to_remove() const
List of particles to be removed from the simulation.
void set_4position(const FourVector &pos)
Set the particle's 4-position directly.
double e_
Energy density in the rest frame.
The ThermLatticeNode class is intended to compute thermodynamical quantities in a cell given a set of...
HadronClass
Specifier to classify the different hadron species according to their quantum numbers.
std::vector< int > mult_int_
Integer multiplicity for each particle type.
double mus_
Net strangeness chemical potential.
std::array< double, 7 > mult_classes_
The different hadron species according to the enum defined in.
A pointer-like interface to global references to ParticleType objects.
double nb() const
Get net baryon density of the cell in the computational frame.
Angles provides a common interface for generating directions: i.e., two angles that should be interpr...
std::unique_ptr< RectangularLattice< ThermLatticeNode > > lat_
The lattice on which the thermodynamic quantities are calculated.
void distribute_isotropically()
Populate the object with a new direction.
The Particles class abstracts the storage and manipulation of particles.
ParticleTypePtrList list_eos_particles() const
Extracts the particles in the hadron gas equation of state from the complete list of particle types i...
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
double mult_class(const HadronClass cl) const
bool is_time_to_thermalize(std::unique_ptr< Clock > &clock) const
Check that the clock is close to n * period of thermalization, since the thermalization only happens ...
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
double cell_volume_
Volume of a single cell, necessary to convert thermal densities to actual particle numbers...
ParticleData contains the dynamic information of a certain particle.
void boost_momentum(const ThreeVector &v)
Apply a Lorentz-boost to only the momentum.
const ParticleTypePtrList eos_typelist_
List of particle types from which equation of state is computed.
const double t_start_
Starting time of the simulation.
FourVector Tmu0() const
Get Four-momentum flow of the cell.