Version: SMASH-3.2
smash::BoxModus Class Reference

#include <boxmodus.h>

BoxModus: Provides a modus for infinite matter calculations.

Matter is confined in a cubical box. Depending on the initial condition, particles are either reflected on the boundaries (not implemented now) or inserted on opposite positions.

To use this modus, choose Modus: Box

General:
Modus: Box

in the configuration file.

Options for BoxModus go in the "Modi"→"Box" section of the configuration:

Modi:
Box:
# definitions here

The following configuration options are understood: Box

Definition at line 46 of file boxmodus.h.

Inheritance diagram for smash::BoxModus:
smash::ModusDefault

Public Member Functions

 BoxModus (Configuration modus_config, const ExperimentParameters &parameters)
 Constructor. More...
 
double initial_conditions (Particles *particles, const ExperimentParameters &parameters)
 Generates initial state of the particles in the system according to specified parameters: number of particles of each species, momentum and coordinate space distributions. More...
 
int impose_boundary_conditions (Particles *particles, const OutputsList &output_list={})
 Enforces that all particles are inside the box at the beginning of an event. More...
 
Grid< GridOptions::PeriodicBoundariescreate_grid (const Particles &particles, double min_cell_length, double timestep_duration, CollisionCriterion crit, const bool include_unformed_particles, CellSizeStrategy strategy=CellSizeStrategy::Optimal) const
 Creates the Grid with normal boundary conditions. More...
 
std::unique_ptr< GrandCanThermalizercreate_grandcan_thermalizer (Configuration &conf) const
 Creates GrandCanThermalizer. More...
 
double max_timestep (double max_transverse_distance_sqr) const
 
double equilibration_time () const
 
bool is_box () const
 
double length () const
 
- Public Member Functions inherited from smash::ModusDefault
int impose_boundary_conditions (Particles *, const OutputsList &={})
 Enforces sensible positions for the particles. More...
 
bool is_collider () const
 
bool is_box () const
 
bool is_list () const
 
bool is_sphere () const
 
double sqrt_s_NN () const
 
double impact_parameter () const
 
void sample_impact () const
 sample impact parameter for collider modus More...
 
double velocity_projectile () const
 
double velocity_target () const
 
FermiMotion fermi_motion () const
 
double max_timestep (double) const
 
double equilibration_time () const
 
double length () const
 
double radius () const
 
bool calculation_frame_is_fixed_target () const
 
double nuclei_passing_time () const
 Get the passing time of the two nuclei in a collision. More...
 
bool is_IC_for_hybrid () const
 
const InitialConditionParametersIC_parameters () const
 
const std::map< int32_t, double > & fluid_background ()
 
const RectangularLattice< EnergyMomentumTensor > & fluid_lattice ()
 
void build_fluidization_lattice ([[maybe_unused]] const double t, [[maybe_unused]] const std::vector< Particles > &ensembles, [[maybe_unused]] const DensityParameters &dens_par)
 Build lattice of energy momentum tensor. More...
 
Grid< GridOptions::Normalcreate_grid (const Particles &particles, double min_cell_length, double timestep_duration, CollisionCriterion crit, const bool include_unformed_particles, CellSizeStrategy strategy=CellSizeStrategy::Optimal) const
 Creates the Grid with normal boundary conditions. More...
 
std::unique_ptr< GrandCanThermalizercreate_grandcan_thermalizer (Configuration &conf) const
 Creates GrandCanThermalizer. More...
 

Private Attributes

const BoxInitialCondition initial_condition_
 Initial momenta distribution: thermal or peaked momenta. More...
 
const double length_
 Length of the cube's edge in fm. More...
 
const double equilibration_time_
 time after which output is written More...
 
const double temperature_
 Temperature of the Box in GeV. More...
 
const double start_time_ = 0.
 Initial time of the box. More...
 
const bool use_thermal_ = false
 Whether to use a thermal initialization for all particles instead of specific numbers. More...
 
const double mub_
 Baryon chemical potential for thermal initialization; only used if use_thermal_ is true. More...
 
const double mus_
 Strange chemical potential for thermal initialization; only used if use_thermal_ is true. More...
 
const double muq_
 Charge chemical potential for thermal initialization; only used if use_thermal_ is true. More...
 
const bool account_for_resonance_widths_
 In case of thermal initialization: true – account for resonance spectral functions, while computing multiplicities and sampling masses, false – simply use pole masses. More...
 
const std::map< PdgCode, int > init_multipl_
 Particle multiplicities at initialization; required if use_thermal_ is false. More...
 
std::map< PdgCode, double > average_multipl_
 Average multiplicities in case of thermal initialization. More...
 
const std::optional< PdgCodejet_pdg_
 Optional PDG code of the particle to use as a jet. More...
 
const double jet_mom_
 Initial momentum of the jet particle; only used if insert_jet_ is true. More...
 

Friends

std::ostream & operator<< (std::ostream &out, const BoxModus &m)
 Console output on startup of box specific parameters; writes the initial state for the box to the output stream. More...
 

Constructor & Destructor Documentation

◆ BoxModus()

smash::BoxModus::BoxModus ( Configuration  modus_config,
const ExperimentParameters parameters 
)
explicit

Constructor.

Gathers all configuration variables for the Box.

Parameters
[in]modus_configThe configuration object that sets all initial conditions of the experiment.
[in]parametersUnused, but necessary because of templated initialization

Definition at line 66 of file boxmodus.cc.

69  modus_config.take(InputKeys::modi_box_initialCondition)),
70  length_(modus_config.take(InputKeys::modi_box_length)),
72  modus_config.take(InputKeys::modi_box_equilibrationTime)),
74  start_time_(modus_config.take(InputKeys::modi_box_startTime)),
84  ? std::map<PdgCode, int>()
85  : modus_config.take(InputKeys::modi_box_initialMultiplicities)),
86  /* Note that it is crucial not to take other keys from the Jet section
87  * before Jet_PDG, since we want here the take to throw in case the user
88  * had a Jet section without the mandatory Jet_PDG key. If all other keys
89  * are taken first, the section is removed from the config because empty,
90  * and has_section(InputSections::m_b_jet) method would return false.
91  */
92  jet_pdg_(modus_config.has_section(InputSections::m_b_jet)
93  ? make_optional<PdgCode>(
94  modus_config.take(InputKeys::modi_box_jet_jetPdg))
95  : std::nullopt),
96 
97  jet_mom_(modus_config.take(InputKeys::modi_box_jet_jetMomentum)) {
98  if (parameters.res_lifetime_factor < 0.) {
99  throw std::invalid_argument(
100  "Resonance lifetime modifier cannot be negative!");
101  }
102  // Check consistency, just in case
103  if (std::abs(length_ - parameters.box_length) > really_small) {
104  throw std::runtime_error("Box length inconsistency");
105  }
106 }
const std::optional< PdgCode > jet_pdg_
Optional PDG code of the particle to use as a jet.
Definition: boxmodus.h:190
const double muq_
Charge chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: boxmodus.h:169
const double jet_mom_
Initial momentum of the jet particle; only used if insert_jet_ is true.
Definition: boxmodus.h:194
const bool use_thermal_
Whether to use a thermal initialization for all particles instead of specific numbers.
Definition: boxmodus.h:154
const double mub_
Baryon chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: boxmodus.h:159
const BoxInitialCondition initial_condition_
Initial momenta distribution: thermal or peaked momenta.
Definition: boxmodus.h:141
const std::map< PdgCode, int > init_multipl_
Particle multiplicities at initialization; required if use_thermal_ is false.
Definition: boxmodus.h:180
const bool account_for_resonance_widths_
In case of thermal initialization: true – account for resonance spectral functions,...
Definition: boxmodus.h:175
const double temperature_
Temperature of the Box in GeV.
Definition: boxmodus.h:147
const double mus_
Strange chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: boxmodus.h:164
const double start_time_
Initial time of the box.
Definition: boxmodus.h:149
const double equilibration_time_
time after which output is written
Definition: boxmodus.h:145
const double length_
Length of the cube's edge in fm.
Definition: boxmodus.h:143
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
static const Key< double > modi_box_temperature
See user guide description for more information.
Definition: input_keys.h:4290
static const Key< bool > modi_box_useThermalMultiplicities
See user guide description for more information.
Definition: input_keys.h:4389
static const Key< double > modi_box_baryonChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4332
static const Key< double > modi_box_length
See user guide description for more information.
Definition: input_keys.h:4265
static const Key< bool > modi_box_accountResonanceWidths
See user guide description for more information.
Definition: input_keys.h:4319
static const Key< double > modi_box_equilibrationTime
See user guide description for more information.
Definition: input_keys.h:4361
static const Key< double > modi_box_jet_jetMomentum
See user guide description for more information.
Definition: input_keys.h:4414
static const Key< double > modi_box_startTime
See user guide description for more information.
Definition: input_keys.h:4278
static const Key< double > modi_box_strangeChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4375
static const Key< double > modi_box_chargeChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4345
static const Key< PdgCode > modi_box_jet_jetPdg
See user guide description for more information.
Definition: input_keys.h:4427
static const Key< BoxInitialCondition > modi_box_initialCondition
See user guide description for more information.
Definition: input_keys.h:4253
static const Key< std::map< PdgCode, int > > modi_box_initialMultiplicities
See user guide description for more information.
Definition: input_keys.h:4233
static const Section m_b_jet
Subsection for the jet in box modus.
Definition: input_keys.h:93

Member Function Documentation

◆ initial_conditions()

double smash::BoxModus::initial_conditions ( Particles particles,
const ExperimentParameters parameters 
)

Generates initial state of the particles in the system according to specified parameters: number of particles of each species, momentum and coordinate space distributions.

Subsequently makes the total 3-momentum 0.

Parameters
[out]particlesAn empty list that gets filled up by this function
[in]parametersThe initialization parameters of the box
Returns
The starting time of the simulation

Initialize formation time

Definition at line 108 of file boxmodus.cc.

109  {
110  double momentum_radial = 0.0, mass = 0.0;
111  Angles phitheta;
112  FourVector momentum_total(0, 0, 0, 0);
113  auto uniform_length = random::make_uniform_distribution(0.0, this->length_);
114  const double T = this->temperature_;
115  const double V = length_ * length_ * length_;
116  /* Create NUMBER OF PARTICLES according to configuration, or thermal case */
117  if (use_thermal_) {
118  if (average_multipl_.empty()) {
119  for (const ParticleType &ptype : ParticleType::list_all()) {
120  if (HadronGasEos::is_eos_particle(ptype)) {
121  const double lifetime_factor =
122  ptype.is_stable() ? 1. : parameters.res_lifetime_factor;
123  const double n = lifetime_factor * HadronGasEos::partial_density(
124  ptype, T, mub_, mus_, muq_,
126  average_multipl_[ptype.pdgcode()] = n * V * parameters.testparticles;
127  }
128  }
129  }
130  double nb_init = 0.0, ns_init = 0.0, nq_init = 0.0;
131  for (const auto &mult : average_multipl_) {
132  const int thermal_mult_int = random::poisson(mult.second);
133  particles->create(thermal_mult_int, mult.first);
134  nb_init += mult.second * mult.first.baryon_number();
135  ns_init += mult.second * mult.first.strangeness();
136  nq_init += mult.second * mult.first.charge();
137  logg[LBox].debug(mult.first, " initial multiplicity ", thermal_mult_int);
138  }
139  logg[LBox].info("Initial hadron gas baryon density ", nb_init);
140  logg[LBox].info("Initial hadron gas strange density ", ns_init);
141  logg[LBox].info("Initial hadron gas charge density ", nq_init);
142  } else {
143  for (const auto &p : init_multipl_) {
144  particles->create(p.second * parameters.testparticles, p.first);
145  logg[LBox].debug("Particle ", p.first, " initial multiplicity ",
146  p.second);
147  }
148  }
149  std::unique_ptr<QuantumSampling> quantum_sampling;
151  quantum_sampling = std::make_unique<QuantumSampling>(init_multipl_, V, T);
152  }
153  for (ParticleData &data : *particles) {
154  /* Set MOMENTUM SPACE distribution */
156  /* initial thermal momentum is the average 3T */
157  momentum_radial = 3.0 * T;
158  mass = data.pole_mass();
159  } else {
160  if (this->initial_condition_ ==
162  /* thermal momentum according Maxwell-Boltzmann distribution */
164  ? data.type().mass()
165  : HadronGasEos::sample_mass_thermal(data.type(), 1.0 / T);
166  momentum_radial = sample_momenta_from_thermal(T, mass);
167  } else if (this->initial_condition_ ==
169  /*
170  * Sampling the thermal momentum according Bose/Fermi/Boltzmann
171  * distribution.
172  * We take the pole mass as the mass.
173  */
174  mass = data.type().mass();
175  momentum_radial = quantum_sampling->sample(data.pdgcode());
176  }
177  }
178  phitheta.distribute_isotropically();
179  logg[LBox].debug(data.type().name(), "(id ", data.id(),
180  ") radial momentum ", momentum_radial, ", direction",
181  phitheta);
182  data.set_4momentum(mass, phitheta.threevec() * momentum_radial);
183  momentum_total += data.momentum();
184 
185  /* Set COORDINATE SPACE distribution */
186  ThreeVector pos{uniform_length(), uniform_length(), uniform_length()};
187  data.set_4position(FourVector(start_time_, pos));
189  data.set_formation_time(start_time_);
190  }
191 
192  /* Make total 3-momentum 0 */
193  for (ParticleData &data : *particles) {
194  data.set_4momentum(data.momentum().abs(),
195  data.momentum().threevec() -
196  momentum_total.threevec() / particles->size());
197  }
198 
199  /* Add a single highly energetic particle in the center of the box (jet) */
200  if (jet_pdg_) {
201  auto &jet_particle = particles->create(jet_pdg_.value());
202  jet_particle.set_formation_time(start_time_);
203  jet_particle.set_4position(FourVector(start_time_, 0., 0., 0.));
204  jet_particle.set_4momentum(ParticleType::find(jet_pdg_.value()).mass(),
205  ThreeVector(jet_mom_, 0., 0.));
206  }
207 
208  /* Recalculate total momentum */
209  momentum_total = FourVector(0, 0, 0, 0);
210  for (ParticleData &data : *particles) {
211  momentum_total += data.momentum();
212  /* IC: debug checks */
213  logg[LBox].debug() << data;
214  }
215  /* allows to check energy conservation */
216  logg[LBox].debug() << "Initial total 4-momentum [GeV]: " << momentum_total;
217  return start_time_;
218 }
std::map< PdgCode, double > average_multipl_
Average multiplicities in case of thermal initialization.
Definition: boxmodus.h:185
static double partial_density(const ParticleType &ptype, double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute partial density of one hadron sort.
Definition: hadgas_eos.cc:270
static double sample_mass_thermal(const ParticleType &ptype, double beta)
Sample resonance mass in a thermal medium.
Definition: hadgas_eos.cc:385
static bool is_eos_particle(const ParticleType &ptype)
Check if a particle belongs to the EoS.
Definition: hadgas_eos.h:355
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
Definition: particletype.cc:99
static const ParticleTypeList & list_all()
Definition: particletype.cc:51
double mass() const
Definition: particletype.h:145
@ ThermalMomentaBoltzmann
A thermalized ensemble is generated, with momenta sampled from a Maxwell-Boltzmann distribution.
@ ThermalMomentaQuantum
A thermalized ensemble is generated, with momenta of baryons(mesons) sampled from a Fermi(Bose) distr...
@ PeakedMomenta
All particles have the same momentum with T being the temperature.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:40
constexpr int p
Proton.
constexpr int n
Neutron.
int poisson(const T &lam)
Returns a Poisson distributed random number.
Definition: random.h:226
uniform_dist< T > make_uniform_distribution(T min, T max)
Definition: random.h:135
double sample_momenta_from_thermal(const double temperature, const double mass)
Samples a momentum from the Maxwell-Boltzmann (thermal) distribution in a faster way,...
static constexpr int LBox
Definition: boxmodus.cc:30

◆ impose_boundary_conditions()

int smash::BoxModus::impose_boundary_conditions ( Particles particles,
const OutputsList &  output_list = {} 
)

Enforces that all particles are inside the box at the beginning of an event.

It checks if the particles were placed correctly inside the box at initialization and places them inside if they are not.

Parameters
[in]particlesparticles to check their position and possibly move it
[in]output_listoutput objects
Returns
The number of particles that were put back into the box

In BoxModus if a particle crosses the wall of the box, it is inserted from the opposite side. However these wall crossings are not performed by this function but in the Experiment constructor when the WallCrossActionsFinder are created. Wall crossings are written to collision output: this is where OutputsList is used.

Definition at line 220 of file boxmodus.cc.

221  {
222  int wraps = 0;
223 
224  for (ParticleData &data : *particles) {
225  FourVector position = data.position();
226  bool wall_hit = enforce_periodic_boundaries(position.begin() + 1,
227  position.end(), length_);
228  if (wall_hit) {
229  const ParticleData incoming_particle(data);
230  data.set_4position(position);
231  ++wraps;
232  ActionPtr action =
233  std::make_unique<WallcrossingAction>(incoming_particle, data);
234  for (const auto &output : output_list) {
235  if (!output->is_dilepton_output() && !output->is_photon_output()) {
236  output->at_interaction(*action, 0.);
237  }
238  }
239  }
240  }
241  logg[LBox].debug("Moved ", wraps, " particles back into the box.");
242  return wraps;
243 }
static bool enforce_periodic_boundaries(Iterator begin, const Iterator &end, typename std::iterator_traits< Iterator >::value_type length)
Enforces periodic boundaries on the given collection of values.
Definition: algorithms.h:53

◆ create_grid()

Grid<GridOptions::PeriodicBoundaries> smash::BoxModus::create_grid ( const Particles particles,
double  min_cell_length,
double  timestep_duration,
CollisionCriterion  crit,
const bool  include_unformed_particles,
CellSizeStrategy  strategy = CellSizeStrategy::Optimal 
) const
inline

Creates the Grid with normal boundary conditions.

Parameters
[in]particlesThe Particles object containing all particles of the currently running Experiment.
[in]min_cell_lengthThe minimal length of the grid cells.
[in]timestep_durationDuration of the timestep. It is necessary for formation times treatment: if particle is fully or partially formed before the end of the timestep, it has to be on the grid.
[in]critCollision criterion (decides if cell number can be limited)
[in]include_unformed_particlesinclude unformed particles from the grid (worsens runtime, necessary for IC output)
[in]strategyThe strategy to determine the cell size
Returns
the Grid object
See also
Grid::Grid

Definition at line 94 of file boxmodus.h.

98  {
100  if (crit == CollisionCriterion::Stochastic) {
102  }
103  return {{{0, 0, 0}, {length_, length_, length_}},
104  particles,
105  min_cell_length,
106  timestep_duration,
107  limit,
108  include_unformed_particles,
109  strategy};
110  }
@ Stochastic
Stochastic Criteiron.
CellNumberLimitation
Identifies whether the number of cells should be limited.
Definition: grid.h:55
@ ParticleNumber
Limit the number of cells to the number of particles.
@ None
No cell number limitation.

◆ create_grandcan_thermalizer()

std::unique_ptr<GrandCanThermalizer> smash::BoxModus::create_grandcan_thermalizer ( Configuration conf) const
inline

Creates GrandCanThermalizer.

(Special Box implementation.)

Parameters
[in]confconfiguration object
Returns
unique pointer to created thermalizer class

Definition at line 118 of file boxmodus.h.

119  {
120  const std::array<double, 3> lat_size = {length_, length_, length_};
121  const std::array<double, 3> origin = {0., 0., 0.};
122  const bool periodicity = true;
123  return std::make_unique<GrandCanThermalizer>(conf, lat_size, origin,
124  periodicity);
125  }

◆ max_timestep()

double smash::BoxModus::max_timestep ( double  max_transverse_distance_sqr) const
inline

Returns
Maximal timestep accepted by this modus. Negative means infinity.

Definition at line 128 of file boxmodus.h.

128  {
129  return 0.5 * std::sqrt(length_ * length_ - max_transverse_distance_sqr);
130  }

◆ equilibration_time()

double smash::BoxModus::equilibration_time ( ) const
inline
Returns
equilibration time of the box

Definition at line 133 of file boxmodus.h.

133 { return equilibration_time_; }

◆ is_box()

bool smash::BoxModus::is_box ( ) const
inline
Returns
whether the modus is box (also, trivially true)

Definition at line 135 of file boxmodus.h.

135 { return true; }

◆ length()

double smash::BoxModus::length ( ) const
inline
Returns
length of the box

Definition at line 137 of file boxmodus.h.

137 { return length_; }

Member Data Documentation

◆ initial_condition_

const BoxInitialCondition smash::BoxModus::initial_condition_
private

Initial momenta distribution: thermal or peaked momenta.

Definition at line 141 of file boxmodus.h.

◆ length_

const double smash::BoxModus::length_
private

Length of the cube's edge in fm.

Definition at line 143 of file boxmodus.h.

◆ equilibration_time_

const double smash::BoxModus::equilibration_time_
private

time after which output is written

Definition at line 145 of file boxmodus.h.

◆ temperature_

const double smash::BoxModus::temperature_
private

Temperature of the Box in GeV.

Definition at line 147 of file boxmodus.h.

◆ start_time_

const double smash::BoxModus::start_time_ = 0.
private

Initial time of the box.

Definition at line 149 of file boxmodus.h.

◆ use_thermal_

const bool smash::BoxModus::use_thermal_ = false
private

Whether to use a thermal initialization for all particles instead of specific numbers.

Definition at line 154 of file boxmodus.h.

◆ mub_

const double smash::BoxModus::mub_
private

Baryon chemical potential for thermal initialization; only used if use_thermal_ is true.

Definition at line 159 of file boxmodus.h.

◆ mus_

const double smash::BoxModus::mus_
private

Strange chemical potential for thermal initialization; only used if use_thermal_ is true.

Definition at line 164 of file boxmodus.h.

◆ muq_

const double smash::BoxModus::muq_
private

Charge chemical potential for thermal initialization; only used if use_thermal_ is true.

Definition at line 169 of file boxmodus.h.

◆ account_for_resonance_widths_

const bool smash::BoxModus::account_for_resonance_widths_
private

In case of thermal initialization: true – account for resonance spectral functions, while computing multiplicities and sampling masses, false – simply use pole masses.

Definition at line 175 of file boxmodus.h.

◆ init_multipl_

const std::map<PdgCode, int> smash::BoxModus::init_multipl_
private

Particle multiplicities at initialization; required if use_thermal_ is false.

Definition at line 180 of file boxmodus.h.

◆ average_multipl_

std::map<PdgCode, double> smash::BoxModus::average_multipl_
private

Average multiplicities in case of thermal initialization.

Saved to avoid recalculating at every event

Definition at line 185 of file boxmodus.h.

◆ jet_pdg_

const std::optional<PdgCode> smash::BoxModus::jet_pdg_
private

Optional PDG code of the particle to use as a jet.

Same setup as in the sphere modus case.

Definition at line 190 of file boxmodus.h.

◆ jet_mom_

const double smash::BoxModus::jet_mom_
private

Initial momentum of the jet particle; only used if insert_jet_ is true.

Definition at line 194 of file boxmodus.h.


The documentation for this class was generated from the following files: