Version: SMASH-1.5
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. More...
 
Grid< GridOptions::PeriodicBoundariescreate_grid (const Particles &particles, double min_cell_length, 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
 
- 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
 
double impact_parameter () const
 
double velocity_projectile () const
 
double velocity_target () const
 
FermiMotion fermi_motion () const
 
double max_timestep (double) const
 
double length () const
 
Grid< GridOptions::Normalcreate_grid (const Particles &particles, double min_cell_length, 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 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 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...
 

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

162  : initial_condition_(modus_config.take({"Box", "Initial_Condition"})),
163  length_(modus_config.take({"Box", "Length"})),
164  temperature_(modus_config.take({"Box", "Temperature"})),
165  start_time_(modus_config.take({"Box", "Start_Time"}, 0.)),
166  use_thermal_(
167  modus_config.take({"Box", "Use_Thermal_Multiplicities"}, false)),
168  mub_(modus_config.take({"Box", "Baryon_Chemical_Potential"}, 0.)),
169  mus_(modus_config.take({"Box", "Strange_Chemical_Potential"}, 0.)),
171  ? std::map<PdgCode, int>()
172  : modus_config.take({"Box", "Init_Multiplicities"})
173  .convert_for(init_multipl_)) {}
const double start_time_
Initial time of the box.
Definition: boxmodus.h:130
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126
const std::map< PdgCode, int > init_multipl_
Particle multiplicities at initialization; required if use_thermal_ is false.
Definition: boxmodus.h:150
const bool use_thermal_
Whether to use a thermal initialization for all particles instead of specific numbers.
Definition: boxmodus.h:135
const double temperature_
Temperature of the Box in GeV.
Definition: boxmodus.h:128
const BoxInitialCondition initial_condition_
Initial momenta distribution: thermal or peaked momenta.
Definition: boxmodus.h:124
const double mus_
Strange chemical potential for thermal initialization; only used if use_thermal_ is true...
Definition: boxmodus.h:145
const double mub_
Baryon chemical potential for thermal initialization; only used if use_thermal_ is true...
Definition: boxmodus.h:140

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

176  {
177  const auto &log = logger<LogArea::Box>();
178  double momentum_radial = 0, mass;
179  Angles phitheta;
180  FourVector momentum_total(0, 0, 0, 0);
181  auto uniform_length = random::make_uniform_distribution(0.0, this->length_);
182  const double T = this->temperature_;
183  /* Create NUMBER OF PARTICLES according to configuration, or thermal case */
184  if (use_thermal_) {
185  const double V = length_ * length_ * length_;
186  if (average_multipl_.empty()) {
187  for (const ParticleType &ptype : ParticleType::list_all()) {
188  if (HadronGasEos::is_eos_particle(ptype)) {
189  const double n = HadronGasEos::partial_density(ptype, T, mub_, mus_);
190  average_multipl_[ptype.pdgcode()] = n * V * parameters.testparticles;
191  }
192  }
193  }
194  double nb_init = 0.0, ns_init = 0.0;
195  for (const auto &mult : average_multipl_) {
196  const int thermal_mult_int = random::poisson(mult.second);
197  particles->create(thermal_mult_int, mult.first);
198  nb_init += mult.second * mult.first.baryon_number();
199  ns_init += mult.second * mult.first.strangeness();
200  log.debug(mult.first, " initial multiplicity ", thermal_mult_int);
201  }
202  log.info("Initial hadron gas baryon density ", nb_init);
203  log.info("Initial hadron gas strange density ", ns_init);
204  } else {
205  for (const auto &p : init_multipl_) {
206  particles->create(p.second * parameters.testparticles, p.first);
207  log.debug("Particle ", p.first, " initial multiplicity ", p.second);
208  }
209  }
210 
211  for (ParticleData &data : *particles) {
212  /* Set MOMENTUM SPACE distribution */
214  /* initial thermal momentum is the average 3T */
215  momentum_radial = 3.0 * T;
216  mass = data.pole_mass();
217  } else {
218  /* thermal momentum according Maxwell-Boltzmann distribution */
219  mass = HadronGasEos::sample_mass_thermal(data.type(), 1.0 / T);
220  momentum_radial = sample_momenta_from_thermal(T, mass);
221  }
222  phitheta.distribute_isotropically();
223  log.debug(data.type().name(), "(id ", data.id(), ") radial momentum ",
224  momentum_radial, ", direction", phitheta);
225  data.set_4momentum(mass, phitheta.threevec() * momentum_radial);
226  momentum_total += data.momentum();
227 
228  /* Set COORDINATE SPACE distribution */
229  ThreeVector pos{uniform_length(), uniform_length(), uniform_length()};
230  data.set_4position(FourVector(start_time_, pos));
232  data.set_formation_time(start_time_);
233  }
234 
235  /* Make total 3-momentum 0 */
236  for (ParticleData &data : *particles) {
237  data.set_4momentum(data.momentum().abs(),
238  data.momentum().threevec() -
239  momentum_total.threevec() / particles->size());
240  }
241 
242  /* Recalculate total momentum */
243  momentum_total = FourVector(0, 0, 0, 0);
244  for (ParticleData &data : *particles) {
245  momentum_total += data.momentum();
246  /* IC: debug checks */
247  log.debug() << data;
248  }
249  /* allows to check energy conservation */
250  log.debug() << "Initial total 4-momentum [GeV]: " << momentum_total;
251  return start_time_;
252 }
const double start_time_
Initial time of the box.
Definition: boxmodus.h:130
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126
const std::map< PdgCode, int > init_multipl_
Particle multiplicities at initialization; required if use_thermal_ is false.
Definition: boxmodus.h:150
const bool use_thermal_
Whether to use a thermal initialization for all particles instead of specific numbers.
Definition: boxmodus.h:135
const double temperature_
Temperature of the Box in GeV.
Definition: boxmodus.h:128
const BoxInitialCondition initial_condition_
Initial momenta distribution: thermal or peaked momenta.
Definition: boxmodus.h:124
const double mus_
Strange chemical potential for thermal initialization; only used if use_thermal_ is true...
Definition: boxmodus.h:145
static bool is_eos_particle(const ParticleType &ptype)
Check if a particle belongs to the EoS.
Definition: hadgas_eos.h:277
static const ParticleTypeList & list_all()
Definition: particletype.cc:55
std::map< PdgCode, double > average_multipl_
Average multiplicities in case of thermal initialization.
Definition: boxmodus.h:155
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.
static double sample_mass_thermal(const ParticleType &ptype, double beta)
Sample resonance mass in a thermal medium.
Definition: hadgas_eos.cc:305
uniform_dist< T > make_uniform_distribution(T min, T max)
Definition: random.h:132
constexpr int p
Proton.
static double partial_density(const ParticleType &ptype, double T, double mub, double mus)
Compute partial density of one hadron sort.
Definition: hadgas_eos.cc:219
int poisson(const T &lam)
Returns a Poisson distributed random number.
Definition: random.h:223
constexpr int n
Neutron.
const double mub_
Baryon chemical potential for thermal initialization; only used if use_thermal_ is true...
Definition: boxmodus.h:140
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.

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. Wall crossings are written to collision output: this is where OutputsList is used.

Definition at line 254 of file boxmodus.cc.

255  {
256  const auto &log = logger<LogArea::Box>();
257  int wraps = 0;
258 
259  for (ParticleData &data : *particles) {
260  FourVector position = data.position();
261  bool wall_hit = enforce_periodic_boundaries(position.begin() + 1,
262  position.end(), length_);
263  if (wall_hit) {
264  const ParticleData incoming_particle(data);
265  data.set_4position(position);
266  ++wraps;
267  ActionPtr action =
268  make_unique<WallcrossingAction>(incoming_particle, data);
269  for (const auto &output : output_list) {
270  if (!output->is_dilepton_output() && !output->is_photon_output()) {
271  output->at_interaction(*action, 0.);
272  }
273  }
274  }
275  }
276  log.debug("Moved ", wraps, " particles back into the box.");
277  return wraps;
278 }
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126
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
Here is the call graph for this function:

◆ create_grid()

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

Definition at line 90 of file boxmodus.h.

92  {
93  return {{{0, 0, 0}, {length_, length_, length_}},
94  particles,
95  min_cell_length,
96  strategy};
97  }
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126

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

106  {
107  const std::array<double, 3> lat_size = {length_, length_, length_};
108  const std::array<double, 3> origin = {0., 0., 0.};
109  const bool periodicity = true;
110  return make_unique<GrandCanThermalizer>(conf, lat_size, origin,
111  periodicity);
112  }
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126

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

115  {
116  return 0.5 * std::sqrt(length_ * length_ - max_transverse_distance_sqr);
117  }
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126

◆ length()

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

Definition at line 120 of file boxmodus.h.

120 { return length_; }
const double length_
Length of the cube&#39;s edge in fm/c.
Definition: boxmodus.h:126

Member Data Documentation

◆ initial_condition_

const BoxInitialCondition smash::BoxModus::initial_condition_
private

Initial momenta distribution: thermal or peaked momenta.

Definition at line 124 of file boxmodus.h.

◆ length_

const double smash::BoxModus::length_
private

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

Definition at line 126 of file boxmodus.h.

◆ temperature_

const double smash::BoxModus::temperature_
private

Temperature of the Box in GeV.

Definition at line 128 of file boxmodus.h.

◆ start_time_

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

Initial time of the box.

Definition at line 130 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 135 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 140 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 145 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 150 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 155 of file boxmodus.h.


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