Version: SMASH-3.3
decayaction.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2013-2021,2024-2025
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #include "smash/decayaction.h"
11 
12 #include "smash/decaymodes.h"
13 #include "smash/logging.h"
14 #include "smash/pdgcode.h"
16 
17 namespace smash {
18 
19 static constexpr int LDecayModes = LogArea::DecayModes::id;
20 
22  SpinInteractionType spin_interaction_type)
23  : Action({p}, time),
24  total_width_(0.),
25  spin_interaction_type_(spin_interaction_type) {}
26 
27 void DecayAction::add_decays(DecayBranchList pv) {
28  add_processes<DecayBranch>(std::move(pv), decay_channels_, total_width_);
29 }
30 
31 void DecayAction::add_decay(DecayBranchPtr p) {
32  add_process<DecayBranch>(p, decay_channels_, total_width_);
33 }
34 
36  logg[LDecayModes].debug("Process: Resonance decay. ");
37  /* Execute a decay process for the selected particle.
38  *
39  * randomly select one of the decay modes of the particle
40  * according to their relative weights. Then decay the particle
41  * by calling function sample_2body_phasespace or sample_manybody_phasespace.
42  */
43  const DecayBranch *proc =
44  choose_channel<DecayBranch>(decay_channels_, total_width_);
46  // set positions of the outgoing particles
47  for (auto &p : outgoing_particles_) {
48  p.set_4position(incoming_particles_[0].position());
49  }
50  process_type_ = proc->get_type();
51  L_ = proc->angular_momentum();
52  partial_width_ = proc->weight();
53 
54  switch (outgoing_particles_.size()) {
55  case 2:
57  if (pot_pointer) {
59  } else {
60  return true;
61  }
62  case 3:
64  return true;
65  default:
66  throw InvalidDecay(
67  "DecayAction::perform: Only 1->2 or 1->3 processes are supported. "
68  "Decay from 1->" +
70  " was requested. (PDGcode=" +
71  incoming_particles_[0].pdgcode().string() + ", mass=" +
72  std::to_string(incoming_particles_[0].effective_mass()) + ")");
73  }
74 }
75 
77  int n_try = 1000;
78  while (n_try--) {
80  break;
81  }
82 
83  const bool core_in_incoming = incoming_particles_[0].is_core();
84  // Set formation time.
85  for (auto &p : outgoing_particles_) {
86  logg[LDecayModes].debug("particle momenta in lrf ", p);
87  // assuming decaying particles are always fully formed
88  p.set_formation_time(time_of_execution_);
89  // Boost to the computational frame
90  p.boost_momentum(-total_momentum_of_outgoing_particles().velocity());
91  logg[LDecayModes].debug("particle momenta in comp ", p);
92  if (core_in_incoming) {
93  p.fluidize();
94  }
95  if (p.type().pdgcode().is_heavy_flavor()) {
96  p.set_perturbative_weight(incoming_particles_[0].perturbative_weight());
97  }
98  }
99 
100  /*
101  * @brief Σ* → Λ + π decay: propagate Λ polarization from the intermediate
102  * resonance.
103  *
104  * During Λ+π → Σ* formation we stored the incoming Λ polarization by
105  * writing it into the Σ* spin 4-vector (optionally applying a Λ spin-flip
106  * probability, cf. arXiv:2404.15890v2). At decay, we must hand this
107  * polarization back to the outgoing Λ to transport Λ polarization through the
108  * resonance stage.
109  */
111  outgoing_particles_.size() == 2) {
112  // Check for Σ* → Λ + π decay channel
113  int lambda_idx = -1;
114  int pion_idx = -1;
115  const bool is_sigmastar_decay = is_sigmastar_to_lambda_pion_decay(
117  lambda_idx, pion_idx);
118 
119  if (is_sigmastar_decay) {
120  auto &sigma_star = incoming_particles_[0];
121  auto &lambda = outgoing_particles_[lambda_idx];
122  auto &pion = outgoing_particles_[pion_idx];
123  // Copy spin vector from Σ* to Λ and boost it to the Λ frame
124  FourVector final_spin_vector =
125  sigma_star.spin_vector().lorentz_boost(lambda.velocity());
126  lambda.set_spin_vector(final_spin_vector);
127  pion.set_spin_vector(FourVector{0., 0., 0., 0.});
128  } else {
129  // Set unpolarized spin vectors
131  }
133  outgoing_particles_.size() != 2) {
134  // Set unpolarized spin vectors
136  }
137 }
138 
140  assert(outgoing_particles_.size() == 2);
142  const double cm_kin_energy = p_tot.abs();
143  const std::pair<double, double> masses = sample_masses(cm_kin_energy);
144 
145  const bool is_valid = !std::isnan(masses.first) && !std::isnan(masses.second);
146 
147  if (pot_pointer) {
149  }
150 
151  sample_angles(masses, cm_kin_energy);
152 }
153 
154 /* This is overridden from the Action class in order to
155  * take care of the angular momentum L_. */
156 std::pair<double, double> DecayAction::sample_masses(
157  double kinetic_energy_cm) const {
158  const ParticleType &t_a = outgoing_particles_[0].type();
159  const ParticleType &t_b = outgoing_particles_[1].type();
160 
161  // start with pole masses
162  std::pair<double, double> masses = {t_a.mass(), t_b.mass()};
163 
164  const bool below_threshold_energy =
165  kinetic_energy_cm < t_a.min_mass_kinematic() + t_b.min_mass_kinematic();
166 
167  const bool return_nan_on_failure = pot_pointer != nullptr;
168 
169  if (below_threshold_energy) {
170  if (return_nan_on_failure) {
171  return {smash_NaN<double>, smash_NaN<double>};
172  } else {
173  const std::string reaction =
174  incoming_particles_[0].type().name() + "→" + t_a.name() + t_b.name();
176  reaction + ": not enough energy, " +
177  std::to_string(kinetic_energy_cm) + " < " +
178  std::to_string(t_a.min_mass_kinematic()) + " + " +
180  }
181  }
182  // If one of the particles is a resonance, sample its mass.
183  if (!t_a.is_stable() && t_b.is_stable()) {
184  masses.first = t_a.sample_resonance_mass(t_b.mass(), kinetic_energy_cm, L_);
185  } else if (!t_b.is_stable() && t_a.is_stable()) {
186  masses.second =
187  t_b.sample_resonance_mass(t_a.mass(), kinetic_energy_cm, L_);
188  } else if (!t_a.is_stable() && !t_b.is_stable()) {
189  // two resonances in final state
190  masses = t_a.sample_resonance_masses(t_b, kinetic_energy_cm, L_);
191  }
192 
193  return masses;
194 }
195 
196 void DecayAction::format_debug_output(std::ostream &out) const {
197  out << "Decay of " << incoming_particles_ << " to " << outgoing_particles_
198  << ", sqrt(s)=" << format(sqrt_s(), "GeV", 11, 9);
199 }
200 
201 } // namespace smash
Thrown for example when ScatterAction is called to perform with a wrong number of final-state particl...
Definition: action.h:330
Action is the base class for a generic process that takes a number of incoming particles and transfor...
Definition: action.h:35
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
Definition: action.cc:160
virtual void sample_angles(std::pair< double, double > masses, double kinetic_energy_cm)
Sample final-state momenta in general X->2 processes (here: using an isotropical angular distribution...
Definition: action.cc:282
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
Definition: action.h:372
virtual void sample_manybody_phasespace()
Sample the full n-body phase-space (masses, momenta, angles) in the center-of-mass frame for the fina...
Definition: action.cc:449
const double time_of_execution_
Time at which the action is supposed to be performed (absolute time in the lab frame in fm).
Definition: action.h:378
void assign_unpolarized_spin_vector_to_outgoing_particles()
Assign an unpolarized spin vector to all outgoing particles.
Definition: action.cc:478
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:271
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:364
ProcessType process_type_
type of process
Definition: action.h:381
bool is_valid(const Particles &particles) const
Check whether the action still applies.
Definition: action.cc:29
Thrown when DecayAction is called to perform with 0 or more than 2 entries in outgoing_particles.
Definition: decayaction.h:106
double partial_width_
partial decay width to the chosen outgoing channel
Definition: decayaction.h:184
double total_width_
total decay width
Definition: decayaction.h:181
bool sample_outgoing_particles()
Sample outgoing particle types, masses and angles return success of sampling.
Definition: decayaction.cc:35
SpinInteractionType spin_interaction_type_
Spin interaction type.
Definition: decayaction.h:191
void add_decays(DecayBranchList pv)
Add several new decays at once.
Definition: decayaction.cc:27
DecayAction(const ParticleData &p, double time, SpinInteractionType spin_interaction_type=SpinInteractionType::Off)
Construct a DecayAction from a particle p.
Definition: decayaction.cc:21
std::optional< bool > was_2body_phase_space_sampled_with_potentials_as_valid_
Optional success flag for sampling outgoing particles.
Definition: decayaction.h:125
void generate_final_state() override
Generate the final state of the decay process.
Definition: decayaction.cc:76
void add_decay(DecayBranchPtr p)
Add one new decay.
Definition: decayaction.cc:31
DecayBranchList decay_channels_
List of possible decays.
Definition: decayaction.h:178
int L_
Angular momentum of the decay.
Definition: decayaction.h:187
void sample_2body_phasespace() override
Sample the full 2-body phase space (masses, momenta, angles) in the center-of-mass frame for the fina...
Definition: decayaction.cc:139
static bool is_sigmastar_to_lambda_pion_decay(const ParticleData &parent, const ParticleData &daughter0, const ParticleData &daughter1, int &lambda_idx, int &pion_idx)
Check for the decay Σ*→Λπ.
Definition: decayaction.h:138
std::pair< double, double > sample_masses(double kinetic_energy_cm) const override
Sample the masses of the final particles.
Definition: decayaction.cc:156
DecayBranch is a derivative of ProcessBranch, which is used to represent decay channels.
ProcessType get_type() const override
int angular_momentum() const
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:33
double abs() const
calculate the lorentz invariant absolute value
Definition: fourvector.h:464
FourVector lorentz_boost(const ThreeVector &v) const
Returns the FourVector boosted with velocity v.
Definition: fourvector.cc:17
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:59
Particle type contains the static properties of a particle species.
Definition: particletype.h:98
double sample_resonance_mass(const double mass_stable, const double cms_energy, int L=0) const
Resonance mass sampling for 2-particle final state with one resonance (type given by 'this') and one ...
double min_mass_kinematic() const
The minimum mass of the resonance that is kinematically allowed.
const std::string & name() const
Definition: particletype.h:142
bool is_stable() const
Definition: particletype.h:249
double mass() const
Definition: particletype.h:145
std::pair< double, double > sample_resonance_masses(const ParticleType &t2, const double cms_energy, int L=0) const
Resonance mass sampling for 2-particle final state with two resonances.
ParticleList particle_list() const
double weight() const
SpinInteractionType
Possible spin interaction types.
@ Off
No spin interactions.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:40
void format_debug_output(std::ostream &out) const override
Writes information about this decay action to the out stream.
Definition: decayaction.cc:196
FormattingHelper< T > format(const T &value, const char *unit, int width=-1, int precision=-1)
Acts as a stream modifier for std::ostream to output an object with an optional suffix string and wit...
Definition: logging.h:217
constexpr int p
Proton.
Definition: action.h:24
std::string to_string(ThermodynamicQuantity quantity)
Convert a ThermodynamicQuantity enum value to its corresponding string.
Definition: stringify.cc:26
Potentials * pot_pointer
Pointer to a Potential class.
static constexpr int LDecayModes
Definition: decayaction.cc:19