Version: SMASH-1.8
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:
[legend]
Collaboration diagram for smash::BoxModus:
[legend]

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, 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 length () const
 
double equilibration_time () const
 
bool is_box () const
 
- Public Member Functions inherited from smash::ModusDefault
int impose_boundary_conditions (Particles *, const OutputsList &={})
 Enforces sensible positions for the particles. More...
 
int total_N_number () const
 
int proj_N_number () const
 
bool cll_in_nucleus () const
 
bool is_collider () const
 
bool is_box () const
 
bool is_list () const
 
double impact_parameter () const
 
double velocity_projectile () const
 
double velocity_target () const
 
FermiMotion fermi_motion () const
 
double max_timestep (double) const
 
double length () const
 
double equilibration_time () const
 
double nuclei_passing_time () const
 Get the passing time of the two nuclei in a collision. More...
 
Grid< GridOptions::Normalcreate_grid (const Particles &particles, double min_cell_length, double timestep_duration, CellSizeStrategy strategy=CellSizeStrategy::Optimal) const
 Creates the Grid with normal boundary conditions. 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/c. 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 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 bool insert_jet_ = false
 Whether to insert a single high energy particle at the center of the box (0,0,0). More...
 
const PdgCode jet_pdg_
 Pdg of the particle to use as a jet; necessary if insert_jet_ is true, unused otherwise. 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)
 

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 214 of file boxmodus.cc.

216  : initial_condition_(modus_config.take({"Box", "Initial_Condition"})),
217  length_(modus_config.take({"Box", "Length"})),
219  modus_config.take({"Box", "Equilibration_Time"}, -1.)),
220  temperature_(modus_config.take({"Box", "Temperature"})),
221  start_time_(modus_config.take({"Box", "Start_Time"}, 0.)),
222  use_thermal_(
223  modus_config.take({"Box", "Use_Thermal_Multiplicities"}, false)),
224  mub_(modus_config.take({"Box", "Baryon_Chemical_Potential"}, 0.)),
225  mus_(modus_config.take({"Box", "Strange_Chemical_Potential"}, 0.)),
227  modus_config.take({"Box", "Account_Resonance_Widths"}, true)),
229  ? std::map<PdgCode, int>()
230  : modus_config.take({"Box", "Init_Multiplicities"})
231  .convert_for(init_multipl_)),
232  insert_jet_(modus_config.has_value({"Box", "Jet", "Jet_PDG"})),
233  jet_pdg_(insert_jet_ ? modus_config.take({"Box", "Jet", "Jet_PDG"})
234  .convert_for(jet_pdg_)
235  : pdg::p), // dummy default; never used
236  jet_mom_(modus_config.take({"Box", "Jet", "Jet_Momentum"}, 20.)) {
237  if (parameters.res_lifetime_factor < 0.) {
238  throw std::invalid_argument(
239  "Resonance lifetime modifier cannot be negative!");
240  }
241 }

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 243 of file boxmodus.cc.

244  {
245  double momentum_radial = 0, mass;
246  Angles phitheta;
247  FourVector momentum_total(0, 0, 0, 0);
248  auto uniform_length = random::make_uniform_distribution(0.0, this->length_);
249  const double T = this->temperature_;
250  /* Create NUMBER OF PARTICLES according to configuration, or thermal case */
251  if (use_thermal_) {
252  const double V = length_ * length_ * length_;
253  if (average_multipl_.empty()) {
254  for (const ParticleType &ptype : ParticleType::list_all()) {
255  if (HadronGasEos::is_eos_particle(ptype)) {
256  const double lifetime_factor =
257  ptype.is_stable() ? 1. : parameters.res_lifetime_factor;
258  const double n = lifetime_factor * HadronGasEos::partial_density(
259  ptype, T, mub_, mus_,
261  average_multipl_[ptype.pdgcode()] = n * V * parameters.testparticles;
262  }
263  }
264  }
265  double nb_init = 0.0, ns_init = 0.0;
266  for (const auto &mult : average_multipl_) {
267  const int thermal_mult_int = random::poisson(mult.second);
268  particles->create(thermal_mult_int, mult.first);
269  nb_init += mult.second * mult.first.baryon_number();
270  ns_init += mult.second * mult.first.strangeness();
271  logg[LBox].debug(mult.first, " initial multiplicity ", thermal_mult_int);
272  }
273  logg[LBox].info("Initial hadron gas baryon density ", nb_init);
274  logg[LBox].info("Initial hadron gas strange density ", ns_init);
275  } else {
276  for (const auto &p : init_multipl_) {
277  particles->create(p.second * parameters.testparticles, p.first);
278  logg[LBox].debug("Particle ", p.first, " initial multiplicity ",
279  p.second);
280  }
281  }
282 
283  for (ParticleData &data : *particles) {
284  /* Set MOMENTUM SPACE distribution */
286  /* initial thermal momentum is the average 3T */
287  momentum_radial = 3.0 * T;
288  mass = data.pole_mass();
289  } else {
290  /* thermal momentum according Maxwell-Boltzmann distribution */
292  ? data.type().mass()
293  : HadronGasEos::sample_mass_thermal(data.type(), 1.0 / T);
294  momentum_radial = sample_momenta_from_thermal(T, mass);
295  }
296  phitheta.distribute_isotropically();
297  logg[LBox].debug(data.type().name(), "(id ", data.id(),
298  ") radial momentum ", momentum_radial, ", direction",
299  phitheta);
300  data.set_4momentum(mass, phitheta.threevec() * momentum_radial);
301  momentum_total += data.momentum();
302 
303  /* Set COORDINATE SPACE distribution */
304  ThreeVector pos{uniform_length(), uniform_length(), uniform_length()};
305  data.set_4position(FourVector(start_time_, pos));
307  data.set_formation_time(start_time_);
308  }
309 
310  /* Make total 3-momentum 0 */
311  for (ParticleData &data : *particles) {
312  data.set_4momentum(data.momentum().abs(),
313  data.momentum().threevec() -
314  momentum_total.threevec() / particles->size());
315  }
316 
317  /* Add a single highly energetic particle in the center of the box (jet) */
318  if (insert_jet_) {
319  auto &jet_particle = particles->create(jet_pdg_);
320  jet_particle.set_formation_time(start_time_);
321  jet_particle.set_4position(FourVector(start_time_, 0., 0., 0.));
322  jet_particle.set_4momentum(ParticleType::find(jet_pdg_).mass(),
323  ThreeVector(jet_mom_, 0., 0.));
324  }
325 
326  /* Recalculate total momentum */
327  momentum_total = FourVector(0, 0, 0, 0);
328  for (ParticleData &data : *particles) {
329  momentum_total += data.momentum();
330  /* IC: debug checks */
331  logg[LBox].debug() << data;
332  }
333  /* allows to check energy conservation */
334  logg[LBox].debug() << "Initial total 4-momentum [GeV]: " << momentum_total;
335  return start_time_;
336 }
Here is the call graph for this function:

◆ 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 338 of file boxmodus.cc.

339  {
340  int wraps = 0;
341 
342  for (ParticleData &data : *particles) {
343  FourVector position = data.position();
344  bool wall_hit = enforce_periodic_boundaries(position.begin() + 1,
345  position.end(), length_);
346  if (wall_hit) {
347  const ParticleData incoming_particle(data);
348  data.set_4position(position);
349  ++wraps;
350  ActionPtr action =
351  make_unique<WallcrossingAction>(incoming_particle, data);
352  for (const auto &output : output_list) {
353  if (!output->is_dilepton_output() && !output->is_photon_output()) {
354  output->at_interaction(*action, 0.);
355  }
356  }
357  }
358  }
359  logg[LBox].debug("Moved ", wraps, " particles back into the box.");
360  return wraps;
361 }
Here is the call graph for this function:

◆ create_grid()

Grid<GridOptions::PeriodicBoundaries> smash::BoxModus::create_grid ( const Particles particles,
double  min_cell_length,
double  timestep_duration,
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]strategyThe strategy to determine the cell size
Returns
the Grid object
See also
Grid::Grid

Definition at line 94 of file boxmodus.h.

97  {
98  return {{{0, 0, 0}, {length_, length_, length_}},
99  particles,
100  min_cell_length,
101  timestep_duration,
102  strategy};
103  }

◆ 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 111 of file boxmodus.h.

112  {
113  const std::array<double, 3> lat_size = {length_, length_, length_};
114  const std::array<double, 3> origin = {0., 0., 0.};
115  const bool periodicity = true;
116  return make_unique<GrandCanThermalizer>(conf, lat_size, origin,
117  periodicity);
118  }

◆ 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 121 of file boxmodus.h.

121  {
122  return 0.5 * std::sqrt(length_ * length_ - max_transverse_distance_sqr);
123  }

◆ length()

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

Definition at line 126 of file boxmodus.h.

126 { return length_; }

◆ equilibration_time()

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

Definition at line 128 of file boxmodus.h.

128 { return equilibration_time_; }

◆ is_box()

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

Definition at line 130 of file boxmodus.h.

130 { return true; }

Member Data Documentation

◆ initial_condition_

const BoxInitialCondition smash::BoxModus::initial_condition_
private

Initial momenta distribution: thermal or peaked momenta.

Definition at line 134 of file boxmodus.h.

◆ length_

const double smash::BoxModus::length_
private

Length of the cube's edge in fm/c.

Definition at line 136 of file boxmodus.h.

◆ equilibration_time_

const double smash::BoxModus::equilibration_time_
private

time after which output is written

Definition at line 138 of file boxmodus.h.

◆ temperature_

const double smash::BoxModus::temperature_
private

Temperature of the Box in GeV.

Definition at line 140 of file boxmodus.h.

◆ start_time_

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

Initial time of the box.

Definition at line 142 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 147 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 152 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 157 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 163 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 168 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 173 of file boxmodus.h.

◆ insert_jet_

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

Whether to insert a single high energy particle at the center of the box (0,0,0).

This particle will initially be moving along the x-axis.

Definition at line 179 of file boxmodus.h.

◆ jet_pdg_

const PdgCode smash::BoxModus::jet_pdg_
private

Pdg of the particle to use as a jet; necessary if insert_jet_ is true, unused otherwise.

Definition at line 184 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 188 of file boxmodus.h.


The documentation for this class was generated from the following files:
smash::BoxModus::account_for_resonance_widths_
const bool account_for_resonance_widths_
In case of thermal initialization: true – account for resonance spectral functions,...
Definition: boxmodus.h:163
smash::BoxModus::mub_
const double mub_
Baryon chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: boxmodus.h:152
smash::HadronGasEos::partial_density
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.
Definition: hadgas_eos.cc:234
smash::BoxModus::init_multipl_
const std::map< PdgCode, int > init_multipl_
Particle multiplicities at initialization; required if use_thermal_ is false.
Definition: boxmodus.h:168
smash::BoxModus::length_
const double length_
Length of the cube's edge in fm/c.
Definition: boxmodus.h:136
smash::BoxModus::use_thermal_
const bool use_thermal_
Whether to use a thermal initialization for all particles instead of specific numbers.
Definition: boxmodus.h:147
smash::BoxModus::jet_pdg_
const PdgCode jet_pdg_
Pdg of the particle to use as a jet; necessary if insert_jet_ is true, unused otherwise.
Definition: boxmodus.h:184
smash::BoxModus::temperature_
const double temperature_
Temperature of the Box in GeV.
Definition: boxmodus.h:140
smash::enforce_periodic_boundaries
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
smash::BoxModus::insert_jet_
const bool insert_jet_
Whether to insert a single high energy particle at the center of the box (0,0,0).
Definition: boxmodus.h:179
smash::BoxModus::initial_condition_
const BoxInitialCondition initial_condition_
Initial momenta distribution: thermal or peaked momenta.
Definition: boxmodus.h:134
smash::logg
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
smash::ParticleType::find
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
Definition: particletype.cc:105
smash::BoxModus::mus_
const double mus_
Strange chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: boxmodus.h:157
smash::HadronGasEos::sample_mass_thermal
static double sample_mass_thermal(const ParticleType &ptype, double beta)
Sample resonance mass in a thermal medium.
Definition: hadgas_eos.cc:325
smash::BoxModus::average_multipl_
std::map< PdgCode, double > average_multipl_
Average multiplicities in case of thermal initialization.
Definition: boxmodus.h:173
smash::random::poisson
int poisson(const T &lam)
Returns a Poisson distributed random number.
Definition: random.h:226
smash::random::make_uniform_distribution
uniform_dist< T > make_uniform_distribution(T min, T max)
Definition: random.h:135
smash::LBox
static constexpr int LBox
Definition: boxmodus.cc:35
smash::HadronGasEos::is_eos_particle
static bool is_eos_particle(const ParticleType &ptype)
Check if a particle belongs to the EoS.
Definition: hadgas_eos.h:308
BoxInitialCondition::PeakedMomenta
smash::pdg::p
constexpr int p
Proton.
Definition: pdgcode_constants.h:28
smash::pdg::n
constexpr int n
Neutron.
Definition: pdgcode_constants.h:30
smash::sample_momenta_from_thermal
double sample_momenta_from_thermal(const double temperature, const double mass)
Samples a momentum from the Maxwell-Boltzmann (thermal) distribution in a faster way,...
Definition: distributions.cc:190
smash::BoxModus::equilibration_time_
const double equilibration_time_
time after which output is written
Definition: boxmodus.h:138
smash::BoxModus::start_time_
const double start_time_
Initial time of the box.
Definition: boxmodus.h:142
smash::ParticleType::list_all
static const ParticleTypeList & list_all()
Definition: particletype.cc:57
smash::BoxModus::jet_mom_
const double jet_mom_
Initial momentum of the jet particle; only used if insert_jet_ is true.
Definition: boxmodus.h:188