Version: SMASH-3.0
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)
 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...
 
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. 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 total_cross_section_
 Total hadronic cross section. 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...
 

Private Attributes

StringProcessstring_process_ = nullptr
 Pointer to interface class for strings. 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 
)

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.

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  box_length_ = box_length;
37 }
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:369
bool isotropic_
Do this collision isotropically?
double string_formation_time_
Time fragments take to be fully formed in hard string excitation.
double total_cross_section_
Total hadronic cross section.

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 39 of file scatteraction.cc.

39  {
40  add_process<CollisionBranch>(p, collision_channels_, total_cross_section_);
41 }
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 43 of file scatteraction.cc.

43  {
44  add_processes<CollisionBranch>(std::move(pv), collision_channels_,
46 }

◆ 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 180 of file scatteraction.cc.

180  {
181  // local copy of particles (since we need to boost them)
182  ParticleData p_a = incoming_particles_[0];
183  ParticleData p_b = incoming_particles_[1];
184  /* Boost particles to center-of-momentum frame. */
185  const ThreeVector velocity = beta_cm();
186  p_a.boost(velocity);
187  p_b.boost(velocity);
188  const ThreeVector pos_diff =
189  p_a.position().threevec() - p_b.position().threevec();
190  const ThreeVector mom_diff =
191  p_a.momentum().threevec() - p_b.momentum().threevec();
192 
193  logg[LScatterAction].debug("Particle ", incoming_particles_,
194  " position difference [fm]: ", pos_diff,
195  ", momentum difference [GeV]: ", mom_diff);
196 
197  const double dp2 = mom_diff.sqr();
198  const double dr2 = pos_diff.sqr();
199  /* Zero momentum leads to infite distance. */
200  if (dp2 < really_small) {
201  return dr2;
202  }
203  const double dpdr = pos_diff * mom_diff;
204 
205  /* UrQMD squared distance criterion, in the center of momentum frame:
206  * position of particle a: x_a
207  * position of particle b: x_b
208  * momentum of particle a: p_a
209  * momentum of particle b: p_b
210  * d^2_{coll} = (x_a - x_b)^2 - ((x_a - x_b) . (p_a - p_b))^2 / (p_a - p_b)^2
211  */
212  const double result = dr2 - dpdr * dpdr / dp2;
213  return result > 0.0 ? result : 0.0;
214 }
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:345
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 216 of file scatteraction.cc.

216  {
217  // local copy of particles (since we need to boost them)
218  ParticleData p_a = incoming_particles_[0];
219  ParticleData p_b = incoming_particles_[1];
220 
221  const FourVector delta_x = p_a.position() - p_b.position();
222  const double mom_diff_sqr =
223  (p_a.momentum().threevec() - p_b.momentum().threevec()).sqr();
224  const double x_sqr = delta_x.sqr();
225 
226  if (mom_diff_sqr < really_small) {
227  return -x_sqr;
228  }
229 
230  const double p_a_sqr = p_a.momentum().sqr();
231  const double p_b_sqr = p_b.momentum().sqr();
232  const double p_a_dot_x = p_a.momentum().Dot(delta_x);
233  const double p_b_dot_x = p_b.momentum().Dot(delta_x);
234  const double p_a_dot_p_b = p_a.momentum().Dot(p_b.momentum());
235 
236  const double b_sqr =
237  -x_sqr -
238  (p_a_sqr * std::pow(p_b_dot_x, 2) + p_b_sqr * std::pow(p_a_dot_x, 2) -
239  2 * p_a_dot_p_b * p_a_dot_x * p_b_dot_x) /
240  (std::pow(p_a_dot_p_b, 2) - p_a_sqr * p_b_sqr);
241  return b_sqr > 0.0 ? b_sqr : 0.0;
242 }

◆ 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 157 of file scatteraction.cc.

157 { return total_momentum().sqr(); }
FourVector total_momentum() const
Sum of 4-momenta of incoming particles.
Definition: action.h:379
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 171 of file scatteraction.cc.

171  {
172  const double m1 = incoming_particles()[0].effective_mass();
173  const double m2 = incoming_particles()[1].effective_mass();
174  const double m_s = mandelstam_s();
175  const double lamb = lambda_tilde(m_s, m1 * m1, m2 * m2);
176  return std::sqrt(lamb) / (2. * incoming_particles()[0].momentum().x0() *
177  incoming_particles()[1].momentum().x0());
178 }
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 48 of file scatteraction.cc.

48  {
49  logg[LScatterAction].debug("Incoming particles: ", incoming_particles_);
50 
51  /* Decide for a particular final state. */
52  const CollisionBranch *proc = choose_channel<CollisionBranch>(
54  process_type_ = proc->get_type();
55  outgoing_particles_ = proc->particle_list();
56  partial_cross_section_ = proc->weight();
57 
58  logg[LScatterAction].debug("Chosen channel: ", process_type_,
60 
61  /* The production point of the new particles. */
62  FourVector middle_point = get_interaction_point();
63 
64  switch (process_type_) {
66  /* 2->2 elastic scattering */
68  break;
70  /* resonance formation */
72  break;
74  /* 2->2 inelastic scattering */
75  /* Sample the particle momenta in CM system. */
77  break;
81  /* 2->m scattering */
83  break;
91  break;
92  default:
93  throw InvalidScatterAction(
94  "ScatterAction::generate_final_state: Invalid process type " +
95  std::to_string(static_cast<int>(process_type_)) + " was requested. " +
96  "(PDGcode1=" + incoming_particles_[0].pdgcode().string() +
97  ", PDGcode2=" + incoming_particles_[1].pdgcode().string() + ")");
98  }
99 
100  for (ParticleData &new_particle : outgoing_particles_) {
101  // Boost to the computational frame
102  new_particle.boost_momentum(
104  /* Set positions of the outgoing particles */
105  if (proc->get_type() != ProcessType::Elastic) {
106  new_particle.set_4position(middle_point);
107  }
108  }
109 }
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:353
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:68
ProcessType process_type_
type of process
Definition: action.h:362
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 139 of file scatteraction.cc.

139  {
140  return total_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
141  incoming_particles_[1].xsec_scaling_factor();
142 }

◆ 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 144 of file scatteraction.cc.

144  {
145  return partial_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
146  incoming_particles_[1].xsec_scaling_factor();
147 }

◆ 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 300 of file scatteraction.cc.

301  {
304  // We potentially have more than two particles, so the following angular
305  // distributions don't work. Instead we just keep the angular
306  // distributions generated by string fragmentation.
307  return;
308  }
309  assert(outgoing_particles_.size() == 2);
310 
311  // NN scattering is anisotropic currently
312  const bool nn_scattering = incoming_particles_[0].type().is_nucleon() &&
313  incoming_particles_[1].type().is_nucleon();
314  /* Elastic process is anisotropic and
315  * the angular distribution is based on the NN elastic scattering. */
316  const bool el_scattering = process_type_ == ProcessType::Elastic;
317 
318  const double mass_in_a = incoming_particles_[0].effective_mass();
319  const double mass_in_b = incoming_particles_[1].effective_mass();
320 
321  ParticleData *p_a = &outgoing_particles_[0];
322  ParticleData *p_b = &outgoing_particles_[1];
323 
324  const double mass_a = masses.first;
325  const double mass_b = masses.second;
326 
327  const std::array<double, 2> t_range = get_t_range<double>(
328  kinetic_energy_cm, mass_in_a, mass_in_b, mass_a, mass_b);
329  Angles phitheta;
330  if (el_scattering && !isotropic_) {
334  double mandelstam_s_new = 0.;
335  if (nn_scattering) {
336  mandelstam_s_new = mandelstam_s();
337  } else {
338  /* In the case of elastic collisions other than NN collisions,
339  * there is an ambiguity on how to get the lab-frame momentum (plab),
340  * since the incoming particles can have different masses.
341  * Right now, we first obtain the center-of-mass momentum
342  * of the collision (pcom_now).
343  * Then, the lab-frame momentum is evaluated from the mandelstam s,
344  * which yields the original center-of-mass momentum
345  * when nucleon mass is assumed. */
346  const double pcm_now = pCM_from_s(mandelstam_s(), mass_in_a, mass_in_b);
347  mandelstam_s_new =
348  4. * std::sqrt(pcm_now * pcm_now + nucleon_mass * nucleon_mass);
349  }
350  double bb, a, plab = plab_from_s(mandelstam_s_new);
351  if (nn_scattering &&
352  p_a->pdgcode().antiparticle_sign() ==
353  p_b->pdgcode().antiparticle_sign() &&
354  std::abs(p_a->type().charge() + p_b->type().charge()) == 1) {
355  // proton-neutron and antiproton-antineutron
356  bb = std::max(Cugnon_bnp(plab), really_small);
357  a = (plab < 0.8) ? 1. : 0.64 / (plab * plab);
358  } else {
359  /* all others including pp, nn and AQM elastic processes
360  * This is applied for all particle pairs, which are allowed to
361  * interact elastically. */
362  bb = std::max(Cugnon_bpp(plab), really_small);
363  a = 1.;
364  }
365  double t = random::expo(bb, t_range[0], t_range[1]);
366  if (random::canonical() > 1. / (1. + a)) {
367  t = t_range[0] + t_range[1] - t;
368  }
369  // determine scattering angles in center-of-mass frame
370  phitheta = Angles(2. * M_PI * random::canonical(),
371  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
372  } else if (nn_scattering && p_a->pdgcode().is_Delta() &&
373  p_b->pdgcode().is_nucleon() &&
374  p_a->pdgcode().antiparticle_sign() ==
375  p_b->pdgcode().antiparticle_sign() &&
376  !isotropic_) {
380  const double plab = plab_from_s(mandelstam_s());
381  const double bb = std::max(Cugnon_bpp(plab), really_small);
382  double t = random::expo(bb, t_range[0], t_range[1]);
383  if (random::canonical() > 0.5) {
384  t = t_range[0] + t_range[1] - t; // symmetrize
385  }
386  phitheta = Angles(2. * M_PI * random::canonical(),
387  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
388  } else if (nn_scattering && p_b->pdgcode().is_nucleon() && !isotropic_ &&
389  (p_a->type().is_Nstar() || p_a->type().is_Deltastar())) {
391  const std::array<double, 4> p{1.46434, 5.80311, -6.89358, 1.94302};
392  const double a = p[0] + mass_a * (p[1] + mass_a * (p[2] + mass_a * p[3]));
393  /* If the resonance is so heavy that the index "a" exceeds 30,
394  * the power function turns out to be too sharp. Take t directly to be
395  * t_0 in such a case. */
396  double t = t_range[0];
397  if (a < 30) {
398  t = random::power(-a, t_range[0], t_range[1]);
399  }
400  if (random::canonical() > 0.5) {
401  t = t_range[0] + t_range[1] - t; // symmetrize
402  }
403  phitheta = Angles(2. * M_PI * random::canonical(),
404  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
405  } else {
406  /* isotropic angular distribution */
407  phitheta.distribute_isotropically();
408  }
409 
410  ThreeVector pscatt = phitheta.threevec();
411  // 3-momentum of first incoming particle in center-of-mass frame
412  ThreeVector pcm =
413  incoming_particles_[0].momentum().lorentz_boost(beta_cm()).threevec();
414  pscatt.rotate_z_axis_to(pcm);
415 
416  // final-state CM momentum
417  const double p_f = pCM(kinetic_energy_cm, mass_a, mass_b);
418  if (!(p_f > 0.0)) {
419  logg[LScatterAction].warn("Particle: ", p_a->pdgcode(),
420  " radial momentum: ", p_f);
421  logg[LScatterAction].warn("Etot: ", kinetic_energy_cm, " m_a: ", mass_a,
422  " m_b: ", mass_b);
423  }
424  p_a->set_4momentum(mass_a, pscatt * p_f);
425  p_b->set_4momentum(mass_b, -pscatt * p_f);
426 
427  /* Debug message is printed before boost, so that p_a and p_b are
428  * the momenta in the center of mass frame and thus opposite to
429  * each other.*/
430  logg[LScatterAction].debug("p_a: ", *p_a, "\np_b: ", *p_b);
431 }
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.

Parameters
[in]finder_parametersparameters for collision finding.

Definition at line 111 of file scatteraction.cc.

112  {
113  CrossSections xs(incoming_particles_, sqrt_s(),
115  CollisionBranchList processes =
116  xs.generate_collision_list(finder_parameters, string_process_);
117 
118  /* Add various subprocesses.*/
119  add_collisions(std::move(processes));
120 
121  /* If the string processes are not triggered by a probability, then they
122  * always happen as long as the parametrized total cross section is larger
123  * than the sum of the cross sections of the non-string processes, and the
124  * square root s exceeds the threshold by at least 0.9 GeV. The cross section
125  * of the string processes are counted by taking the difference between the
126  * parametrized total and the sum of the non-strings. */
127  if (!finder_parameters.strings_with_probability &&
128  xs.string_probability(finder_parameters)) {
129  const double xs_diff =
130  xs.high_energy(finder_parameters.transition_high_energy) -
131  cross_section();
132  if (xs_diff > 0.) {
133  add_collisions(xs.string_excitation(xs_diff, string_process_,
134  finder_parameters.use_AQM));
135  }
136  }
137 }
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
StringProcess * string_process_
Pointer to interface class for strings.
void add_collisions(CollisionBranchList pv)
Add several new collision channels at once.
virtual double cross_section() const
Get the total cross section of the scattering particles.

◆ collision_channels()

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

Get list of possible collision channels.

Returns
list of possible collision channels.

Definition at line 152 of file scatteraction.h.

152  {
153  return collision_channels_;
154  }

◆ 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 172 of file scatteraction.h.

172  {
173  string_process_ = str_proc;
174  }

◆ cross_section()

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

Get the total cross section of the scattering particles.

Returns
total cross section.

Definition at line 181 of file scatteraction.h.

181 { return total_cross_section_; }

◆ 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 159 of file scatteraction.cc.

159  {
160  const double m1 = incoming_particles_[0].effective_mass();
161  const double m2 = incoming_particles_[1].effective_mass();
162  return pCM(sqrt_s(), m1, m2);
163 }

◆ 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 165 of file scatteraction.cc.

165  {
166  const double m1 = incoming_particles_[0].effective_mass();
167  const double m2 = incoming_particles_[1].effective_mass();
168  return pCM_sqr(sqrt_s(), m1, m2);
169 }
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 149 of file scatteraction.cc.

149  {
150  return total_momentum().velocity();
151 }
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 153 of file scatteraction.cc.

153  {
154  return (1. / std::sqrt(1.0 - beta_cm().sqr()));
155 }

◆ elastic_scattering()

void smash::ScatterAction::elastic_scattering ( )
protected

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

Definition at line 433 of file scatteraction.cc.

433  {
434  // copy initial particles into final state
437  // resample momenta
438  sample_angles({outgoing_particles_[0].effective_mass(),
439  outgoing_particles_[1].effective_mass()},
440  sqrt_s());
441 }
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 443 of file scatteraction.cc.

443  {
444  // create new particles
447 }
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 449 of file scatteraction.cc.

449  {
452  logg[LScatterAction].debug("2->", outgoing_particles_.size(),
453  " scattering:", incoming_particles_, " -> ",
455 }
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:442

◆ 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 479 of file scatteraction.cc.

479  {
482  /* Check momentum difference for debugging */
483  FourVector out_mom;
484  for (ParticleData data : outgoing_particles_) {
485  out_mom += data.momentum();
486  }
487  logg[LPythia].debug("Incoming momenta string:", total_momentum());
488  logg[LPythia].debug("Outgoing momenta string:", out_mom);
489 }
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 495 of file scatteraction.cc.

495  {
496  assert(incoming_particles_.size() == 2);
497  // Disable floating point exception trap for Pythia
498  {
499  DisableFloatTraps guard;
500  /* initialize the string_process_ object for this particular collision */
502  /* implement collision */
503  bool success = false;
504  int ntry = 0;
505  const int ntry_max = 10000;
506  while (!success && ntry < ntry_max) {
507  ntry++;
508  switch (process_type_) {
510  /* single diffractive to A+X */
511  success = string_process_->next_SDiff(true);
512  break;
514  /* single diffractive to X+B */
515  success = string_process_->next_SDiff(false);
516  break;
518  /* double diffractive */
519  success = string_process_->next_DDiff();
520  break;
522  /* soft non-diffractive */
523  success = string_process_->next_NDiffSoft();
524  break;
526  /* soft BBbar 2 mesonic annihilation */
527  success = string_process_->next_BBbarAnn();
528  break;
530  success = string_process_->next_NDiffHard();
531  break;
532  default:
533  logg[LPythia].error("Unknown string process required.");
534  success = false;
535  }
536  }
537  if (ntry == ntry_max) {
538  /* If pythia fails to form a string, it is usually because the energy
539  * is not large enough. In this case, annihilation is then enforced. If
540  * this process still does not not produce any results, it defaults to
541  * an elastic collision. */
542  bool success_newtry = false;
543 
544  /* Check if the initial state is a baryon-antibaryon state.*/
545  PdgCode part1 = incoming_particles_[0].pdgcode(),
546  part2 = incoming_particles_[1].pdgcode();
547  bool is_BBbar_Pair = (part1.baryon_number() != 0) &&
548  (part1.baryon_number() == -part2.baryon_number());
549 
550  /* Decide on the new process .*/
551  if (is_BBbar_Pair) {
553  } else {
555  }
556  /* Perform the new process*/
557  int ntry_new = 0;
558  while (!success_newtry && ntry_new < ntry_max) {
559  ntry_new++;
560  if (is_BBbar_Pair) {
561  success_newtry = string_process_->next_BBbarAnn();
562  } else {
563  success_newtry = string_process_->next_DDiff();
564  }
565  }
566 
567  if (success_newtry) {
569  }
570 
571  if (!success_newtry) {
572  /* If annihilation fails:
573  * Particles are normally added after process selection for
574  * strings, outgoing_particles is still uninitialized, and memory
575  * needs to be allocated. We also shift the process_type_ to elastic
576  * so that sample_angles does a proper treatment. */
577  outgoing_particles_.reserve(2);
578  outgoing_particles_.push_back(ParticleData{incoming_particles_[0]});
579  outgoing_particles_.push_back(ParticleData{incoming_particles_[1]});
582  }
583  } else {
585  }
586  }
587 }
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:359
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 457 of file scatteraction.cc.

457  {
458  if (outgoing_particles_.size() != 1) {
459  std::string s =
460  "resonance_formation: "
461  "Incorrect number of particles in final state: ";
462  s += std::to_string(outgoing_particles_.size()) + " (";
463  s += incoming_particles_[0].pdgcode().string() + " + ";
464  s += incoming_particles_[1].pdgcode().string() + ")";
465  throw InvalidResonanceFormation(s);
466  }
467  // Set the momentum of the formed resonance in its rest frame.
468  outgoing_particles_[0].set_4momentum(
469  total_momentum_of_outgoing_particles().abs(), 0., 0., 0.);
471  /* this momentum is evaluated in the computational frame. */
472  logg[LScatterAction].debug("Momentum of the new particle: ",
473  outgoing_particles_[0].momentum());
474 }

Member Data Documentation

◆ collision_channels_

CollisionBranchList smash::ScatterAction::collision_channels_
protected

List of possible collisions.

Definition at line 242 of file scatteraction.h.

◆ total_cross_section_

double smash::ScatterAction::total_cross_section_
protected

Total hadronic cross section.

Definition at line 245 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 248 of file scatteraction.h.

◆ isotropic_

bool smash::ScatterAction::isotropic_ = false
protected

Do this collision isotropically?

Definition at line 251 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 254 of file scatteraction.h.

◆ string_process_

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

Pointer to interface class for strings.

Definition at line 271 of file scatteraction.h.


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