Version: SMASH-3.3
smash::SphereModus Class Reference

#include <spheremodus.h>

SphereModus: Provides a modus for expanding matter calculations.

Matter is put in a sphere of radius R with uniform density; isotropic thermal momenta are typically used for initialization, although other initial momentum states are also included, see Bazow:2016oky [8] and Tindall:2016try [63]

To use this modus, choose

General:
Modus: Sphere

in the configuration file.

Options for SphereModus go in the "Modi"→"Sphere" section of the configuration:

Modi:
Sphere:
# definitions here

The following configuration options are understood: Sphere

Definition at line 49 of file spheremodus.h.

Inheritance diagram for smash::SphereModus:
smash::ModusDefault

Public Member Functions

 SphereModus (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...
 
bool is_sphere () const
 
double radius () 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

double radius_
 Sphere radius (in fm) More...
 
double sphere_temperature_
 Temperature for momentum distribution (in GeV) More...
 
const double start_time_ = 0.
 Starting time for the Sphere. 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 double hf_multiplier_
 Multiplicative factor for thermal multiplicity of heavy flavored hadrons; only used if use_thermal_ is true. More...
 
const bool account_for_resonance_widths_
 In case of thermal initialization: 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 SphereInitialCondition init_distr_
 Initialization scheme for momenta in the sphere; used for expanding metric setup. More...
 
const double radial_velocity_
 Parameter \( u_0\) in the initial flow velocity profile of particles in the sphere, which has the form \( u = u_0 (r / R)^n\). More...
 
const double radial_velocity_exponent_
 Parameter \( n\) in the initial flow velocity profile of particles in the sphere, which has the form \( u = u_0 (r / R)^n\). More...
 
const std::optional< PdgCodejet_pdg_
 Optional PDG code of the particle to use as a jet, i.e. More...
 
const double jet_mom_
 Initial momentum of the jet particle; only used if jet_pdg_ is not nullopt. More...
 
const ThreeVector jet_pos_
 Initial position of the jet particle; only used if jet_pdg_ is not nullopt. More...
 
const bool jet_back_
 Create the back to back jet with the corresponding antiparticle; only used if jet_pdg_ is not nullopt. More...
 
const double jet_back_separation_ = smash_NaN<double>
 Initial separation between the back to back jets; can only be set by the user if jet_back_ is true. More...
 
const SpinInteractionType spin_interaction_type_
 Spin interaction type. More...
 

Friends

std::ostream & operator<< (std::ostream &out, const SphereModus &m)
 Writes the initial state for the Sphere to the output stream. More...
 

Constructor & Destructor Documentation

◆ SphereModus()

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

Constructor.

Takes all there is to take from the (truncated!) configuration object (only contains configuration for this modus).

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 37 of file spheremodus.cc.

39  : radius_(modus_config.take(InputKeys::modi_sphere_radius)),
41  modus_config.take(InputKeys::modi_sphere_temperature)),
53  ? std::map<PdgCode, int>()
54  : modus_config.take(
58  modus_config.take(InputKeys::modi_sphere_addRadialVelocity)),
61  /* Note that it is crucial not to take other keys from the Jet section
62  * before Jet_PDG, since we want here the take to throw in case the user
63  * had a Jet section without the mandatory Jet_PDG key. If all other keys
64  * are taken first, the section is removed from the config because empty,
65  * and has_section(InputSections::m_s_jet) method would return false.
66  */
67  jet_pdg_(modus_config.has_section(InputSections::m_s_jet)
68  ? make_optional<PdgCode>(
69  modus_config.take(InputKeys::modi_sphere_jet_jetPdg))
70  : std::nullopt),
75  jet_back_ ? modus_config.take(
77  : 0),
79  modus_config.take(InputKeys::collTerm_spinInteractions)) {
80  if (!jet_back_ &&
81  modus_config.has_value(InputKeys::modi_sphere_jet_backToBackSeparation)) {
82  throw std::invalid_argument(
83  "In order to specify 'Back_To_Back_Separation', 'Back_To_Back' must be "
84  "true.");
85  }
86  if (radial_velocity_ > 1.0) {
87  throw std::invalid_argument(
88  "Additional velocity cannot be greater than 1!");
89  }
90  if (sphere_temperature_ <= 0.0) {
91  throw std::invalid_argument("Temperature must be positive!");
92  }
93  if (radial_velocity_exponent_ < 0.0) {
94  throw std::invalid_argument("Flow velocity exponent cannot be negative!");
95  }
96 }
const bool account_for_resonance_widths_
In case of thermal initialization:
Definition: spheremodus.h:121
const bool use_thermal_
Whether to use a thermal initialization for all particles instead of specific numbers.
Definition: spheremodus.h:94
const double muq_
Charge chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: spheremodus.h:109
const ThreeVector jet_pos_
Initial position of the jet particle; only used if jet_pdg_ is not nullopt.
Definition: spheremodus.h:162
double sphere_temperature_
Temperature for momentum distribution (in GeV)
Definition: spheremodus.h:87
const double start_time_
Starting time for the Sphere.
Definition: spheremodus.h:89
const double radial_velocity_exponent_
Parameter in the initial flow velocity profile of particles in the sphere, which has the form .
Definition: spheremodus.h:146
const SphereInitialCondition init_distr_
Initialization scheme for momenta in the sphere; used for expanding metric setup.
Definition: spheremodus.h:136
const std::optional< PdgCode > jet_pdg_
Optional PDG code of the particle to use as a jet, i.e.
Definition: spheremodus.h:154
const bool jet_back_
Create the back to back jet with the corresponding antiparticle; only used if jet_pdg_ is not nullopt...
Definition: spheremodus.h:167
const SpinInteractionType spin_interaction_type_
Spin interaction type.
Definition: spheremodus.h:174
const double hf_multiplier_
Multiplicative factor for thermal multiplicity of heavy flavored hadrons; only used if use_thermal_ i...
Definition: spheremodus.h:114
const double mub_
Baryon chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: spheremodus.h:99
const double jet_mom_
Initial momentum of the jet particle; only used if jet_pdg_ is not nullopt.
Definition: spheremodus.h:158
const std::map< PdgCode, int > init_multipl_
Particle multiplicities at initialization; required if use_thermal_ is false.
Definition: spheremodus.h:126
const double jet_back_separation_
Initial separation between the back to back jets; can only be set by the user if jet_back_ is true.
Definition: spheremodus.h:172
const double mus_
Strange chemical potential for thermal initialization; only used if use_thermal_ is true.
Definition: spheremodus.h:104
const double radial_velocity_
Parameter in the initial flow velocity profile of particles in the sphere, which has the form .
Definition: spheremodus.h:141
double radius_
Sphere radius (in fm)
Definition: spheremodus.h:85
static const Key< double > modi_sphere_addRadialVelocity
See user guide description for more information.
Definition: input_keys.h:4248
static const Key< double > modi_sphere_addRadialVelocityExponent
See user guide description for more information.
Definition: input_keys.h:4263
static const Key< std::map< PdgCode, int > > modi_sphere_initialMultiplicities
See user guide description for more information.
Definition: input_keys.h:4170
static const Key< double > modi_sphere_chargeChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4293
static const Key< double > modi_sphere_temperature
See user guide description for more information.
Definition: input_keys.h:4206
static const Key< std::array< double, 3 > > modi_sphere_jet_jetPosition
See user guide description for more information.
Definition: input_keys.h:4427
static const Key< PdgCode > modi_sphere_jet_jetPdg
See user guide description for more information.
Definition: input_keys.h:4401
static const Key< double > modi_sphere_jet_jetMomentum
See user guide description for more information.
Definition: input_keys.h:4413
static const Key< bool > modi_sphere_accountResonanceWidths
See user guide description for more information.
Definition: input_keys.h:4228
static const Key< double > modi_sphere_startTime
See user guide description for more information.
Definition: input_keys.h:4194
static const Key< double > modi_sphere_heavyFlavorMultiplier
See user guide description for more information.
Definition: input_keys.h:4357
static const Key< SphereInitialCondition > modi_sphere_initialCondition
See user guide description for more information.
Definition: input_keys.h:4317
static const Key< double > modi_sphere_radius
See user guide description for more information.
Definition: input_keys.h:4182
static const Key< bool > modi_sphere_jet_backToBack
See user guide description for more information.
Definition: input_keys.h:4443
static const Key< double > modi_sphere_baryonChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4278
static const Key< SpinInteractionType > collTerm_spinInteractions
See user guide description for more information.
Definition: input_keys.h:2549
static const Key< double > modi_sphere_jet_backToBackSeparation
See user guide description for more information.
Definition: input_keys.h:4459
static const Key< bool > modi_sphere_useThermalMultiplicities
See user guide description for more information.
Definition: input_keys.h:4378
static const Key< double > modi_sphere_strangeChemicalPotential
See user guide description for more information.
Definition: input_keys.h:4334
static const Section m_s_jet
Subsection for the jet in sphere modus.
Definition: input_keys.h:137

Member Function Documentation

◆ initial_conditions()

double smash::SphereModus::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.

Susbsequently 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

Definition at line 153 of file spheremodus.cc.

154  {
155  FourVector momentum_total(0, 0, 0, 0);
156  const double T = this->sphere_temperature_;
157  const double V = 4.0 / 3.0 * M_PI * radius_ * radius_ * radius_;
158  /* Create NUMBER OF PARTICLES according to configuration */
159  if (use_thermal_) {
160  if (average_multipl_.empty()) {
161  for (const ParticleType &ptype : ParticleType::list_all()) {
162  const bool is_eos_particle = HadronGasEos::is_eos_particle(ptype);
163  const bool use_heavy_flavor = ptype.pdgcode().is_heavy_flavor() &&
165  if (is_eos_particle || use_heavy_flavor) {
166  const double n = HadronGasEos::partial_density(
168  average_multipl_[ptype.pdgcode()] = n * V * parameters.testparticles;
169  if (ptype.pdgcode().is_heavy_flavor()) {
170  average_multipl_[ptype.pdgcode()] *= hf_multiplier_;
171  }
172  }
173  }
174  }
175  double nb_init = 0.0, ns_init = 0.0, nq_init = 0.0;
176  for (const auto &mult : average_multipl_) {
177  const int thermal_mult_int = random::poisson(mult.second);
178  particles->create(thermal_mult_int, mult.first);
179  nb_init += mult.second * mult.first.baryon_number();
180  ns_init += mult.second * mult.first.strangeness();
181  nq_init += mult.second * mult.first.charge();
182  logg[LSphere].debug(mult.first, " initial multiplicity ",
183  thermal_mult_int);
184  }
185  logg[LSphere].info("Initial hadron gas baryon density ", nb_init);
186  logg[LSphere].info("Initial hadron gas strange density ", ns_init);
187  logg[LSphere].info("Initial hadron gas charge density ", nq_init);
188  } else {
189  for (const auto &p : init_multipl_) {
190  particles->create(p.second * parameters.testparticles, p.first);
191  logg[LSphere].debug("Particle ", p.first, " initial multiplicity ",
192  p.second);
193  }
194  }
195  std::unique_ptr<QuantumSampling> quantum_sampling;
197  quantum_sampling = std::make_unique<QuantumSampling>(init_multipl_, V, T);
198  }
199  /* loop over particle data to fill in momentum and position information */
200  for (ParticleData &data : *particles) {
201  Angles phitheta;
202  /* thermal momentum according Maxwell-Boltzmann distribution */
203  double momentum_radial = 0.0, mass = data.pole_mass();
204  /* assign momentum_radial according to requested distribution */
205  switch (init_distr_) {
207  momentum_radial = sample_momenta_IC_ES(T);
208  break;
210  momentum_radial = sample_momenta_1M_IC(T, mass);
211  break;
213  momentum_radial = sample_momenta_2M_IC(T, mass);
214  break;
216  momentum_radial = sample_momenta_non_eq_mass(T, mass);
217  break;
219  default:
221  ? data.type().mass()
222  : HadronGasEos::sample_mass_thermal(data.type(), 1.0 / T);
223  momentum_radial = sample_momenta_from_thermal(T, mass);
224  break;
226  /*
227  * **********************************************************************
228  * Sampling the thermal momentum according Bose/Fermi/Boltzmann
229  * distribution.
230  * We take the pole mass as the mass.
231  * **********************************************************************
232  */
233  mass = data.type().mass();
234  momentum_radial = quantum_sampling->sample(data.pdgcode());
235  break;
236  }
237  phitheta.distribute_isotropically();
238  logg[LSphere].debug(data.type().name(), "(id ", data.id(),
239  ") radial momentum ", momentum_radial, ", direction",
240  phitheta);
241  data.set_4momentum(mass, phitheta.threevec() * momentum_radial);
242  momentum_total += data.momentum();
243  /* uniform sampling in a sphere with radius r */
244  double position_radial;
245  position_radial = std::cbrt(random::canonical()) * radius_;
246  Angles pos_phitheta;
247  pos_phitheta.distribute_isotropically();
248  data.set_4position(
249  FourVector(start_time_, pos_phitheta.threevec() * position_radial));
250  data.set_formation_time(start_time_);
251  }
252 
253  /* Boost in radial direction with an underlying velocity field of the form
254  * u_r = u_0 * (r / R)^n
255  */
256  if (radial_velocity_ > 0.0) {
257  for (ParticleData &data : *particles) {
258  double particle_radius = std::sqrt(data.position().sqr3());
259  auto e_r = data.position().threevec() / particle_radius;
260  auto radial_velocity =
261  -1.0 * radial_velocity_ * e_r *
262  std::pow(particle_radius / radius_, radial_velocity_exponent_);
263  data.set_4momentum(data.momentum().lorentz_boost(radial_velocity));
264  momentum_total += data.momentum();
265  }
266  }
267 
268  /* Make total 3-momentum in sphere 0 and set unpolarized spin vector if spin
269  * interactions are enabled */
270  for (ParticleData &data : *particles) {
271  data.set_4momentum(data.momentum().abs(),
272  data.momentum().threevec() -
273  momentum_total.threevec() / particles->size());
275  data.set_unpolarized_spin_vector();
276  }
277  }
278 
279  /* Add a single highly energetic particle in the center of the sphere (jet) */
280  if (jet_pdg_) {
281  auto &pdg = jet_pdg_.value();
282  auto &jet_particle = particles->create(pdg);
283  auto displacement = ThreeVector(jet_back_separation_ / 2., 0., 0.);
284  jet_particle.set_formation_time(start_time_);
285  jet_particle.set_4position(
286  FourVector(start_time_, jet_pos_ + displacement));
287  jet_particle.set_4momentum(ParticleType::find(pdg).mass(),
288  ThreeVector(jet_mom_, 0., 0.));
289  if (jet_back_) {
290  auto &anti_pdg = pdg.has_antiparticle() ? pdg.get_antiparticle() : pdg;
291  auto &jet_antiparticle = particles->create(anti_pdg);
292  jet_antiparticle.set_formation_time(start_time_);
293  jet_antiparticle.set_4position(
294  FourVector(start_time_, jet_pos_ - displacement));
295  jet_antiparticle.set_4momentum(ParticleType::find(anti_pdg).mass(),
296  ThreeVector(-jet_mom_, 0., 0.));
297  }
299  jet_particle.set_unpolarized_spin_vector();
300  }
301  }
302 
303  /* Recalculate total momentum */
304  momentum_total = FourVector(0, 0, 0, 0);
305  for (ParticleData &data : *particles) {
306  momentum_total += data.momentum();
307  /* IC: debug checks */
308  logg[LSphere].debug() << data;
309  }
310  /* allows to check energy conservation */
311  logg[LSphere].debug() << "Sphere initial total 4-momentum [GeV]: "
312  << momentum_total;
313  return start_time_;
314 }
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
std::map< PdgCode, double > average_multipl_
Average multiplicities in case of thermal initialization.
Definition: spheremodus.h:131
@ ThermalMomentaBoltzmann
A thermalized ensemble is generated, with momenta sampled from a Maxwell-Boltzmann distribution.
@ IC_ES
Off-equilibrium distribution used in massless comparisons of SMASH to the extended universe metric.
@ ThermalMomentaQuantum
A thermalized ensemble is generated, with momenta of baryons(mesons) sampled from a Fermi(Bose) distr...
@ IC_Massive
A generalization of IC_ES for the non-zero mass case; note that there is currently no analytical comp...
@ IC_2M
Off-equilibrium distribution used in massless comparisons of SMASH to the extended universe metric.
@ IC_1M
Off-equilibrium distribution used in massless comparisons of SMASH to the extended universe metric.
@ Off
No spin interactions.
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
T canonical()
Definition: random.h:113
static constexpr int LSphere
Definition: spheremodus.cc:35
double sample_momenta_from_thermal(const double temperature, const double mass)
Samples a momentum from the Maxwell-Boltzmann (thermal) distribution in a faster way,...
double sample_momenta_IC_ES(const double temperature)
Sample momenta according to the momentum distribution in Bazow:2016oky .
double sample_momenta_non_eq_mass(const double temperature, const double mass)
Samples a momentum via rejection method from the non-equilibrium distribution.
double sample_momenta_1M_IC(const double temperature, const double mass)
Samples a momentum from the non-equilibrium distribution 1M_IC from Bazow:2016oky .
double sample_momenta_2M_IC(const double temperature, const double mass)
Samples a momentum from the non-equilibrium distribution 2M_IC from Bazow:2016oky .
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:41

◆ is_sphere()

bool smash::SphereModus::is_sphere ( ) const
inline
Returns
If the modus is sphere modus, which is always true

Definition at line 79 of file spheremodus.h.

79 { return true; }

◆ radius()

double smash::SphereModus::radius ( ) const
inline
Returns
radius

Definition at line 81 of file spheremodus.h.

81 { return radius_; }

Member Data Documentation

◆ radius_

double smash::SphereModus::radius_
private

Sphere radius (in fm)

Definition at line 85 of file spheremodus.h.

◆ sphere_temperature_

double smash::SphereModus::sphere_temperature_
private

Temperature for momentum distribution (in GeV)

Definition at line 87 of file spheremodus.h.

◆ start_time_

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

Starting time for the Sphere.

Definition at line 89 of file spheremodus.h.

◆ use_thermal_

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

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

Definition at line 94 of file spheremodus.h.

◆ mub_

const double smash::SphereModus::mub_
private

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

Definition at line 99 of file spheremodus.h.

◆ mus_

const double smash::SphereModus::mus_
private

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

Definition at line 104 of file spheremodus.h.

◆ muq_

const double smash::SphereModus::muq_
private

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

Definition at line 109 of file spheremodus.h.

◆ hf_multiplier_

const double smash::SphereModus::hf_multiplier_
private

Multiplicative factor for thermal multiplicity of heavy flavored hadrons; only used if use_thermal_ is true.

Definition at line 114 of file spheremodus.h.

◆ account_for_resonance_widths_

const bool smash::SphereModus::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 121 of file spheremodus.h.

◆ init_multipl_

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

Particle multiplicities at initialization; required if use_thermal_ is false.

Definition at line 126 of file spheremodus.h.

◆ average_multipl_

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

Average multiplicities in case of thermal initialization.

Saved to avoid recalculating at every event

Definition at line 131 of file spheremodus.h.

◆ init_distr_

const SphereInitialCondition smash::SphereModus::init_distr_
private

Initialization scheme for momenta in the sphere; used for expanding metric setup.

Definition at line 136 of file spheremodus.h.

◆ radial_velocity_

const double smash::SphereModus::radial_velocity_
private

Parameter \( u_0\) in the initial flow velocity profile of particles in the sphere, which has the form \( u = u_0 (r / R)^n\).

Definition at line 141 of file spheremodus.h.

◆ radial_velocity_exponent_

const double smash::SphereModus::radial_velocity_exponent_
private

Parameter \( n\) in the initial flow velocity profile of particles in the sphere, which has the form \( u = u_0 (r / R)^n\).

Definition at line 146 of file spheremodus.h.

◆ jet_pdg_

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

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

a single high energy particle at the center (0,0,0) of the expanding sphere. This particle will be placed at t=0 and will initially be moving along the x axis, in the positive or negative direction depending on the sign of its momentum.

Definition at line 154 of file spheremodus.h.

◆ jet_mom_

const double smash::SphereModus::jet_mom_
private

Initial momentum of the jet particle; only used if jet_pdg_ is not nullopt.

Definition at line 158 of file spheremodus.h.

◆ jet_pos_

const ThreeVector smash::SphereModus::jet_pos_
private

Initial position of the jet particle; only used if jet_pdg_ is not nullopt.

Definition at line 162 of file spheremodus.h.

◆ jet_back_

const bool smash::SphereModus::jet_back_
private

Create the back to back jet with the corresponding antiparticle; only used if jet_pdg_ is not nullopt.

Definition at line 167 of file spheremodus.h.

◆ jet_back_separation_

const double smash::SphereModus::jet_back_separation_ = smash_NaN<double>
private

Initial separation between the back to back jets; can only be set by the user if jet_back_ is true.

Definition at line 172 of file spheremodus.h.

◆ spin_interaction_type_

const SpinInteractionType smash::SphereModus::spin_interaction_type_
private

Spin interaction type.

Definition at line 174 of file spheremodus.h.


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