Version: SMASH-3.1
smash::ScatterAction Class Reference

#include <scatteraction.h>

ScatterAction is a special action which takes two incoming particles and performs a scattering, producing one or more final-state particles.

Definition at line 30 of file scatteraction.h.

Inheritance diagram for smash::ScatterAction:
smash::Action smash::BremsstrahlungAction smash::ScatterActionPhoton

Classes

class  InvalidScatterAction
 Thrown when ScatterAction is called to perform with unknown ProcessType. More...
 

Public Member Functions

 ScatterAction (const ParticleData &in_part1, const ParticleData &in_part2, double time, bool isotropic=false, double string_formation_time=1.0, double box_length=-1.0, bool is_total_parametrized=false)
 Construct a ScatterAction object. More...
 
void add_collision (CollisionBranchPtr p)
 Add a new collision channel. More...
 
void add_collisions (CollisionBranchList pv)
 Add several new collision channels at once. More...
 
double transverse_distance_sqr () const
 Calculate the transverse distance of the two incoming particles in their local rest frame. More...
 
double cov_transverse_distance_sqr () const
 Calculate the transverse distance of the two incoming particles in their local rest frame written in a covariant form. More...
 
double mandelstam_s () const
 Determine the Mandelstam s variable,. More...
 
double relative_velocity () const
 Get the relative velocity of the two incoming particles. More...
 
void generate_final_state () override
 Generate the final-state of the scattering process. More...
 
double get_total_weight () const override
 Get the total cross section of scattering particles. More...
 
double get_partial_weight () const override
 Get the partial cross section of the chosen channel. More...
 
void sample_angles (std::pair< double, double > masses, double kinetic_energy_cm) override
 Sample final-state angles in a 2->2 collision (possibly anisotropic). More...
 
void add_all_scatterings (const ScatterActionsFinderParameters &finder_parameters)
 Add all possible scattering subprocesses for this action object. More...
 
void set_parametrized_total_cross_section (const ScatterActionsFinderParameters &finder_parameters)
 Given the incoming particles, assigns the correct parametrization of the total cross section. More...
 
const CollisionBranchList & collision_channels ()
 Get list of possible collision channels. More...
 
void set_string_interface (StringProcess *str_proc)
 Set the StringProcess object to be used. More...
 
virtual double cross_section () const
 Get the total cross section of the scattering particles, either from a parametrization, or from the sum of partials. More...
 
- Public Member Functions inherited from smash::Action
 Action (const ParticleList &in_part, double time)
 Construct an action object with incoming particles and relative time. More...
 
 Action (const ParticleData &in_part, const ParticleData &out_part, double time, ProcessType type)
 Construct an action object with the incoming particles, relative time, and the already known outgoing particles and type of the process. More...
 
 Action (const ParticleList &in_part, const ParticleList &out_part, double absolute_execution_time, ProcessType type)
 Construct an action object with the incoming particles, absolute time, and the already known outgoing particles and type of the process. More...
 
 Action (const Action &)=delete
 Copying is disabled. Use pointers or create a new Action. More...
 
virtual ~Action ()
 Virtual Destructor. More...
 
bool operator< (const Action &rhs) const
 Determine whether one action takes place before another in time. More...
 
virtual ProcessType get_type () const
 Get the process type. More...
 
template<typename Branch >
void add_process (ProcessBranchPtr< Branch > &p, ProcessBranchList< Branch > &subprocesses, double &total_weight)
 Add a new subprocess. More...
 
template<typename Branch >
void add_processes (ProcessBranchList< Branch > pv, ProcessBranchList< Branch > &subprocesses, double &total_weight)
 Add several new subprocesses at once. More...
 
virtual double perform (Particles *particles, uint32_t id_process)
 Actually perform the action, e.g. More...
 
bool is_valid (const Particles &particles) const
 Check whether the action still applies. More...
 
bool is_pauli_blocked (const std::vector< Particles > &ensembles, const PauliBlocker &p_bl) const
 Check if the action is Pauli-blocked. More...
 
const ParticleList & incoming_particles () const
 Get the list of particles that go into the action. More...
 
void update_incoming (const Particles &particles)
 Update the incoming particles that are stored in this action to the state they have in the global particle list. More...
 
const ParticleList & outgoing_particles () const
 Get the list of particles that resulted from the action. More...
 
double time_of_execution () const
 Get the time at which the action is supposed to be performed. More...
 
virtual double check_conservation (const uint32_t id_process) const
 Check various conservation laws. More...
 
double sqrt_s () const
 Determine the total energy in the center-of-mass frame [GeV]. More...
 
FourVector total_momentum_of_outgoing_particles () const
 Calculate the total kinetic momentum of the outgoing particles. More...
 
FourVector get_interaction_point () const
 Get the interaction point. More...
 
std::pair< FourVector, FourVectorget_potential_at_interaction_point () const
 Get the skyrme and asymmetry potential at the interaction point. More...
 
void set_stochastic_pos_idx ()
 Setter function that stores a random incoming particle index latter used to determine the interaction point. More...
 

Protected Member Functions

double cm_momentum () const
 Get the momentum of the center of mass of the incoming particles in the calculation frame. More...
 
double cm_momentum_squared () const
 Get the squared momentum of the center of mass of the incoming particles in the calculation frame. More...
 
ThreeVector beta_cm () const
 Get the velocity of the center of mass of the scattering/incoming particles in the calculation frame. More...
 
double gamma_cm () const
 Get the gamma factor corresponding to a boost to the center of mass frame of the colliding particles. More...
 
void elastic_scattering ()
 Perform an elastic two-body scattering, i.e. just exchange momentum. More...
 
void inelastic_scattering ()
 Perform an inelastic two-body scattering, i.e. new particles are formed. More...
 
void two_to_many_scattering ()
 Perform an inelastic two-to-many-body scattering (more than 2) More...
 
void create_string_final_state ()
 Creates the final states for string-processes after they are performed. More...
 
void string_excitation ()
 Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string excitation and decay. More...
 
void format_debug_output (std::ostream &out) const override
 Writes information about this scatter action to the out stream. More...
 
- Protected Member Functions inherited from smash::Action
FourVector total_momentum () const
 Sum of 4-momenta of incoming particles. More...
 
template<typename Branch >
const Branch * choose_channel (const ProcessBranchList< Branch > &subprocesses, double total_weight)
 Decide for a particular final-state channel via Monte-Carlo and return it as a ProcessBranch. More...
 
virtual std::pair< double, double > sample_masses (double kinetic_energy_cm) const
 Sample final-state masses in general X->2 processes (thus also fixing the absolute c.o.m. More...
 
void sample_2body_phasespace ()
 Sample the full 2-body phase-space (masses, momenta, angles) in the center-of-mass frame for the final state particles. More...
 
virtual void sample_manybody_phasespace ()
 Sample the full n-body phase-space (masses, momenta, angles) in the center-of-mass frame for the final state particles. More...
 
void assign_formation_time_to_outgoing_particles ()
 Assign the formation time to the outgoing particles. More...
 

Protected Attributes

CollisionBranchList collision_channels_
 List of possible collisions. More...
 
double sum_of_partial_cross_sections_
 Current sum of partial hadronic cross sections. More...
 
double partial_cross_section_
 Partial cross-section to the chosen outgoing channel. More...
 
bool isotropic_ = false
 Do this collision isotropically? More...
 
double string_formation_time_ = 1.0
 Time fragments take to be fully formed in hard string excitation. More...
 
- Protected Attributes inherited from smash::Action
ParticleList incoming_particles_
 List with data of incoming particles. More...
 
ParticleList outgoing_particles_
 Initially this stores only the PDG codes of final-state particles. More...
 
const double time_of_execution_
 Time at which the action is supposed to be performed (absolute time in the lab frame in fm). More...
 
ProcessType process_type_
 type of process More...
 
double box_length_ = -1.0
 Box length: needed to determine coordinates of collision correctly in case of collision through the wall. More...
 
int stochastic_position_idx_ = -1
 This stores a randomly-chosen index to an incoming particle. More...
 

Private Member Functions

bool is_elastic () const
 Check if the scattering is elastic. More...
 
void resonance_formation ()
 Perform a 2->1 resonance-formation process. More...
 
void rescale_outgoing_branches ()
 Loop over the possible branches and rescales their weight according to the desired total cross section. More...
 
ParticleTypePtr try_find_pseudoresonance (const PseudoResonance method, const StringTransitionParameters &transition) const
 Try to find a pseudo-resonance that can be created from the incoming particles using a given method. More...
 

Private Attributes

StringProcessstring_process_ = nullptr
 Pointer to interface class for strings. More...
 
bool is_total_parametrized_ = false
 Whether the total cross section is parametrized. More...
 
std::optional< double > parametrized_total_cross_section_ = std::nullopt
 If cross section is parametrized, store the value. More...
 
bool were_processes_added_ = false
 Lock for calling add_all_scatterings only once. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from smash::Action
static double lambda_tilde (double a, double b, double c)
 Little helper function that calculates the lambda function (sometimes written with a tilde to better distinguish it) that appears e.g. More...
 
static void sample_manybody_phasespace_impl (double sqrts, const std::vector< double > &m, std::vector< FourVector > &sampled_momenta)
 Implementation of the full n-body phase-space sampling (masses, momenta, angles) in the center-of-mass frame for the final state particles. More...
 

Constructor & Destructor Documentation

◆ ScatterAction()

smash::ScatterAction::ScatterAction ( const ParticleData in_part1,
const ParticleData in_part2,
double  time,
bool  isotropic = false,
double  string_formation_time = 1.0,
double  box_length = -1.0,
bool  is_total_parametrized = false 
)

Construct a ScatterAction object.

Parameters
[in]in_part1first scattering partner
[in]in_part2second scattering partner
[in]timeTime at which the action is supposed to take place
[in]isotropicif true, do the collision isotropically
[in]string_formation_timeTime string fragments take to form
[in]box_lengthPassing box length to determine coordinate of the collision, in case it happened through the wall in a box. If negative, then there is no wrapping.
[in]is_total_parametrizedWhether the total cross section used for collision finding is parametrized

Definition at line 28 of file scatteraction.cc.

32  : Action({in_part_a, in_part_b}, time),
34  isotropic_(isotropic),
35  string_formation_time_(string_formation_time),
36  is_total_parametrized_(is_total_parametrized) {
37  box_length_ = box_length;
40  }
41 }
Action(const ParticleList &in_part, double time)
Construct an action object with incoming particles and relative time.
Definition: action.h:44
double box_length_
Box length: needed to determine coordinates of collision correctly in case of collision through the w...
Definition: action.h:379
bool isotropic_
Do this collision isotropically?
double string_formation_time_
Time fragments take to be fully formed in hard string excitation.
std::optional< double > parametrized_total_cross_section_
If cross section is parametrized, store the value.
bool is_total_parametrized_
Whether the total cross section is parametrized.
double sum_of_partial_cross_sections_
Current sum of partial hadronic cross sections.

Member Function Documentation

◆ add_collision()

void smash::ScatterAction::add_collision ( CollisionBranchPtr  p)

Add a new collision channel.

Parameters
[in]pChannel to be added.

Definition at line 43 of file scatteraction.cc.

43  {
44  add_process<CollisionBranch>(p, collision_channels_,
46 }
CollisionBranchList collision_channels_
List of possible collisions.
constexpr int p
Proton.

◆ add_collisions()

void smash::ScatterAction::add_collisions ( CollisionBranchList  pv)

Add several new collision channels at once.

Parameters
[in]pvlist of channels to be added.

Definition at line 48 of file scatteraction.cc.

48  {
49  add_processes<CollisionBranch>(std::move(pv), collision_channels_,
51 }

◆ transverse_distance_sqr()

double smash::ScatterAction::transverse_distance_sqr ( ) const

Calculate the transverse distance of the two incoming particles in their local rest frame.

According to UrQMD criterion, Bass:1998ca [5] eq. (3.27):

  • position of particle a: \(\mathbf{x}_a\)
  • position of particle b: \(\mathbf{x}_b\)
  • momentum of particle a: \(\mathbf{p}_a\)
  • momentum of particle b: \(\mathbf{p}_b\)

\[ d^2_\mathrm{coll} = (\mathbf{x}_a - \mathbf{x}_b)^2 - \frac{\bigl[(\mathbf{x}_a - \mathbf{x}_b) \cdot (\mathbf{p}_a - \mathbf{p}_b)\bigr]^2 } {(\mathbf{p}_a - \mathbf{p}_b)^2} \]

Returns
squared distance \(d^2_\mathrm{coll}\).

Definition at line 336 of file scatteraction.cc.

336  {
337  // local copy of particles (since we need to boost them)
338  ParticleData p_a = incoming_particles_[0];
339  ParticleData p_b = incoming_particles_[1];
340  /* Boost particles to center-of-momentum frame. */
341  const ThreeVector velocity = beta_cm();
342  p_a.boost(velocity);
343  p_b.boost(velocity);
344  const ThreeVector pos_diff =
345  p_a.position().threevec() - p_b.position().threevec();
346  const ThreeVector mom_diff =
347  p_a.momentum().threevec() - p_b.momentum().threevec();
348 
349  logg[LScatterAction].debug("Particle ", incoming_particles_,
350  " position difference [fm]: ", pos_diff,
351  ", momentum difference [GeV]: ", mom_diff);
352 
353  const double dp2 = mom_diff.sqr();
354  const double dr2 = pos_diff.sqr();
355  /* Zero momentum leads to infite distance. */
356  if (dp2 < really_small) {
357  return dr2;
358  }
359  const double dpdr = pos_diff * mom_diff;
360 
361  /* UrQMD squared distance criterion, in the center of momentum frame:
362  * position of particle a: x_a
363  * position of particle b: x_b
364  * momentum of particle a: p_a
365  * momentum of particle b: p_b
366  * d^2_{coll} = (x_a - x_b)^2 - ((x_a - x_b) . (p_a - p_b))^2 / (p_a - p_b)^2
367  */
368  const double result = dr2 - dpdr * dpdr / dp2;
369  return result > 0.0 ? result : 0.0;
370 }
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:355
ThreeVector beta_cm() const
Get the velocity of the center of mass of the scattering/incoming particles in the calculation frame.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
static constexpr int LScatterAction

◆ cov_transverse_distance_sqr()

double smash::ScatterAction::cov_transverse_distance_sqr ( ) const

Calculate the transverse distance of the two incoming particles in their local rest frame written in a covariant form.

Equivalent to the UrQMD transverse distance. See Hirano:2012yy [26] (5.6)-(5.11).

Returns
squared distance \(d^2_\mathrm{coll}\).

Definition at line 372 of file scatteraction.cc.

372  {
373  // local copy of particles (since we need to boost them)
374  ParticleData p_a = incoming_particles_[0];
375  ParticleData p_b = incoming_particles_[1];
376 
377  const FourVector delta_x = p_a.position() - p_b.position();
378  const double mom_diff_sqr =
379  (p_a.momentum().threevec() - p_b.momentum().threevec()).sqr();
380  const double x_sqr = delta_x.sqr();
381 
382  if (mom_diff_sqr < really_small) {
383  return -x_sqr;
384  }
385 
386  const double p_a_sqr = p_a.momentum().sqr();
387  const double p_b_sqr = p_b.momentum().sqr();
388  const double p_a_dot_x = p_a.momentum().Dot(delta_x);
389  const double p_b_dot_x = p_b.momentum().Dot(delta_x);
390  const double p_a_dot_p_b = p_a.momentum().Dot(p_b.momentum());
391 
392  const double b_sqr =
393  -x_sqr -
394  (p_a_sqr * std::pow(p_b_dot_x, 2) + p_b_sqr * std::pow(p_a_dot_x, 2) -
395  2 * p_a_dot_p_b * p_a_dot_x * p_b_dot_x) /
396  (std::pow(p_a_dot_p_b, 2) - p_a_sqr * p_b_sqr);
397  return b_sqr > 0.0 ? b_sqr : 0.0;
398 }

◆ mandelstam_s()

double smash::ScatterAction::mandelstam_s ( ) const

Determine the Mandelstam s variable,.

\[s = (p_a + p_b)^2\]

Equal to the square of CMS energy.

Returns
Mandelstam s

Definition at line 313 of file scatteraction.cc.

313 { return total_momentum().sqr(); }
FourVector total_momentum() const
Sum of 4-momenta of incoming particles.
Definition: action.h:389
double sqr() const
calculate the square of the vector (which is a scalar)
Definition: fourvector.h:460

◆ relative_velocity()

double smash::ScatterAction::relative_velocity ( ) const

Get the relative velocity of the two incoming particles.

For a defintion see e.g. Seifert:2017oyb [50], eq. (5)

Returns
relative velocity.

Definition at line 327 of file scatteraction.cc.

327  {
328  const double m1 = incoming_particles()[0].effective_mass();
329  const double m2 = incoming_particles()[1].effective_mass();
330  const double m_s = mandelstam_s();
331  const double lamb = lambda_tilde(m_s, m1 * m1, m2 * m2);
332  return std::sqrt(lamb) / (2. * incoming_particles()[0].momentum().x0() *
333  incoming_particles()[1].momentum().x0());
334 }
const ParticleList & incoming_particles() const
Get the list of particles that go into the action.
Definition: action.cc:58
static double lambda_tilde(double a, double b, double c)
Little helper function that calculates the lambda function (sometimes written with a tilde to better ...
Definition: action.h:315
double mandelstam_s() const
Determine the Mandelstam s variable,.

◆ generate_final_state()

void smash::ScatterAction::generate_final_state ( )
overridevirtual

Generate the final-state of the scattering process.

Performs either elastic or inelastic scattering.

Exceptions
InvalidScatterAction

Implements smash::Action.

Reimplemented in smash::ScatterActionPhoton.

Definition at line 53 of file scatteraction.cc.

53  {
54  logg[LScatterAction].debug("Incoming particles: ", incoming_particles_);
55 
56  const CollisionBranch *proc = choose_channel<CollisionBranch>(
60  process_type_ = proc->get_type();
61  outgoing_particles_ = proc->particle_list();
62  partial_cross_section_ = proc->weight();
63 
64  logg[LScatterAction].debug("Chosen channel: ", process_type_,
66 
67  /* The production point of the new particles. */
68  FourVector middle_point = get_interaction_point();
69 
70  switch (process_type_) {
72  /* 2->2 elastic scattering */
74  break;
76  /* resonance formation */
78  break;
80  /* 2->2 inelastic scattering */
81  /* Sample the particle momenta in CM system. */
83  break;
87  /* 2->m scattering */
89  break;
97  break;
98  default:
99  throw InvalidScatterAction(
100  "ScatterAction::generate_final_state: Invalid process type " +
101  std::to_string(static_cast<int>(process_type_)) + " was requested. " +
102  "(PDGcode1=" + incoming_particles_[0].pdgcode().string() +
103  ", PDGcode2=" + incoming_particles_[1].pdgcode().string() + ")");
104  }
105 
106  for (ParticleData &new_particle : outgoing_particles_) {
107  // Boost to the computational frame
108  new_particle.boost_momentum(
110  /* Set positions of the outgoing particles */
111  if (proc->get_type() != ProcessType::Elastic) {
112  new_particle.set_4position(middle_point);
113  }
114  }
115 }
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
Definition: action.cc:157
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
Definition: action.h:363
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:68
ProcessType process_type_
type of process
Definition: action.h:372
void resonance_formation()
Perform a 2->1 resonance-formation process.
double partial_cross_section_
Partial cross-section to the chosen outgoing channel.
void two_to_many_scattering()
Perform an inelastic two-to-many-body scattering (more than 2)
void string_excitation()
Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string e...
void elastic_scattering()
Perform an elastic two-body scattering, i.e. just exchange momentum.
void inelastic_scattering()
Perform an inelastic two-body scattering, i.e. new particles are formed.
@ TwoToOne
See here for a short description.
@ StringSoftDoubleDiffractive
See here for a short description.
@ TwoToFive
See here for a short description.
@ StringSoftSingleDiffractiveXB
See here for a short description.
@ TwoToTwo
See here for a short description.
@ Elastic
See here for a short description.
@ TwoToFour
See here for a short description.
@ StringSoftAnnihilation
See here for a short description.
@ StringSoftNonDiffractive
See here for a short description.
@ StringSoftSingleDiffractiveAX
See here for a short description.
@ StringHard
See here for a short description.
@ TwoToThree
See here for a short description.

◆ get_total_weight()

double smash::ScatterAction::get_total_weight ( ) const
overridevirtual

Get the total cross section of scattering particles.

Returns
total cross section.

Implements smash::Action.

Reimplemented in smash::ScatterActionPhoton.

Definition at line 294 of file scatteraction.cc.

294  {
296  incoming_particles_[0].xsec_scaling_factor() *
297  incoming_particles_[1].xsec_scaling_factor();
298 }

◆ get_partial_weight()

double smash::ScatterAction::get_partial_weight ( ) const
overridevirtual

Get the partial cross section of the chosen channel.

Returns
partial cross section.

Implements smash::Action.

Definition at line 300 of file scatteraction.cc.

300  {
301  return partial_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
302  incoming_particles_[1].xsec_scaling_factor();
303 }

◆ sample_angles()

void smash::ScatterAction::sample_angles ( std::pair< double, double >  masses,
double  kinetic_energy_cm 
)
overridevirtual

Sample final-state angles in a 2->2 collision (possibly anisotropic).

NN → NN: Choose angular distribution according to Cugnon parametrization, see Cugnon:1996kh [18].

NN → NΔ: Sample scattering angles in center-of-mass frame from an anisotropic angular distribution, using the same distribution as for elastic pp scattering, as suggested in Cugnon:1996kh [18].

NN → NR: Fit to HADES data, see Agakishiev:2014wqa [2].

Reimplemented from smash::Action.

Definition at line 456 of file scatteraction.cc.

457  {
460  // We potentially have more than two particles, so the following angular
461  // distributions don't work. Instead we just keep the angular
462  // distributions generated by string fragmentation.
463  return;
464  }
465  assert(outgoing_particles_.size() == 2);
466 
467  // NN scattering is anisotropic currently
468  const bool nn_scattering = incoming_particles_[0].type().is_nucleon() &&
469  incoming_particles_[1].type().is_nucleon();
470  /* Elastic process is anisotropic and
471  * the angular distribution is based on the NN elastic scattering. */
472  const bool el_scattering = process_type_ == ProcessType::Elastic;
473 
474  const double mass_in_a = incoming_particles_[0].effective_mass();
475  const double mass_in_b = incoming_particles_[1].effective_mass();
476 
477  ParticleData *p_a = &outgoing_particles_[0];
478  ParticleData *p_b = &outgoing_particles_[1];
479 
480  const double mass_a = masses.first;
481  const double mass_b = masses.second;
482 
483  const std::array<double, 2> t_range = get_t_range<double>(
484  kinetic_energy_cm, mass_in_a, mass_in_b, mass_a, mass_b);
485  Angles phitheta;
486  if (el_scattering && !isotropic_) {
490  double mandelstam_s_new = 0.;
491  if (nn_scattering) {
492  mandelstam_s_new = mandelstam_s();
493  } else {
494  /* In the case of elastic collisions other than NN collisions,
495  * there is an ambiguity on how to get the lab-frame momentum (plab),
496  * since the incoming particles can have different masses.
497  * Right now, we first obtain the center-of-mass momentum
498  * of the collision (pcom_now).
499  * Then, the lab-frame momentum is evaluated from the mandelstam s,
500  * which yields the original center-of-mass momentum
501  * when nucleon mass is assumed. */
502  const double pcm_now = pCM_from_s(mandelstam_s(), mass_in_a, mass_in_b);
503  mandelstam_s_new =
504  4. * std::sqrt(pcm_now * pcm_now + nucleon_mass * nucleon_mass);
505  }
506  double bb, a, plab = plab_from_s(mandelstam_s_new);
507  if (nn_scattering &&
508  p_a->pdgcode().antiparticle_sign() ==
509  p_b->pdgcode().antiparticle_sign() &&
510  std::abs(p_a->type().charge() + p_b->type().charge()) == 1) {
511  // proton-neutron and antiproton-antineutron
512  bb = std::max(Cugnon_bnp(plab), really_small);
513  a = (plab < 0.8) ? 1. : 0.64 / (plab * plab);
514  } else {
515  /* all others including pp, nn and AQM elastic processes
516  * This is applied for all particle pairs, which are allowed to
517  * interact elastically. */
518  bb = std::max(Cugnon_bpp(plab), really_small);
519  a = 1.;
520  }
521  double t = random::expo(bb, t_range[0], t_range[1]);
522  if (random::canonical() > 1. / (1. + a)) {
523  t = t_range[0] + t_range[1] - t;
524  }
525  // determine scattering angles in center-of-mass frame
526  phitheta = Angles(2. * M_PI * random::canonical(),
527  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
528  } else if (nn_scattering && p_a->pdgcode().is_Delta() &&
529  p_b->pdgcode().is_nucleon() &&
530  p_a->pdgcode().antiparticle_sign() ==
531  p_b->pdgcode().antiparticle_sign() &&
532  !isotropic_) {
536  const double plab = plab_from_s(mandelstam_s());
537  const double bb = std::max(Cugnon_bpp(plab), really_small);
538  double t = random::expo(bb, t_range[0], t_range[1]);
539  if (random::canonical() > 0.5) {
540  t = t_range[0] + t_range[1] - t; // symmetrize
541  }
542  phitheta = Angles(2. * M_PI * random::canonical(),
543  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
544  } else if (nn_scattering && p_b->pdgcode().is_nucleon() && !isotropic_ &&
545  (p_a->type().is_Nstar() || p_a->type().is_Deltastar())) {
547  const std::array<double, 4> p{1.46434, 5.80311, -6.89358, 1.94302};
548  const double a = p[0] + mass_a * (p[1] + mass_a * (p[2] + mass_a * p[3]));
549  /* If the resonance is so heavy that the index "a" exceeds 30,
550  * the power function turns out to be too sharp. Take t directly to be
551  * t_0 in such a case. */
552  double t = t_range[0];
553  if (a < 30) {
554  t = random::power(-a, t_range[0], t_range[1]);
555  }
556  if (random::canonical() > 0.5) {
557  t = t_range[0] + t_range[1] - t; // symmetrize
558  }
559  phitheta = Angles(2. * M_PI * random::canonical(),
560  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
561  } else {
562  /* isotropic angular distribution */
563  phitheta.distribute_isotropically();
564  }
565 
566  ThreeVector pscatt = phitheta.threevec();
567  // 3-momentum of first incoming particle in center-of-mass frame
568  ThreeVector pcm =
569  incoming_particles_[0].momentum().lorentz_boost(beta_cm()).threevec();
570  pscatt.rotate_z_axis_to(pcm);
571 
572  // final-state CM momentum
573  const double p_f = pCM(kinetic_energy_cm, mass_a, mass_b);
574  if (!(p_f > 0.0)) {
575  logg[LScatterAction].warn("Particle: ", p_a->pdgcode(),
576  " radial momentum: ", p_f);
577  logg[LScatterAction].warn("Etot: ", kinetic_energy_cm, " m_a: ", mass_a,
578  " m_b: ", mass_b);
579  }
580  p_a->set_4momentum(mass_a, pscatt * p_f);
581  p_b->set_4momentum(mass_b, -pscatt * p_f);
582 
583  /* Debug message is printed before boost, so that p_a and p_b are
584  * the momenta in the center of mass frame and thus opposite to
585  * each other.*/
586  logg[LScatterAction].debug("p_a: ", *p_a, "\np_b: ", *p_b);
587 }
T power(T n, T xMin, T xMax)
Draws a random number according to a power-law distribution ~ x^n.
Definition: random.h:203
T expo(T A, T x1, T x2)
Draws a random number x from an exponential distribution exp(A*x), where A is assumed to be positive,...
Definition: random.h:166
T canonical()
Definition: random.h:113
double plab_from_s(double mandelstam_s, double mass)
Convert Mandelstam-s to p_lab in a fixed-target collision.
Definition: kinematics.h:157
static double Cugnon_bnp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic np...
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:79
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:58
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:66
bool is_string_soft_process(ProcessType p)
Check if a given process type is a soft string excitation.
static double Cugnon_bpp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic pp...

◆ add_all_scatterings()

void smash::ScatterAction::add_all_scatterings ( const ScatterActionsFinderParameters finder_parameters)

Add all possible scattering subprocesses for this action object.

This can only be called once per ScatterAction instance.

Parameters
[in]finder_parametersparameters for collision finding.

Definition at line 117 of file scatteraction.cc.

118  {
119  if (were_processes_added_) {
120  logg[LScatterAction].fatal() << "Trying to add processes again.";
121  throw std::logic_error(
122  "add_all_scatterings should be called only once per ScatterAction "
123  "instance");
124  } else {
125  were_processes_added_ = true;
126  }
127  CrossSections xs(incoming_particles_, sqrt_s(),
129  CollisionBranchList processes =
130  xs.generate_collision_list(finder_parameters, string_process_);
131 
132  /* Add various subprocesses.*/
133  add_collisions(std::move(processes));
134 
135  /* If the string processes are not triggered by a probability, then they
136  * always happen as long as the parametrized total cross section is larger
137  * than the sum of the cross sections of the non-string processes, and the
138  * square root s exceeds the threshold by at least 0.9 GeV. The cross section
139  * of the string processes are counted by taking the difference between the
140  * parametrized total and the sum of the non-strings. */
141  if (!finder_parameters.strings_with_probability &&
142  xs.string_probability(finder_parameters)) {
143  const double xs_diff =
144  xs.high_energy(finder_parameters.transition_high_energy) -
146  if (xs_diff > 0.) {
147  add_collisions(xs.string_excitation(xs_diff, string_process_,
148  finder_parameters.use_AQM));
149  }
150  }
151 
152  ParticleTypePtr pseudoresonance =
153  try_find_pseudoresonance(finder_parameters.pseudoresonance_method,
154  finder_parameters.transition_high_energy);
155  if (pseudoresonance && finder_parameters.two_to_one) {
156  const double xs_total =
159  : xs.high_energy(finder_parameters.transition_high_energy);
160  const double xs_gap = xs_total - sum_of_partial_cross_sections_;
161  // The pseudo-resonance is only created if there is a (positive) cross
162  // section gap
163  if (xs_gap > really_small) {
164  auto pseudoresonance_branch = std::make_unique<CollisionBranch>(
165  *pseudoresonance, xs_gap, ProcessType::TwoToOne);
166  add_collision(std::move(pseudoresonance_branch));
167  logg[LScatterAction].debug()
168  << "Pseudoresonance between " << incoming_particles_[0].type().name()
169  << " and " << incoming_particles_[1].type().name() << "is"
170  << pseudoresonance->name() << " with cross section " << xs_gap
171  << "mb.";
172  }
173  }
174  were_processes_added_ = true;
175  // Rescale the branches so that their sum matches the parametrization
178  }
179 }
std::pair< FourVector, FourVector > get_potential_at_interaction_point() const
Get the skyrme and asymmetry potential at the interaction point.
Definition: action.cc:112
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:271
void rescale_outgoing_branches()
Loop over the possible branches and rescales their weight according to the desired total cross sectio...
void add_collision(CollisionBranchPtr p)
Add a new collision channel.
StringProcess * string_process_
Pointer to interface class for strings.
void add_collisions(CollisionBranchList pv)
Add several new collision channels at once.
bool were_processes_added_
Lock for calling add_all_scatterings only once.
ParticleTypePtr try_find_pseudoresonance(const PseudoResonance method, const StringTransitionParameters &transition) const
Try to find a pseudo-resonance that can be created from the incoming particles using a given method.

◆ set_parametrized_total_cross_section()

void smash::ScatterAction::set_parametrized_total_cross_section ( const ScatterActionsFinderParameters finder_parameters)

Given the incoming particles, assigns the correct parametrization of the total cross section.

Parameters
[in]finder_parametersParameters for collision finding.

Definition at line 277 of file scatteraction.cc.

278  {
279  CrossSections xs(incoming_particles_, sqrt_s(),
281 
284  xs.parametrized_total(finder_parameters);
285  } else {
286  logg[LScatterAction].fatal()
287  << "Trying to parametrize total cross section when it shouldn't be.";
288  throw std::logic_error(
289  "This function can only be called on ScatterAction objects with "
290  "parametrized cross section.");
291  }
292 }

◆ collision_channels()

const CollisionBranchList& smash::ScatterAction::collision_channels ( )
inline

Get list of possible collision channels.

Returns
list of possible collision channels.

Definition at line 165 of file scatteraction.h.

165  {
166  return collision_channels_;
167  }

◆ set_string_interface()

void smash::ScatterAction::set_string_interface ( StringProcess str_proc)
inline

Set the StringProcess object to be used.

The StringProcess object is used to handle string excitation and to generate final state particles.

Parameters
[in]str_procString process object to be used.

Definition at line 185 of file scatteraction.h.

185  {
186  string_process_ = str_proc;
187  }

◆ cross_section()

virtual double smash::ScatterAction::cross_section ( ) const
inlinevirtual

Get the total cross section of the scattering particles, either from a parametrization, or from the sum of partials.

Returns
total cross section.

Definition at line 195 of file scatteraction.h.

195  {
198  }
200  }

◆ cm_momentum()

double smash::ScatterAction::cm_momentum ( ) const
protected

Get the momentum of the center of mass of the incoming particles in the calculation frame.

Returns
center of mass momentum.

Definition at line 315 of file scatteraction.cc.

315  {
316  const double m1 = incoming_particles_[0].effective_mass();
317  const double m2 = incoming_particles_[1].effective_mass();
318  return pCM(sqrt_s(), m1, m2);
319 }

◆ cm_momentum_squared()

double smash::ScatterAction::cm_momentum_squared ( ) const
protected

Get the squared momentum of the center of mass of the incoming particles in the calculation frame.

Returns
center of mass momentum squared.

Definition at line 321 of file scatteraction.cc.

321  {
322  const double m1 = incoming_particles_[0].effective_mass();
323  const double m2 = incoming_particles_[1].effective_mass();
324  return pCM_sqr(sqrt_s(), m1, m2);
325 }
T pCM_sqr(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:91

◆ beta_cm()

ThreeVector smash::ScatterAction::beta_cm ( ) const
protected

Get the velocity of the center of mass of the scattering/incoming particles in the calculation frame.

Note: Do not use this function to boost the outgoing particles. Use total_momentum_of_outgoing_particles(), which corrects for the effect of potentials on intial and final state.

Returns
boost velocity between center of mass and calculation frame.

Definition at line 305 of file scatteraction.cc.

305  {
306  return total_momentum().velocity();
307 }
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:333

◆ gamma_cm()

double smash::ScatterAction::gamma_cm ( ) const
protected

Get the gamma factor corresponding to a boost to the center of mass frame of the colliding particles.

Returns
gamma factor.

Definition at line 309 of file scatteraction.cc.

309  {
310  return (1. / std::sqrt(1.0 - beta_cm().sqr()));
311 }

◆ elastic_scattering()

void smash::ScatterAction::elastic_scattering ( )
protected

Perform an elastic two-body scattering, i.e. just exchange momentum.

Definition at line 589 of file scatteraction.cc.

589  {
590  // copy initial particles into final state
593  // resample momenta
594  sample_angles({outgoing_particles_[0].effective_mass(),
595  outgoing_particles_[1].effective_mass()},
596  sqrt_s());
597 }
void sample_angles(std::pair< double, double > masses, double kinetic_energy_cm) override
Sample final-state angles in a 2->2 collision (possibly anisotropic).

◆ inelastic_scattering()

void smash::ScatterAction::inelastic_scattering ( )
protected

Perform an inelastic two-body scattering, i.e. new particles are formed.

Definition at line 599 of file scatteraction.cc.

599  {
600  // create new particles
603 }
void sample_2body_phasespace()
Sample the full 2-body phase-space (masses, momenta, angles) in the center-of-mass frame for the fina...
Definition: action.cc:302
void assign_formation_time_to_outgoing_particles()
Assign the formation time to the outgoing particles.
Definition: action.cc:188

◆ two_to_many_scattering()

void smash::ScatterAction::two_to_many_scattering ( )
protected

Perform an inelastic two-to-many-body scattering (more than 2)

Definition at line 605 of file scatteraction.cc.

605  {
608  logg[LScatterAction].debug("2->", outgoing_particles_.size(),
609  " scattering:", incoming_particles_, " -> ",
611 }
virtual void sample_manybody_phasespace()
Sample the full n-body phase-space (masses, momenta, angles) in the center-of-mass frame for the fina...
Definition: action.cc:446

◆ create_string_final_state()

void smash::ScatterAction::create_string_final_state ( )
protected

Creates the final states for string-processes after they are performed.

Definition at line 635 of file scatteraction.cc.

635  {
638  /* Check momentum difference for debugging */
639  FourVector out_mom;
640  for (ParticleData data : outgoing_particles_) {
641  out_mom += data.momentum();
642  }
643  logg[LPythia].debug("Incoming momenta string:", total_momentum());
644  logg[LPythia].debug("Outgoing momenta string:", out_mom);
645 }
ParticleList get_final_state()
a function to get the final state particle list which is called after the collision
static constexpr int LPythia
Definition: stringprocess.h:26

◆ string_excitation()

void smash::ScatterAction::string_excitation ( )
protected

Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string excitation and decay.

Definition at line 651 of file scatteraction.cc.

651  {
652  assert(incoming_particles_.size() == 2);
653  // Disable floating point exception trap for Pythia
654  {
655  DisableFloatTraps guard;
656  /* initialize the string_process_ object for this particular collision */
658  /* implement collision */
659  bool success = false;
660  int ntry = 0;
661  const int ntry_max = 10000;
662  while (!success && ntry < ntry_max) {
663  ntry++;
664  switch (process_type_) {
666  /* single diffractive to A+X */
667  success = string_process_->next_SDiff(true);
668  break;
670  /* single diffractive to X+B */
671  success = string_process_->next_SDiff(false);
672  break;
674  /* double diffractive */
675  success = string_process_->next_DDiff();
676  break;
678  /* soft non-diffractive */
679  success = string_process_->next_NDiffSoft();
680  break;
682  /* soft BBbar 2 mesonic annihilation */
683  success = string_process_->next_BBbarAnn();
684  break;
686  success = string_process_->next_NDiffHard();
687  break;
688  default:
689  logg[LPythia].error("Unknown string process required.");
690  success = false;
691  }
692  }
693  if (ntry == ntry_max) {
694  /* If pythia fails to form a string, it is usually because the energy
695  * is not large enough. In this case, annihilation is then enforced. If
696  * this process still does not not produce any results, it defaults to
697  * an elastic collision. */
698  bool success_newtry = false;
699 
700  /* Check if the initial state is a baryon-antibaryon state.*/
701  PdgCode part1 = incoming_particles_[0].pdgcode(),
702  part2 = incoming_particles_[1].pdgcode();
703  bool is_BBbar_Pair = (part1.baryon_number() != 0) &&
704  (part1.baryon_number() == -part2.baryon_number());
705 
706  /* Decide on the new process .*/
707  if (is_BBbar_Pair) {
709  } else {
711  }
712  /* Perform the new process*/
713  int ntry_new = 0;
714  while (!success_newtry && ntry_new < ntry_max) {
715  ntry_new++;
716  if (is_BBbar_Pair) {
717  success_newtry = string_process_->next_BBbarAnn();
718  } else {
719  success_newtry = string_process_->next_DDiff();
720  }
721  }
722 
723  if (success_newtry) {
725  }
726 
727  if (!success_newtry) {
728  /* If annihilation fails:
729  * Particles are normally added after process selection for
730  * strings, outgoing_particles is still uninitialized, and memory
731  * needs to be allocated. We also shift the process_type_ to elastic
732  * so that sample_angles does a proper treatment. */
733  outgoing_particles_.reserve(2);
734  outgoing_particles_.push_back(ParticleData{incoming_particles_[0]});
735  outgoing_particles_.push_back(ParticleData{incoming_particles_[1]});
738  }
739  } else {
741  }
742  }
743 }
const double time_of_execution_
Time at which the action is supposed to be performed (absolute time in the lab frame in fm).
Definition: action.h:369
void create_string_final_state()
Creates the final states for string-processes after they are performed.
bool next_SDiff(bool is_AB_to_AX)
Single-diffractive process is based on single pomeron exchange described in Ingelman:1984ns .
bool next_NDiffSoft()
Soft Non-diffractive process is modelled in accordance with dual-topological approach Capella:1978ig ...
bool next_DDiff()
Double-diffractive process ( A + B -> X + X ) is similar to the single-diffractive process,...
bool next_BBbarAnn()
Baryon-antibaryon annihilation process Based on what UrQMD Bass:1998ca , Bleicher:1999xi does,...
void init(const ParticleList &incoming, double tcoll)
initialization feed intial particles, time of collision and gamma factor of the center of mass.
bool next_NDiffHard()
Hard Non-diffractive process is based on PYTHIA 8 with partonic showers and interactions.
@ FailedString
See here for a short description.

◆ is_elastic()

bool smash::ScatterAction::is_elastic ( ) const
private

Check if the scattering is elastic.

Returns
whether the scattering is elastic.

◆ resonance_formation()

void smash::ScatterAction::resonance_formation ( )
private

Perform a 2->1 resonance-formation process.

Exceptions
InvalidResonanceFormation

Definition at line 613 of file scatteraction.cc.

613  {
614  if (outgoing_particles_.size() != 1) {
615  std::string s =
616  "resonance_formation: "
617  "Incorrect number of particles in final state: ";
618  s += std::to_string(outgoing_particles_.size()) + " (";
619  s += incoming_particles_[0].pdgcode().string() + " + ";
620  s += incoming_particles_[1].pdgcode().string() + ")";
621  throw InvalidResonanceFormation(s);
622  }
623  // Set the momentum of the formed resonance in its rest frame.
624  outgoing_particles_[0].set_4momentum(
625  total_momentum_of_outgoing_particles().abs(), 0., 0., 0.);
627  /* this momentum is evaluated in the computational frame. */
628  logg[LScatterAction].debug("Momentum of the new particle: ",
629  outgoing_particles_[0].momentum());
630 }

◆ rescale_outgoing_branches()

void smash::ScatterAction::rescale_outgoing_branches ( )
private

Loop over the possible branches and rescales their weight according to the desired total cross section.

In case the current sum of partials is close to 0, a warning is issued as this would not happen in an usual run, and an elastic process is added to match the total.

Definition at line 181 of file scatteraction.cc.

181  {
182  if (!were_processes_added_) {
183  logg[LScatterAction].fatal()
184  << "Trying to rescale branches before adding processes.";
185  throw std::logic_error(
186  "This function can only be called after having added processes.");
187  }
189  logg[LScatterAction].warn() << "Current total cross section is roughly "
190  "zero and no rescaling to match "
191  "the parametrized one will be done.\nAn "
192  "elastic process will be added,"
193  "instead, to match the total cross section.";
194  auto elastic_branch = std::make_unique<CollisionBranch>(
195  incoming_particles_[0].type(), incoming_particles_[1].type(),
197  add_collision(std::move(elastic_branch));
198  } else {
199  const double reweight =
201  logg[LScatterAction].debug("Reweighting ", sum_of_partial_cross_sections_,
203  for (auto &proc : collision_channels_) {
204  proc->set_weight(proc->weight() * reweight);
205  }
206  }
207 }

◆ try_find_pseudoresonance()

ParticleTypePtr smash::ScatterAction::try_find_pseudoresonance ( const PseudoResonance  method,
const StringTransitionParameters transition 
) const
private

Try to find a pseudo-resonance that can be created from the incoming particles using a given method.

Parameters
[in]methodused to select the pseudo-resonance among possible candidates. See user guide description for more information.
[in]transitionparameters for the string transition region, which are also used to determine when a pseudo-resonance can be created.
Returns
the appropriate pseudo-resonance, if there is any, or an invalid pointer otherwise.

Definition at line 209 of file scatteraction.cc.

211  {
212  const ParticleTypePtr type_a = &incoming_particles_[0].type();
213  const ParticleTypePtr type_b = &incoming_particles_[1].type();
214  const double desired_mass = sqrt_s();
215 
216  double string_offset = 0.5 * transition.sqrts_add_lower;
217  const bool nucleon_and_pion = (type_a->is_nucleon() && type_b->is_pion()) ||
218  (type_a->is_pion() && type_b->is_nucleon());
219  const bool two_nucleons = type_a->is_nucleon() && type_b->is_nucleon();
220  const bool two_pions = type_a->is_pion() && type_b->is_pion();
221  const bool nucleon_and_kaon = (type_a->is_nucleon() && type_b->is_kaon()) ||
222  (type_a->is_kaon() && type_b->is_nucleon());
223  if (nucleon_and_pion) {
224  string_offset = transition.sqrts_range_Npi.first - pion_mass - nucleon_mass;
225  } else if (two_nucleons) {
226  string_offset = transition.sqrts_range_NN.first - 2 * nucleon_mass;
227  } else if (two_pions) {
228  string_offset = transition.pipi_offset;
229  } else if (nucleon_and_kaon) {
230  string_offset = transition.KN_offset;
231  }
232  /*
233  * Artificial cutoff to create a pseudo-resonance only close to the string
234  * transition, where data or first-principle models are unhelpful
235  */
236  if (desired_mass < incoming_particles_[0].effective_mass() +
237  incoming_particles_[1].effective_mass() +
238  string_offset) {
239  return {};
240  }
241 
242  if (method == PseudoResonance::None) {
243  return {};
244  } else if (method == PseudoResonance::LargestFromUnstable ||
246  if (type_a->is_stable() && type_b->is_stable()) {
247  return {};
248  }
249  }
250 
251  // If this list is empty, there are no possible pseudo-resonances.
252  ParticleTypePtrList list = list_possible_resonances(type_a, type_b);
253  if (std::empty(list)) {
254  return {};
255  }
256 
257  if (method == PseudoResonance::Largest ||
259  auto largest = *std::max_element(list.begin(), list.end(),
260  [](ParticleTypePtr a, ParticleTypePtr b) {
261  return a->mass() < b->mass();
262  });
263  return largest;
264  } else if (method == PseudoResonance::Closest ||
266  auto comparison = [&desired_mass](ParticleTypePtr a, ParticleTypePtr b) {
267  return std::abs(a->mass() - desired_mass) <
268  std::abs(b->mass() - desired_mass);
269  };
270  auto closest = *std::min_element(list.begin(), list.end(), comparison);
271  return closest;
272  } else {
273  throw std::logic_error("Unknown method for selecting pseudoresonance.");
274  }
275 }
@ Closest
Resonance with the pole mass closest from the invariant mass of incoming particles for all processes.
@ ClosestFromUnstable
Closest resonance for a given mass from processes with at least one resonance in the incoming particl...
@ None
No pseudo-resonance is created.
@ LargestFromUnstable
Heaviest possible resonance from processes with at least one resonance in the incoming particles.
@ Largest
Resonance of largest mass for all processes.
ParticleTypePtrList list_possible_resonances(const ParticleTypePtr type_a, const ParticleTypePtr type_b)
Lists the possible resonances that decay into two particles.
constexpr double pion_mass
Pion mass in GeV.
Definition: constants.h:65

Member Data Documentation

◆ collision_channels_

CollisionBranchList smash::ScatterAction::collision_channels_
protected

List of possible collisions.

Definition at line 261 of file scatteraction.h.

◆ sum_of_partial_cross_sections_

double smash::ScatterAction::sum_of_partial_cross_sections_
protected

Current sum of partial hadronic cross sections.

Definition at line 264 of file scatteraction.h.

◆ partial_cross_section_

double smash::ScatterAction::partial_cross_section_
protected

Partial cross-section to the chosen outgoing channel.

Definition at line 267 of file scatteraction.h.

◆ isotropic_

bool smash::ScatterAction::isotropic_ = false
protected

Do this collision isotropically?

Definition at line 270 of file scatteraction.h.

◆ string_formation_time_

double smash::ScatterAction::string_formation_time_ = 1.0
protected

Time fragments take to be fully formed in hard string excitation.

Definition at line 273 of file scatteraction.h.

◆ string_process_

StringProcess* smash::ScatterAction::string_process_ = nullptr
private

Pointer to interface class for strings.

Definition at line 313 of file scatteraction.h.

◆ is_total_parametrized_

bool smash::ScatterAction::is_total_parametrized_ = false
private

Whether the total cross section is parametrized.

Definition at line 316 of file scatteraction.h.

◆ parametrized_total_cross_section_

std::optional<double> smash::ScatterAction::parametrized_total_cross_section_ = std::nullopt
private

If cross section is parametrized, store the value.

Definition at line 319 of file scatteraction.h.

◆ were_processes_added_

bool smash::ScatterAction::were_processes_added_ = false
private

Lock for calling add_all_scatterings only once.

Definition at line 322 of file scatteraction.h.


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