Version: SMASH-3.0
scatteractionsfinder.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2023
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/scatteraction.h"
22 #include "smash/stringfunctions.h"
23 
24 namespace smash {
25 static constexpr int LFindScatter = LogArea::FindScatter::id;
34  Configuration& config, const ExperimentParameters& parameters)
35  : finder_parameters_(create_finder_parameters(config, parameters)),
36  isotropic_(config.take({"Collision_Term", "Isotropic"}, false)),
37  box_length_(parameters.box_length),
38  string_formation_time_(config.take(
39  {"Collision_Term", "String_Parameters", "Formation_Time"}, 1.)) {
40  if (is_constant_elastic_isotropic()) {
41  logg[LFindScatter].info(
42  "Constant elastic isotropic cross-section mode:", " using ",
43  finder_parameters_.elastic_parameter, " mb as maximal cross-section.");
44  }
45  if (finder_parameters_.included_multi.any() &&
46  finder_parameters_.coll_crit != CollisionCriterion::Stochastic) {
47  throw std::invalid_argument(
48  "Multi-body reactions (like e.g. 3->1 or 3->2) are only possible with "
49  "the stochastic "
50  "collision "
51  "criterion. Change your config accordingly.");
52  }
53 
54  if (finder_parameters_
56  1 &&
57  (finder_parameters_
58  .included_2to2[IncludedReactions::PiDeuteron_to_pidprime] == 1 ||
59  finder_parameters_
60  .included_2to2[IncludedReactions::NDeuteron_to_Ndprime] == 1)) {
61  throw std::invalid_argument(
62  "To prevent double counting it is not possible to enable deuteron 3->2 "
63  "reactions\nand reactions involving the d' at the same time\ni.e. to "
64  "include \"Deuteron_3to2\" in `Multi_Particle_Reactions` and\n "
65  "\"PiDeuteron_to_pidprime\" "
66  "or \"NDeuteron_to_Ndprime\" in `Included_2to2` at the same time.\n"
67  "Change your config accordingly.");
68  }
69 
70  if (finder_parameters_
72  1 &&
74  throw std::invalid_argument(
75  "Do not use the d' resonance and enable \"Deuteron_3to2\" "
76  "`Multi_Particle_Reactions` at the same time. Either use the direct "
77  "3-to-2 reactions or the d' together with \"PiDeuteron_to_pidprime\" "
78  "and \"NDeuteron_to_Ndprime\" in `Included_2to2`. Otherwise the "
79  "deuteron 3-to-2 reactions would be double counted.");
80  }
81 
82  if ((finder_parameters_.nnbar_treatment == NNbarTreatment::TwoToFive &&
83  finder_parameters_
85  1) ||
86  (finder_parameters_
88  1 &&
89  finder_parameters_.nnbar_treatment != NNbarTreatment::TwoToFive)) {
90  throw std::invalid_argument(
91  "In order to conserve detailed balance, when \"NNbar_5to2\" is "
92  "included in\n`Multi_Particle_Reactions`, the `NNbarTreatment` has to "
93  "be set to \"two to five\" and vice versa.");
94  }
95 
96  if (finder_parameters_.nnbar_treatment == NNbarTreatment::Resonances &&
97  finder_parameters_.included_2to2[IncludedReactions::NNbar] != 1) {
98  throw std::invalid_argument(
99  "'NNbar' has to be in the list of allowed 2 to 2 processes "
100  "to enable annihilation to go through resonances");
101  }
102 
103  if (finder_parameters_.strings_switch) {
104  auto subconfig = config.extract_sub_configuration(
105  {"Collision_Term", "String_Parameters"}, Configuration::GetEmpty::Yes);
106  string_process_interface_ = std::make_unique<StringProcess>(
107  subconfig.take({"String_Tension"}, 1.0), string_formation_time_,
108  subconfig.take({"Gluon_Beta"}, 0.5),
109  subconfig.take({"Gluon_Pmin"}, 0.001),
110  subconfig.take({"Quark_Alpha"}, 2.0),
111  subconfig.take({"Quark_Beta"}, 7.0),
112  subconfig.take({"Strange_Supp"}, 0.16),
113  subconfig.take({"Diquark_Supp"}, 0.036),
114  subconfig.take({"Sigma_Perp"}, 0.42),
115  subconfig.take({"StringZ_A_Leading"}, 0.2),
116  subconfig.take({"StringZ_B_Leading"}, 2.0),
117  subconfig.take({"StringZ_A"}, 2.0), subconfig.take({"StringZ_B"}, 0.55),
118  subconfig.take({"String_Sigma_T"}, 0.5),
119  subconfig.take({"Form_Time_Factor"}, 1.0),
120  subconfig.take({"Mass_Dependent_Formation_Times"}, false),
121  subconfig.take({"Prob_proton_to_d_uu"}, 1. / 3.),
122  subconfig.take({"Separate_Fragment_Baryon"}, true),
123  subconfig.take({"Popcorn_Rate"}, 0.15),
124  subconfig.take({"Use_Monash_Tune"},
125  parameters.use_monash_tune_default.value()));
126  }
127 }
128 
130  Configuration& config, const ExperimentParameters& parameters) {
131  std::pair<double, double> sqrts_range_Npi =
132  config.take({"Collision_Term", "String_Transition", "Sqrts_Range_Npi"},
134  std::pair<double, double> sqrts_range_NN =
135  config.take({"Collision_Term", "String_Transition", "Sqrts_Range_NN"},
137  if (sqrts_range_Npi.first < nucleon_mass + pion_mass) {
138  sqrts_range_Npi.first = nucleon_mass + pion_mass;
139  if (sqrts_range_Npi.second < sqrts_range_Npi.first)
140  sqrts_range_Npi.second = sqrts_range_Npi.first;
141  logg[LFindScatter].warn(
142  "Lower bound of Sqrts_Range_Npi too small, setting it to mass "
143  "threshold. New range is [",
144  sqrts_range_Npi.first, ',', sqrts_range_Npi.second, "] GeV");
145  }
146  if (sqrts_range_NN.first < 2 * nucleon_mass) {
147  sqrts_range_NN.first = 2 * nucleon_mass;
148  if (sqrts_range_NN.second < sqrts_range_NN.first)
149  sqrts_range_NN.second = sqrts_range_NN.first;
150  logg[LFindScatter].warn(
151  "Lower bound of Sqrts_Range_NN too small, setting it to mass "
152  "threshold. New range is [",
153  sqrts_range_NN.first, ',', sqrts_range_NN.second, "] GeV.");
154  }
155  return {
156  config.take({"Collision_Term", "Elastic_Cross_Section"}, -1.),
157  parameters.low_snn_cut,
158  parameters.scale_xs,
159  config.take({"Collision_Term", "Additional_Elastic_Cross_Section"}, 0.0),
160  parameters.maximum_cross_section,
161  parameters.coll_crit,
162  parameters.nnbar_treatment,
163  parameters.included_2to2,
164  parameters.included_multi,
165  parameters.testparticles,
166  parameters.two_to_one,
167  config.take({"Modi", "Collider", "Collisions_Within_Nucleus"}, false),
168  parameters.strings_switch,
169  config.take({"Collision_Term", "Use_AQM"}, true),
170  config.take({"Collision_Term", "Strings_with_Probability"}, true),
171  config.take({"Collision_Term", "Only_Warn_For_High_Probability"}, false),
173  sqrts_range_Npi, sqrts_range_NN,
174  config.take({"Collision_Term", "String_Transition", "Sqrts_Lower"},
176  config.take(
177  {"Collision_Term", "String_Transition", "Sqrts_Range_Width"},
179  config.take(
180  {"Collision_Term", "String_Transition", "PiPi_Offset"},
182  config.take(
183  {"Collision_Term", "String_Transition", "KN_Offset"},
185 }
186 
188  const ParticleData& data_a, const ParticleData& data_b, double dt,
189  const std::vector<FourVector>& beam_momentum,
190  const double gcell_vol) const {
191  /* If the two particles
192  * 1) belong to one of the two colliding nuclei, and
193  * 2) both of them have never experienced any collisions,
194  * then the collisions between them are banned. */
196  assert(data_a.id() >= 0);
197  assert(data_b.id() >= 0);
198  bool in_same_nucleus = (data_a.belongs_to() == BelongsTo::Projectile &&
199  data_b.belongs_to() == BelongsTo::Projectile) ||
200  (data_a.belongs_to() == BelongsTo::Target &&
201  data_b.belongs_to() == BelongsTo::Target);
202  bool never_interacted_before =
203  data_a.get_history().collisions_per_particle == 0 &&
204  data_b.get_history().collisions_per_particle == 0;
205  if (in_same_nucleus && never_interacted_before) {
206  return nullptr;
207  }
208  }
209 
210  // No grid or search in cell means no collision for stochastic criterion
212  gcell_vol < really_small) {
213  return nullptr;
214  }
215 
216  // Determine time of collision.
217  const double time_until_collision =
218  collision_time(data_a, data_b, dt, beam_momentum);
219 
220  // Check that collision happens in this timestep.
221  if (time_until_collision < 0. || time_until_collision >= dt) {
222  return nullptr;
223  }
224 
225  // Create ScatterAction object.
226  ScatterActionPtr act = std::make_unique<ScatterAction>(
227  data_a, data_b, time_until_collision, isotropic_, string_formation_time_,
228  box_length_);
229 
231  act->set_stochastic_pos_idx();
232  }
233 
235  act->set_string_interface(string_process_interface_.get());
236  }
237 
238  // Distance squared calculation not needed for stochastic criterion
239  const double distance_squared =
241  ? act->transverse_distance_sqr()
243  ? act->cov_transverse_distance_sqr()
244  : 0.0;
245 
246  // Don't calculate cross section if the particles are very far apart.
247  // Not needed for stochastic criterion because of cell structure.
249  distance_squared >=
251  return nullptr;
252  }
253 
254  // Add various subprocesses.
255  act->add_all_scatterings(finder_parameters_);
256  double xs = act->cross_section() * fm2_mb /
257  static_cast<double>(finder_parameters_.testparticles);
258 
259  // Take cross section scaling factors into account
260  xs *= data_a.xsec_scaling_factor(time_until_collision);
261  xs *= data_b.xsec_scaling_factor(time_until_collision);
262 
264  const double v_rel = act->relative_velocity();
265  /* Collision probability for 2-particle scattering, see
266  * \iref{Staudenmaier:2021lrg}. */
267  const double prob = xs * v_rel * dt / gcell_vol;
268 
269  logg[LFindScatter].debug(
270  "Stochastic collison criterion parameters (2-particles):\nprob = ",
271  prob, ", xs = ", xs, ", v_rel = ", v_rel, ", dt = ", dt,
272  ", gcell_vol = ", gcell_vol,
273  ", testparticles = ", finder_parameters_.testparticles);
274 
275  if (prob > 1.) {
276  std::stringstream err;
277  err << "Probability larger than 1 for stochastic rates. ( P_22 = " << prob
278  << " )\n"
279  << data_a.type().name() << data_b.type().name() << " with masses "
280  << data_a.momentum().abs() << " and " << data_b.momentum().abs()
281  << " at sqrts[GeV] = " << act->sqrt_s()
282  << " with xs[fm^2]/Ntest = " << xs
283  << "\nConsider using smaller timesteps.";
285  logg[LFindScatter].warn(err.str());
286  } else {
287  throw std::runtime_error(err.str());
288  }
289  }
290 
291  // probability criterion
292  double random_no = random::uniform(0., 1.);
293  if (random_no > prob) {
294  return nullptr;
295  }
296 
299  // just collided with this particle
300  if (data_a.id_process() > 0 && data_a.id_process() == data_b.id_process()) {
301  logg[LFindScatter].debug("Skipping collided particles at time ",
302  data_a.position().x0(), " due to process ",
303  data_a.id_process(), "\n ", data_a, "\n<-> ",
304  data_b);
305 
306  return nullptr;
307  }
308 
309  // Cross section for collision criterion
310  const double cross_section_criterion = xs * M_1_PI;
311 
312  // distance criterion according to cross_section
313  if (distance_squared >= cross_section_criterion) {
314  return nullptr;
315  }
316 
317  logg[LFindScatter].debug("particle distance squared: ", distance_squared,
318  "\n ", data_a, "\n<-> ", data_b);
319  }
320 
321  return act;
322 }
323 
325  const ParticleList& plist, double dt, const double gcell_vol) const {
326  /* If all particles
327  * 1) belong to the two colliding nuclei
328  * 2) are within the same nucleus
329  * 3) have never experienced any collisons,
330  * then the collision between them are banned also for multi-particle
331  * interactions. */
333  bool all_projectile =
334  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
335  return data.belongs_to() == BelongsTo::Projectile;
336  });
337  bool all_target =
338  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
339  return data.belongs_to() == BelongsTo::Target;
340  });
341  bool none_collided =
342  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
343  return data.get_history().collisions_per_particle == 0;
344  });
345  if ((all_projectile || all_target) && none_collided) {
346  return nullptr;
347  }
348  }
349  // No grid or search in cell
350  if (gcell_vol < really_small) {
351  return nullptr;
352  }
353 
354  /* Optimisation for later: Already check here at the beginning
355  * if collision with plist is possible before constructing actions. */
356 
357  // 1. Determine time of collision.
358  const double time_until_collision = dt * random::uniform(0., 1.);
359 
360  // 2. Create ScatterAction object.
361  ScatterActionMultiPtr act =
362  std::make_unique<ScatterActionMulti>(plist, time_until_collision);
363 
364  act->set_stochastic_pos_idx();
365 
366  // 3. Add possible final states (dt and gcell_vol for probability calculation)
367  act->add_possible_reactions(dt, gcell_vol, finder_parameters_.included_multi);
368 
369  /* 4. Return total collision probability
370  * Scales with 1 over the number of testpartciles to the power of the
371  * number of incoming particles - 1 */
372  const double prob =
373  act->get_total_weight() /
374  std::pow(finder_parameters_.testparticles, plist.size() - 1);
375 
376  // 5. Check that probability is smaller than one
377  if (prob > 1.) {
378  std::stringstream err;
379  err << "Probability " << prob << " larger than 1 for stochastic rates for ";
380  for (const ParticleData& data : plist) {
381  err << data.type().name();
382  }
383  err << " at sqrts[GeV] = " << act->sqrt_s()
384  << "\nConsider using smaller timesteps.";
386  logg[LFindScatter].warn(err.str());
387  } else {
388  throw std::runtime_error(err.str());
389  }
390  }
391 
392  // 6. Perform probability decisions
393  double random_no = random::uniform(0., 1.);
394  if (random_no > prob) {
395  return nullptr;
396  }
397 
398  return act;
399 }
400 
402  const ParticleList& search_list, double dt, const double gcell_vol,
403  const std::vector<FourVector>& beam_momentum) const {
404  std::vector<ActionPtr> actions;
405  for (const ParticleData& p1 : search_list) {
406  for (const ParticleData& p2 : search_list) {
407  // Check for 2 particle scattering
408  if (p1.id() < p2.id()) {
409  ActionPtr act =
410  check_collision_two_part(p1, p2, dt, beam_momentum, gcell_vol);
411  if (act) {
412  actions.push_back(std::move(act));
413  }
414  }
416  // Also, check for 3 particle scatterings with stochastic criterion
417  for (const ParticleData& p3 : search_list) {
422  if (p1.id() < p2.id() && p2.id() < p3.id()) {
423  ActionPtr act =
424  check_collision_multi_part({p1, p2, p3}, dt, gcell_vol);
425  if (act) {
426  actions.push_back(std::move(act));
427  }
428  }
429  }
430  for (const ParticleData& p4 : search_list) {
433  if (p1.id() < p2.id() && p2.id() < p3.id() && p3.id() < p4.id()) {
434  ActionPtr act =
435  check_collision_multi_part({p1, p2, p3, p4}, dt, gcell_vol);
436  if (act) {
437  actions.push_back(std::move(act));
438  }
439  }
440  }
443  search_list.size() >= 5) {
444  for (const ParticleData& p5 : search_list) {
445  if ((p1.id() < p2.id() && p2.id() < p3.id() &&
446  p3.id() < p4.id() && p4.id() < p5.id()) &&
447  (p1.is_pion() && p2.is_pion() && p3.is_pion() &&
448  p4.is_pion() && p5.is_pion())) {
449  // at the moment only pure pion 5-body reactions
450  ActionPtr act = check_collision_multi_part(
451  {p1, p2, p3, p4, p5}, dt, gcell_vol);
452  if (act) {
453  actions.push_back(std::move(act));
454  }
455  }
456  }
457  }
458  }
459  }
460  }
461  }
462  }
463  return actions;
464 }
465 
467  const ParticleList& search_list, const ParticleList& neighbors_list,
468  double dt, const std::vector<FourVector>& beam_momentum) const {
469  std::vector<ActionPtr> actions;
471  // Only search in cells
472  return actions;
473  }
474  for (const ParticleData& p1 : search_list) {
475  for (const ParticleData& p2 : neighbors_list) {
476  assert(p1.id() != p2.id());
477  // Check if a collision is possible.
478  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
479  if (act) {
480  actions.push_back(std::move(act));
481  }
482  }
483  }
484  return actions;
485 }
486 
488  const ParticleList& search_list, const Particles& surrounding_list,
489  double dt, const std::vector<FourVector>& beam_momentum) const {
490  std::vector<ActionPtr> actions;
492  // Only search in cells
493  return actions;
494  }
495  for (const ParticleData& p2 : surrounding_list) {
496  /* don't look for collisions if the particle from the surrounding list is
497  * also in the search list */
498  auto result = std::find_if(
499  search_list.begin(), search_list.end(),
500  [&p2](const ParticleData& p) { return p.id() == p2.id(); });
501  if (result != search_list.end()) {
502  continue;
503  }
504  for (const ParticleData& p1 : search_list) {
505  // Check if a collision is possible.
506  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
507  if (act) {
508  actions.push_back(std::move(act));
509  }
510  }
511  }
512  return actions;
513 }
514 
516  constexpr double time = 0.0;
517 
518  const size_t N_isotypes = IsoParticleType::list_all().size();
519  const size_t N_pairs = N_isotypes * (N_isotypes - 1) / 2;
520 
521  std::cout << N_isotypes << " iso-particle types." << std::endl;
522  std::cout << "They can make " << N_pairs << " pairs." << std::endl;
523  std::vector<double> momentum_scan_list = {0.1, 0.3, 0.5, 1.0,
524  2.0, 3.0, 5.0, 10.0};
525  for (const IsoParticleType& A_isotype : IsoParticleType::list_all()) {
526  for (const IsoParticleType& B_isotype : IsoParticleType::list_all()) {
527  if (&A_isotype > &B_isotype) {
528  continue;
529  }
530  bool any_nonzero_cs = false;
531  std::vector<std::string> r_list;
532  for (const ParticleTypePtr A_type : A_isotype.get_states()) {
533  for (const ParticleTypePtr B_type : B_isotype.get_states()) {
534  if (A_type > B_type) {
535  continue;
536  }
537  ParticleData A(*A_type), B(*B_type);
538  for (auto mom : momentum_scan_list) {
539  A.set_4momentum(A.pole_mass(), mom, 0.0, 0.0);
540  B.set_4momentum(B.pole_mass(), -mom, 0.0, 0.0);
541  ScatterActionPtr act = std::make_unique<ScatterAction>(
542  A, B, time, isotropic_, string_formation_time_);
544  act->set_string_interface(string_process_interface_.get());
545  }
546  act->add_all_scatterings(finder_parameters_);
547  const double total_cs = act->cross_section();
548  if (total_cs <= 0.0) {
549  continue;
550  }
551  any_nonzero_cs = true;
552  for (const auto& channel : act->collision_channels()) {
553  const auto type = channel->get_type();
554  std::string r;
555  if (is_string_soft_process(type) ||
556  type == ProcessType::StringHard) {
557  r = A_type->name() + B_type->name() + std::string(" → strings");
558  } else {
559  std::string r_type =
560  (type == ProcessType::Elastic) ? std::string(" (el)")
561  : (channel->get_type() == ProcessType::TwoToTwo)
562  ? std::string(" (inel)")
563  : std::string(" (?)");
564  r = A_type->name() + B_type->name() + std::string(" → ") +
565  channel->particle_types()[0]->name() +
566  channel->particle_types()[1]->name() + r_type;
567  }
568  isoclean(r);
569  r_list.push_back(r);
570  }
571  }
572  }
573  }
574  std::sort(r_list.begin(), r_list.end());
575  r_list.erase(std::unique(r_list.begin(), r_list.end()), r_list.end());
576  if (any_nonzero_cs) {
577  for (auto r : r_list) {
578  std::cout << r;
579  if (r_list.back() != r) {
580  std::cout << ", ";
581  }
582  }
583  std::cout << std::endl;
584  }
585  }
586  }
587 }
588 
592  std::string name_;
593 
596 
598  double mass_;
599 
608  FinalStateCrossSection(const std::string& name, double cross_section,
609  double mass)
610  : name_(name), cross_section_(cross_section), mass_(mass) {}
611 };
612 
613 namespace decaytree {
614 
626 struct Node {
627  public:
629  std::string name_;
630 
632  double weight_;
633 
635  ParticleTypePtrList initial_particles_;
636 
638  ParticleTypePtrList final_particles_;
639 
641  ParticleTypePtrList state_;
642 
644  std::vector<Node> children_;
645 
647  Node(const Node&) = delete;
649  Node(Node&&) = default;
650 
661  Node(const std::string& name, double weight,
662  ParticleTypePtrList&& initial_particles,
663  ParticleTypePtrList&& final_particles, ParticleTypePtrList&& state,
664  std::vector<Node>&& children)
665  : name_(name),
666  weight_(weight),
667  initial_particles_(std::move(initial_particles)),
668  final_particles_(std::move(final_particles)),
669  state_(std::move(state)),
670  children_(std::move(children)) {}
671 
683  Node& add_action(const std::string& name, double weight,
684  ParticleTypePtrList&& initial_particles,
685  ParticleTypePtrList&& final_particles) {
686  // Copy parent state and update it.
687  ParticleTypePtrList state(state_);
688  for (const auto& p : initial_particles) {
689  state.erase(std::find(state.begin(), state.end(), p));
690  }
691  for (const auto& p : final_particles) {
692  state.push_back(p);
693  }
694  // Sort the state to normalize the output.
695  std::sort(state.begin(), state.end(),
697  return a->name() < b->name();
698  });
699  // Push new node to children.
700  Node new_node(name, weight, std::move(initial_particles),
701  std::move(final_particles), std::move(state), {});
702  children_.emplace_back(std::move(new_node));
703  return children_.back();
704  }
705 
707  void print() const { print_helper(0); }
708 
712  std::vector<FinalStateCrossSection> final_state_cross_sections() const {
713  std::vector<FinalStateCrossSection> result;
714  final_state_cross_sections_helper(0, result, "", 1.);
715  return result;
716  }
717 
718  private:
725  void print_helper(uint64_t depth) const {
726  for (uint64_t i = 0; i < depth; i++) {
727  std::cout << " ";
728  }
729  std::cout << name_ << " " << weight_ << std::endl;
730  for (const auto& child : children_) {
731  child.print_helper(depth + 1);
732  }
733  }
734 
747  uint64_t depth, std::vector<FinalStateCrossSection>& result,
748  const std::string& name, double weight,
749  bool show_intermediate_states = false) const {
750  // The first node corresponds to the total cross section and has to be
751  // ignored. The second node corresponds to the partial cross section. All
752  // further nodes correspond to branching ratios.
753  if (depth > 0) {
754  weight *= weight_;
755  }
756 
757  std::string new_name;
758  double mass = 0.;
759 
760  if (show_intermediate_states) {
761  new_name = name;
762  if (!new_name.empty()) {
763  new_name += "->";
764  }
765  new_name += name_;
766  new_name += "{";
767  } else {
768  new_name = "";
769  }
770  for (const auto& s : state_) {
771  new_name += s->name();
772  mass += s->mass();
773  }
774  if (show_intermediate_states) {
775  new_name += "}";
776  }
777 
778  if (children_.empty()) {
779  result.emplace_back(FinalStateCrossSection(new_name, weight, mass));
780  return;
781  }
782  for (const auto& child : children_) {
783  child.final_state_cross_sections_helper(depth + 1, result, new_name,
784  weight, show_intermediate_states);
785  }
786  }
787 };
788 
797 static std::string make_decay_name(const std::string& res_name,
798  const DecayBranchPtr& decay,
799  ParticleTypePtrList& final_state) {
800  std::stringstream name;
801  name << "[" << res_name << "->";
802  for (const auto& p : decay->particle_types()) {
803  name << p->name();
804  final_state.push_back(p);
805  }
806  name << "]";
807  return name.str();
808 }
809 
817 static void add_decays(Node& node, double sqrts) {
818  // If there is more than one unstable particle in the current state, then
819  // there will be redundant paths in the decay tree, corresponding to
820  // reorderings of the decays. To avoid double counting, we normalize by the
821  // number of possible decay orderings. Normalizing by the number of unstable
822  // particles recursively corresponds to normalizing by the factorial that
823  // gives the number of reorderings.
824  //
825  // Ideally, the redundant paths should never be added to the decay tree, but
826  // we never have more than two redundant paths, so it probably does not
827  // matter much.
828  uint32_t n_unstable = 0;
829  double sqrts_minus_masses = sqrts;
830  for (const ParticleTypePtr ptype : node.state_) {
831  if (!ptype->is_stable()) {
832  n_unstable += 1;
833  }
834  sqrts_minus_masses -= ptype->mass();
835  }
836  const double norm =
837  n_unstable != 0 ? 1. / static_cast<double>(n_unstable) : 1.;
838 
839  for (const ParticleTypePtr ptype : node.state_) {
840  if (!ptype->is_stable()) {
841  const double sqrts_decay = sqrts_minus_masses + ptype->mass();
842  bool can_decay = false;
843  for (const auto& decay : ptype->decay_modes().decay_mode_list()) {
844  // Make sure to skip kinematically impossible decays.
845  // In principle, we would have to integrate over the mass of the
846  // resonance, but as an approximation we just assume it at its pole.
847  double final_state_mass = 0.;
848  for (const auto& p : decay->particle_types()) {
849  final_state_mass += p->mass();
850  }
851  if (final_state_mass > sqrts_decay) {
852  continue;
853  }
854  can_decay = true;
855 
856  ParticleTypePtrList parts;
857  const auto name = make_decay_name(ptype->name(), decay, parts);
858  auto& new_node = node.add_action(name, norm * decay->weight(), {ptype},
859  std::move(parts));
860  add_decays(new_node, sqrts_decay);
861  }
862  if (!can_decay) {
863  // Remove final-state cross sections with resonances that cannot
864  // decay due to our "mass = pole mass" approximation.
865  node.weight_ = 0;
866  return;
867  }
868  }
869  }
870 }
871 
872 } // namespace decaytree
873 
879 static void deduplicate(std::vector<FinalStateCrossSection>& final_state_xs) {
880  std::sort(final_state_xs.begin(), final_state_xs.end(),
881  [](const FinalStateCrossSection& a,
882  const FinalStateCrossSection& b) { return a.name_ < b.name_; });
883  auto current = final_state_xs.begin();
884  while (current != final_state_xs.end()) {
885  auto adjacent = std::adjacent_find(
886  current, final_state_xs.end(),
887  [](const FinalStateCrossSection& a, const FinalStateCrossSection& b) {
888  return a.name_ == b.name_;
889  });
890  current = adjacent;
891  if (adjacent != final_state_xs.end()) {
892  adjacent->cross_section_ += (adjacent + 1)->cross_section_;
893  final_state_xs.erase(adjacent + 1);
894  }
895  }
896 }
897 
899  const ParticleType& a, const ParticleType& b, double m_a, double m_b,
900  bool final_state, std::vector<double>& plab) const {
901  typedef std::vector<std::pair<double, double>> xs_saver;
902  std::map<std::string, xs_saver> xs_dump;
903  std::map<std::string, double> outgoing_total_mass;
904 
905  ParticleData a_data(a), b_data(b);
906  int n_momentum_points = 200;
907  constexpr double momentum_step = 0.02;
908  if (plab.size() > 0) {
909  n_momentum_points = plab.size();
910  // Remove duplicates.
911  std::sort(plab.begin(), plab.end());
912  plab.erase(std::unique(plab.begin(), plab.end()), plab.end());
913  }
914  for (int i = 0; i < n_momentum_points; i++) {
915  double momentum;
916  if (plab.size() > 0) {
917  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
918  } else {
919  momentum = momentum_step * (i + 1);
920  }
921  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
922  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
923  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
924  const ParticleList incoming = {a_data, b_data};
925  ScatterActionPtr act = std::make_unique<ScatterAction>(
926  a_data, b_data, 0., isotropic_, string_formation_time_);
928  act->set_string_interface(string_process_interface_.get());
929  }
930  act->add_all_scatterings(finder_parameters_);
931  decaytree::Node tree(a.name() + b.name(), act->cross_section(), {&a, &b},
932  {&a, &b}, {&a, &b}, {});
933  const CollisionBranchList& processes = act->collision_channels();
934  for (const auto& process : processes) {
935  const double xs = process->weight();
936  if (xs <= 0.0) {
937  continue;
938  }
939  if (!final_state) {
940  std::stringstream process_description_stream;
941  process_description_stream << *process;
942  const std::string& description = process_description_stream.str();
943  double m_tot = 0.0;
944  for (const auto& ptype : process->particle_types()) {
945  m_tot += ptype->mass();
946  }
947  outgoing_total_mass[description] = m_tot;
948  if (!xs_dump[description].empty() &&
949  std::abs(xs_dump[description].back().first - sqrts) <
950  really_small) {
951  xs_dump[description].back().second += xs;
952  } else {
953  xs_dump[description].push_back(std::make_pair(sqrts, xs));
954  }
955  } else {
956  std::stringstream process_description_stream;
957  process_description_stream << *process;
958  const std::string& description = process_description_stream.str();
959  ParticleTypePtrList initial_particles = {&a, &b};
960  ParticleTypePtrList final_particles = process->particle_types();
961  auto& process_node =
962  tree.add_action(description, xs, std::move(initial_particles),
963  std::move(final_particles));
964  decaytree::add_decays(process_node, sqrts);
965  }
966  }
967  xs_dump["total"].push_back(std::make_pair(sqrts, act->cross_section()));
968  // Total cross-section should be the first in the list -> negative mass
969  outgoing_total_mass["total"] = -1.0;
970  if (final_state) {
971  // tree.print();
972  auto final_state_xs = tree.final_state_cross_sections();
973  deduplicate(final_state_xs);
974  for (const auto& p : final_state_xs) {
975  // Don't print empty columns.
976  //
977  // FIXME(steinberg): The better fix would be to not have them in the
978  // first place.
979  if (p.name_ == "") {
980  continue;
981  }
982  outgoing_total_mass[p.name_] = p.mass_;
983  xs_dump[p.name_].push_back(std::make_pair(sqrts, p.cross_section_));
984  }
985  }
986  }
987  // Get rid of cross sections that are zero.
988  // (This only happens if their is a resonance in the final state that cannot
989  // decay with our simplified assumptions.)
990  for (auto it = begin(xs_dump); it != end(xs_dump);) {
991  // Sum cross section over all energies.
992  const xs_saver& xs = (*it).second;
993  double sum = 0;
994  for (const auto& p : xs) {
995  sum += p.second;
996  }
997  if (sum == 0.) {
998  it = xs_dump.erase(it);
999  } else {
1000  ++it;
1001  }
1002  }
1003 
1004  // Nice ordering of channels by summed pole mass of products
1005  std::vector<std::string> all_channels;
1006  for (const auto& channel : xs_dump) {
1007  all_channels.push_back(channel.first);
1008  }
1009  std::sort(all_channels.begin(), all_channels.end(),
1010  [&](const std::string& str_a, const std::string& str_b) {
1011  return outgoing_total_mass[str_a] < outgoing_total_mass[str_b];
1012  });
1013 
1014  // Print header
1015  std::cout << "# Dumping partial " << a.name() << b.name()
1016  << " cross-sections in mb, energies in GeV" << std::endl;
1017  std::cout << " sqrt_s";
1018  // Align everything to 16 unicode characters.
1019  // This should be enough for the longest channel name (7 final-state
1020  // particles).
1021  for (const auto& channel : all_channels) {
1022  std::cout << utf8::fill_left(channel, 16, ' ');
1023  }
1024  std::cout << std::endl;
1025 
1026  // Print out all partial cross-sections in mb
1027  for (int i = 0; i < n_momentum_points; i++) {
1028  double momentum;
1029  if (plab.size() > 0) {
1030  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
1031  } else {
1032  momentum = momentum_step * (i + 1);
1033  }
1034  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
1035  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
1036  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
1037  std::printf("%9.6f", sqrts);
1038  for (const auto& channel : all_channels) {
1039  const xs_saver energy_and_xs = xs_dump[channel];
1040  size_t j = 0;
1041  for (; j < energy_and_xs.size() && energy_and_xs[j].first < sqrts; j++) {
1042  }
1043  double xs = 0.0;
1044  if (j < energy_and_xs.size() &&
1045  std::abs(energy_and_xs[j].first - sqrts) < really_small) {
1046  xs = energy_and_xs[j].second;
1047  }
1048  std::printf("%16.6f", xs); // Same alignment as in the header.
1049  }
1050  std::printf("\n");
1051  }
1052 }
1053 
1054 } // namespace smash
Interface to the SMASH configuration files.
Value take(std::initializer_list< const char * > keys)
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.
default_type default_value() const
Get the default value of the key.
Definition: input_keys.h:132
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:85
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:339
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:672
Particle type contains the static properties of a particle species.
Definition: particletype.h:97
const DecayModes & decay_modes() const
static const ParticleTypePtr try_find(PdgCode pdgcode)
Returns the ParticleTypePtr for the given pdgcode.
Definition: particletype.cc:89
const std::string & name() const
Definition: particletype.h:141
bool is_stable() const
Definition: particletype.h:242
double mass() const
Definition: particletype.h:144
The Particles class abstracts the storage and manipulation of particles.
Definition: particles.h:33
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...
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.
@ NNbar_5to2
@ A3_Nuclei_4to2
@ Deuteron_3to2
@ Meson_3to1
@ Stochastic
Stochastic Criteiron.
@ Geometric
(Default) 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:39
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
ScatterActionsFinderParameters create_finder_parameters(Configuration &config, const ExperimentParameters &parameters)
Gather all relevant parameters for a ScatterActionsFinder either getting them from an ExperimentParam...
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
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.
const CollisionCriterion coll_crit
Employed collision criterion.
NNbarTreatment nnbar_treatment
This indicates how NN̅ annihilation should be treated; options are to neglect it, make it conserve de...
double maximum_cross_section
The maximal cross section (in mb) for which it is guaranteed that all collisions with this cross sect...
bool strings_switch
This indicates whether string fragmentation is switched on.
const ReactionsBitSet included_2to2
This indicates which two to two reactions are switched off.
double scale_xs
Global factor which all cross sections are scaled with.
bool two_to_one
This indicates whether string fragmentation is switched on.
int testparticles
Number of test-particles.
double low_snn_cut
Elastic collisions between the nucleons with the square root s below low_snn_cut are excluded.
const MultiParticleReactionsBitSet included_multi
This indicates which multi-particle reactions are switched on.
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
static const Key< double > collTerm_stringTrans_lower
See user guide description for more information.
Definition: input_keys.h:2340
static const Key< std::pair< double, double > > collTerm_stringTrans_rangeNpi
See user guide description for more information.
Definition: input_keys.h:2375
static const Key< double > collTerm_stringTrans_range_width
See user guide description for more information.
Definition: input_keys.h:2391
static const Key< double > collTerm_stringTrans_pipiOffset
See user guide description for more information.
Definition: input_keys.h:2327
static const Key< double > collTerm_stringTrans_KNOffset
See user guide description for more information.
Definition: input_keys.h:2312
static const Key< std::pair< double, double > > collTerm_stringTrans_rangeNN
See user guide description for more information.
Definition: input_keys.h:2358
Helper structure for ScatterActionsFinder.
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.
const bool strings_switch
Indicates whether string fragmentation is switched on.
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 CollisionCriterion coll_crit
Specifies which collision criterion is used.
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.