Belle II Software  release-08-01-10
DecayDescriptor.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <analysis/DecayDescriptor/DecayString.h>
12 #include <analysis/DecayDescriptor/DecayDescriptorParticle.h>
13 #include <analysis/dataobjects/Particle.h>
14 
15 #include <vector>
16 #include <string>
17 
18 namespace Belle2 {
23  class MCParticle;
24  class Particle;
25 
34  private:
41  std::vector<DecayDescriptor> m_daughters;
45  bool m_isNULL;
47  template <class T>
48  int match(const T* p, int iDaughter_p);
63  std::vector<std::vector<std::pair<int, std::string>>> m_hierarchy;
64 
66  bool m_isInitOK;
67 
68  public:
70  const static DecayDescriptor& s_NULL;
71 
73  operator DecayDescriptor* ()
74  {
75  return m_isNULL ? nullptr : this;
76  }
79 
81  DecayDescriptor(const DecayDescriptor&) = default;
82 
85 
87  std::vector<std::vector<std::pair<int, std::string>>> getHierarchyOfSelected();
88 
90  std::vector<std::vector<std::pair<int, std::string>>> getHierarchyOfSelected(const std::vector<std::pair<int, std::string>>&
91  currentPath);
92 
95  bool init(const std::string& str);
96 
100  bool init(const DecayString& s);
101 
108  int match(const Particle* p) {return match<Particle>(p, -1);}
109 
111  int match(const MCParticle* p) {return match<MCParticle>(p, -1);}
112 
115 
117  void resetMatch();
118 
120  std::vector<const Particle*> getSelectionParticles(const Particle* particle);
124  std::vector<std::string> getSelectionNames();
127  std::vector<int> getSelectionPDGCodes();
130  {
131  return &m_mother;
132  }
134  int getNDaughters() const
135  {
136  return m_daughters.size();
137  }
139  const DecayDescriptor* getDaughter(int i) const
140  {
141  return (i < getNDaughters()) ? &(m_daughters[i]) : nullptr;
142  }
144  int getProperty() const
145  {
146  return m_properties;
147  }
150  {
151  return (m_properties & Particle::PropertyFlags::c_IsIgnoreRadiatedPhotons) > 0;
152  }
154  bool isIgnoreIntermediate() const
155  {
156  return (m_properties & Particle::PropertyFlags::c_IsIgnoreIntermediate) > 0;
157  }
159  bool isIgnoreMassive() const
160  {
161  return (m_properties & Particle::PropertyFlags::c_IsIgnoreMassive) > 0;
162  }
164  bool isIgnoreNeutrino() const
165  {
166  return (m_properties & Particle::PropertyFlags::c_IsIgnoreNeutrino) > 0;
167  }
169  bool isIgnoreGamma() const
170  {
171  return (m_properties & Particle::PropertyFlags::c_IsIgnoreGamma) > 0;
172  }
174  bool isIgnoreBrems() const
175  {
176  return (m_properties & Particle::PropertyFlags::c_IsIgnoreBrems) > 0;
177  }
178 
180  bool isSelfConjugated() const;
181 
183  bool isInitOK() const
184  {
185  return m_isInitOK;
186  }
187 
191  std::vector<const Particle*>& selparticles,
192  std::vector<std::string>& selnames);
193  };
194 
196 }
197 
Represents a particle in the DecayDescriptor.
The DecayDescriptor stores information about a decay tree or parts of a decay tree.
DecayDescriptor()
Default ctor.
bool isIgnoreBrems() const
Check if added Brems gammas shall be ignored.
bool isIgnoreRadiatedPhotons() const
Check if additional radiated photons shall be ignored.
bool init(const std::string &str)
Initialise the DecayDescriptor from given string.
bool isSelfConjugated() const
Is the decay or the particle self conjugated.
const DecayDescriptorParticle * getMother() const
return mother.
DecayDescriptorParticle m_mother
Mother of the decay ('left side').
bool m_isInitOK
Is this object initialized correctly?
bool getSelectionParticlesAndNames(const Particle *particle, std::vector< const Particle * > &selparticles, std::vector< std::string > &selnames)
Takes as input argument a (reconstructed) Particle, tries to match with this DecayDescriptorElement a...
bool isIgnoreNeutrino() const
Check if missing neutrinos shall be ignored.
std::vector< int > getSelectionPDGCodes()
Return list of PDG codes of selected particles.
bool isInitOK() const
Check if the object initialized correctly.
int match(const MCParticle *p)
See match(const Particle* p).
int match(const Particle *p)
Check if the DecayDescriptor matches with the given Particle.
int getMatchedDaughter()
Particle daughter ID set by previous call of match(const Particle*) function.
int getNDaughters() const
return number of direct daughters.
int getProperty() const
return property of the particle.
void resetMatch()
Reset results from previous call of the match() function.
bool m_isNULL
Is this the NULL object?
bool isIgnoreMassive() const
Check if missing massive final state particles shall be ignored.
const DecayDescriptor * getDaughter(int i) const
return i-th daughter (0 based index).
std::vector< std::vector< std::pair< int, std::string > > > getHierarchyOfSelected()
Function to get hierarchy of selected particles and their names (for python use)
bool isIgnoreGamma() const
Check if missing gammas shall be ignored.
static const DecayDescriptor & s_NULL
Singleton object representing NULL.
std::vector< std::string > getSelectionNames()
Return list of human readable names of selected particles.
std::vector< DecayDescriptor > m_daughters
Direct daughters of the decaying particle.
std::vector< std::vector< std::pair< int, std::string > > > m_hierarchy
Collection of hierarchy paths of selected particles.
std::vector< const Particle * > getSelectionParticles(const Particle *particle)
Get a vector of pointers with selected daughters in the decay tree.
int m_iDaughter_p
ID of the Daughter Particle* matched to this DecayDescriptor.
int match(const T *p, int iDaughter_p)
Internally called by match(Particle*) and match(MCParticle*) function.
bool isIgnoreIntermediate() const
Check if intermediate resonances/particles shall be ignored.
int m_properties
Particle property.
DecayDescriptor(const DecayDescriptor &)=default
Want the default copy ctor.
DecayDescriptor & operator=(const DecayDescriptor &)=default
Want the default assignment operator.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Class to store reconstructed particles.
Definition: Particle.h:75
boost::variant< boost::recursive_wrapper< DecayStringDecay >, DecayStringParticle > DecayString
The DecayStringElement can be either a DecayStringDecay or a vector of mother particles.
Definition: DecayString.h:18
Abstract base class for different kinds of events.