Version: SMASH-2.2
scatteractionsfinder.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2022
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/cxx14compat.h"
18 #include "smash/decaymodes.h"
19 #include "smash/logging.h"
20 #include "smash/scatteraction.h"
23 #include "smash/stringfunctions.h"
24 
25 namespace smash {
26 static constexpr int LFindScatter = LogArea::FindScatter::id;
309  Configuration config, const ExperimentParameters& parameters)
310  : coll_crit_(parameters.coll_crit),
311  elastic_parameter_(
312  config.take({"Collision_Term", "Elastic_Cross_Section"}, -1.)),
313  testparticles_(parameters.testparticles),
314  isotropic_(config.take({"Collision_Term", "Isotropic"}, false)),
315  two_to_one_(parameters.two_to_one),
316  incl_set_(parameters.included_2to2),
317  incl_multi_set_(parameters.included_multi),
318  scale_xs_(parameters.scale_xs),
319  additional_el_xs_(parameters.additional_el_xs),
320  low_snn_cut_(parameters.low_snn_cut),
321  strings_switch_(parameters.strings_switch),
322  use_AQM_(parameters.use_AQM),
323  strings_with_probability_(parameters.strings_with_probability),
324  nnbar_treatment_(parameters.nnbar_treatment),
325  box_length_(parameters.box_length),
326  string_formation_time_(config.take(
327  {"Collision_Term", "String_Parameters", "Formation_Time"}, 1.)),
328  maximum_cross_section_(parameters.maximum_cross_section),
329  allow_first_collisions_within_nucleus_(
330  parameters.allow_collisions_within_nucleus),
331  only_warn_for_high_prob_(config.take(
332  {"Collision_Term", "Only_Warn_For_High_Probability"}, false)) {
333  if (is_constant_elastic_isotropic()) {
334  logg[LFindScatter].info(
335  "Constant elastic isotropic cross-section mode:", " using ",
336  elastic_parameter_, " mb as maximal cross-section.");
337  }
338  if (incl_multi_set_.any() && coll_crit_ != CollisionCriterion::Stochastic) {
339  throw std::invalid_argument(
340  "Multi-body reactions (like e.g. 3->1 or 3->2) are only possible with "
341  "the stochastic "
342  "collision "
343  "criterion. Change your config accordingly.");
344  }
345 
346  if (incl_multi_set_[IncludedMultiParticleReactions::Deuteron_3to2] == 1 &&
348  incl_set_[IncludedReactions::NDeuteron_to_Ndprime] == 1)) {
349  throw std::invalid_argument(
350  "To prevent double counting it is not possible to enable deuteron 3->2 "
351  "reactions\nand reactions involving the d' at the same time\ni.e. to "
352  "include \"Deuteron_3to2\" in `Multi_Particle_Reactions` and\n "
353  "\"PiDeuteron_to_pidprime\" "
354  "or \"NDeuteron_to_Ndprime\" in `Included_2to2` at the same time.\n"
355  "Change your config accordingly.");
356  }
357 
358  if (incl_multi_set_[IncludedMultiParticleReactions::Deuteron_3to2] == 1 &&
360  throw std::invalid_argument(
361  "Do not use the d' resonance and enable \"Deuteron_3to2\" "
362  "`Multi_Particle_Reactions` at the same time. Either use the direct "
363  "3-to-2 reactions or the d' together with \"PiDeuteron_to_pidprime\" "
364  "and \"NDeuteron_to_Ndprime\" in `Included_2to2`. Otherwise the "
365  "deuteron 3-to-2 reactions would be double counted.");
366  }
367 
368  if ((nnbar_treatment_ == NNbarTreatment::TwoToFive &&
369  incl_multi_set_[IncludedMultiParticleReactions::NNbar_5to2] != 1) ||
370  (incl_multi_set_[IncludedMultiParticleReactions::NNbar_5to2] == 1 &&
371  nnbar_treatment_ != NNbarTreatment::TwoToFive)) {
372  throw std::invalid_argument(
373  "In order to conserve detailed balance, when \"NNbar_5to2\" is "
374  "included in\n`Multi_Particle_Reactions`, the `NNbarTreatment` has to "
375  "be set to \"two to five\" and vice versa.");
376  }
377 
378  if (nnbar_treatment_ == NNbarTreatment::Resonances &&
379  incl_set_[IncludedReactions::NNbar] != 1) {
380  throw std::invalid_argument(
381  "'NNbar' has to be in the list of allowed 2 to 2 processes "
382  "to enable annihilation to go through resonances");
383  }
384 
385  if (strings_switch_) {
386  auto subconfig = config["Collision_Term"]["String_Parameters"];
387  string_process_interface_ = make_unique<StringProcess>(
388  subconfig.take({"String_Tension"}, 1.0), string_formation_time_,
389  subconfig.take({"Gluon_Beta"}, 0.5),
390  subconfig.take({"Gluon_Pmin"}, 0.001),
391  subconfig.take({"Quark_Alpha"}, 2.0),
392  subconfig.take({"Quark_Beta"}, 7.0),
393  subconfig.take({"Strange_Supp"}, 0.16),
394  subconfig.take({"Diquark_Supp"}, 0.036),
395  subconfig.take({"Sigma_Perp"}, 0.42),
396  subconfig.take({"StringZ_A_Leading"}, 0.2),
397  subconfig.take({"StringZ_B_Leading"}, 2.0),
398  subconfig.take({"StringZ_A"}, 2.0), subconfig.take({"StringZ_B"}, 0.55),
399  subconfig.take({"String_Sigma_T"}, 0.5),
400  subconfig.take({"Form_Time_Factor"}, 1.0),
401  subconfig.take({"Mass_Dependent_Formation_Times"}, false),
402  subconfig.take({"Prob_proton_to_d_uu"}, 1. / 3.),
403  subconfig.take({"Separate_Fragment_Baryon"}, true),
404  subconfig.take({"Popcorn_Rate"}, 0.15));
405  }
406 }
407 
409  const ParticleData& data_a, const ParticleData& data_b, double dt,
410  const std::vector<FourVector>& beam_momentum,
411  const double gcell_vol) const {
412  /* If the two particles
413  * 1) belong to one of the two colliding nuclei, and
414  * 2) both of them have never experienced any collisions,
415  * then the collisions between them are banned. */
417  assert(data_a.id() >= 0);
418  assert(data_b.id() >= 0);
419  bool in_same_nucleus = (data_a.belongs_to() == BelongsTo::Projectile &&
420  data_b.belongs_to() == BelongsTo::Projectile) ||
421  (data_a.belongs_to() == BelongsTo::Target &&
422  data_b.belongs_to() == BelongsTo::Target);
423  bool never_interacted_before =
424  data_a.get_history().collisions_per_particle == 0 &&
425  data_b.get_history().collisions_per_particle == 0;
426  if (in_same_nucleus && never_interacted_before) {
427  return nullptr;
428  }
429  }
430 
431  // No grid or search in cell means no collision for stochastic criterion
433  gcell_vol < really_small) {
434  return nullptr;
435  }
436 
437  // Determine time of collision.
438  const double time_until_collision =
439  collision_time(data_a, data_b, dt, beam_momentum);
440 
441  // Check that collision happens in this timestep.
442  if (time_until_collision < 0. || time_until_collision >= dt) {
443  return nullptr;
444  }
445 
446  // Create ScatterAction object.
447  ScatterActionPtr act = make_unique<ScatterAction>(
448  data_a, data_b, time_until_collision, isotropic_, string_formation_time_,
449  box_length_);
450 
452  act->set_stochastic_pos_idx();
453  }
454 
455  if (strings_switch_) {
456  act->set_string_interface(string_process_interface_.get());
457  }
458 
459  // Distance squared calculation not needed for stochastic criterion
460  const double distance_squared =
462  ? act->transverse_distance_sqr()
464  ? act->cov_transverse_distance_sqr()
465  : 0.0;
466 
467  // Don't calculate cross section if the particles are very far apart.
468  // Not needed for stochastic criterion because of cell structure.
470  distance_squared >= max_transverse_distance_sqr(testparticles_)) {
471  return nullptr;
472  }
473 
474  // Add various subprocesses.
475  act->add_all_scatterings(elastic_parameter_, two_to_one_, incl_set_,
479 
480  double xs =
481  act->cross_section() * fm2_mb / static_cast<double>(testparticles_);
482 
483  // Take cross section scaling factors into account
484  xs *= data_a.xsec_scaling_factor(time_until_collision);
485  xs *= data_b.xsec_scaling_factor(time_until_collision);
486 
488  const double v_rel = act->relative_velocity();
489  /* Collision probability for 2-particle scattering, see
490  * \iref{Staudenmaier:2021lrg}. */
491  const double prob = xs * v_rel * dt / gcell_vol;
492 
493  logg[LFindScatter].debug(
494  "Stochastic collison criterion parameters (2-particles):\nprob = ",
495  prob, ", xs = ", xs, ", v_rel = ", v_rel, ", dt = ", dt,
496  ", gcell_vol = ", gcell_vol, ", testparticles = ", testparticles_);
497 
498  if (prob > 1.) {
499  std::stringstream err;
500  err << "Probability larger than 1 for stochastic rates. ( P_22 = " << prob
501  << " )\nConsider using smaller timesteps.";
503  logg[LFindScatter].warn(err.str());
504  } else {
505  throw std::runtime_error(err.str());
506  }
507  }
508 
509  // probability criterion
510  double random_no = random::uniform(0., 1.);
511  if (random_no > prob) {
512  return nullptr;
513  }
514 
517  // just collided with this particle
518  if (data_a.id_process() > 0 && data_a.id_process() == data_b.id_process()) {
519  logg[LFindScatter].debug("Skipping collided particles at time ",
520  data_a.position().x0(), " due to process ",
521  data_a.id_process(), "\n ", data_a, "\n<-> ",
522  data_b);
523 
524  return nullptr;
525  }
526 
527  // Cross section for collision criterion
528  const double cross_section_criterion = xs * M_1_PI;
529 
530  // distance criterion according to cross_section
531  if (distance_squared >= cross_section_criterion) {
532  return nullptr;
533  }
534 
535  logg[LFindScatter].debug("particle distance squared: ", distance_squared,
536  "\n ", data_a, "\n<-> ", data_b);
537  }
538 
539  // Using std::move here is redundant with newer compilers, but required for
540  // supporting GCC 4.8. Once we drop this support, std::move should be removed.
541  return std::move(act);
542 }
543 
545  const ParticleList& plist, double dt, const double gcell_vol) const {
546  /* If all particles
547  * 1) belong to the two colliding nuclei
548  * 2) are within the same nucleus
549  * 3) have never experienced any collisons,
550  * then the collision between them are banned also for multi-particle
551  * interactions. */
553  bool all_projectile =
554  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
555  return data.belongs_to() == BelongsTo::Projectile;
556  });
557  bool all_target =
558  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
559  return data.belongs_to() == BelongsTo::Target;
560  });
561  bool none_collided =
562  std::all_of(plist.begin(), plist.end(), [&](const ParticleData& data) {
563  return data.get_history().collisions_per_particle == 0;
564  });
565  if ((all_projectile || all_target) && none_collided) {
566  return nullptr;
567  }
568  }
569  // No grid or search in cell
570  if (gcell_vol < really_small) {
571  return nullptr;
572  }
573 
574  /* Optimisation for later: Already check here at the beginning
575  * if collision with plist is possible before constructing actions. */
576 
577  // 1. Determine time of collision.
578  const double time_until_collision = dt * random::uniform(0., 1.);
579 
580  // 2. Create ScatterAction object.
581  ScatterActionMultiPtr act =
582  make_unique<ScatterActionMulti>(plist, time_until_collision);
583 
584  act->set_stochastic_pos_idx();
585 
586  // 3. Add possible final states (dt and gcell_vol for probability calculation)
587  act->add_possible_reactions(dt, gcell_vol, incl_multi_set_);
588 
589  /* 4. Return total collision probability
590  * Scales with 1 over the number of testpartciles to the power of the
591  * number of incoming particles - 1 */
592  const double prob =
593  act->get_total_weight() / std::pow(testparticles_, plist.size() - 1);
594 
595  // 5. Check that probability is smaller than one
596  if (prob > 1.) {
597  std::stringstream err;
598  err << "Probability larger than 1 for stochastic rates. ( P_nm = " << prob
599  << " )\nConsider using smaller timesteps.";
601  logg[LFindScatter].warn(err.str());
602  } else {
603  throw std::runtime_error(err.str());
604  }
605  }
606 
607  // 6. Perform probability decisions
608  double random_no = random::uniform(0., 1.);
609  if (random_no > prob) {
610  return nullptr;
611  }
612 
613  return std::move(act);
614 }
615 
617  const ParticleList& search_list, double dt, const double gcell_vol,
618  const std::vector<FourVector>& beam_momentum) const {
619  std::vector<ActionPtr> actions;
620  for (const ParticleData& p1 : search_list) {
621  for (const ParticleData& p2 : search_list) {
622  // Check for 2 particle scattering
623  if (p1.id() < p2.id()) {
624  ActionPtr act =
625  check_collision_two_part(p1, p2, dt, beam_momentum, gcell_vol);
626  if (act) {
627  actions.push_back(std::move(act));
628  }
629  }
630  if (incl_multi_set_.any()) {
631  // Also, check for 3 particle scatterings with stochastic criterion
632  for (const ParticleData& p3 : search_list) {
634  1 ||
636  1) {
637  if (p1.id() < p2.id() && p2.id() < p3.id()) {
638  ActionPtr act =
639  check_collision_multi_part({p1, p2, p3}, dt, gcell_vol);
640  if (act) {
641  actions.push_back(std::move(act));
642  }
643  }
644  }
645  for (const ParticleData& p4 : search_list) {
646  if (incl_multi_set_
648  if (p1.id() < p2.id() && p2.id() < p3.id() && p3.id() < p4.id()) {
649  ActionPtr act =
650  check_collision_multi_part({p1, p2, p3, p4}, dt, gcell_vol);
651  if (act) {
652  actions.push_back(std::move(act));
653  }
654  }
655  }
657  1 &&
658  search_list.size() >= 5) {
659  for (const ParticleData& p5 : search_list) {
660  if ((p1.id() < p2.id() && p2.id() < p3.id() &&
661  p3.id() < p4.id() && p4.id() < p5.id()) &&
662  (p1.is_pion() && p2.is_pion() && p3.is_pion() &&
663  p4.is_pion() && p5.is_pion())) {
664  // at the moment only pure pion 5-body reactions
665  ActionPtr act = check_collision_multi_part(
666  {p1, p2, p3, p4, p5}, dt, gcell_vol);
667  if (act) {
668  actions.push_back(std::move(act));
669  }
670  }
671  }
672  }
673  }
674  }
675  }
676  }
677  }
678  return actions;
679 }
680 
682  const ParticleList& search_list, const ParticleList& neighbors_list,
683  double dt, const std::vector<FourVector>& beam_momentum) const {
684  std::vector<ActionPtr> actions;
686  // Only search in cells
687  return actions;
688  }
689  for (const ParticleData& p1 : search_list) {
690  for (const ParticleData& p2 : neighbors_list) {
691  assert(p1.id() != p2.id());
692  // Check if a collision is possible.
693  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
694  if (act) {
695  actions.push_back(std::move(act));
696  }
697  }
698  }
699  return actions;
700 }
701 
703  const ParticleList& search_list, const Particles& surrounding_list,
704  double dt, const std::vector<FourVector>& beam_momentum) const {
705  std::vector<ActionPtr> actions;
707  // Only search in cells
708  return actions;
709  }
710  for (const ParticleData& p2 : surrounding_list) {
711  /* don't look for collisions if the particle from the surrounding list is
712  * also in the search list */
713  auto result = std::find_if(
714  search_list.begin(), search_list.end(),
715  [&p2](const ParticleData& p) { return p.id() == p2.id(); });
716  if (result != search_list.end()) {
717  continue;
718  }
719  for (const ParticleData& p1 : search_list) {
720  // Check if a collision is possible.
721  ActionPtr act = check_collision_two_part(p1, p2, dt, beam_momentum);
722  if (act) {
723  actions.push_back(std::move(act));
724  }
725  }
726  }
727  return actions;
728 }
729 
731  constexpr double time = 0.0;
732 
733  const size_t N_isotypes = IsoParticleType::list_all().size();
734  const size_t N_pairs = N_isotypes * (N_isotypes - 1) / 2;
735 
736  std::cout << N_isotypes << " iso-particle types." << std::endl;
737  std::cout << "They can make " << N_pairs << " pairs." << std::endl;
738  std::vector<double> momentum_scan_list = {0.1, 0.3, 0.5, 1.0,
739  2.0, 3.0, 5.0, 10.0};
740  for (const IsoParticleType& A_isotype : IsoParticleType::list_all()) {
741  for (const IsoParticleType& B_isotype : IsoParticleType::list_all()) {
742  if (&A_isotype > &B_isotype) {
743  continue;
744  }
745  bool any_nonzero_cs = false;
746  std::vector<std::string> r_list;
747  for (const ParticleTypePtr A_type : A_isotype.get_states()) {
748  for (const ParticleTypePtr B_type : B_isotype.get_states()) {
749  if (A_type > B_type) {
750  continue;
751  }
752  ParticleData A(*A_type), B(*B_type);
753  for (auto mom : momentum_scan_list) {
754  A.set_4momentum(A.pole_mass(), mom, 0.0, 0.0);
755  B.set_4momentum(B.pole_mass(), -mom, 0.0, 0.0);
756  ScatterActionPtr act = make_unique<ScatterAction>(
757  A, B, time, isotropic_, string_formation_time_);
758  if (strings_switch_) {
759  act->set_string_interface(string_process_interface_.get());
760  }
761  act->add_all_scatterings(
766  const double total_cs = act->cross_section();
767  if (total_cs <= 0.0) {
768  continue;
769  }
770  any_nonzero_cs = true;
771  for (const auto& channel : act->collision_channels()) {
772  const auto type = channel->get_type();
773  std::string r;
774  if (is_string_soft_process(type) ||
775  type == ProcessType::StringHard) {
776  r = A_type->name() + B_type->name() + std::string(" → strings");
777  } else {
778  std::string r_type =
779  (type == ProcessType::Elastic)
780  ? std::string(" (el)")
781  : (channel->get_type() == ProcessType::TwoToTwo)
782  ? std::string(" (inel)")
783  : std::string(" (?)");
784  r = A_type->name() + B_type->name() + std::string(" → ") +
785  channel->particle_types()[0]->name() +
786  channel->particle_types()[1]->name() + r_type;
787  }
788  isoclean(r);
789  r_list.push_back(r);
790  }
791  }
792  }
793  }
794  std::sort(r_list.begin(), r_list.end());
795  r_list.erase(std::unique(r_list.begin(), r_list.end()), r_list.end());
796  if (any_nonzero_cs) {
797  for (auto r : r_list) {
798  std::cout << r;
799  if (r_list.back() != r) {
800  std::cout << ", ";
801  }
802  }
803  std::cout << std::endl;
804  }
805  }
806  }
807 }
808 
812  std::string name_;
813 
816 
818  double mass_;
819 
828  FinalStateCrossSection(const std::string& name, double cross_section,
829  double mass)
830  : name_(name), cross_section_(cross_section), mass_(mass) {}
831 };
832 
833 namespace decaytree {
834 
846 struct Node {
847  public:
849  std::string name_;
850 
852  double weight_;
853 
855  ParticleTypePtrList initial_particles_;
856 
858  ParticleTypePtrList final_particles_;
859 
861  ParticleTypePtrList state_;
862 
864  std::vector<Node> children_;
865 
867  Node(const Node&) = delete;
869  Node(Node&&) = default;
870 
881  Node(const std::string& name, double weight,
882  ParticleTypePtrList&& initial_particles,
883  ParticleTypePtrList&& final_particles, ParticleTypePtrList&& state,
884  std::vector<Node>&& children)
885  : name_(name),
886  weight_(weight),
887  initial_particles_(std::move(initial_particles)),
888  final_particles_(std::move(final_particles)),
889  state_(std::move(state)),
890  children_(std::move(children)) {}
891 
903  Node& add_action(const std::string& name, double weight,
904  ParticleTypePtrList&& initial_particles,
905  ParticleTypePtrList&& final_particles) {
906  // Copy parent state and update it.
907  ParticleTypePtrList state(state_);
908  for (const auto& p : initial_particles) {
909  state.erase(std::find(state.begin(), state.end(), p));
910  }
911  for (const auto& p : final_particles) {
912  state.push_back(p);
913  }
914  // Sort the state to normalize the output.
915  std::sort(state.begin(), state.end(),
917  return a->name() < b->name();
918  });
919  // Push new node to children.
920  Node new_node(name, weight, std::move(initial_particles),
921  std::move(final_particles), std::move(state), {});
922  children_.emplace_back(std::move(new_node));
923  return children_.back();
924  }
925 
927  void print() const { print_helper(0); }
928 
932  std::vector<FinalStateCrossSection> final_state_cross_sections() const {
933  std::vector<FinalStateCrossSection> result;
934  final_state_cross_sections_helper(0, result, "", 1.);
935  return result;
936  }
937 
938  private:
945  void print_helper(uint64_t depth) const {
946  for (uint64_t i = 0; i < depth; i++) {
947  std::cout << " ";
948  }
949  std::cout << name_ << " " << weight_ << std::endl;
950  for (const auto& child : children_) {
951  child.print_helper(depth + 1);
952  }
953  }
954 
967  uint64_t depth, std::vector<FinalStateCrossSection>& result,
968  const std::string& name, double weight,
969  bool show_intermediate_states = false) const {
970  // The first node corresponds to the total cross section and has to be
971  // ignored. The second node corresponds to the partial cross section. All
972  // further nodes correspond to branching ratios.
973  if (depth > 0) {
974  weight *= weight_;
975  }
976 
977  std::string new_name;
978  double mass = 0.;
979 
980  if (show_intermediate_states) {
981  new_name = name;
982  if (!new_name.empty()) {
983  new_name += "->";
984  }
985  new_name += name_;
986  new_name += "{";
987  } else {
988  new_name = "";
989  }
990  for (const auto& s : state_) {
991  new_name += s->name();
992  mass += s->mass();
993  }
994  if (show_intermediate_states) {
995  new_name += "}";
996  }
997 
998  if (children_.empty()) {
999  result.emplace_back(FinalStateCrossSection(new_name, weight, mass));
1000  return;
1001  }
1002  for (const auto& child : children_) {
1003  child.final_state_cross_sections_helper(depth + 1, result, new_name,
1004  weight, show_intermediate_states);
1005  }
1006  }
1007 };
1008 
1017 static std::string make_decay_name(const std::string& res_name,
1018  const DecayBranchPtr& decay,
1019  ParticleTypePtrList& final_state) {
1020  std::stringstream name;
1021  name << "[" << res_name << "->";
1022  for (const auto& p : decay->particle_types()) {
1023  name << p->name();
1024  final_state.push_back(p);
1025  }
1026  name << "]";
1027  return name.str();
1028 }
1029 
1037 static void add_decays(Node& node, double sqrts) {
1038  // If there is more than one unstable particle in the current state, then
1039  // there will be redundant paths in the decay tree, corresponding to
1040  // reorderings of the decays. To avoid double counting, we normalize by the
1041  // number of possible decay orderings. Normalizing by the number of unstable
1042  // particles recursively corresponds to normalizing by the factorial that
1043  // gives the number of reorderings.
1044  //
1045  // Ideally, the redundant paths should never be added to the decay tree, but
1046  // we never have more than two redundant paths, so it probably does not
1047  // matter much.
1048  uint32_t n_unstable = 0;
1049  double sqrts_minus_masses = sqrts;
1050  for (const ParticleTypePtr ptype : node.state_) {
1051  if (!ptype->is_stable()) {
1052  n_unstable += 1;
1053  }
1054  sqrts_minus_masses -= ptype->mass();
1055  }
1056  const double norm =
1057  n_unstable != 0 ? 1. / static_cast<double>(n_unstable) : 1.;
1058 
1059  for (const ParticleTypePtr ptype : node.state_) {
1060  if (!ptype->is_stable()) {
1061  const double sqrts_decay = sqrts_minus_masses + ptype->mass();
1062  bool can_decay = false;
1063  for (const auto& decay : ptype->decay_modes().decay_mode_list()) {
1064  // Make sure to skip kinematically impossible decays.
1065  // In principle, we would have to integrate over the mass of the
1066  // resonance, but as an approximation we just assume it at its pole.
1067  double final_state_mass = 0.;
1068  for (const auto& p : decay->particle_types()) {
1069  final_state_mass += p->mass();
1070  }
1071  if (final_state_mass > sqrts_decay) {
1072  continue;
1073  }
1074  can_decay = true;
1075 
1076  ParticleTypePtrList parts;
1077  const auto name = make_decay_name(ptype->name(), decay, parts);
1078  auto& new_node = node.add_action(name, norm * decay->weight(), {ptype},
1079  std::move(parts));
1080  add_decays(new_node, sqrts_decay);
1081  }
1082  if (!can_decay) {
1083  // Remove final-state cross sections with resonances that cannot
1084  // decay due to our "mass = pole mass" approximation.
1085  node.weight_ = 0;
1086  return;
1087  }
1088  }
1089  }
1090 }
1091 
1092 } // namespace decaytree
1093 
1099 static void deduplicate(std::vector<FinalStateCrossSection>& final_state_xs) {
1100  std::sort(final_state_xs.begin(), final_state_xs.end(),
1101  [](const FinalStateCrossSection& a,
1102  const FinalStateCrossSection& b) { return a.name_ < b.name_; });
1103  auto current = final_state_xs.begin();
1104  while (current != final_state_xs.end()) {
1105  auto adjacent = std::adjacent_find(
1106  current, final_state_xs.end(),
1107  [](const FinalStateCrossSection& a, const FinalStateCrossSection& b) {
1108  return a.name_ == b.name_;
1109  });
1110  current = adjacent;
1111  if (adjacent != final_state_xs.end()) {
1112  adjacent->cross_section_ += (adjacent + 1)->cross_section_;
1113  final_state_xs.erase(adjacent + 1);
1114  }
1115  }
1116 }
1117 
1119  const ParticleType& a, const ParticleType& b, double m_a, double m_b,
1120  bool final_state, std::vector<double>& plab) const {
1121  typedef std::vector<std::pair<double, double>> xs_saver;
1122  std::map<std::string, xs_saver> xs_dump;
1123  std::map<std::string, double> outgoing_total_mass;
1124 
1125  ParticleData a_data(a), b_data(b);
1126  int n_momentum_points = 200;
1127  constexpr double momentum_step = 0.02;
1128  /*
1129  // Round to output precision.
1130  for (auto& p : plab) {
1131  p = std::floor((p * 100000) + 0.5) / 100000;
1132  }
1133  */
1134  if (plab.size() > 0) {
1135  n_momentum_points = plab.size();
1136  // Remove duplicates.
1137  std::sort(plab.begin(), plab.end());
1138  plab.erase(std::unique(plab.begin(), plab.end()), plab.end());
1139  }
1140  for (int i = 0; i < n_momentum_points; i++) {
1141  double momentum;
1142  if (plab.size() > 0) {
1143  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
1144  } else {
1145  momentum = momentum_step * (i + 1);
1146  }
1147  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
1148  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
1149  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
1150  const ParticleList incoming = {a_data, b_data};
1151  ScatterActionPtr act = make_unique<ScatterAction>(
1152  a_data, b_data, 0., isotropic_, string_formation_time_);
1153  if (strings_switch_) {
1154  act->set_string_interface(string_process_interface_.get());
1155  }
1156  act->add_all_scatterings(elastic_parameter_, two_to_one_, incl_set_,
1160  decaytree::Node tree(a.name() + b.name(), act->cross_section(), {&a, &b},
1161  {&a, &b}, {&a, &b}, {});
1162  const CollisionBranchList& processes = act->collision_channels();
1163  for (const auto& process : processes) {
1164  const double xs = process->weight();
1165  if (xs <= 0.0) {
1166  continue;
1167  }
1168  if (!final_state) {
1169  std::stringstream process_description_stream;
1170  process_description_stream << *process;
1171  const std::string& description = process_description_stream.str();
1172  double m_tot = 0.0;
1173  for (const auto& ptype : process->particle_types()) {
1174  m_tot += ptype->mass();
1175  }
1176  outgoing_total_mass[description] = m_tot;
1177  if (!xs_dump[description].empty() &&
1178  std::abs(xs_dump[description].back().first - sqrts) <
1179  really_small) {
1180  xs_dump[description].back().second += xs;
1181  } else {
1182  xs_dump[description].push_back(std::make_pair(sqrts, xs));
1183  }
1184  } else {
1185  std::stringstream process_description_stream;
1186  process_description_stream << *process;
1187  const std::string& description = process_description_stream.str();
1188  ParticleTypePtrList initial_particles = {&a, &b};
1189  ParticleTypePtrList final_particles = process->particle_types();
1190  auto& process_node =
1191  tree.add_action(description, xs, std::move(initial_particles),
1192  std::move(final_particles));
1193  decaytree::add_decays(process_node, sqrts);
1194  }
1195  }
1196  xs_dump["total"].push_back(std::make_pair(sqrts, act->cross_section()));
1197  // Total cross-section should be the first in the list -> negative mass
1198  outgoing_total_mass["total"] = -1.0;
1199  if (final_state) {
1200  // tree.print();
1201  auto final_state_xs = tree.final_state_cross_sections();
1202  deduplicate(final_state_xs);
1203  for (const auto& p : final_state_xs) {
1204  // Don't print empty columns.
1205  //
1206  // FIXME(steinberg): The better fix would be to not have them in the
1207  // first place.
1208  if (p.name_ == "") {
1209  continue;
1210  }
1211  outgoing_total_mass[p.name_] = p.mass_;
1212  xs_dump[p.name_].push_back(std::make_pair(sqrts, p.cross_section_));
1213  }
1214  }
1215  }
1216  // Get rid of cross sections that are zero.
1217  // (This only happens if their is a resonance in the final state that cannot
1218  // decay with our simplified assumptions.)
1219  for (auto it = begin(xs_dump); it != end(xs_dump);) {
1220  // Sum cross section over all energies.
1221  const xs_saver& xs = (*it).second;
1222  double sum = 0;
1223  for (const auto& p : xs) {
1224  sum += p.second;
1225  }
1226  if (sum == 0.) {
1227  it = xs_dump.erase(it);
1228  } else {
1229  ++it;
1230  }
1231  }
1232 
1233  // Nice ordering of channels by summed pole mass of products
1234  std::vector<std::string> all_channels;
1235  for (const auto& channel : xs_dump) {
1236  all_channels.push_back(channel.first);
1237  }
1238  std::sort(all_channels.begin(), all_channels.end(),
1239  [&](const std::string& str_a, const std::string& str_b) {
1240  return outgoing_total_mass[str_a] < outgoing_total_mass[str_b];
1241  });
1242 
1243  // Print header
1244  std::cout << "# Dumping partial " << a.name() << b.name()
1245  << " cross-sections in mb, energies in GeV" << std::endl;
1246  std::cout << " sqrt_s";
1247  // Align everything to 16 unicode characters.
1248  // This should be enough for the longest channel name (7 final-state
1249  // particles).
1250  for (const auto& channel : all_channels) {
1251  std::cout << utf8::fill_left(channel, 16, ' ');
1252  }
1253  std::cout << std::endl;
1254 
1255  // Print out all partial cross-sections in mb
1256  for (int i = 0; i < n_momentum_points; i++) {
1257  double momentum;
1258  if (plab.size() > 0) {
1259  momentum = pCM_from_s(s_from_plab(plab.at(i), m_a, m_b), m_a, m_b);
1260  } else {
1261  momentum = momentum_step * (i + 1);
1262  }
1263  a_data.set_4momentum(m_a, momentum, 0.0, 0.0);
1264  b_data.set_4momentum(m_b, -momentum, 0.0, 0.0);
1265  const double sqrts = (a_data.momentum() + b_data.momentum()).abs();
1266  std::printf("%9.6f", sqrts);
1267  for (const auto& channel : all_channels) {
1268  const xs_saver energy_and_xs = xs_dump[channel];
1269  size_t j = 0;
1270  for (; j < energy_and_xs.size() && energy_and_xs[j].first < sqrts; j++) {
1271  }
1272  double xs = 0.0;
1273  if (j < energy_and_xs.size() &&
1274  std::abs(energy_and_xs[j].first - sqrts) < really_small) {
1275  xs = energy_and_xs[j].second;
1276  }
1277  std::printf("%16.6f", xs); // Same alignment as in the header.
1278  }
1279  std::printf("\n");
1280  }
1281 }
1282 
1283 } // namespace smash
Interface to the SMASH configuration files.
const DecayBranchList & decay_mode_list() const
Definition: decaymodes.h:63
double x0() const
Definition: fourvector.h:308
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
double xsec_scaling_factor(double delta_time=0.) const
Return the cross section scaling factor at a given time.
Definition: particledata.cc:82
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:671
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
const int testparticles_
Number of test particles.
const MultiParticleReactionsBitSet incl_multi_set_
List of included multi-particle reactions.
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.
const double elastic_parameter_
Elastic cross section parameter (in mb).
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 strings_switch_
Switch to turn off string excitation.
const double additional_el_xs_
Additional constant elastic cross section.
const double scale_xs_
Factor by which all (partial) cross sections are scaled.
const bool allow_first_collisions_within_nucleus_
If particles within nucleus are allowed to collide for their first time.
const NNbarTreatment nnbar_treatment_
Switch for NNbar reactions.
const double low_snn_cut_
Elastic collsions between two nucleons with sqrt_s below low_snn_cut_ are excluded.
const bool use_AQM_
Switch to control whether to use AQM or not.
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,...
const bool only_warn_for_high_prob_
Switch to turn off throwing an exception for collision probabilities larger than 1.
const CollisionCriterion coll_crit_
Specifies which collision criterion is used.
ScatterActionsFinder(Configuration config, const ExperimentParameters &parameters)
Constructor of the finder with the given parameters.
const ReactionsBitSet incl_set_
List of included 2<->2 reactions.
const bool two_to_one_
Enable 2->1 processes.
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.
const bool strings_with_probability_
Decide whether to implement string fragmentation based on a probability.
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.
const PdgCode dprime(PdgCode::from_decimal(1000010021))
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
static void deduplicate(std::vector< FinalStateCrossSection > &final_state_xs)
Deduplicate the final-state cross sections by summing.
@ TwoToTwo
2->2 inelastic scattering
@ Elastic
elastic scattering: particles remain the same, only momenta change
@ StringHard
hard string process involving 2->2 QCD process by PYTHIA.
bool all_of(Container &&c, UnaryPredicate &&p)
Convenience wrapper for std::all_of that operates on a complete container.
Definition: algorithms.h:80
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.
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
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.