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