Version: SMASH-1.8
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 31 of file scatteraction.h.

Inheritance diagram for smash::ScatterAction:
[legend]
Collaboration diagram for smash::ScatterAction:
[legend]

Classes

class  InvalidScatterAction
 

Public Member Functions

 ScatterAction (const ParticleData &in_part1, const ParticleData &in_part2, double time, bool isotropic=false, double string_formation_time=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 mandelstam_s () const
 Determine the Mandelstam s variable,. 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 (double elastic_parameter, bool two_to_one, ReactionsBitSet included_2to2, double low_snn_cut, bool strings_switch, bool use_AQM, bool strings_with_probability, NNbarTreatment nnbar_treatment)
 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 void 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 Particles &particles, 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 void 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...
 

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 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 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
 
- 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_3body_phasespace ()
 Sample the full 3-body phase-space (masses, momenta, angles) in the center-of-mass frame for the final state 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/c). More...
 
ProcessType process_type_
 type of process 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...
 

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 
)

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

Definition at line 33 of file scatteraction.cc.

36  : Action({in_part_a, in_part_b}, time),
38  isotropic_(isotropic),
39  string_formation_time_(string_formation_time) {}

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

41  {
42  add_process<CollisionBranch>(p, collision_channels_, total_cross_section_);
43 }
Here is the caller graph for this function:

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

45  {
46  add_processes<CollisionBranch>(std::move(pv), collision_channels_,
48 }
Here is the caller graph for this function:

◆ 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 position of particle a: x_a
position of particle b: x_b
momentum of particle a: p_a
momentum of particle b: p_b

\[d^2_\mathrm{coll} = (\vec{x_a} - \vec{x_b})^2 - \frac{((\vec{x_a} - \vec{x_b}) \cdot (\vec{p_a} - \vec{p_b}))^2 } {(\vec{p_a} - \vec{p_b})^2}\]

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

UrQMD squared distance criterion: Bass:1998ca [4] (3.27): in center of momentum frame position of particle a: x_a position of particle b: x_b momentum of particle a: p_a momentum of particle b: p_b d^2_{coll} = (x_a - x_b)^2 - ((x_a - x_b) . (p_a - p_b))^2 / (p_a - p_b)^2

Definition at line 168 of file scatteraction.cc.

168  {
169  // local copy of particles (since we need to boost them)
170  ParticleData p_a = incoming_particles_[0];
171  ParticleData p_b = incoming_particles_[1];
172  /* Boost particles to center-of-momentum frame. */
173  const ThreeVector velocity = beta_cm();
174  p_a.boost(velocity);
175  p_b.boost(velocity);
176  const ThreeVector pos_diff =
177  p_a.position().threevec() - p_b.position().threevec();
178  const ThreeVector mom_diff =
179  p_a.momentum().threevec() - p_b.momentum().threevec();
180 
181  logg[LScatterAction].debug("Particle ", incoming_particles_,
182  " position difference [fm]: ", pos_diff,
183  ", momentum difference [GeV]: ", mom_diff);
184 
185  const double dp2 = mom_diff.sqr();
186  const double dr2 = pos_diff.sqr();
187  /* Zero momentum leads to infite distance. */
188  if (dp2 < really_small) {
189  return dr2;
190  }
191  const double dpdr = pos_diff * mom_diff;
192 
201  const double result = dr2 - dpdr * dpdr / dp2;
202  return result > 0.0 ? result : 0.0;
203 }
Here is the call graph for this function:

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

154 { return total_momentum().sqr(); }
Here is the call graph for this function:
Here is the caller graph for this function:

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

50  {
51  logg[LScatterAction].debug("Incoming particles: ", incoming_particles_);
52 
53  /* Decide for a particular final state. */
54  const CollisionBranch *proc = choose_channel<CollisionBranch>(
56  process_type_ = proc->get_type();
57  outgoing_particles_ = proc->particle_list();
58  partial_cross_section_ = proc->weight();
59 
60  logg[LScatterAction].debug("Chosen channel: ", process_type_,
62 
63  /* The production point of the new particles. */
64  FourVector middle_point = get_interaction_point();
65 
66  switch (process_type_) {
68  /* 2->2 elastic scattering */
70  break;
72  /* resonance formation */
74  break;
76  /* 2->2 inelastic scattering */
77  /* Sample the particle momenta in CM system. */
79  break;
87  break;
88  default:
89  throw InvalidScatterAction(
90  "ScatterAction::generate_final_state: Invalid process type " +
91  std::to_string(static_cast<int>(process_type_)) + " was requested. " +
92  "(PDGcode1=" + incoming_particles_[0].pdgcode().string() +
93  ", PDGcode2=" + incoming_particles_[1].pdgcode().string() + ")");
94  }
95 
96  for (ParticleData &new_particle : outgoing_particles_) {
97  // Boost to the computational frame
98  new_particle.boost_momentum(
100  /* Set positions of the outgoing particles */
101  if (proc->get_type() != ProcessType::Elastic) {
102  new_particle.set_4position(middle_point);
103  }
104  }
105 }
Here is the call graph for this function:

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

136  {
137  return total_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
138  incoming_particles_[1].xsec_scaling_factor();
139 }

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

141  {
142  return partial_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
143  incoming_particles_[1].xsec_scaling_factor();
144 }

◆ 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 [14].

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 [14].

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

Reimplemented from smash::Action.

Definition at line 249 of file scatteraction.cc.

250  {
253  // We potentially have more than two particles, so the following angular
254  // distributions don't work. Instead we just keep the angular
255  // distributions generated by string fragmentation.
256  return;
257  }
258  assert(outgoing_particles_.size() == 2);
259 
260  // NN scattering is anisotropic currently
261  const bool nn_scattering = incoming_particles_[0].type().is_nucleon() &&
262  incoming_particles_[1].type().is_nucleon();
263  /* Elastic process is anisotropic and
264  * the angular distribution is based on the NN elastic scattering. */
265  const bool el_scattering = process_type_ == ProcessType::Elastic;
266 
267  const double mass_in_a = incoming_particles_[0].effective_mass();
268  const double mass_in_b = incoming_particles_[1].effective_mass();
269 
270  ParticleData *p_a = &outgoing_particles_[0];
271  ParticleData *p_b = &outgoing_particles_[1];
272 
273  const double mass_a = masses.first;
274  const double mass_b = masses.second;
275 
276  const std::array<double, 2> t_range = get_t_range<double>(
277  kinetic_energy_cm, mass_in_a, mass_in_b, mass_a, mass_b);
278  Angles phitheta;
279  if (el_scattering && !isotropic_) {
283  double mandelstam_s_new = 0.;
284  if (nn_scattering) {
285  mandelstam_s_new = mandelstam_s();
286  } else {
287  /* In the case of elastic collisions other than NN collisions,
288  * there is an ambiguity on how to get the lab-frame momentum (plab),
289  * since the incoming particles can have different masses.
290  * Right now, we first obtain the center-of-mass momentum
291  * of the collision (pcom_now).
292  * Then, the lab-frame momentum is evaluated from the mandelstam s,
293  * which yields the original center-of-mass momentum
294  * when nucleon mass is assumed. */
295  const double pcm_now = pCM_from_s(mandelstam_s(), mass_in_a, mass_in_b);
296  mandelstam_s_new =
297  4. * std::sqrt(pcm_now * pcm_now + nucleon_mass * nucleon_mass);
298  }
299  double bb, a, plab = plab_from_s(mandelstam_s_new);
300  if (nn_scattering &&
301  p_a->pdgcode().antiparticle_sign() ==
302  p_b->pdgcode().antiparticle_sign() &&
303  std::abs(p_a->type().charge() + p_b->type().charge()) == 1) {
304  // proton-neutron and antiproton-antineutron
305  bb = std::max(Cugnon_bnp(plab), really_small);
306  a = (plab < 0.8) ? 1. : 0.64 / (plab * plab);
307  } else {
308  /* all others including pp, nn and AQM elastic processes
309  * This is applied for all particle pairs, which are allowed to
310  * interact elastically. */
311  bb = std::max(Cugnon_bpp(plab), really_small);
312  a = 1.;
313  }
314  double t = random::expo(bb, t_range[0], t_range[1]);
315  if (random::canonical() > 1. / (1. + a)) {
316  t = t_range[0] + t_range[1] - t;
317  }
318  // determine scattering angles in center-of-mass frame
319  phitheta = Angles(2. * M_PI * random::canonical(),
320  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
321  } else if (nn_scattering && p_a->pdgcode().is_Delta() &&
322  p_b->pdgcode().is_nucleon() &&
323  p_a->pdgcode().antiparticle_sign() ==
324  p_b->pdgcode().antiparticle_sign() &&
325  !isotropic_) {
329  const double plab = plab_from_s(mandelstam_s());
330  const double bb = std::max(Cugnon_bpp(plab), really_small);
331  double t = random::expo(bb, t_range[0], t_range[1]);
332  if (random::canonical() > 0.5) {
333  t = t_range[0] + t_range[1] - t; // symmetrize
334  }
335  phitheta = Angles(2. * M_PI * random::canonical(),
336  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
337  } else if (nn_scattering && p_b->pdgcode().is_nucleon() && !isotropic_ &&
338  (p_a->type().is_Nstar() || p_a->type().is_Deltastar())) {
340  const std::array<double, 4> p{1.46434, 5.80311, -6.89358, 1.94302};
341  const double a = p[0] + mass_a * (p[1] + mass_a * (p[2] + mass_a * p[3]));
342  /* If the resonance is so heavy that the index "a" exceeds 30,
343  * the power function turns out to be too sharp. Take t directly to be
344  * t_0 in such a case. */
345  double t = t_range[0];
346  if (a < 30) {
347  t = random::power(-a, t_range[0], t_range[1]);
348  }
349  if (random::canonical() > 0.5) {
350  t = t_range[0] + t_range[1] - t; // symmetrize
351  }
352  phitheta = Angles(2. * M_PI * random::canonical(),
353  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
354  } else {
355  /* isotropic angular distribution */
356  phitheta.distribute_isotropically();
357  }
358 
359  ThreeVector pscatt = phitheta.threevec();
360  // 3-momentum of first incoming particle in center-of-mass frame
361  ThreeVector pcm =
362  incoming_particles_[0].momentum().lorentz_boost(beta_cm()).threevec();
363  pscatt.rotate_z_axis_to(pcm);
364 
365  // final-state CM momentum
366  const double p_f = pCM(kinetic_energy_cm, mass_a, mass_b);
367  if (!(p_f > 0.0)) {
368  logg[LScatterAction].warn("Particle: ", p_a->pdgcode(),
369  " radial momentum: ", p_f);
370  logg[LScatterAction].warn("Etot: ", kinetic_energy_cm, " m_a: ", mass_a,
371  " m_b: ", mass_b);
372  }
373  p_a->set_4momentum(mass_a, pscatt * p_f);
374  p_b->set_4momentum(mass_b, -pscatt * p_f);
375 
376  /* Debug message is printed before boost, so that p_a and p_b are
377  * the momenta in the center of mass frame and thus opposite to
378  * each other.*/
379  logg[LScatterAction].debug("p_a: ", *p_a, "\np_b: ", *p_b);
380 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ add_all_scatterings()

void smash::ScatterAction::add_all_scatterings ( double  elastic_parameter,
bool  two_to_one,
ReactionsBitSet  included_2to2,
double  low_snn_cut,
bool  strings_switch,
bool  use_AQM,
bool  strings_with_probability,
NNbarTreatment  nnbar_treatment 
)

Add all possible scattering subprocesses for this action object.

Parameters
[in]elastic_parameterIf non-zero, given global elastic cross section.
[in]two_to_one2->1 reactions enabled?
[in]included_2to2Which 2->2 reactions are enabled?
[in]low_snn_cutElastic collisions with CME below are forbidden.
[in]strings_switchAre string processes enabled?
[in]use_AQMuse elastic cross sections via AQM?
[in]strings_with_probabilityAre string processes triggered according to a probability?
[in]nnbar_treatmentNNbar treatment through resonance, strings or none

Definition at line 107 of file scatteraction.cc.

110  {
111  CrossSections xs(incoming_particles_, sqrt_s(),
113  CollisionBranchList processes = xs.generate_collision_list(
114  elastic_parameter, two_to_one, included_2to2, low_snn_cut, strings_switch,
115  use_AQM, strings_with_probability, nnbar_treatment, string_process_);
116 
117  /* Add various subprocesses.*/
118  add_collisions(std::move(processes));
119 
120  /* If the string processes are not triggered by a probability, then they
121  * always happen as long as the parametrized total cross section is larger
122  * than the sum of the cross sections of the non-string processes, and the
123  * square root s exceeds the threshold by at least 0.9 GeV. The cross section
124  * of the string processes are counted by taking the difference between the
125  * parametrized total and the sum of the non-strings. */
126  if (!strings_with_probability &&
127  xs.string_probability(strings_switch, strings_with_probability, use_AQM,
128  nnbar_treatment == NNbarTreatment::Strings) == 1.) {
129  const double xs_diff = xs.high_energy() - cross_section();
130  if (xs_diff > 0.) {
131  add_collisions(xs.string_excitation(xs_diff, string_process_, use_AQM));
132  }
133  }
134 }
Here is the call graph for this function:

◆ collision_channels()

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

Get list of possible collision channels.

Returns
list of possible collision channels.

Definition at line 142 of file scatteraction.h.

142  {
143  return collision_channels_;
144  }

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

162  {
163  string_process_ = str_proc;
164  }

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

171 { return total_cross_section_; }
Here is the caller graph for this function:

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

156  {
157  const double m1 = incoming_particles_[0].effective_mass();
158  const double m2 = incoming_particles_[1].effective_mass();
159  return pCM(sqrt_s(), m1, m2);
160 }
Here is the call graph for this function:
Here is the caller graph for this function:

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

162  {
163  const double m1 = incoming_particles_[0].effective_mass();
164  const double m2 = incoming_particles_[1].effective_mass();
165  return pCM_sqr(sqrt_s(), m1, m2);
166 }
Here is the call graph for this function:

◆ beta_cm()

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

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

Returns
boost velocity between center of mass and calculation frame.

Definition at line 146 of file scatteraction.cc.

146  {
147  return total_momentum().velocity();
148 }
Here is the call graph for this function:
Here is the caller graph for this function:

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

150  {
151  return (1. / std::sqrt(1.0 - beta_cm().sqr()));
152 }
Here is the call graph for this function:

◆ elastic_scattering()

void smash::ScatterAction::elastic_scattering ( )
protected

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

Definition at line 382 of file scatteraction.cc.

382  {
383  // copy initial particles into final state
386  // resample momenta
387  sample_angles({outgoing_particles_[0].effective_mass(),
388  outgoing_particles_[1].effective_mass()},
389  sqrt_s());
390 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ inelastic_scattering()

void smash::ScatterAction::inelastic_scattering ( )
protected

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

Definition at line 392 of file scatteraction.cc.

392  {
393  // create new particles
395  /* Set the formation time of the 2 particles to the larger formation time of
396  * the incoming particles, if it is larger than the execution time;
397  * execution time is otherwise taken to be the formation time */
398  const double t0 = incoming_particles_[0].formation_time();
399  const double t1 = incoming_particles_[1].formation_time();
400 
401  const size_t index_tmax = (t0 > t1) ? 0 : 1;
402  const double form_time_begin =
403  incoming_particles_[index_tmax].begin_formation_time();
404  const double sc =
405  incoming_particles_[index_tmax].initial_xsec_scaling_factor();
406  if (t0 > time_of_execution_ || t1 > time_of_execution_) {
407  /* The newly produced particles are supposed to continue forming exactly
408  * like the latest forming ingoing particle. Therefore the details on the
409  * formation are adopted. The initial cross section scaling factor of the
410  * incoming particles is considered to also be the scaling factor of the
411  * newly produced outgoing particles.
412  */
413  outgoing_particles_[0].set_slow_formation_times(form_time_begin,
414  std::max(t0, t1));
415  outgoing_particles_[1].set_slow_formation_times(form_time_begin,
416  std::max(t0, t1));
417  outgoing_particles_[0].set_cross_section_scaling_factor(sc);
418  outgoing_particles_[1].set_cross_section_scaling_factor(sc);
419  } else {
420  outgoing_particles_[0].set_formation_time(time_of_execution_);
421  outgoing_particles_[1].set_formation_time(time_of_execution_);
422  }
423 }
Here is the call graph for this function:
Here is the caller graph for this function:

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

464  {
465  assert(incoming_particles_.size() == 2);
466  // Disable floating point exception trap for Pythia
467  {
468  DisableFloatTraps guard;
469  /* initialize the string_process_ object for this particular collision */
471  /* implement collision */
472  bool success = false;
473  int ntry = 0;
474  const int ntry_max = 10000;
475  while (!success && ntry < ntry_max) {
476  ntry++;
477  switch (process_type_) {
479  /* single diffractive to A+X */
480  success = string_process_->next_SDiff(true);
481  break;
483  /* single diffractive to X+B */
484  success = string_process_->next_SDiff(false);
485  break;
487  /* double diffractive */
488  success = string_process_->next_DDiff();
489  break;
491  /* soft non-diffractive */
492  success = string_process_->next_NDiffSoft();
493  break;
495  /* soft BBbar 2 mesonic annihilation */
496  success = string_process_->next_BBbarAnn();
497  break;
499  success = string_process_->next_NDiffHard();
500  break;
501  default:
502  logg[LPythia].error("Unknown string process required.");
503  success = false;
504  }
505  }
506  if (ntry == ntry_max) {
507  /* If pythia fails to form a string, it is usually because the energy
508  * is not large enough. In this case, an elastic scattering happens.
509  *
510  * Since particles are normally added after process selection for
511  * strings, outgoing_particles is still uninitialized, and memory
512  * needs to be allocated. We also shift the process_type_ to elastic
513  * so that sample_angles does a proper treatment. */
514  outgoing_particles_.reserve(2);
515  outgoing_particles_.push_back(ParticleData{incoming_particles_[0]});
516  outgoing_particles_.push_back(ParticleData{incoming_particles_[1]});
519  } else {
521  /* If the incoming particles already were unformed, the formation
522  * times and cross section scaling factors need to be adjusted */
523  const double tform_in = std::max(incoming_particles_[0].formation_time(),
524  incoming_particles_[1].formation_time());
525  if (tform_in > time_of_execution_) {
526  const double fin =
527  (incoming_particles_[0].formation_time() >
528  incoming_particles_[1].formation_time())
529  ? incoming_particles_[0].initial_xsec_scaling_factor()
530  : incoming_particles_[1].initial_xsec_scaling_factor();
531  for (size_t i = 0; i < outgoing_particles_.size(); i++) {
532  const double tform_out = outgoing_particles_[i].formation_time();
533  const double fout =
534  outgoing_particles_[i].initial_xsec_scaling_factor();
535  /* The new cross section scaling factor will be the product of the
536  * cross section scaling factor of the ingoing particles and of the
537  * outgoing ones (since the outgoing ones are also string fragments
538  * and thus take time to form). */
539  outgoing_particles_[i].set_cross_section_scaling_factor(fin * fout);
540  /* If the unformed incoming particles' formation time is larger than
541  * the current outgoing particle's formation time, then the latter
542  * is overwritten by the former*/
543  if (tform_in > tform_out) {
544  outgoing_particles_[i].set_slow_formation_times(time_of_execution_,
545  tform_in);
546  }
547  }
548  }
549  /* Check momentum difference for debugging */
550  FourVector out_mom;
551  for (ParticleData data : outgoing_particles_) {
552  out_mom += data.momentum();
553  }
554  logg[LPythia].debug("Incoming momenta string:", total_momentum());
555  logg[LPythia].debug("Outgoing momenta string:", out_mom);
556  }
557  }
558 }
Here is the call graph for this function:
Here is the caller graph for this function:

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

425  {
426  if (outgoing_particles_.size() != 1) {
427  std::string s =
428  "resonance_formation: "
429  "Incorrect number of particles in final state: ";
430  s += std::to_string(outgoing_particles_.size()) + " (";
431  s += incoming_particles_[0].pdgcode().string() + " + ";
432  s += incoming_particles_[1].pdgcode().string() + ")";
433  throw InvalidResonanceFormation(s);
434  }
435  // Set the momentum of the formed resonance in its rest frame.
436  outgoing_particles_[0].set_4momentum(
437  total_momentum_of_outgoing_particles().abs(), 0., 0., 0.);
438  /* Set the formation time of the resonance to the larger formation time of the
439  * incoming particles, if it is larger than the execution time; execution time
440  * is otherwise taken to be the formation time */
441  const double t0 = incoming_particles_[0].formation_time();
442  const double t1 = incoming_particles_[1].formation_time();
443 
444  const size_t index_tmax = (t0 > t1) ? 0 : 1;
445  const double begin_form_time =
446  incoming_particles_[index_tmax].begin_formation_time();
447  const double sc =
448  incoming_particles_[index_tmax].initial_xsec_scaling_factor();
449  if (t0 > time_of_execution_ || t1 > time_of_execution_) {
450  outgoing_particles_[0].set_slow_formation_times(begin_form_time,
451  std::max(t0, t1));
452  outgoing_particles_[0].set_cross_section_scaling_factor(sc);
453  } else {
454  outgoing_particles_[0].set_formation_time(time_of_execution_);
455  }
456  /* this momentum is evaluated in the computational frame. */
457  logg[LScatterAction].debug("Momentum of the new particle: ",
458  outgoing_particles_[0].momentum());
459 }
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ collision_channels_

CollisionBranchList smash::ScatterAction::collision_channels_
protected

List of possible collisions.

Definition at line 222 of file scatteraction.h.

◆ total_cross_section_

double smash::ScatterAction::total_cross_section_
protected

Total hadronic cross section.

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

◆ isotropic_

bool smash::ScatterAction::isotropic_ = false
protected

Do this collision isotropically?

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

◆ string_process_

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

Pointer to interface class for strings.

Definition at line 251 of file scatteraction.h.


The documentation for this class was generated from the following files:
smash::Action::incoming_particles_
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:304
smash::ProcessType::StringHard
hard string process involving 2->2 QCD process by PYTHIA.
smash::Action::total_momentum
FourVector total_momentum() const
Sum of 4-momenta of incoming particles.
Definition: action.h:324
smash::Action::time_of_execution_
const double time_of_execution_
Time at which the action is supposed to be performed (absolute time in the lab frame in fm/c).
Definition: action.h:318
smash::StringProcess::get_final_state
ParticleList get_final_state()
a function to get the final state particle list which is called after the collision
Definition: processstring.h:766
smash::plab_from_s
double plab_from_s(double mandelstam_s, double mass)
Convert Mandelstam-s to p_lab in a fixed-target collision.
Definition: kinematics.h:157
smash::ScatterAction::resonance_formation
void resonance_formation()
Perform a 2->1 resonance-formation process.
Definition: scatteraction.cc:425
smash::StringProcess::init
void init(const ParticleList &incoming, double tcoll)
initialization feed intial particles, time of collision and gamma factor of the center of mass.
Definition: processstring.cc:200
smash::random::expo
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
smash::Action::sample_2body_phasespace
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:206
smash::StringProcess::next_NDiffHard
bool next_NDiffHard()
Hard Non-diffractive process is based on PYTHIA 8 with partonic showers and interactions.
Definition: processstring.cc:513
smash::ProcessType::StringSoftDoubleDiffractive
double diffractive. Two strings are formed, one from A and one from B.
NNbarTreatment::Strings
Use string fragmentation.
smash::ProcessType::StringSoftNonDiffractive
non-diffractive. Two strings are formed both have ends in A and B.
smash::is_string_soft_process
bool is_string_soft_process(ProcessType p)
Check if a given process type is a soft string excitation.
Definition: processbranch.cc:18
smash::FourVector::sqr
double sqr() const
calculate the square of the vector (which is a scalar)
Definition: fourvector.h:450
smash::ScatterAction::collision_channels_
CollisionBranchList collision_channels_
List of possible collisions.
Definition: scatteraction.h:222
smash::ScatterAction::elastic_scattering
void elastic_scattering()
Perform an elastic two-body scattering, i.e. just exchange momentum.
Definition: scatteraction.cc:382
smash::nucleon_mass
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:55
smash::ProcessType::TwoToTwo
2->2 inelastic scattering
smash::ScatterAction::beta_cm
ThreeVector beta_cm() const
Get the velocity of the center of mass of the scattering particles in the calculation frame.
Definition: scatteraction.cc:146
smash::ProcessType::StringSoftSingleDiffractiveXB
single diffractive AB->XB.
smash::logg
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
smash::ScatterAction::string_process_
StringProcess * string_process_
Pointer to interface class for strings.
Definition: scatteraction.h:251
smash::pCM
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:79
smash::really_small
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
smash::Cugnon_bpp
static double Cugnon_bpp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic pp...
Definition: scatteraction.cc:218
smash::ScatterAction::string_excitation
void string_excitation()
Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string e...
Definition: scatteraction.cc:464
smash::StringProcess::next_SDiff
bool next_SDiff(bool is_AB_to_AX)
Single-diffractive process is based on single pomeron exchange described in Ingelman:1984ns .
Definition: processstring.cc:229
smash::ScatterAction::string_formation_time_
double string_formation_time_
Time fragments take to be fully formed in hard string excitation.
Definition: scatteraction.h:234
smash::StringProcess::next_DDiff
bool next_DDiff()
Double-diffractive process ( A + B -> X + X ) is similar to the single-diffractive process,...
Definition: processstring.cc:381
smash::Action::process_type_
ProcessType process_type_
type of process
Definition: action.h:321
smash::ScatterAction::total_cross_section_
double total_cross_section_
Total hadronic cross section.
Definition: scatteraction.h:225
smash::LPythia
static constexpr int LPythia
Definition: processstring.h:25
smash::LScatterAction
static constexpr int LScatterAction
Definition: bremsstrahlungaction.cc:17
smash::ScatterAction::cross_section
virtual double cross_section() const
Get the total cross section of the scattering particles.
Definition: scatteraction.h:171
smash::Action::get_potential_at_interaction_point
std::pair< FourVector, FourVector > get_potential_at_interaction_point() const
Get the skyrme and asymmetry potential at the interaction point.
Definition: action.cc:79
smash::Action::total_momentum_of_outgoing_particles
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
Definition: action.cc:123
smash::FourVector::velocity
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:323
smash::ScatterAction::sample_angles
void sample_angles(std::pair< double, double > masses, double kinetic_energy_cm) override
Sample final-state angles in a 2->2 collision (possibly anisotropic).
Definition: scatteraction.cc:249
smash::ScatterAction::isotropic_
bool isotropic_
Do this collision isotropically?
Definition: scatteraction.h:231
smash::Action::outgoing_particles_
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
Definition: action.h:312
smash::ProcessType::StringSoftSingleDiffractiveAX
(41-45) soft string excitations.
smash::Action::sqrt_s
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:266
smash::ProcessType::StringSoftAnnihilation
a special case of baryon-antibaryon annihilation.
smash::random::power
T power(T n, T xMin, T xMax)
Draws a random number according to a power-law distribution ~ x^n.
Definition: random.h:203
smash::ProcessType::Elastic
elastic scattering: particles remain the same, only momenta change
smash::ScatterAction::inelastic_scattering
void inelastic_scattering()
Perform an inelastic two-body scattering, i.e. new particles are formed.
Definition: scatteraction.cc:392
smash::Action::Action
Action(const ParticleList &in_part, double time)
Construct an action object with incoming particles and relative time.
Definition: action.h:44
smash::Action::get_interaction_point
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:69
smash::StringProcess::next_BBbarAnn
bool next_BBbarAnn()
Baryon-antibaryon annihilation process Based on what UrQMD Bass:1998ca , Bleicher:1999xi does,...
Definition: processstring.cc:1504
smash::pdg::p
constexpr int p
Proton.
Definition: pdgcode_constants.h:28
smash::ScatterAction::partial_cross_section_
double partial_cross_section_
Partial cross-section to the chosen outgoing channel.
Definition: scatteraction.h:228
smash::ProcessType::TwoToOne
resonance formation (2->1)
smash::random::canonical
T canonical()
Definition: random.h:113
smash::Cugnon_bnp
static double Cugnon_bnp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic np...
Definition: scatteraction.cc:237
smash::ScatterAction::add_collisions
void add_collisions(CollisionBranchList pv)
Add several new collision channels at once.
Definition: scatteraction.cc:45
smash::StringProcess::next_NDiffSoft
bool next_NDiffSoft()
Soft Non-diffractive process is modelled in accordance with dual-topological approach Capella:1978ig ...
Definition: processstring.cc:436
smash::ScatterAction::mandelstam_s
double mandelstam_s() const
Determine the Mandelstam s variable,.
Definition: scatteraction.cc:154
smash::pCM_sqr
T pCM_sqr(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:91
smash::pCM_from_s
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:66