30                           ParticleTypePtrList particle_types) {
 
   34     if (type == &mode->type()) {
 
   35       mode->set_weight(mode->weight() + ratio);
 
   40   decay_modes_.push_back(std::make_unique<DecayBranch>(*type, ratio));
 
   44                                       ParticleTypePtrList particle_types,
 
   50     if (type->has_mother(mother) && type->has_particles(particle_types) &&
 
   51         type->angular_momentum() == L) {
 
   57   switch (particle_types.size()) {
 
   60                       particle_types[1]->pdgcode())) {
 
   62             std::make_unique<TwoBodyDecayDilepton>(particle_types, L));
 
   63       } 
else if (particle_types[0]->is_stable() &&
 
   64                  particle_types[1]->is_stable()) {
 
   66             std::make_unique<TwoBodyDecayStable>(particle_types, L));
 
   67       } 
else if (particle_types[0]->is_stable() ||
 
   68                  particle_types[1]->is_stable()) {
 
   70             std::make_unique<TwoBodyDecaySemistable>(particle_types, L));
 
   73             std::make_unique<TwoBodyDecayUnstable>(particle_types, L));
 
   78                           particle_types[1]->pdgcode(),
 
   79                           particle_types[2]->pdgcode())) {
 
   81             mother, particle_types, L));
 
   84             std::make_unique<ThreeBodyDecay>(particle_types, L));
 
   89           "DecayModes::get_decay_type was instructed to add a decay mode " 
   92           " particles. This is an invalid input.");
 
  100   bool is_large_renormalization = 
false;
 
  102     sum += mode->weight();
 
  106                             ": Extremely small renormalization constant: ", sum,
 
  107                             "\n=> Skipping the renormalization.");
 
  109     is_large_renormalization = (std::abs(sum - 1.) > 0.01);
 
  111                             ": Renormalizing decay modes with ", sum);
 
  112     double new_sum = 0.0;
 
  114       mode->set_weight(mode->weight() / sum);
 
  115       new_sum += mode->weight();
 
  117     logg[
LDecayModes].debug(
"After renormalization sum of ratios is ", new_sum);
 
  119   return is_large_renormalization;
 
  139   int min_L = std::min(std::abs(s0 - s1 - s2), std::abs(s0 - s1 + s2));
 
  140   min_L = std::min(min_L, std::abs(s0 + s1 - s2));
 
  141   if (min_L % 2 != 0) {
 
  142     throw std::runtime_error(
 
  143         "min_angular_momentum: sum of spins should be integer");
 
  150       std::min(std::abs(s0 - s1 + s2 + s3), std::abs(s0 + s1 - s2 + s3));
 
  151   min_L = std::min(min_L, std::abs(s0 + s1 + s2 - s3));
 
  152   min_L = std::min(min_L, std::abs(s0 - s1 - s2 + s3));
 
  153   min_L = std::min(min_L, std::abs(s0 - s1 + s2 - s3));
 
  154   min_L = std::min(min_L, std::abs(s0 + s1 - s2 - s3));
 
  155   min_L = std::min(min_L, std::abs(s0 - s1 - s2 - s3));
 
  156   if (min_L % 2 != 0) {
 
  157     throw std::runtime_error(
 
  158         "min_angular_momentum: sum of spins should be integer");
 
  166   static std::vector<DecayTypePtr> decaytypes;
 
  172   static std::vector<DecayModes> decaymodes;
 
  178   ParticleTypePtrList mother_states;
 
  179   std::vector<DecayModes> decay_modes_to_add;  
 
  180   int total_large_renormalized = 0;
 
  182   const auto end_of_decaymodes = [&]() {
 
  183     if (isotype_mother == 
nullptr) {  
 
  187     for (
size_t m = 0; m < mother_states.size(); m++) {
 
  188       if (decay_modes_to_add[m].
is_empty() && !mother_states[m]->is_stable()) {
 
  190                             mother_states[m]->name());
 
  192       bool is_large_renorm =
 
  193           decay_modes_to_add[m].renormalize(mother_states[m]->name());
 
  194       total_large_renormalized += is_large_renorm;
 
  195       PdgCode pdgcode = mother_states[m]->pdgcode();
 
  197       decaymodes[
find_offset(pdgcode)] = std::move(decay_modes_to_add[m]);
 
  201       logg[
LDecayModes].debug(
"generating decay modes for anti-multiplet: " +
 
  202                               isotype_mother->
name());
 
  203       for (
const auto &state : mother_states) {
 
  204         PdgCode pdg = state->pdgcode();
 
  210           ParticleTypePtrList list = mode->particle_types();
 
  211           for (
auto &type : list) {
 
  212             if (type->has_antiparticle()) {
 
  213               type = type->get_antiparticle();
 
  216           decay_modes_anti.
add_mode(&type_anti, mode->weight(),
 
  217                                     mode->angular_momentum(), list);
 
  225   uint64_t linenumber = 1;
 
  227     const auto trimmed = 
trim(line.text);
 
  228     assert(!trimmed.empty());  
 
  230     if (trimmed.find_first_of(
" \t") == std::string::npos) {
 
  233       std::string name = 
trim(line.text);
 
  236       decay_modes_to_add.clear();
 
  237       decay_modes_to_add.resize(mother_states.size());
 
  240       for (
size_t m = 0; m < mother_states.size(); m++) {
 
  241         PdgCode pdgcode = mother_states[m]->pdgcode();
 
  248       std::istringstream lineinput(line.text);
 
  249       std::vector<std::string> decay_particles;
 
  250       decay_particles.reserve(3);
 
  259                           ": '" + line.text + 
"'");
 
  266         decay_particles.emplace_back(name);
 
  268         const bool is_multiplet = isotype;
 
  270         if (!is_multiplet && !is_state) {
 
  273               " is neither an isospin multiplet nor a particle." + 
" (line " +
 
  276         const bool is_hadronic_multiplet =
 
  277             is_multiplet && isotype->get_states()[0]->is_hadron();
 
  278         multi &= is_hadronic_multiplet;
 
  282       const int s0 = isotype_mother->
spin();
 
  288         switch (decay_particles.size()) {
 
  294             parity = isotype_daughter_1.
parity() * isotype_daughter_2.
parity();
 
  295             const int s1 = isotype_daughter_1.
spin();
 
  296             const int s2 = isotype_daughter_2.
spin();
 
  298             max_L = (s0 + s1 + s2) / 2;
 
  300             bool forbidden_by_isospin = 
true;
 
  301             for (
size_t m = 0; m < mother_states.size(); m++) {
 
  302               for (
const auto &daughter1 : isotype_daughter_1.
get_states()) {
 
  303                 for (
const auto &daughter2 : isotype_daughter_2.
get_states()) {
 
  306                       *daughter1, *daughter2, *mother_states[m]);
 
  310                         "decay mode generated: " + mother_states[m]->name() +
 
  311                         " -> " + daughter1->name() + 
" " + daughter2->name() +
 
  313                     decay_modes_to_add[m].add_mode(mother_states[m],
 
  315                                                    {daughter1, daughter2});
 
  316                     forbidden_by_isospin = 
false;
 
  321             if (forbidden_by_isospin) {
 
  323               s << 
",\nwhere isospin mother: " << isotype_mother->
isospin()
 
  324                 << 
", daughters: " << isotype_daughter_1.
isospin() << 
" " 
  325                 << isotype_daughter_2.
isospin();
 
  327                                  " decay mode is forbidden by isospin: \"" +
 
  328                                  line.text + 
"\"" + s.str());
 
  339             parity = isotype_daughter_1.
parity() * isotype_daughter_2.
parity() *
 
  340                      isotype_daughter_3.
parity();
 
  341             const int s1 = isotype_daughter_1.
spin();
 
  342             const int s2 = isotype_daughter_2.
spin();
 
  343             const int s3 = isotype_daughter_2.
spin();
 
  345             max_L = (s0 + s1 + s2 + s3) / 2;
 
  347             for (
size_t m = 0; m < mother_states.size(); m++) {
 
  348               for (
const auto &daughter1 : isotype_daughter_1.
get_states()) {
 
  349                 for (
const auto &daughter2 : isotype_daughter_2.
get_states()) {
 
  350                   for (
const auto &daughter3 :
 
  353                         *daughter1, *daughter2, *daughter3, *mother_states[m]);
 
  357                           "decay mode generated: " + mother_states[m]->name() +
 
  358                           " -> " + daughter1->name() + 
" " + daughter2->name() +
 
  359                           " " + daughter3->name() + 
" (" +
 
  361                       decay_modes_to_add[m].add_mode(
 
  362                           mother_states[m], ratio * cg_sqr, L,
 
  363                           {daughter1, daughter2, daughter3});
 
  372             throw std::runtime_error(
 
  373                 "References to isospin multiplets only " 
  374                 "allowed in two-body or three-body decays: " +
 
  381         ParticleTypePtrList types;
 
  384         for (
auto part : decay_particles) {
 
  387           } 
catch (std::runtime_error &e) {
 
  388             throw std::runtime_error(std::string() + e.what() + 
" (line " +
 
  392           charge += types.back()->charge();
 
  393           parity *= types.back()->parity();
 
  395         if (types.size() == 2) {
 
  396           const int s1 = types[0]->spin();
 
  397           const int s2 = types[1]->spin();
 
  399           max_L = (s0 + s1 + s2) / 2;
 
  400         } 
else if (types.size() == 3) {
 
  401           const int s1 = types[0]->spin();
 
  402           const int s2 = types[1]->spin();
 
  403           const int s3 = types[2]->spin();
 
  405           max_L = (s0 + s1 + s2 + s3) / 2;
 
  408                              " decay mode has an invalid number of particles" 
  409                              " in the final state " +
 
  413         bool no_decays = 
true;
 
  414         for (
size_t m = 0; m < mother_states.size(); m++) {
 
  415           if (mother_states[m]->charge() == charge) {
 
  417                 "decay mode found: ", mother_states[m]->name(), 
" -> ",
 
  418                 decay_particles.size());
 
  419             decay_modes_to_add[m].add_mode(mother_states[m], ratio, L, types);
 
  425                              " decay mode violates charge conservation " +
 
  437       if (decay_particles.size() == 2 && parity != mother_states[0]->parity()) {
 
  439                            " decay mode violates parity conservation " +
 
  444       if (L < min_L || L > max_L) {
 
  446             mother_states[0]->name() +
 
  447             " decay mode violates angular momentum conservation: " +
 
  450             ": \"" + trimmed + 
"\")");
 
  460   for (
const auto &mother : particles) {
 
  461     if (mother.is_stable()) {
 
  464     const auto &decays = mother.decay_modes().decay_mode_list();
 
  465     for (
const auto &decay : decays) {
 
  466       if (mother.mass() <= decay->threshold()) {
 
  468         s << mother.name() << 
" →  ";
 
  469         for (
const auto &
p : decay->particle_types()) {
 
  472         s << 
" with " << mother.mass() << 
" ≤ " << decay->threshold();
 
  474             "For all decays, the minimum mass of daughters" 
  475             "must be smaller\nthan the mother's pole mass " 
  476             "(Manley-Saleski Ansatz)\n" 
  477             "Violated by the following decay: " +
 
  482   if (total_large_renormalized > 0) {
 
  484         "Branching ratios of ", total_large_renormalized,
 
  485         " hadrons were renormalized by more than 1% to have sum 1.");
 
The DecayModes class is used to store and update information about decay branches (i....
 
const DecayBranchList & decay_mode_list() const
 
void add_mode(ParticleTypePtr mother, double ratio, int L, ParticleTypePtrList particle_types)
Add a decay mode using all necessary information.
 
bool renormalize(const std::string &name)
Renormalize the branching ratios to add up to 1.
 
static std::vector< DecayModes > * all_decay_modes
A list of all DecayModes objects using the same indexing as all_particle_types.
 
DecayBranchList decay_modes_
Vector of decay modes.
 
static DecayType * get_decay_type(ParticleTypePtr mother, ParticleTypePtrList particle_types, int L)
Retrieve a decay type.
 
static void load_decaymodes(const std::string &input)
Loads the DecayModes map as described in the input string.
 
DecayType is the abstract base class for all decay types.
 
IsoParticleType is a class to represent isospin multiplets.
 
static const IsoParticleType * try_find(const std::string &name)
Returns the IsoParticleType pointer for the given name.
 
static const IsoParticleType & find(const std::string &name)
Returns the IsoParticleType object for the given name.
 
int isospin() const
Returns twice the total isospin of the multiplet.
 
bool has_anti_multiplet() const
Check if there is a multiplet of antiparticles, which is different from the original multiplet.
 
const std::string & name() const
Returns the name of the multiplet.
 
unsigned int spin() const
Returns twice the spin of the multiplet.
 
ParticleTypePtrList get_states() const
Returns list of states that form part of the multiplet.
 
static const ParticleTypePtr find_state(const std::string &name)
Returns the ParticleType object for the given name, by first finding the correct multiplet and then l...
 
A pointer-like interface to global references to ParticleType objects.
 
Particle type contains the static properties of a particle species.
 
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
 
static bool exists(PdgCode pdgcode)
 
static const ParticleTypeList & list_all()
 
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
 
PdgCode get_antiparticle() const
Construct the antiparticle to a given PDG code.
 
Collection of useful constants that are known at compile time.
 
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
 
std::size_t find_offset(PdgCode pdg)
Passes back the address offset of a particletype in the list of all particles.
 
bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2)
 
static int min_angular_momentum(int s0, int s1, int s2)
 
double isospin_clebsch_gordan_sqr_3to1(const ParticleType &p_a, const ParticleType &p_b, const ParticleType &p_c, const ParticleType &Res)
Calculate the squared isospin Clebsch-Gordan coefficient for three particles p_a, p_b and p_c couplin...
 
build_vector_< Line > line_parser(const std::string &input)
Helper function for parsing particles.txt and decaymodes.txt.
 
std::string trim(const std::string &s)
Strip leading and trailing whitespaces.
 
Parity
Represent the parity of a particle type.
 
bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2, const PdgCode pdg3)
 
std::string to_string(ThermodynamicQuantity quantity)
Convert a ThermodynamicQuantity enum value to its corresponding string.
 
std::vector< DecayTypePtr > * all_decay_types
Global pointer to the decay types list.
 
constexpr double really_small
Numerical error tolerance.
 
double isospin_clebsch_gordan_sqr_2to1(const ParticleType &p_a, const ParticleType &p_b, const ParticleType &Res)
Calculate the squared isospin Clebsch-Gordan coefficient for two particles p_a and p_b coupling to a ...
 
static constexpr int LDecayModes
 
Line consists of a line number and the contents of that line.