Version: SMASH-3.2
scatteractionsfinder.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2025
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
11 
12 #include <algorithm>
13 #include <map>
14 #include <vector>
15 
16 #include "smash/constants.h"
17 #include "smash/decaymodes.h"
18 #include "smash/logging.h"
19 #include "smash/parametrizations.h"
20 #include "smash/scatteraction.h"
23 #include "smash/stringfunctions.h"
24 
25 namespace smash {
26 static constexpr int LFindScatter = LogArea::FindScatter::id;
35  Configuration& config, const ExperimentParameters& parameters)
36  : finder_parameters_(config, parameters),
37  isotropic_(config.take(InputKeys::collTerm_isotropic)),
38  box_length_(parameters.box_length),
39  string_formation_time_(
40  config.take(InputKeys::collTerm_stringParam_formationTime)) {
42  logg[LFindScatter].info(
43  "Constant elastic isotropic cross-section mode:", " using ",
44  finder_parameters_.elastic_parameter, " mb as maximal cross-section.");
45  }
48  throw std::invalid_argument(
49  "Multi-body reactions (like e.g. 3->1 or 3->2) are only possible with "
50  "the stochastic "
51  "collision "
52  "criterion. Change your config accordingly.");
53  }
54 
57  1 &&
59  .included_2to2[IncludedReactions::PiDeuteron_to_pidprime] == 1 ||
61  .included_2to2[IncludedReactions::NDeuteron_to_Ndprime] == 1)) {
62  throw std::invalid_argument(
63  "To prevent double counting it is not possible to enable deuteron 3->2 "
64  "reactions\nand reactions involving the d' at the same time\ni.e. to "
65  "include \"Deuteron_3to2\" in `Multi_Particle_Reactions` and\n "
66  "\"PiDeuteron_to_pidprime\" "
67  "or \"NDeuteron_to_Ndprime\" in `Included_2to2` at the same time.\n"
68  "Change your config accordingly.");
69  }
70 
73  1 &&
75  throw std::invalid_argument(
76  "Do not use the d' resonance and enable \"Deuteron_3to2\" "
77  "`Multi_Particle_Reactions` at the same time. Either use the direct "
78  "3-to-2 reactions or the d' together with \"PiDeuteron_to_pidprime\" "
79  "and \"NDeuteron_to_Ndprime\" in `Included_2to2`. Otherwise the "
80  "deuteron 3-to-2 reactions would be double counted.");
81  }
82 
86  1) ||
89  1 &&
91  throw std::invalid_argument(
92  "In order to conserve detailed balance, when \"NNbar_5to2\" is "
93  "included in\n`Multi_Particle_Reactions`, the `NNbarTreatment` has to "
94  "be set to \"two to five\" and vice versa.");
95  }
96 
99  throw std::invalid_argument(
100  "'NNbar' has to be in the list of allowed 2 to 2 processes "
101  "to enable annihilation to go through resonances");
102  }
103 
105  string_process_interface_ = std::make_unique<StringProcess>(
126  parameters.use_monash_tune_default.value()));
127  }
128 }
129 
131  Configuration& config) {
132  auto sqrts_range_Npi = config.take(InputKeys::collTerm_stringTrans_rangeNpi);
133  auto sqrts_range_NN = config.take(InputKeys::collTerm_stringTrans_rangeNN);
134 
135  if (sqrts_range_Npi.first < nucleon_mass + pion_mass) {
136  sqrts_range_Npi.first = nucleon_mass + pion_mass;
137  if (sqrts_range_Npi.second < sqrts_range_Npi.first)
138  sqrts_range_Npi.second = sqrts_range_Npi.first;
139  logg[LFindScatter].warn(
140  "Lower bound of Sqrts_Range_Npi too small, setting it to mass "
141  "threshold. New range is [",
142  sqrts_range_Npi.first, ',', sqrts_range_Npi.second, "] GeV");
143  }
144  if (sqrts_range_NN.first < 2 * nucleon_mass) {
145  sqrts_range_NN.first = 2 * nucleon_mass;
146  if (sqrts_range_NN.second < sqrts_range_NN.first)
147  sqrts_range_NN.second = sqrts_range_NN.first;
148  logg[LFindScatter].warn(
149  "Lower bound of Sqrts_Range_NN too small, setting it to mass "
150  "threshold. New range is [",
151  sqrts_range_NN.first, ',', sqrts_range_NN.second, "] GeV.");
152  }
153 
154  return {sqrts_range_Npi,
155  sqrts_range_NN,
160 }
161 
163  Configuration& config, const ExperimentParameters& parameters)
164  : elastic_parameter(config.take(InputKeys::collTerm_elasticCrossSection)),
165  low_snn_cut(parameters.low_snn_cut),
166  scale_xs(parameters.scale_xs),
167  additional_el_xs(
168  config.take(InputKeys::collTerm_additionalElasticCrossSection)),
169  maximum_cross_section(parameters.maximum_cross_section),
170  coll_crit(parameters.coll_crit),
171  nnbar_treatment(parameters.nnbar_treatment),
172  included_2to2(parameters.included_2to2),
173  included_multi(parameters.included_multi),
174  testparticles(parameters.testparticles),
175  two_to_one(parameters.two_to_one),
176  allow_collisions_within_nucleus(
177  config.take(InputKeys::modi_collider_collisionWithinNucleus)),
178  strings_switch(parameters.strings_switch),
179  use_AQM(config.take(InputKeys::collTerm_useAQM)),
180  strings_with_probability(
181  config.take(InputKeys::collTerm_stringsWithProbability)),
182  only_warn_for_high_prob(
183  config.take(InputKeys::collTerm_onlyWarnForHighProbability)),
184  transition_high_energy{create_string_transition_parameters(config)},
185  total_xs_strategy(config.take(InputKeys::collTerm_totXsStrategy)),
186  pseudoresonance_method(config.take(InputKeys::collTerm_pseudoresonance)),
187  AQM_charm_suppression(
188  config.take(InputKeys::collTerm_HF_AQMcSuppression)),
189  AQM_bottom_suppression(
190  config.take(InputKeys::collTerm_HF_AQMbSuppression)) {
192  logg[LFindScatter].info(
193  "Evaluating total cross sections from partial processes.");
194  } else if (parameters.included_2to2[IncludedReactions::Elastic] == 1 &&
195  parameters.included_2to2.count() == 1) {
196  throw std::invalid_argument(
197  "The BottomUp strategy for total cross section evaluation is needed to "
198  "have only elastic interactions, please change the configuration "
199  "accordingly.");
201  logg[LFindScatter].info(
202  "Evaluating total cross sections from parametrizations.");
204  logg[LFindScatter].info(
205  "Evaluating total cross sections from parametrizations only for "
206  "measured processes.");
207  }
208 
211  throw std::invalid_argument(
212  "Suppression factors for AQM should be between 0 and 1.");
213  }
214 }
215 
217  const ParticleData& data_a, const ParticleData& data_b, double dt,
218  const std::vector<FourVector>& beam_momentum,
219  const double gcell_vol) const {
220  /* If the two particles
221  * 1) belong to one of the two colliding nuclei, and
222  * 2) both of them have never experienced any collisions,
223  * then the collisions between them are banned. */
225  assert(data_a.id() >= 0);
226  assert(data_b.id() >= 0);
227  bool in_same_nucleus = (data_a.belongs_to() == BelongsTo::Projectile &&
228  data_b.belongs_to() == BelongsTo::Projectile) ||
229  (data_a.belongs_to() == BelongsTo::Target &&
230  data_b.belongs_to() == BelongsTo::Target);
231  bool never_interacted_before =
232  data_a.get_history().collisions_per_particle == 0 &&
233  data_b.get_history().collisions_per_particle == 0;
234  if (in_same_nucleus && never_interacted_before) {
235  return nullptr;
236  }
237  }
238 
239  // No grid or search in cell means no collision for stochastic criterion
241  gcell_vol < really_small) {
242  return nullptr;
243  }
244 
245  // Determine time of collision.
246  const double time_until_collision =
247  collision_time(data_a, data_b, dt, beam_momentum);
248 
249  // Check that collision happens in this timestep.
250  if (time_until_collision < 0. || time_until_collision >= dt) {
251  return nullptr;
252  }
253 
254  // Determine which total cross section to use
255  bool incoming_parametrized = (finder_parameters_.total_xs_strategy ==
259  const PdgCode& pdg_a = data_a.type().pdgcode();
260  const PdgCode& pdg_b = data_b.type().pdgcode();
261  incoming_parametrized = parametrization_exists(pdg_a, pdg_b);
262  }
263 
264  // Create ScatterAction object.
265  ScatterActionPtr act = std::make_unique<ScatterAction>(
266  data_a, data_b, time_until_collision, isotropic_, string_formation_time_,
267  box_length_, incoming_parametrized);
268 
270  act->set_stochastic_pos_idx();
271  }
272 
274  act->set_string_interface(string_process_interface_.get());
275  }
276 
277  // Distance squared calculation not needed for stochastic criterion
278  const double distance_squared =
280  ? act->transverse_distance_sqr()
282  ? act->cov_transverse_distance_sqr()
283  : 0.0;
284 
285  // Don't calculate cross section if the particles are very far apart.
286  // Not needed for stochastic criterion because of cell structure.
288  distance_squared >=
290  return nullptr;
291  }
292 
293  if (incoming_parametrized) {
294  act->set_parametrized_total_cross_section(finder_parameters_);
295  } else {
296  // Add various subprocesses.
297  act->add_all_scatterings(finder_parameters_);
298  }
299 
300  double xs = act->cross_section() * fm2_mb /
301  static_cast<double>(finder_parameters_.testparticles);
302 
303  // Take cross section scaling factors into account
304  xs *= data_a.xsec_scaling_factor(time_until_collision);
305  xs *= data_b.xsec_scaling_factor(time_until_collision);
306 
308  const double v_rel = act->relative_velocity();
309  /* Collision probability for 2-particle scattering, see
310  * \iref{Staudenmaier:2021lrg}. */
311  const double prob = xs * v_rel * dt / gcell_vol;
312 
313  logg[LFindScatter].debug(
314  "Stochastic collison criterion parameters (2-particles):\nprob = ",
315  prob, ", xs = ", xs, ", v_rel = ", v_rel, ", dt = ", dt,
316  ", gcell_vol = ", gcell_vol,
317  ", testparticles = ", finder_parameters_.testparticles);
318 
319  if (prob > 1.) {
320  std::stringstream err;
321  err << "Probability larger than 1 for stochastic rates. ( P_22 = " << prob
322  << " )\n"
323  << data_a.type().name() << data_b.type().name() << " with masses "
324  << data_a.momentum().abs() << " and " << data_b.momentum().abs()
325  << " at sqrts[GeV] = " << act->sqrt_s()
326  << " with xs[fm^2]/Ntest = " << xs
327  << "\nConsider using smaller timesteps.";
329  logg[LFindScatter].warn(err.str());
330  } else {
331  throw std::runtime_error(err.str());
332  }
333  }
334 
335  // probability criterion
336  double random_no = random::uniform(0., 1.);
337  if (random_no > prob) {
338  return nullptr;
339  }
340 
343  // just collided with this particle
344  if (data_a.id_process() > 0 && data_a.id_process() == data_b.id_process()) {
345  logg[LFindScatter].debug("Skipping collided particles at time ",
346  data_a.position().x0(), " due to process ",
347  data_a.id_process(), "\n ", data_a, "\n<-> ",
348  data_b);
349 
350  return nullptr;
351  }
352 
353  // Cross section for collision criterion
354  const double cross_section_criterion = xs * M_1_PI;
355 
356  // distance criterion according to cross_section
357  if (distance_squared >= cross_section_criterion) {
358  return nullptr;
359  }
360 
361  logg[LFindScatter].debug("particle distance squared: ", distance_squared,
362  "\n ", data_a, "\n<-> ", data_b);
363  }
364 
365  // Include possible outgoing branches
366  if (incoming_parametrized) {
367  act->add_all_scatterings(finder_parameters_);
368  }
369 
370  return act;
371 }
372 
374  const ParticleList& plist, double dt, const double gcell_vol) const {
375  /* If all particles
376  * 1) belong to the two colliding nuclei
377  * 2) are within the same nucleus
378  * 3) have never experienced any collisons,
379  * then the collision between them are banned also for multi-particle
380  * interactions. */
382  bool all_projectile =
383  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
384  return data.belongs_to() == BelongsTo::Projectile;
385  });
386  bool all_target =
387  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
388  return data.belongs_to() == BelongsTo::Target;
389  });
390  bool none_collided =
391  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
392  return data.get_history().collisions_per_particle == 0;
393  });
394  if ((all_projectile || all_target) && none_collided) {
395  return nullptr;
396  }
397  }
398  // No grid or search in cell
399  if (gcell_vol < really_small) {
400  return nullptr;
401  }
402 
403  /* Optimisation for later: Already check here at the beginning
404  * if collision with plist is possible before constructing actions. */
405 
406  // 1. Determine time of collision.
407  const double time_until_collision = dt * random::uniform(0., 1.);
408 
409  // 2. Create ScatterAction object.
410  ScatterActionMultiPtr act =
411  std::make_unique<ScatterActionMulti>(plist, time_until_collision);
412 
413  act->set_stochastic_pos_idx();
414 
415  // 3. Add possible final states (dt and gcell_vol for probability calculation)
416  act->add_possible_reactions(dt, gcell_vol, finder_parameters_.included_multi);
417 
418  /* 4. Return total collision probability
419  * Scales with 1 over the number of testpartciles to the power of the
420  * number of incoming particles - 1 */
421  const double prob =
422  act->get_total_weight() /
423  std::pow(finder_parameters_.testparticles, plist.size() - 1);
424 
425  // 5. Check that probability is smaller than one
426  if (prob > 1.) {
427  std::stringstream err;
428  err << "Probability " << prob << " larger than 1 for stochastic rates for ";
429  for (const ParticleData& data : plist) {
430  err << data.type().name();
431  }
432  err << " at sqrts[GeV] = " << act->sqrt_s()
433  << "\nConsider using smaller timesteps.";
435  logg[LFindScatter].warn(err.str());
436  } else {
437  throw std::runtime_error(err.str());
438  }
439  }
440 
441  // 6. Perform probability decisions
442  double random_no = random::uniform(0., 1.);
443  if (random_no > prob) {
444  return nullptr;
445  }
446 
447  return act;
448 }
449 
451  const ParticleList& search_list, double dt, const double gcell_vol,
452  const std::vector<FourVector>& beam_momentum) const {
453  std::vector<ActionPtr> actions;
454  for (const ParticleData& p1 : search_list) {
455  for (const ParticleData& p2 : search_list) {
456  // Check for 2 particle scattering
457  if (p1.id() < p2.id()) {
458  ActionPtr act =
459  check_collision_two_part(p1, p2, dt, beam_momentum, gcell_vol);
460  if (act) {
461  actions.push_back(std::move(act));
462  }
463  }
465  // Also, check for 3 particle scatterings with stochastic criterion
466  for (const ParticleData& p3 : search_list) {
471  if (p1.id() < p2.id() && p2.id() < p3.id()) {
472  ActionPtr act =
473  check_collision_multi_part({p1, p2, p3}, dt, gcell_vol);
474  if (act) {
475  actions.push_back(std::move(act));
476  }
477  }
478  }
479  for (const ParticleData& p4 : search_list) {
482  if (p1.id() < p2.id() && p2.id() < p3.id() && p3.id() < p4.id()) {
483  ActionPtr act =
484  check_collision_multi_part({p1, p2, p3, p4}, dt, gcell_vol);
485  if (act) {
486  actions.push_back(std::move(act));
487  }
488  }
489  }
492  search_list.size() >= 5) {
493  for (const ParticleData& p5 : search_list) {
494  if ((p1.id() < p2.id() && p2.id() < p3.id() &&
495  p3.id() < p4.id() && p4.id() < p5.id()) &&
496  (p1.is_pion() && p2.is_pion() && p3.is_pion() &&
497  p4.is_pion() && p5.is_pion())) {
498  // at the moment only pure pion 5-body reactions
499  ActionPtr act = check_collision_multi_part(
500  {p1, p2, p3, p4, p5}, dt, gcell_vol);
501  if (act) {
502  actions.push_back(std::move(act));
503  }
504  }
505  }
506  }
507  }
508  }
509  }
510  }
511  }
512  return actions;
513 }
514 
516  const ParticleList& search_list, const ParticleList& neighbors_list,
517  double dt, const std::vector<FourVector>& beam_momentum) const {
518  std::vector<ActionPtr> actions;
520  // Only search in cells
521  return actions;
522  }
523  for (const ParticleData& p1 : search_list) {
524  for (const ParticleData& p2 : neighbors_list) {
525  assert(p1.id() != p2.id());
526  // Check if a collision is possible.
527  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
528  if (act) {
529  actions.push_back(std::move(act));
530  }
531  }
532  }
533  return actions;
534 }
535 
537  const ParticleList& search_list, const Particles& surrounding_list,
538  double dt, const std::vector<FourVector>& beam_momentum) const {
539  std::vector<ActionPtr> actions;
541  // Only search in cells
542  return actions;
543  }
544  for (const ParticleData& p2 : surrounding_list) {
545  /* don't look for collisions if the particle from the surrounding list is
546  * also in the search list */
547  auto result = std::find_if(
548  search_list.begin(), search_list.end(),
549  [&p2](const ParticleData& p) { return p.id() == p2.id(); });
550  if (result != search_list.end()) {
551  continue;
552  }
553  for (const ParticleData& p1 : search_list) {
554  // Check if a collision is possible.
555  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
556  if (act) {
557  actions.push_back(std::move(act));
558  }
559  }
560  }
561  return actions;
562 }
563 
565  constexpr double time = 0.0;
566 
567  const size_t N_isotypes = IsoParticleType::list_all().size();
568  const size_t N_pairs = N_isotypes * (N_isotypes - 1) / 2;
569 
570  std::cout << N_isotypes << " iso-particle types." << std::endl;
571  std::cout << "They can make " << N_pairs << " pairs." << std::endl;
572  std::vector<double> momentum_scan_list = {0.1, 0.3, 0.5, 1.0,
573  2.0, 3.0, 5.0, 10.0};
574  for (const IsoParticleType& A_isotype : IsoParticleType::list_all()) {
575  for (const IsoParticleType& B_isotype : IsoParticleType::list_all()) {
576  if (&A_isotype > &B_isotype) {
577  continue;
578  }
579  bool any_nonzero_cs = false;
580  std::vector<std::string> r_list;
581  for (const ParticleTypePtr A_type : A_isotype.get_states()) {
582  for (const ParticleTypePtr B_type : B_isotype.get_states()) {
583  if (A_type > B_type) {
584  continue;
585  }
586  ParticleData A(*A_type), B(*B_type);
587  for (auto mom : momentum_scan_list) {
588  A.set_4momentum(A.pole_mass(), mom, 0.0, 0.0);
589  B.set_4momentum(B.pole_mass(), -mom, 0.0, 0.0);
590  ScatterActionPtr act = std::make_unique<ScatterAction>(
591  A, B, time, isotropic_, string_formation_time_, -1, false);
593  act->set_string_interface(string_process_interface_.get());
594  }
595  act->add_all_scatterings(finder_parameters_);
596  const double total_cs = act->cross_section();
597  if (total_cs <= 0.0) {
598  continue;
599  }
600  any_nonzero_cs = true;
601  for (const auto& channel : act->collision_channels()) {
602  const auto type = channel->get_type();
603  std::string r;
604  if (is_string_soft_process(type) ||
605  type == ProcessType::StringHard) {
606  r = A_type->name() + B_type->name() + std::string(" → strings");
607  } else {
608  std::string r_type =
609  (type == ProcessType::Elastic) ? std::string(" (el)")
610  : (channel->get_type() == ProcessType::TwoToTwo)
611  ? std::string(" (inel)")
612  : std::string(" (?)");
613  r = A_type->name() + B_type->name() + std::string(" → ") +
614  channel->particle_types()[0]->name() +
615  channel->particle_types()[1]->name() + r_type;
616  }
617  isoclean(r);
618  r_list.push_back(r);
619  }
620  }
621  }
622  }
623  std::sort(r_list.begin(), r_list.end());
624  r_list.erase(std::unique(r_list.begin(), r_list.end()), r_list.end());
625  if (any_nonzero_cs) {
626  for (auto r : r_list) {
627  std::cout << r;
628  if (r_list.back() != r) {
629  std::cout << ", ";
630  }
631  }
632  std::cout << std::endl;
633  }
634  }
635  }
636 }
637 
641  std::string name_;
642 
645 
647  double mass_;
648 
657  FinalStateCrossSection(const std::string& name, double cross_section,
658  double mass)
659  : name_(name), cross_section_(cross_section), mass_(mass) {}
660 };
661 
662 namespace decaytree {
663 
675 struct Node {
676  public:
678  std::string name_;
679 
681  double weight_;
682 
684  ParticleTypePtrList initial_particles_;
685 
687  ParticleTypePtrList final_particles_;
688 
690  ParticleTypePtrList state_;
691 
693  std::vector<Node> children_;
694 
696  Node(const Node&) = delete;
698  Node(Node&&) = default;
699 
710  Node(const std::string& name, double weight,
711  ParticleTypePtrList&& initial_particles,
712  ParticleTypePtrList&& final_particles, ParticleTypePtrList&& state,
713  std::vector<Node>&& children)
714  : name_(name),
715  weight_(weight),
716  initial_particles_(std::move(initial_particles)),
717  final_particles_(std::move(final_particles)),
718  state_(std::move(state)),
719  children_(std::move(children)) {}
720 
732  Node& add_action(const std::string& name, double weight,
733  ParticleTypePtrList&& initial_particles,
734  ParticleTypePtrList&& final_particles) {
735  // Copy parent state and update it.
736  ParticleTypePtrList state(state_);
737  for (const auto& p : initial_particles) {
738  state.erase(std::find(state.begin(), state.end(), p));
739  }
740  for (const auto& p : final_particles) {
741  state.push_back(p);
742  }
743  // Sort the state to normalize the output.
744  std::sort(state.begin(), state.end(),
746  return a->name() < b->name();
747  });
748  // Push new node to children.
749  Node new_node(name, weight, std::move(initial_particles),
750  std::move(final_particles), std::move(state), {});
751  children_.emplace_back(std::move(new_node));
752  return children_.back();
753  }
754 
756  void print() const { print_helper(0); }
757 
761  std::vector<FinalStateCrossSection> final_state_cross_sections() const {
762  std::vector<FinalStateCrossSection> result;
763  final_state_cross_sections_helper(0, result, "", 1.);
764  return result;
765  }
766 
767  private:
774  void print_helper(uint64_t depth) const {
775  for (uint64_t i = 0; i < depth; i++) {
776  std::cout << " ";
777  }
778  std::cout << name_ << " " << weight_ << std::endl;
779  for (const auto& child : children_) {
780  child.print_helper(depth + 1);
781  }
782  }
783 
796  uint64_t depth, std::vector<FinalStateCrossSection>& result,
797  const std::string& name, double weight,
798  bool show_intermediate_states = false) const {
799  // The first node corresponds to the total cross section and has to be
800  // ignored. The second node corresponds to the partial cross section. All
801  // further nodes correspond to branching ratios.
802  if (depth > 0) {
803  weight *= weight_;
804  }
805 
806  std::string new_name;
807  double mass = 0.;
808 
809  if (show_intermediate_states) {
810  new_name = name;
811  if (!new_name.empty()) {
812  new_name += "->";
813  }
814  new_name += name_;
815  new_name += "{";
816  } else {
817  new_name = "";
818  }
819  for (const auto& s : state_) {
820  new_name += s->name();
821  mass += s->mass();
822  }
823  if (show_intermediate_states) {
824  new_name += "}";
825  }
826 
827  if (children_.empty()) {
828  result.emplace_back(FinalStateCrossSection(new_name, weight, mass));
829  return;
830  }
831  for (const auto& child : children_) {
832  child.final_state_cross_sections_helper(depth + 1, result, new_name,
833  weight, show_intermediate_states);
834  }
835  }
836 };
837 
846 static std::string make_decay_name(const std::string& res_name,
847  const DecayBranchPtr& decay,
848  ParticleTypePtrList& final_state) {
849  std::stringstream name;
850  name << "[" << res_name << "->";
851  for (const auto& p : decay->particle_types()) {
852  name << p->name();
853  final_state.push_back(p);
854  }
855  name << "]";
856  return name.str();
857 }
858 
866 static void add_decays(Node& node, double sqrts) {
867  // If there is more than one unstable particle in the current state, then
868  // there will be redundant paths in the decay tree, corresponding to
869  // reorderings of the decays. To avoid double counting, we normalize by the
870  // number of possible decay orderings. Normalizing by the number of unstable
871  // particles recursively corresponds to normalizing by the factorial that
872  // gives the number of reorderings.
873  //
874  // Ideally, the redundant paths should never be added to the decay tree, but
875  // we never have more than two redundant paths, so it probably does not
876  // matter much.
877  uint32_t n_unstable = 0;
878  double sqrts_minus_masses = sqrts;
879  for (const ParticleTypePtr ptype : node.state_) {
880  if (!ptype->is_stable()) {
881  n_unstable += 1;
882  }
883  sqrts_minus_masses -= ptype->mass();
884  }
885  const double norm =
886  n_unstable != 0 ? 1. / static_cast<double>(n_unstable) : 1.;
887 
888  for (const ParticleTypePtr ptype : node.state_) {
889  if (!ptype->is_stable()) {
890  const double sqrts_decay = sqrts_minus_masses + ptype->mass();
891  bool can_decay = false;
892  for (const auto& decay : ptype->decay_modes().decay_mode_list()) {
893  // Make sure to skip kinematically impossible decays.
894  // In principle, we would have to integrate over the mass of the
895  // resonance, but as an approximation we just assume it at its pole.
896  double final_state_mass = 0.;
897  for (const auto& p : decay->particle_types()) {
898  final_state_mass += p->mass();
899  }
900  if (final_state_mass > sqrts_decay) {
901  continue;
902  }
903  can_decay = true;
904 
905  ParticleTypePtrList parts;
906  const auto name = make_decay_name(ptype->name(), decay, parts);
907  auto& new_node = node.add_action(name, norm * decay->weight(), {ptype},
908  std::move(parts));
909  add_decays(new_node, sqrts_decay);
910  }
911  if (!can_decay) {
912  // Remove final-state cross sections with resonances that cannot
913  // decay due to our "mass = pole mass" approximation.
914  node.weight_ = 0;
915  return;
916  }
917  }
918  }
919 }
920 
921 } // namespace decaytree
922 
928 static void deduplicate(std::vector<FinalStateCrossSection>& final_state_xs) {
929  std::sort(final_state_xs.begin(), final_state_xs.end(),
930  [](const FinalStateCrossSection& a,
931  const FinalStateCrossSection& b) { return a.name_ < b.name_; });
932  auto current = final_state_xs.begin();
933  while (current != final_state_xs.end()) {
934  auto adjacent = std::adjacent_find(
935  current, final_state_xs.end(),
936  [](const FinalStateCrossSection& a, const FinalStateCrossSection& b) {
937  return a.name_ == b.name_;
938  });
939  current = adjacent;
940  if (adjacent != final_state_xs.end()) {
941  adjacent->cross_section_ += (adjacent + 1)->cross_section_;
942  final_state_xs.erase(adjacent + 1);
943  }
944  }
945 }
946 
948  const ParticleType& a, const ParticleType& b, double m_a, double m_b,
949  bool final_state, std::vector<double>& plab) const {
950  typedef std::vector<std::pair<double, double>> xs_saver;
951  std::map<std::string, xs_saver> xs_dump;
952  std::map<std::string, double> outgoing_total_mass;
953 
954  ParticleData a_data(a), b_data(b);
955  int n_momentum_points = 200;
956  constexpr double momentum_step = 0.02;
957  if (plab.size() > 0) {
958  n_momentum_points = plab.size();
959  // Remove duplicates.
960  std::sort(plab.begin(), plab.end());
961  plab.erase(std::unique(plab.begin(), plab.end()), plab.end());
962  }
963  for (int i = 0; i < n_momentum_points; i++) {
964  double momentum;
965  if (plab.size() > 0) {
966  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
967  } else {
968  momentum = momentum_step * (i + 1);
969  }
970  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
971  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
972  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
973  const ParticleList incoming = {a_data, b_data};
974  ScatterActionPtr act = std::make_unique<ScatterAction>(
975  a_data, b_data, 0., isotropic_, string_formation_time_, -1, false);
977  act->set_string_interface(string_process_interface_.get());
978  }
979  act->add_all_scatterings(finder_parameters_);
980  decaytree::Node tree(a.name() + b.name(), act->cross_section(), {&a, &b},
981  {&a, &b}, {&a, &b}, {});
982  const CollisionBranchList& processes = act->collision_channels();
983  for (const auto& process : processes) {
984  const double xs = process->weight();
985  if (xs <= 0.0) {
986  continue;
987  }
988  if (!final_state) {
989  std::stringstream process_description_stream;
990  process_description_stream << *process;
991  const std::string& description = process_description_stream.str();
992  double m_tot = 0.0;
993  for (const auto& ptype : process->particle_types()) {
994  m_tot += ptype->mass();
995  }
996  outgoing_total_mass[description] = m_tot;
997  if (!xs_dump[description].empty() &&
998  std::abs(xs_dump[description].back().first - sqrts) <
999  really_small) {
1000  xs_dump[description].back().second += xs;
1001  } else {
1002  xs_dump[description].push_back(std::make_pair(sqrts, xs));
1003  }
1004  } else {
1005  std::stringstream process_description_stream;
1006  process_description_stream << *process;
1007  const std::string& description = process_description_stream.str();
1008  ParticleTypePtrList initial_particles = {&a, &b};
1009  ParticleTypePtrList final_particles = process->particle_types();
1010  auto& process_node =
1011  tree.add_action(description, xs, std::move(initial_particles),
1012  std::move(final_particles));
1013  decaytree::add_decays(process_node, sqrts);
1014  }
1015  }
1016  xs_dump["total"].push_back(std::make_pair(sqrts, act->cross_section()));
1017  // Total cross-section should be the first in the list -> negative mass
1018  outgoing_total_mass["total"] = -1.0;
1019  if (final_state) {
1020  // tree.print();
1021  auto final_state_xs = tree.final_state_cross_sections();
1022  deduplicate(final_state_xs);
1023  for (const auto& p : final_state_xs) {
1024  // Don't print empty columns.
1025  //
1026  // FIXME(steinberg): The better fix would be to not have them in the
1027  // first place.
1028  if (p.name_ == "") {
1029  continue;
1030  }
1031  outgoing_total_mass[p.name_] = p.mass_;
1032  xs_dump[p.name_].push_back(std::make_pair(sqrts, p.cross_section_));
1033  }
1034  }
1035  }
1036  // Get rid of cross sections that are zero.
1037  // (This only happens if their is a resonance in the final state that cannot
1038  // decay with our simplified assumptions.)
1039  for (auto it = begin(xs_dump); it != end(xs_dump);) {
1040  // Sum cross section over all energies.
1041  const xs_saver& xs = (*it).second;
1042  double sum = 0;
1043  for (const auto& p : xs) {
1044  sum += p.second;
1045  }
1046  if (sum == 0.) {
1047  it = xs_dump.erase(it);
1048  } else {
1049  ++it;
1050  }
1051  }
1052 
1053  // Nice ordering of channels by summed pole mass of products
1054  std::vector<std::string> all_channels;
1055  for (const auto& channel : xs_dump) {
1056  all_channels.push_back(channel.first);
1057  }
1058  std::sort(all_channels.begin(), all_channels.end(),
1059  [&](const std::string& str_a, const std::string& str_b) {
1060  return outgoing_total_mass[str_a] < outgoing_total_mass[str_b];
1061  });
1062 
1063  // Print header
1064  std::cout << "# Dumping partial " << a.name() << b.name()
1065  << " cross-sections in mb, energies in GeV" << std::endl;
1066  std::cout << " sqrt_s";
1067  // Align everything to 24 unicode characters.
1068  // This should be enough for the longest channel name: 7 final-state
1069  // particles, or 2 of the longest named resonances (currently Ds0*(2317)⁺).
1070  for (const auto& channel : all_channels) {
1071  std::cout << utf8::fill_left(channel, 24, ' ');
1072  }
1073  std::cout << std::endl;
1074 
1075  // Print out all partial cross-sections in mb
1076  for (int i = 0; i < n_momentum_points; i++) {
1077  double momentum;
1078  if (plab.size() > 0) {
1079  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
1080  } else {
1081  momentum = momentum_step * (i + 1);
1082  }
1083  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
1084  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
1085  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
1086  std::printf("%9.6f", sqrts);
1087  for (const auto& channel : all_channels) {
1088  const xs_saver energy_and_xs = xs_dump[channel];
1089  size_t j = 0;
1090  for (; j < energy_and_xs.size() && energy_and_xs[j].first < sqrts; j++) {
1091  }
1092  double xs = 0.0;
1093  if (j < energy_and_xs.size() &&
1094  std::abs(energy_and_xs[j].first - sqrts) < really_small) {
1095  xs = energy_and_xs[j].second;
1096  }
1097  std::printf("%16.6f", xs); // Same alignment as in the header.
1098  }
1099  std::printf("\n");
1100  }
1101 }
1102 
1103 } // namespace smash
Interface to the SMASH configuration files.
T take(const Key< T > &key)
The default interface for SMASH to read configuration values.
const DecayBranchList & decay_mode_list() const
Definition: decaymodes.h:63
double abs() const
calculate the lorentz invariant absolute value
Definition: fourvector.h:464
double x0() const
Definition: fourvector.h:313
IsoParticleType is a class to represent isospin multiplets.
static const IsoParticleTypeList & list_all()
Returns a list of all IsoParticleTypes.
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:58
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
Definition: particledata.h:164
const ParticleType & type() const
Get the type of the particle.
Definition: particledata.h:128
double xsec_scaling_factor(double delta_time=0.) const
Return the cross section scaling factor at a given time.
Definition: particledata.cc:86
uint32_t id_process() const
Get the id of the last action.
Definition: particledata.h:134
const FourVector & momentum() const
Get the particle's 4-momentum.
Definition: particledata.h:158
BelongsTo belongs_to() const
Getter for belongs_to label.
Definition: particledata.h:345
int32_t id() const
Get the id of the particle.
Definition: particledata.h:76
HistoryData get_history() const
Get history information.
Definition: particledata.h:139
double pole_mass() const
Get the particle's pole mass ("on-shell").
Definition: particledata.h:115
const FourVector & position() const
Get the particle's position in Minkowski space.
Definition: particledata.h:204
A pointer-like interface to global references to ParticleType objects.
Definition: particletype.h:679
Particle type contains the static properties of a particle species.
Definition: particletype.h:98
const DecayModes & decay_modes() const
static const ParticleTypePtr try_find(PdgCode pdgcode)
Returns the ParticleTypePtr for the given pdgcode.
Definition: particletype.cc:89
PdgCode pdgcode() const
Definition: particletype.h:157
const std::string & name() const
Definition: particletype.h:142
bool is_stable() const
Definition: particletype.h:246
double mass() const
Definition: particletype.h:145
The Particles class abstracts the storage and manipulation of particles.
Definition: particles.h:33
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:111
const bool only_warn_for_high_prob
Switch to turn off throwing an exception for collision probabilities larger than 1.
const int testparticles
Number of test particles.
ScatterActionsFinderParameters(Configuration &config, const ExperimentParameters &parameters)
Class constructor.
const ReactionsBitSet included_2to2
List of included 2<->2 reactions.
const double elastic_parameter
Elastic cross section parameter (in mb).
const bool strings_switch
Indicates whether string fragmentation is switched on.
const TotalCrossSectionStrategy total_xs_strategy
Method used to evaluate total cross sections for collision finding.
const double AQM_charm_suppression
Factor to reduce cross sections for charmed hadrons.
const MultiParticleReactionsBitSet included_multi
List of included multi-particle reactions.
const bool allow_collisions_within_nucleus
If particles within the same nucleus are allowed to collide for their first time.
const NNbarTreatment nnbar_treatment
Switch for NNbar reactions.
const double AQM_bottom_suppression
Factor to reduce cross sections for bottomed hadrons.
const CollisionCriterion coll_crit
Specifies which collision criterion is used.
ScatterActionsFinder(Configuration &config, const ExperimentParameters &parameters)
Constructor of the finder with the given parameters.
ScatterActionsFinderParameters finder_parameters_
Struct collecting several parameters.
ActionList find_actions_in_cell(const ParticleList &search_list, double dt, const double gcell_vol, const std::vector< FourVector > &beam_momentum) const override
Search for all the possible collisions within one cell.
void dump_cross_sections(const ParticleType &a, const ParticleType &b, double m_a, double m_b, bool final_state, std::vector< double > &plab) const
Print out partial cross-sections of all processes that can occur in the collision of a(mass = m_a) an...
const bool isotropic_
Do all collisions isotropically.
double max_transverse_distance_sqr(int testparticles) const
The maximal distance over which particles can interact in case of the geometric criterion,...
ActionPtr check_collision_two_part(const ParticleData &data_a, const ParticleData &data_b, double dt, const std::vector< FourVector > &beam_momentum={}, const double gcell_vol=0.0) const
Check for a single pair of particles (id_a, id_b) if a collision will happen in the next timestep and...
const double string_formation_time_
Parameter for formation time.
std::unique_ptr< StringProcess > string_process_interface_
Class that deals with strings, interfacing Pythia.
void dump_reactions() const
Prints out all the 2-> n (n > 1) reactions with non-zero cross-sections between all possible pairs of...
bool is_constant_elastic_isotropic() const
If there is only one particle sort, no decays (only elastic scatterings are possible),...
ActionPtr check_collision_multi_part(const ParticleList &plist, double dt, const double gcell_vol) const
Check for multiple i.e.
const double box_length_
Box length: needed to determine coordinates of collision correctly in case of collision through the w...
ActionList find_actions_with_surrounding_particles(const ParticleList &search_list, const Particles &surrounding_list, double dt, const std::vector< FourVector > &beam_momentum) const override
Search for all the possible secondary collisions between the outgoing particles and the rest.
ActionList find_actions_with_neighbors(const ParticleList &search_list, const ParticleList &neighbors_list, double dt, const std::vector< FourVector > &beam_momentum) const override
Search for all the possible collisions among the neighboring cells.
double collision_time(const ParticleData &p1, const ParticleData &p2, double dt, const std::vector< FourVector > &beam_momentum) const
Determine the collision time of the two particles.
Collection of useful constants that are known at compile time.
@ TwoToFive
Directly create 5 pions, use with multi-particle reactions.
@ Resonances
Use intermediate Resonances.
@ TopDownMeasured
Mix the two above, using the parametrizations only for measured processes, and summing up partials fo...
@ TopDown
Use parametrizations based on existing data, rescaling with AQM for unmeasured processes.
@ BottomUp
Sum the existing partial contributions.
@ NNbar_5to2
@ A3_Nuclei_4to2
@ Deuteron_3to2
@ Meson_3to1
@ Stochastic
Stochastic Criteiron.
@ Geometric
Geometric criterion.
@ Covariant
Covariant Criterion.
@ PiDeuteron_to_pidprime
@ NDeuteron_to_Ndprime
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:40
static std::string make_decay_name(const std::string &res_name, const DecayBranchPtr &decay, ParticleTypePtrList &final_state)
Generate name for decay and update final state.
static void add_decays(Node &node, double sqrts)
Add nodes for all decays possible from the given node and all of its children.
constexpr int p
Proton.
constexpr int64_t dprime
Deuteron-prime resonance.
T uniform(T min, T max)
Definition: random.h:88
std::string fill_left(const std::string &s, size_t width, char fill=' ')
Fill string with characters to the left until the given width is reached.
Definition: action.h:24
bool parametrization_exists(const PdgCode &pdg_a, const PdgCode &pdg_b)
Checks if supplied codes have existing parametrizations of total cross sections.
static void deduplicate(std::vector< FinalStateCrossSection > &final_state_xs)
Deduplicate the final-state cross sections by summing.
@ TwoToTwo
See here for a short description.
@ Elastic
See here for a short description.
@ StringHard
See here for a short description.
bool all_of(Container &&c, UnaryPredicate &&p)
Convenience wrapper for std::all_of that operates on a complete container.
Definition: algorithms.h:80
static StringTransitionParameters create_string_transition_parameters(Configuration &config)
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:58
constexpr double pion_mass
Pion mass in GeV.
Definition: constants.h:65
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:66
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
void isoclean(std::string &s)
Remove ⁺, ⁻, ⁰ from string.
static constexpr int LFindScatter
bool is_string_soft_process(ProcessType p)
Check if a given process type is a soft string excitation.
double s_from_plab(double plab, double m_P, double m_T)
Convert p_lab to Mandelstam-s for a fixed-target setup, with a projectile of mass m_P and momentum pl...
Definition: kinematics.h:265
constexpr double fm2_mb
mb <-> fm^2 conversion factor.
Definition: constants.h:28
Helper structure for Experiment.
std::optional< bool > use_monash_tune_default
Bool for the default usage of the monash tune in the collider modus.
const ReactionsBitSet included_2to2
This indicates which two to two reactions are switched off.
Represent a final-state cross section.
FinalStateCrossSection(const std::string &name, double cross_section, double mass)
Construct a final-state cross section.
std::string name_
Name of the final state.
double cross_section_
Corresponding cross section in mb.
double mass_
Total mass of final state particles.
int32_t collisions_per_particle
Collision counter per particle, zero only for initially present particles.
Definition: particledata.h:32
A container to keep track of all ever existed input keys.
Definition: input_keys.h:1067
static const Key< double > collTerm_stringParam_probabilityPToDUU
See user guide description for more information.
Definition: input_keys.h:2849
static const Key< double > collTerm_stringTrans_lower
See user guide description for more information.
Definition: input_keys.h:2644
static const Key< std::pair< double, double > > collTerm_stringTrans_rangeNpi
See user guide description for more information.
Definition: input_keys.h:2679
static const Key< double > collTerm_stringTrans_range_width
See user guide description for more information.
Definition: input_keys.h:2695
static const Key< bool > collTerm_stringParam_useMonashTune
See user guide description for more information.
Definition: input_keys.h:3002
static const Key< double > collTerm_stringParam_sigmaPerp
See user guide description for more information.
Definition: input_keys.h:2884
static const Key< double > collTerm_stringParam_gluonPMin
See user guide description for more information.
Definition: input_keys.h:2760
static const Key< double > collTerm_stringParam_stringZALeading
See user guide description for more information.
Definition: input_keys.h:2957
static const Key< bool > collTerm_stringParam_mDependentFormationTimes
See user guide description for more information.
Definition: input_keys.h:2774
static const Key< double > collTerm_stringParam_diquarkSuppression
See user guide description for more information.
Definition: input_keys.h:2708
static const Key< double > collTerm_stringParam_formTimeFactor
See user guide description for more information.
Definition: input_keys.h:2721
static const Key< double > collTerm_stringParam_quarkAlpha
See user guide description for more information.
Definition: input_keys.h:2789
static const Key< double > collTerm_stringParam_stringZA
See user guide description for more information.
Definition: input_keys.h:2943
static const Key< double > collTerm_stringParam_strangeSuppression
See user guide description for more information.
Definition: input_keys.h:2902
static const Key< double > collTerm_stringParam_stringSigmaT
See user guide description for more information.
Definition: input_keys.h:2915
static const Key< double > collTerm_stringParam_stringZB
See user guide description for more information.
Definition: input_keys.h:2970
static const Key< double > collTerm_stringParam_quarkBeta
See user guide description for more information.
Definition: input_keys.h:2802
static const Key< double > collTerm_stringParam_popcornRate
See user guide description for more information.
Definition: input_keys.h:2817
static const Key< double > collTerm_stringParam_stringZBLeading
See user guide description for more information.
Definition: input_keys.h:2985
static const Key< double > collTerm_stringTrans_pipiOffset
See user guide description for more information.
Definition: input_keys.h:2631
static const Key< double > collTerm_stringParam_stringTension
See user guide description for more information.
Definition: input_keys.h:2930
static const Key< double > collTerm_stringTrans_KNOffset
See user guide description for more information.
Definition: input_keys.h:2616
static const Key< std::pair< double, double > > collTerm_stringTrans_rangeNN
See user guide description for more information.
Definition: input_keys.h:2662
static const Key< double > collTerm_stringParam_gluonBeta
See user guide description for more information.
Definition: input_keys.h:2746
static const Key< bool > collTerm_stringParam_separateFragmentBaryon
See user guide description for more information.
Definition: input_keys.h:2864
Constants related to transition between low and high collision energies.
Node of a decay tree, representing a possible action (2-to-2 or 1-to-2).
std::vector< Node > children_
Possible actions after this action.
void final_state_cross_sections_helper(uint64_t depth, std::vector< FinalStateCrossSection > &result, const std::string &name, double weight, bool show_intermediate_states=false) const
Internal helper function for final_state_cross_sections, to be called recursively to calculate all fi...
ParticleTypePtrList final_particles_
Final-state particle types in this action.
std::vector< FinalStateCrossSection > final_state_cross_sections() const
Node & add_action(const std::string &name, double weight, ParticleTypePtrList &&initial_particles, ParticleTypePtrList &&final_particles)
Add an action to the children of this node.
ParticleTypePtrList initial_particles_
Initial-state particle types in this action.
void print() const
Print the decay tree starting with this node.
void print_helper(uint64_t depth) const
Internal helper function for print, to be called recursively to print all nodes.
double weight_
Weight (cross section or branching ratio).
Node(const Node &)=delete
Cannot be copied.
Node(const std::string &name, double weight, ParticleTypePtrList &&initial_particles, ParticleTypePtrList &&final_particles, ParticleTypePtrList &&state, std::vector< Node > &&children)
Node(Node &&)=default
Move constructor.
ParticleTypePtrList state_
Particle types corresponding to the global state after this action.
std::string name_
Name for printing.