Version: SMASH-1.5
scatteraction.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2018
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #include "smash/scatteraction.h"
11 
12 #include <cmath>
13 
14 #include "Pythia8/Pythia.h"
15 
16 #include "smash/angles.h"
17 #include "smash/constants.h"
18 #include "smash/crosssections.h"
19 #include "smash/cxx14compat.h"
20 #include "smash/fpenvironment.h"
21 #include "smash/kinematics.h"
22 #include "smash/logging.h"
23 #include "smash/parametrizations.h"
24 #include "smash/pdgcode.h"
26 #include "smash/pow.h"
27 #include "smash/processstring.h"
28 #include "smash/random.h"
29 
30 namespace smash {
31 
33  const ParticleData &in_part_b, double time,
34  bool isotropic, double string_formation_time)
35  : Action({in_part_a, in_part_b}, time),
36  total_cross_section_(0.),
37  isotropic_(isotropic),
38  string_formation_time_(string_formation_time) {}
39 
40 void ScatterAction::add_collision(CollisionBranchPtr p) {
41  add_process<CollisionBranch>(p, collision_channels_, total_cross_section_);
42 }
43 
44 void ScatterAction::add_collisions(CollisionBranchList pv) {
45  add_processes<CollisionBranch>(std::move(pv), collision_channels_,
47 }
48 
50  const auto &log = logger<LogArea::ScatterAction>();
51 
52  log.debug("Incoming particles: ", incoming_particles_);
53 
54  /* Decide for a particular final state. */
55  const CollisionBranch *proc = choose_channel<CollisionBranch>(
57  process_type_ = proc->get_type();
60 
61  log.debug("Chosen channel: ", process_type_, outgoing_particles_);
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:
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 }
106 
108  double elastic_parameter, bool two_to_one, ReactionsBitSet included_2to2,
109  double low_snn_cut, bool strings_switch, bool use_AQM,
110  bool strings_with_probability, NNbarTreatment nnbar_treatment) {
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 }
135 
137  return total_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
138  incoming_particles_[1].xsec_scaling_factor();
139 }
140 
142  return partial_cross_section_ * incoming_particles_[0].xsec_scaling_factor() *
143  incoming_particles_[1].xsec_scaling_factor();
144 }
145 
147  return total_momentum().velocity();
148 }
149 
150 double ScatterAction::gamma_cm() const {
151  return (1. / std::sqrt(1.0 - beta_cm().sqr()));
152 }
153 
154 double ScatterAction::mandelstam_s() const { return total_momentum().sqr(); }
155 
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 }
161 
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 }
167 
169  // local copy of particles (since we need to boost them)
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  const auto &log = logger<LogArea::ScatterAction>();
182  log.debug("Particle ", incoming_particles_,
183  " position difference [fm]: ", pos_diff,
184  ", momentum difference [GeV]: ", mom_diff);
185 
186  const double dp2 = mom_diff.sqr();
187  const double dr2 = pos_diff.sqr();
188  /* Zero momentum leads to infite distance. */
189  if (dp2 < really_small) {
190  return dr2;
191  }
192  const double dpdr = pos_diff * mom_diff;
193 
202  return dr2 - dpdr * dpdr / dp2;
203 }
204 
218 static double Cugnon_bpp(double plab) {
219  if (plab < 2.) {
220  double p8 = pow_int(plab, 8);
221  return 5.5 * p8 / (7.7 + p8);
222  } else {
223  return std::min(9.0, 5.334 + 0.67 * (plab - 2.));
224  }
225 }
226 
237 static double Cugnon_bnp(double plab) {
238  if (plab < 0.225) {
239  return 0.;
240  } else if (plab < 0.6) {
241  return 16.53 * (plab - 0.225);
242  } else if (plab < 1.6) {
243  return -1.63 * plab + 7.16;
244  } else {
245  return Cugnon_bpp(plab);
246  }
247 }
248 
249 void ScatterAction::sample_angles(std::pair<double, double> masses,
250  double kinetic_energy_cm) {
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  const auto &log = logger<LogArea::ScatterAction>();
260 
261  // NN scattering is anisotropic currently
262  const bool nn_scattering = incoming_particles_[0].type().is_nucleon() &&
263  incoming_particles_[1].type().is_nucleon();
264  /* Elastic process is anisotropic and
265  * the angular distribution is based on the NN elastic scattering. */
266  const bool el_scattering = process_type_ == ProcessType::Elastic;
267 
268  const double mass_in_a = incoming_particles_[0].effective_mass();
269  const double mass_in_b = incoming_particles_[1].effective_mass();
270 
273 
274  const double mass_a = masses.first;
275  const double mass_b = masses.second;
276 
277  const std::array<double, 2> t_range = get_t_range<double>(
278  kinetic_energy_cm, mass_in_a, mass_in_b, mass_a, mass_b);
279  Angles phitheta;
280  if (el_scattering && !isotropic_) {
284  double mandelstam_s_new = 0.;
285  if (nn_scattering) {
286  mandelstam_s_new = mandelstam_s();
287  } else {
288  /* In the case of elastic collisions other than NN collisions,
289  * there is an ambiguity on how to get the lab-frame momentum (plab),
290  * since the incoming particles can have different masses.
291  * Right now, we first obtain the center-of-mass momentum
292  * of the collision (pcom_now).
293  * Then, the lab-frame momentum is evaluated from the mandelstam s,
294  * which yields the original center-of-mass momentum
295  * when nucleon mass is assumed. */
296  const double pcm_now = pCM_from_s(mandelstam_s(), mass_in_a, mass_in_b);
297  mandelstam_s_new =
298  4. * std::sqrt(pcm_now * pcm_now + nucleon_mass * nucleon_mass);
299  }
300  double bb, a, plab = plab_from_s(mandelstam_s_new);
301  if (nn_scattering &&
302  p_a->pdgcode().antiparticle_sign() ==
303  p_b->pdgcode().antiparticle_sign() &&
304  std::abs(p_a->type().charge() + p_b->type().charge()) == 1) {
305  // proton-neutron and antiproton-antineutron
306  bb = std::max(Cugnon_bnp(plab), really_small);
307  a = (plab < 0.8) ? 1. : 0.64 / (plab * plab);
308  } else {
309  /* all others including pp, nn and AQM elastic processes
310  * This is applied for all particle pairs, which are allowed to
311  * interact elastically. */
312  bb = std::max(Cugnon_bpp(plab), really_small);
313  a = 1.;
314  }
315  double t = random::expo(bb, t_range[0], t_range[1]);
316  if (random::canonical() > 1. / (1. + a)) {
317  t = t_range[0] + t_range[1] - t;
318  }
319  // determine scattering angles in center-of-mass frame
320  phitheta = Angles(2. * M_PI * random::canonical(),
321  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
322  } else if (nn_scattering && p_a->pdgcode().is_Delta() &&
323  p_b->pdgcode().is_nucleon() &&
324  p_a->pdgcode().antiparticle_sign() ==
325  p_b->pdgcode().antiparticle_sign() &&
326  !isotropic_) {
330  const double plab = plab_from_s(mandelstam_s());
331  const double bb = std::max(Cugnon_bpp(plab), really_small);
332  double t = random::expo(bb, t_range[0], t_range[1]);
333  if (random::canonical() > 0.5) {
334  t = t_range[0] + t_range[1] - t; // symmetrize
335  }
336  phitheta = Angles(2. * M_PI * random::canonical(),
337  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
338  } else if (nn_scattering && p_b->pdgcode().is_nucleon() && !isotropic_ &&
339  (p_a->type().is_Nstar() || p_a->type().is_Deltastar())) {
341  const std::array<double, 4> p{1.46434, 5.80311, -6.89358, 1.94302};
342  const double a = p[0] + mass_a * (p[1] + mass_a * (p[2] + mass_a * p[3]));
343  /* If the resonance is so heavy that the index "a" exceeds 30,
344  * the power function turns out to be too sharp. Take t directly to be
345  * t_0 in such a case. */
346  double t = t_range[0];
347  if (a < 30) {
348  t = random::power(-a, t_range[0], t_range[1]);
349  }
350  if (random::canonical() > 0.5) {
351  t = t_range[0] + t_range[1] - t; // symmetrize
352  }
353  phitheta = Angles(2. * M_PI * random::canonical(),
354  1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
355  } else {
356  /* isotropic angular distribution */
357  phitheta.distribute_isotropically();
358  }
359 
360  ThreeVector pscatt = phitheta.threevec();
361  // 3-momentum of first incoming particle in center-of-mass frame
362  ThreeVector pcm =
363  incoming_particles_[0].momentum().LorentzBoost(beta_cm()).threevec();
364  pscatt.rotate_z_axis_to(pcm);
365 
366  // final-state CM momentum
367  const double p_f = pCM(kinetic_energy_cm, mass_a, mass_b);
368  if (!(p_f > 0.0)) {
369  log.warn("Particle: ", p_a->pdgcode(), " radial momentum: ", p_f);
370  log.warn("Etot: ", kinetic_energy_cm, " m_a: ", mass_a, " m_b: ", mass_b);
371  }
372  p_a->set_4momentum(mass_a, pscatt * p_f);
373  p_b->set_4momentum(mass_b, -pscatt * p_f);
374 
375  /* Debug message is printed before boost, so that p_a and p_b are
376  * the momenta in the center of mass frame and thus opposite to
377  * each other.*/
378  log.debug("p_a: ", *p_a, "\np_b: ", *p_b);
379 }
380 
382  // copy initial particles into final state
385  // resample momenta
386  sample_angles({outgoing_particles_[0].effective_mass(),
387  outgoing_particles_[1].effective_mass()},
388  sqrt_s());
389 }
390 
392  // create new particles
394  /* Set the formation time of the 2 particles to the larger formation time of
395  * the incoming particles, if it is larger than the execution time;
396  * execution time is otherwise taken to be the formation time */
397  const double t0 = incoming_particles_[0].formation_time();
398  const double t1 = incoming_particles_[1].formation_time();
399 
400  const size_t index_tmax = (t0 > t1) ? 0 : 1;
401  const double form_time_begin =
402  incoming_particles_[index_tmax].begin_formation_time();
403  const double sc =
404  incoming_particles_[index_tmax].initial_xsec_scaling_factor();
405  if (t0 > time_of_execution_ || t1 > time_of_execution_) {
406  /* The newly produced particles are supposed to continue forming exactly
407  * like the latest forming ingoing particle. Therefore the details on the
408  * formation are adopted. The initial cross section scaling factor of the
409  * incoming particles is considered to also be the scaling factor of the
410  * newly produced outgoing particles.
411  */
412  outgoing_particles_[0].set_slow_formation_times(form_time_begin,
413  std::max(t0, t1));
414  outgoing_particles_[1].set_slow_formation_times(form_time_begin,
415  std::max(t0, t1));
416  outgoing_particles_[0].set_cross_section_scaling_factor(sc);
417  outgoing_particles_[1].set_cross_section_scaling_factor(sc);
418  } else {
419  outgoing_particles_[0].set_formation_time(time_of_execution_);
420  outgoing_particles_[1].set_formation_time(time_of_execution_);
421  }
422 }
423 
425  const auto &log = logger<LogArea::ScatterAction>();
426 
427  if (outgoing_particles_.size() != 1) {
428  std::string s =
429  "resonance_formation: "
430  "Incorrect number of particles in final state: ";
431  s += std::to_string(outgoing_particles_.size()) + " (";
432  s += incoming_particles_[0].pdgcode().string() + " + ";
433  s += incoming_particles_[1].pdgcode().string() + ")";
434  throw InvalidResonanceFormation(s);
435  }
436  // Set the momentum of the formed resonance in its rest frame.
437  outgoing_particles_[0].set_4momentum(
438  total_momentum_of_outgoing_particles().abs(), 0., 0., 0.);
439  /* Set the formation time of the resonance to the larger formation time of the
440  * incoming particles, if it is larger than the execution time; execution time
441  * is otherwise taken to be the formation time */
442  const double t0 = incoming_particles_[0].formation_time();
443  const double t1 = incoming_particles_[1].formation_time();
444 
445  const size_t index_tmax = (t0 > t1) ? 0 : 1;
446  const double begin_form_time =
447  incoming_particles_[index_tmax].begin_formation_time();
448  const double sc =
449  incoming_particles_[index_tmax].initial_xsec_scaling_factor();
450  if (t0 > time_of_execution_ || t1 > time_of_execution_) {
451  outgoing_particles_[0].set_slow_formation_times(begin_form_time,
452  std::max(t0, t1));
453  outgoing_particles_[0].set_cross_section_scaling_factor(sc);
454  } else {
455  outgoing_particles_[0].set_formation_time(time_of_execution_);
456  }
457  /* this momentum is evaluated in the computational frame. */
458  log.debug("Momentum of the new particle: ",
459  outgoing_particles_[0].momentum());
460 }
461 
462 /* This function will generate outgoing particles in computational frame
463  * from a hard process.
464  * The way to excite soft strings is based on the UrQMD model */
466  assert(incoming_particles_.size() == 2);
467  const auto &log = logger<LogArea::Pythia>();
468  // Disable floating point exception trap for Pythia
469  {
470  DisableFloatTraps guard;
471  /* initialize the string_process_ object for this particular collision */
473  /* implement collision */
474  bool success = false;
475  int ntry = 0;
476  const int ntry_max = 10000;
477  while (!success && ntry < ntry_max) {
478  ntry++;
479  switch (process_type_) {
481  /* single diffractive to A+X */
482  success = string_process_->next_SDiff(true);
483  break;
485  /* single diffractive to X+B */
486  success = string_process_->next_SDiff(false);
487  break;
489  /* double diffractive */
490  success = string_process_->next_DDiff();
491  break;
493  /* soft non-diffractive */
494  success = string_process_->next_NDiffSoft();
495  break;
497  /* soft BBbar 2 mesonic annihilation */
498  success = string_process_->next_BBbarAnn();
499  break;
501  success = string_process_->next_NDiffHard();
502  break;
503  default:
504  log.error("Unknown string process required.");
505  success = false;
506  }
507  }
508  if (ntry == ntry_max) {
509  /* If pythia fails to form a string, it is usually because the energy
510  * is not large enough. In this case, an elastic scattering happens.
511  *
512  * Since particles are normally added after process selection for
513  * strings, outgoing_particles is still uninitialized, and memory
514  * needs to be allocated. We also shift the process_type_ to elastic
515  * so that sample_angles does a proper treatment. */
516  outgoing_particles_.reserve(2);
521  } else {
523  /* If the incoming particles already were unformed, the formation
524  * times and cross section scaling factors need to be adjusted */
525  const double tform_in = std::max(incoming_particles_[0].formation_time(),
526  incoming_particles_[1].formation_time());
527  if (tform_in > time_of_execution_) {
528  const double fin =
529  (incoming_particles_[0].formation_time() >
530  incoming_particles_[1].formation_time())
531  ? incoming_particles_[0].initial_xsec_scaling_factor()
532  : incoming_particles_[1].initial_xsec_scaling_factor();
533  for (size_t i = 0; i < outgoing_particles_.size(); i++) {
534  const double tform_out = outgoing_particles_[i].formation_time();
535  const double fout =
536  outgoing_particles_[i].initial_xsec_scaling_factor();
537  /* The new cross section scaling factor will be the product of the
538  * cross section scaling factor of the ingoing particles and of the
539  * outgoing ones (since the outgoing ones are also string fragments
540  * and thus take time to form). */
541  outgoing_particles_[i].set_cross_section_scaling_factor(fin * fout);
542  /* If the unformed incoming particles' formation time is larger than
543  * the current outgoing particle's formation time, then the latter
544  * is overwritten by the former*/
545  if (tform_in > tform_out) {
546  outgoing_particles_[i].set_slow_formation_times(time_of_execution_,
547  tform_in);
548  }
549  }
550  }
551  /* Check momentum difference for debugging */
552  FourVector out_mom;
553  for (ParticleData data : outgoing_particles_) {
554  out_mom += data.momentum();
555  }
556  log.debug("Incoming momenta string:", total_momentum());
557  log.debug("Outgoing momenta string:", out_mom);
558  }
559  }
560 }
561 
562 void ScatterAction::format_debug_output(std::ostream &out) const {
563  out << "Scatter of " << incoming_particles_;
564  if (outgoing_particles_.empty()) {
565  out << " (not performed)";
566  } else {
567  out << " to " << outgoing_particles_;
568  }
569 }
570 
571 } // namespace smash
void add_collisions(CollisionBranchList pv)
Add several new collision channels at once.
void add_collision(CollisionBranchPtr p)
Add a new collision channel.
bool is_Delta() const
Definition: pdgcode.h:336
double string_probability(bool strings_switch, bool use_transition_probability, bool use_AQM, bool treat_nnbar_with_strings) const
double high_energy() const
Determine the parametrized total cross section at high energies for the given collision, which is non-zero for Baryon-Baryon and Nucleon-Pion scatterings currently.
bool is_Deltastar() const
Definition: particletype.h:220
double partial_cross_section_
Partial cross-section to the chosen outgoing channel.
bool next_NDiffSoft()
Soft Non-diffractive process is modelled in accordance with dual-topological approach Capella:1978ig...
int antiparticle_sign() const
Definition: pdgcode.h:537
The ThreeVector class represents a physical three-vector with the components .
Definition: threevector.h:30
void generate_final_state() override
Generate the final-state of the scattering process.
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:163
double diffractive. Two strings are formed, one from A and one from B.
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:68
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:34
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:310
bool is_string_soft_process(ProcessType p)
Check if a given process type is a soft string excitation.
double cm_momentum() const
Get the momentum of the center of mass of the incoming particles in the calculation frame...
ProcessType process_type_
type of process
Definition: action.h:320
double gamma_cm() const
Get the gamma factor corresponding to a boost to the center of mass frame of the colliding particles...
Thrown when ScatterAction is called to perform with unknown ProcessType.
ThreeVector beta_cm() const
Get the velocity of the center of mass of the scattering particles in the calculation frame...
Guard type that safely disables floating point traps for the scope in which it is placed...
Definition: fpenvironment.h:79
Collection of useful constants that are known at compile time.
bool is_nucleon() const
Definition: pdgcode.h:324
a special case of baryon-antibaryon annihilation.
T pCM_sqr(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:91
StringProcess * string_process_
Pointer to interface class for strings.
FourVector total_momentum() const
Sum of 4-momenta of incoming particles.
Definition: action.h:323
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:66
constexpr T pow_int(const T base, unsigned const exponent)
Efficient template for calculating integer powers using squaring.
Definition: pow.h:23
2->2 inelastic scattering
const FourVector & momentum() const
Get the particle&#39;s 4-momentum.
Definition: particledata.h:139
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:52
double get_partial_weight() const override
Get the partial cross section of the chosen channel.
CollisionBranchList string_excitation(double total_string_xs, StringProcess *string_process, bool use_AQM) const
Determine the cross section for string excitations, which is given by the difference between the para...
ProcessType get_type() const override
void resonance_formation()
Perform a 2->1 resonance-formation process.
NNbarTreatment
Treatment of N Nbar Annihilation.
void string_excitation()
Todo(ryu): document better - it is not really UrQMD-based, isn&#39;t it? Perform the UrQMD-based string e...
void boost(const ThreeVector &v)
Apply a full Lorentz boost of momentum and position.
Definition: particledata.h:304
T canonical()
Definition: random.h:110
static double Cugnon_bnp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic np...
double plab_from_s(double mandelstam_s, double mass)
Convert Mandelstam-s to p_lab in a fixed-target collision.
Definition: kinematics.h:157
double sqr() const
calculate the square of the vector (which is a scalar)
Definition: fourvector.h:437
ScatterAction(const ParticleData &in_part1, const ParticleData &in_part2, double time, bool isotropic=false, double string_formation_time=1.0)
Construct a ScatterAction object.
double total_cross_section_
Total hadronic cross section.
void rotate_z_axis_to(ThreeVector &r)
Rotate the z-axis onto the vector r.
Definition: threevector.h:309
double mandelstam_s() const
Determine the Mandelstam s variable,.
double weight() const
ThreeVector threevec() const
Definition: fourvector.h:306
ParticleList particle_list() const
elastic scattering: particles remain the same, only momenta change
void init(const ParticleList &incoming, double tcoll)
initialization feed intial particles, time of collision and gamma factor of the center of mass...
int charge() const
The charge of the particle.
Definition: particletype.h:178
Thrown for example when ScatterAction is called to perform with a wrong number of final-state particl...
Definition: action.h:297
CollisionBranch is a derivative of ProcessBranch, which is used to represent particular final-state c...
bool next_SDiff(bool is_AB_to_AX)
Single-diffractive process is based on single pomeron exchange described in Ingelman:1984ns.
const FourVector & position() const
Get the particle&#39;s position in Minkowski space.
Definition: particledata.h:185
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
Definition: action.h:311
bool next_DDiff()
Double-diffractive process ( A + B -> X + X ) is similar to the single-diffractive process...
void set_4momentum(const FourVector &momentum_vector)
Set the particle&#39;s 4-momentum directly.
Definition: particledata.h:145
double sqr() const
Definition: threevector.h:249
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:303
hard string process involving 2->2 QCD process by PYTHIA.
std::pair< FourVector, FourVector > get_potential_at_interaction_point() const
Get the skyrme and asymmetry potential at the interaction point.
Definition: action.cc:78
CollisionBranchList generate_collision_list(double elastic_parameter, bool two_to_one_switch, ReactionsBitSet included_2to2, double low_snn_cut, bool strings_switch, bool use_AQM, bool strings_with_probability, NNbarTreatment nnbar_treatment, StringProcess *string_process) const
Generate a list of all possible collisions between the incoming particles with the given c...
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
Definition: action.cc:123
double cm_momentum_squared() const
Get the squared momentum of the center of mass of the incoming particles in the calculation frame...
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.
void elastic_scattering()
Perform an elastic two-body scattering, i.e. just exchange momentum.
std::bitset< 6 > ReactionsBitSet
Container for the 2 to 2 reactions in the code.
Action is the base class for a generic process that takes a number of incoming particles and transfor...
Definition: action.h:34
bool isotropic_
Do this collision isotropically?
ParticleList get_final_state()
a function to get the final state particle list which is called after the collision ...
Use string fragmentation.
static double Cugnon_bpp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic pp...
T power(T n, T xMin, T xMax)
Draws a random number according to a power-law distribution ~ x^n.
Definition: random.h:200
The cross section class assembels everything that is needed to calculate the cross section and return...
Definition: crosssections.h:58
bool is_Nstar() const
Definition: particletype.h:211
constexpr int p
Proton.
const ParticleType & type() const
Get the type of the particle.
Definition: particledata.h:109
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:208
PdgCode pdgcode() const
Get the pdgcode of the particle.
Definition: particledata.h:81
resonance formation (2->1)
CollisionBranchList collision_channels_
List of possible collisions.
double get_total_weight() const override
Get the total cross section of scattering particles.
Angles provides a common interface for generating directions: i.e., two angles that should be interpr...
Definition: angles.h:59
virtual double cross_section() const
Get the total cross section of the scattering particles.
bool next_BBbarAnn()
Baryon-antibaryon annihilation process Based on what UrQMD Bass:1998ca, Bleicher:1999xi does...
void sample_angles(std::pair< double, double > masses, double kinetic_energy_cm) override
Sample final-state angles in a 2->2 collision (possibly anisotropic).
non-diffractive. Two strings are formed both have ends in A and B.
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:79
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:32
bool next_NDiffHard()
Hard Non-diffractive process is based on PYTHIA 8 with partonic showers and interactions.
double transverse_distance_sqr() const
Calculate the transverse distance of the two incoming particles in their local rest frame...
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:52
(41-45) soft string excitations.
Definition: action.h:24
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:317
void format_debug_output(std::ostream &out) const override
Writes information about this scatter action to the out stream.
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:265
void inelastic_scattering()
Perform an inelastic two-body scattering, i.e. new particles are formed.