Belle II Software  release-08-01-10
Particle.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 <framework/datastore/RelationsObject.h>
12 #include <framework/gearbox/Const.h>
13 #include <mdst/dataobjects/ECLCluster.h>
14 
15 #include <Math/Vector3D.h>
16 #include <Math/Vector4D.h>
17 #include <TMatrixFfwd.h>
18 #include <TMatrixFSymfwd.h>
19 
20 #include <vector>
21 
22 class TClonesArray;
23 
24 namespace Belle2 {
30  // forward declarations
31  class KLMCluster;
32  class Track;
33  class TrackFitResult;
34  class MCParticle;
35  class PIDLikelihood;
36  class V0;
37 
75  class Particle : public RelationsObject {
76 
77  public:
78 
83  c_Undefined = 0,
84  c_Track = 1,
85  c_ECLCluster = 2,
86  c_KLMCluster = 3,
87  c_V0 = 4,
88  c_MCParticle = 5,
89  c_Composite = 6,
90  c_NoMDSTSource = 7
91  };
92 
94  enum EFlavorType {
96  c_Flavored = 1,
97  };
98 
102  enum {c_DimPosition = 3, c_DimMomentum = 4, c_DimMatrix = 7,
103  c_SizeMatrix = c_DimMatrix * (c_DimMatrix + 1) / 2
104  };
105 
110  enum {c_Px, c_Py, c_Pz, c_E, c_X, c_Y, c_Z};
111 
117  c_Ordinary = 0,
127  };
128 
133  Particle();
134 
141  Particle(const ROOT::Math::PxPyPzEVector& momentum, const int pdgCode);
142 
152  Particle(const ROOT::Math::PxPyPzEVector& momentum,
153  const int pdgCode,
154  EFlavorType flavorType,
155  const EParticleSourceObject particleType,
156  const unsigned mdstIndex);
157 
167  Particle(const ROOT::Math::PxPyPzEVector& momentum,
168  const int pdgCode,
169  EFlavorType flavorType,
170  const std::vector<int>& daughterIndices,
171  TClonesArray* arrayPointer = nullptr);
172 
183  Particle(const ROOT::Math::PxPyPzEVector& momentum,
184  const int pdgCode,
185  EFlavorType flavorType,
186  const std::vector<int>& daughterIndices,
187  int properties,
188  TClonesArray* arrayPointer = nullptr);
189 
201  Particle(const ROOT::Math::PxPyPzEVector& momentum,
202  const int pdgCode,
203  EFlavorType flavorType,
204  const std::vector<int>& daughterIndices,
205  int properties,
206  const std::vector<int>& daughterProperties,
207  TClonesArray* arrayPointer = nullptr);
208 
214  Particle(const Track* track,
215  const Const::ChargedStable& chargedStable);
216 
225  Particle(int trackArrayIndex, const TrackFitResult* trackFit,
226  const Const::ChargedStable& chargedStable);
227 
233  explicit Particle(const ECLCluster* eclCluster,
234  const Const::ParticleType& type = Const::photon);
235 
241  explicit Particle(const KLMCluster* klmCluster, const int pdgCode = Const::Klong.getPDGCode());
242 
247  explicit Particle(const MCParticle* MCparticle);
248 
253 
254  public:
255 
256  // setters
257 
262  void setPDGCode(const int pdg)
263  {
264  m_pdgCode = pdg;
265  }
266 
271  void set4Vector(const ROOT::Math::PxPyPzEVector& p4)
272  {
273  m_px = p4.Px();
274  m_py = p4.Py();
275  m_pz = p4.Pz();
276  m_mass = p4.M();
277  }
278 
283  void set4VectorDividingByMomentumScaling(const ROOT::Math::PxPyPzEVector& p4)
284  {
285  m_px = p4.Px() / m_momentumScale;
286  m_py = p4.Py() / m_momentumScale;
287  m_pz = p4.Pz() / m_momentumScale;
288  m_mass = p4.M();
289  }
290 
295  void setVertex(const ROOT::Math::XYZVector& vertex)
296  {
297  m_x = vertex.X();
298  m_y = vertex.Y();
299  m_z = vertex.Z();
300  };
301 
306  void setMomentumScalingFactor(double momentumScalingFactor)
307  {
308  m_momentumScalingFactor = momentumScalingFactor;
310  }
311 
316  void setMomentumSmearingFactor(double momentumSmearingFactor)
317  {
318  m_momentumSmearingFactor = momentumSmearingFactor;
320  }
321 
326  void setJacobiMatrix(const TMatrixF& jacobiMatrix);
327 
332  void setMomentumVertexErrorMatrix(const TMatrixFSym& errMatrix);
333 
338  void setPValue(double pValue)
339  {
340  m_pValue = pValue;
341  }
342 
346  void setProperty(const int properties)
347  {
348  m_properties = properties;
349  }
350 
358  void updateMomentum(const ROOT::Math::PxPyPzEVector& p4,
359  const ROOT::Math::XYZVector& vertex,
360  const TMatrixFSym& errMatrix,
361  double pValue)
362  {
363  set4Vector(p4);
364  setVertex(vertex);
365  setMomentumVertexErrorMatrix(errMatrix);
366  m_pValue = pValue;
367  }
368 
373  void updateMass(const int pdgCode);
374 
381  void appendDaughter(const Particle* daughter, const bool updateType = true, const int daughterProperty = c_Ordinary);
382 
388  void appendDaughter(int particleIndex, const bool updateType = true)
389  {
390  if (updateType) {
391  // is it a composite particle or fsr corrected?
392  m_particleSource = c_Composite;
393  }
394  m_daughterIndices.push_back(particleIndex);
395  m_daughterProperties.push_back(c_Ordinary);
396  }
397 
403  void removeDaughter(const Particle* daughter, const bool updateType = true);
404 
410  bool replaceDaughter(const Particle* oldDaughter, Particle* newDaughter);
411 
418  bool replaceDaughterRecursively(const Particle* oldDaughter, Particle* newDaughter);
419 
420  // getters
421 
426  int getPDGCode(void) const
427  {
428  return m_pdgCode;
429  }
430 
435  double getCharge(void) const;
436 
442  {
443  return m_flavorType;
444  }
445 
451  {
452  return m_particleSource;
453  }
454 
459  unsigned getMdstArrayIndex(void) const
460  {
461  return m_mdstIndex;
462  }
463 
470  int getProperty() const
471  {
472  return m_properties;
473  }
474 
479  double getMass() const
480  {
481  return m_mass;
482  }
483 
489  //double getMassError() const;
490 
495  double getPDGMass(void) const;
496 
501  double getPDGLifetime() const;
502 
507  double getEnergy() const
508  {
511  }
512 
517  ROOT::Math::PxPyPzEVector get4Vector() const
518  {
519  return ROOT::Math::PxPyPzEVector(m_momentumScale * m_px, m_momentumScale * m_py, m_momentumScale * m_pz, getEnergy());
520  }
521 
526  ROOT::Math::XYZVector getMomentum() const
527  {
528  return m_momentumScale * ROOT::Math::XYZVector(m_px, m_py, m_pz);
529  };
530 
535  double getMomentumMagnitude() const
536  {
537  return m_momentumScale * sqrt(m_px * m_px + m_py * m_py + m_pz * m_pz);
538  };
539 
544  double getP() const
545  {
546  return m_momentumScale * sqrt(m_px * m_px + m_py * m_py + m_pz * m_pz);
547  };
548 
553  double getPx() const
554  {
555  return m_momentumScale * m_px;
556  }
557 
562  double getPy() const
563  {
564  return m_momentumScale * m_py;
565  }
566 
571  double getPz() const
572  {
573  return m_momentumScale * m_pz;
574  }
575 
581  {
582  return m_momentumScale;
583  }
584 
589  ROOT::Math::XYZVector getVertex() const
590  {
591  return ROOT::Math::XYZVector(m_x, m_y, m_z);
592  };
593 
598  double getX() const
599  {
600  return m_x;
601  }
602 
607  double getY() const
608  {
609  return m_y;
610  }
611 
616  double getZ() const
617  {
618  return m_z;
619  }
620 
625  double getPValue() const
626  {
627  return m_pValue;
628  }
629 
634  TMatrixFSym getMomentumVertexErrorMatrix() const;
635 
640  TMatrixFSym getMomentumErrorMatrix() const;
641 
646  TMatrixFSym getVertexErrorMatrix() const;
647 
657  double getCosHelicity(const Particle* mother = nullptr) const;
658 
665  double getCosHelicityDaughter(unsigned iDaughter, unsigned iGrandDaughter = 0) const;
666 
672  double getAcoplanarity() const;
673 
674 
679  int getMdstSource() const;
680 
685  unsigned getNDaughters(void) const
686  {
687  return m_daughterIndices.size();
688  }
689 
694  const std::vector<int>& getDaughterIndices() const
695  {
696  return m_daughterIndices;
697  }
698 
703  const std::vector<int>& getDaughterProperties() const
704  {
705  return m_daughterProperties;
706  }
707 
713  const Particle* getDaughter(unsigned i) const;
714 
724  bool forEachDaughter(const std::function<bool(const Particle*)>& function,
725  bool recursive = true, bool includeSelf = true) const;
726 
731  std::vector<Belle2::Particle*> getDaughters() const;
732  //Need namespace qualifier because ROOT CINT has troubles otherwise
733 
738  std::vector<const Belle2::Particle*> getFinalStateDaughters() const;
739  //Need namespace qualifier because ROOT CINT has troubles otherwise
740 
745  std::vector<const Belle2::Particle*> getAllDaughters() const;
746  //Need namespace qualifier because ROOT CINT has troubles otherwise
747 
754  std::vector<int> getMdstArrayIndices(EParticleSourceObject type) const;
755 
761  bool overlapsWith(const Particle* oParticle) const;
762 
784  bool isCopyOf(const Particle* oParticle, bool doDetailedComparison = false) const;
785 
791  const Track* getTrack() const;
792 
798  const TrackFitResult* getTrackFitResult() const;
799 
806  const V0* getV0() const;
807 
816  const PIDLikelihood* getPIDLikelihood() const;
817 
824  const ECLCluster* getECLCluster() const;
825 
831  double getECLClusterEnergy() const;
832 
839  const KLMCluster* getKLMCluster() const;
840 
848  const MCParticle* getMCParticle() const;
849 
851  std::string getName() const override;
852 
854  std::string getInfoHTML() const override;
855 
859  void print() const;
860 
862  std::vector<std::string> getExtraInfoNames() const;
863 
867  void removeExtraInfo();
868 
873  double getExtraInfo(const std::string& name) const;
874 
876  bool hasExtraInfo(const std::string& name) const;
877 
879  int getExtraInfoMap() const
880  {
881  return !m_extraInfo.empty() ? m_extraInfo[0] : -1;
882  }
883 
885  unsigned int getExtraInfoSize() const
886  {
887  return m_extraInfo.size();
888  }
889 
893  void writeExtraInfo(const std::string& name, const double value);
894 
899  void setExtraInfo(const std::string& name, double value);
900 
905  void addExtraInfo(const std::string& name, double value);
906 
911  TClonesArray* getArrayPointer() const
912  {
913  if (!m_arrayPointer)
915  return m_arrayPointer;
916  }
917 
924  {
925  return std::abs(m_pdgCodeUsedForFit);
926  }
927 
934  {
935  return std::abs(m_pdgCodeUsedForFit) == std::abs(m_pdgCode);
936  }
937 
942  bool isMostLikely() const;
943 
948  std::pair<Const::ChargedStable, const TrackFitResult*> getMostLikelyTrackFitResult() const;
949 
954  bool isMostLikelyTrackFitResult() const;
955 
960  {
961  const int pdg = abs(getPDGCode());
962  if ((pdg == Const::photon.getPDGCode())
963  or (pdg == Const::electron.getPDGCode())
964  or (pdg == Const::muon.getPDGCode())
965  or (pdg == Const::pion.getPDGCode())
966  or (pdg == Const::kaon.getPDGCode())
967  or (pdg == Const::proton.getPDGCode())
968  or (pdg == Const::deuteron.getPDGCode())) {
970  } else if ((pdg == Const::Klong.getPDGCode())
971  or (pdg == Const::neutron.getPDGCode())) {
973  } else {
975  }
976  }
977 
985  const Particle* getParticleFromGeneralizedIndexString(const std::string& generalizedIndex) const;
986 
990  void updateJacobiMatrix();
991 
998  void fillFSPDaughters(std::vector<const Belle2::Particle*>& fspDaughters) const;
999 
1006  void fillAllDaughters(std::vector<const Belle2::Particle*>& allDaughters) const;
1007 
1008 
1009  private:
1010 
1011  // persistent data members
1014  double m_mass;
1015  double m_px;
1016  double m_py;
1017  double m_pz;
1018  double m_momentumScale = 1.0;
1021  double m_x;
1022  double m_y;
1023  double m_z;
1024  double m_errMatrix[c_SizeMatrix] = {};
1025  double m_jacobiMatrix[c_SizeMatrix] = {};
1026  double m_pValue;
1027  std::vector<int> m_daughterIndices;
1030  unsigned m_mdstIndex;
1032  std::vector<int> m_daughterProperties;
1041  int m_identifier = -1;
1042 
1048  std::vector<double> m_extraInfo;
1049 
1050  // transient data members
1057  mutable TClonesArray* m_arrayPointer;
1058 
1059  // private methods
1063  void setMomentumPositionErrorMatrix(const TrackFitResult* trackFit);
1064 
1069  void resetErrorMatrix();
1070 
1075  void resetJacobiMatrix();
1076 
1081  void storeErrorMatrix(const TMatrixFSym& errMatrix);
1082 
1087  void storeJacobiMatrix(const TMatrixF& jacobiMatrix);
1088 
1096  // TODO: this can be optimized for speed
1097  void fillDecayChain(std::vector<int>& decayChain) const;
1098 
1102  void setFlavorType();
1103 
1107  void setMdstArrayIndex(const int arrayIndex);
1108 
1114  int generatePDGCodeFromCharge(const int chargedSign, const Const::ChargedStable& chargedStable);
1115 
1117  // v8: added identifier, changed getMdstSource
1118  // v9: added m_pdgCodeUsedForFit
1119  // v10: added m_properties
1120  // v11: added m_daughterProperties
1121  // v12: renamed EParticleType m_particleType to EParticleSourceObject m_particleSource
1122  // v13: added m_momentumScale
1123  // v14: added m_jacobiMatrix
1124  // v15: added m_momentumScalingFactor and m_momentumSmearingFactor
1125  // v16: use double precision for private members
1126 
1127  friend class ParticleSubset;
1128  };
1129 
1131 } // end namespace Belle2
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:580
The ParticleType class for identifying different particle types.
Definition: Const.h:399
static const ParticleType neutron
neutron particle
Definition: Const.h:666
static const ChargedStable muon
muon particle
Definition: Const.h:651
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:669
static const ChargedStable proton
proton particle
Definition: Const.h:654
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:653
static const ParticleType photon
photon particle
Definition: Const.h:664
static const ChargedStable electron
electron particle
Definition: Const.h:650
static const ChargedStable deuteron
deuteron particle
Definition: Const.h:655
ECL cluster data.
Definition: ECLCluster.h:27
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition: ECLCluster.h:31
@ c_nPhotons
CR is split into n photons (N1)
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
@ c_none
None as initializer.
KLM cluster data.
Definition: KLMCluster.h:28
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Class to collect log likelihoods from TOP, ARICH, dEdx, ECL and KLM aimed for output to mdst includes...
Definition: PIDLikelihood.h:26
Specialised SelectSubset<Particle> that also fixes daughter indices and all ParticleLists.
Class to store reconstructed particles.
Definition: Particle.h:75
void setProperty(const int properties)
sets m_properties
Definition: Particle.h:346
void removeExtraInfo()
Remove all stored extra info fields.
Definition: Particle.cc:1285
double getPx() const
Returns x component of momentum.
Definition: Particle.h:553
void updateJacobiMatrix()
Propagate the photon energy scaling to jacobian elements that were calculated using energy.
Definition: Particle.cc:257
TMatrixFSym getVertexErrorMatrix() const
Returns the 3x3 position error sub-matrix.
Definition: Particle.cc:451
bool isMostLikely() const
Returns true if the (track-based) particle is created with its most likely mass hypothesis based on P...
Definition: Particle.cc:1392
const KLMCluster * getKLMCluster() const
Returns the pointer to the KLMCluster object that was used to create this Particle (ParticleType == c...
Definition: Particle.cc:930
double m_momentumSmearingFactor
momentum smearing factor
Definition: Particle.h:1020
bool replaceDaughter(const Particle *oldDaughter, Particle *newDaughter)
Replace index of given daughter with new daughter, return true if a replacement is made.
Definition: Particle.cc:709
const Track * getTrack() const
Returns the pointer to the Track object that was used to create this Particle (ParticleType == c_Trac...
Definition: Particle.cc:849
Particle(const ROOT::Math::PxPyPzEVector &momentum, const int pdgCode, EFlavorType flavorType, const std::vector< int > &daughterIndices, int properties, const std::vector< int > &daughterProperties, TClonesArray *arrayPointer=nullptr)
Constructor for composite particles.
void setMomentumPositionErrorMatrix(const TrackFitResult *trackFit)
Sets the momentum, position and error matrix for this particle (created from charged Track)
Definition: Particle.cc:1005
void storeErrorMatrix(const TMatrixFSym &errMatrix)
Stores 7x7 error matrix into private member m_errMatrix.
Definition: Particle.cc:1095
bool forEachDaughter(const std::function< bool(const Particle *)> &function, bool recursive=true, bool includeSelf=true) const
Apply a function to all daughters of this particle.
Definition: Particle.cc:1359
std::vector< const Belle2::Particle * > getFinalStateDaughters() const
Returns a vector of pointers to Final State daughter particles.
Definition: Particle.cc:653
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
Definition: Particle.cc:946
void appendDaughter(const Particle *daughter, const bool updateType=true, const int daughterProperty=c_Ordinary)
Appends index of daughter to daughters index array.
Definition: Particle.cc:680
const ECLCluster * getECLCluster() const
Returns the pointer to the ECLCluster object that was used to create this Particle (if ParticleType =...
Definition: Particle.cc:895
void setExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1318
EParticleSourceObject
particle source enumerators
Definition: Particle.h:82
double m_pValue
chi^2 probability of the fit.
Definition: Particle.h:1026
double m_pz
momentum component z
Definition: Particle.h:1017
int m_pdgCodeUsedForFit
PDG code used for the track fit.
Definition: Particle.h:1013
double m_momentumScale
effective momentum scale factor
Definition: Particle.h:1018
std::string getName() const override
Return name of this particle.
Definition: Particle.cc:1169
bool overlapsWith(const Particle *oParticle) const
Returns true if final state ancestors of oParticle overlap.
Definition: Particle.cc:741
double getCosHelicityDaughter(unsigned iDaughter, unsigned iGrandDaughter=0) const
Returns cosine of the helicity angle of the given daughter defined by given grand daughter.
Definition: Particle.cc:506
TClonesArray * getArrayPointer() const
Returns the pointer to the store array which holds the daughter particles.
Definition: Particle.h:911
double getPz() const
Returns z component of momentum.
Definition: Particle.h:571
std::vector< std::string > getExtraInfoNames() const
get a list of the extra info names
Definition: Particle.cc:1179
double getX() const
Returns x component of vertex position.
Definition: Particle.h:598
bool isCopyOf(const Particle *oParticle, bool doDetailedComparison=false) const
Returns true if this Particle and oParticle are copies of each other.
Definition: Particle.cc:756
void writeExtraInfo(const std::string &name, const double value)
Sets the user defined extraInfo.
Definition: Particle.cc:1309
void setVertex(const ROOT::Math::XYZVector &vertex)
Sets position (decay vertex)
Definition: Particle.h:295
std::vector< int > getMdstArrayIndices(EParticleSourceObject type) const
Returns a vector of StoreArray indices of given MDST dataobjects.
Definition: Particle.cc:669
double m_jacobiMatrix[c_SizeMatrix]
error matrix (1D representation)
Definition: Particle.h:1025
void resetJacobiMatrix()
Resets 4x6 error matrix All elements are set to 0.0.
Definition: Particle.cc:1089
void set4VectorDividingByMomentumScaling(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector dividing by the momentum scaling factor.
Definition: Particle.h:283
void storeJacobiMatrix(const TMatrixF &jacobiMatrix)
Stores 4x6 Jacobi matrix into private member m_jacobiMatrix.
Definition: Particle.cc:1106
bool isMostLikelyTrackFitResult() const
Returns true if the (track-based) particle is created with its most likely mass hypothesis based on T...
Definition: Particle.cc:1417
std::vector< double > m_extraInfo
Stores associated user defined values.
Definition: Particle.h:1048
double getEffectiveMomentumScale() const
Returns effective momentum scale which is the product of the momentum scaling and smearing factors.
Definition: Particle.h:580
unsigned m_mdstIndex
0-based index of MDST store array object
Definition: Particle.h:1030
double getPValue() const
Returns chi^2 probability of fit if done or -1.
Definition: Particle.h:625
Particle(const ROOT::Math::PxPyPzEVector &momentum, const int pdgCode)
Constructor from a Lorentz vector and PDG code.
void setFlavorType()
sets m_flavorType using m_pdgCode
Definition: Particle.cc:1154
double getEnergy() const
Returns total energy.
Definition: Particle.h:507
void setMomentumScalingFactor(double momentumScalingFactor)
Sets momentum scaling.
Definition: Particle.h:306
int getMdstSource() const
Returns unique identifier of final state particle (needed in particle combiner)
Definition: Particle.cc:371
ROOT::Math::XYZVector getVertex() const
Returns vertex position (POCA for charged, IP for neutral FS particles)
Definition: Particle.h:589
void fillAllDaughters(std::vector< const Belle2::Particle * > &allDaughters) const
Fill all generations' daughters into a vector.
Definition: Particle.cc:1131
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1267
std::string getInfoHTML() const override
Return a short summary of this object's contents in HTML format.
Definition: Particle.cc:1192
void setMdstArrayIndex(const int arrayIndex)
set mdst array index
Definition: Particle.cc:350
const Particle * getParticleFromGeneralizedIndexString(const std::string &generalizedIndex) const
Explores the decay tree of the particle and returns the (grand^n)daughter identified by a generalized...
Definition: Particle.cc:960
unsigned getMdstArrayIndex(void) const
Returns 0-based index of MDST store array object (0 for composite particles)
Definition: Particle.h:459
double m_py
momentum component y
Definition: Particle.h:1016
int getPDGCodeUsedForFit() const
Return the always positive PDG code which was used for the track fit (if there was a track fit) of th...
Definition: Particle.h:923
double m_x
position component x
Definition: Particle.h:1021
double m_px
momentum component x
Definition: Particle.h:1015
std::pair< Const::ChargedStable, const TrackFitResult * > getMostLikelyTrackFitResult() const
For a (track-based) particle, returns the charged stable mass hypothesis associated to the most proba...
Definition: Particle.cc:1399
EFlavorType getFlavorType() const
Returns flavor type of the decay (for FS particles: flavor type of particle)
Definition: Particle.h:441
~Particle()
Destructor.
void fillFSPDaughters(std::vector< const Belle2::Particle * > &fspDaughters) const
Fill final state particle daughters into a vector.
Definition: Particle.cc:1118
const PIDLikelihood * getPIDLikelihood() const
Returns the pointer to the PIDLikelihood object that is related to the Track, which was used to creat...
Definition: Particle.cc:875
const V0 * getV0() const
Returns the pointer to the V0 object that was used to create this Particle (if ParticleType == c_V0).
Definition: Particle.cc:884
TClonesArray * m_arrayPointer
Internal pointer to DataStore array containing the daughters of this particle.
Definition: Particle.h:1057
EParticleSourceObject m_particleSource
(mdst) source of particle
Definition: Particle.h:1029
void fillDecayChain(std::vector< int > &decayChain) const
Fill vector with (PDGCode, MdstSource) pairs for the entire decay chain.
Definition: Particle.cc:1144
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:426
double getPy() const
Returns y component of momentum.
Definition: Particle.h:562
bool wasExactFitHypothesisUsed() const
Returns true if the type represented by this Particle object was used use as a mass hypothesis during...
Definition: Particle.h:933
int getProperty() const
Returns particle property as a bit pattern The values are defined in the PropertyFlags enum and descr...
Definition: Particle.h:470
double getCosHelicity(const Particle *mother=nullptr) const
Returns cosine of the helicity angle The helicity angle is defined in the rest frame of the particle ...
Definition: Particle.cc:463
double m_momentumScalingFactor
momentum scaling factor
Definition: Particle.h:1019
unsigned getNDaughters(void) const
Returns number of daughter particles.
Definition: Particle.h:685
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:641
const std::vector< int > & getDaughterIndices() const
Returns a vector of store array indices of daughter particles.
Definition: Particle.h:694
double getZ() const
Returns z component of vertex position.
Definition: Particle.h:616
Particle(const ROOT::Math::PxPyPzEVector &momentum, const int pdgCode, EFlavorType flavorType, const std::vector< int > &daughterIndices, int properties, TClonesArray *arrayPointer=nullptr)
Constructor for composite particles.
void resetErrorMatrix()
Resets 7x7 error matrix All elements are set to 0.0.
Definition: Particle.cc:1083
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:608
void setMomentumSmearingFactor(double momentumSmearingFactor)
Sets momentum smearing.
Definition: Particle.h:316
EFlavorType m_flavorType
flavor type.
Definition: Particle.h:1028
int generatePDGCodeFromCharge(const int chargedSign, const Const::ChargedStable &chargedStable)
Generate the PDG code with correct sign, using the charge.
Definition: Particle.cc:1383
double getCharge(void) const
Returns particle charge.
Definition: Particle.cc:626
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:517
void appendDaughter(int particleIndex, const bool updateType=true)
Appends index of daughter to daughters index array.
Definition: Particle.h:388
Particle()
Default constructor.
Definition: Particle.cc:45
TMatrixFSym getMomentumErrorMatrix() const
Returns the 4x4 momentum error matrix.
Definition: Particle.cc:439
double m_mass
particle (invariant) mass
Definition: Particle.h:1014
std::vector< const Belle2::Particle * > getAllDaughters() const
Returns a vector of pointers to all generations' daughter particles.
Definition: Particle.cc:661
double getAcoplanarity() const
Returns acoplanarity angle defined as the angle between the decay planes of the grand daughters in th...
Definition: Particle.cc:537
std::vector< int > m_daughterProperties
daughter particle properties
Definition: Particle.h:1032
double m_errMatrix[c_SizeMatrix]
error matrix (1D representation)
Definition: Particle.h:1024
std::vector< int > m_daughterIndices
daughter particle indices
Definition: Particle.h:1027
double getY() const
Returns y component of vertex position.
Definition: Particle.h:607
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1337
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets Lorentz vector.
Definition: Particle.h:271
void setJacobiMatrix(const TMatrixF &jacobiMatrix)
Sets 4x6 jacobi matrix.
Definition: Particle.cc:411
int getExtraInfoMap() const
Return the id of the associated ParticleExtraInfoMap or -1 if no map is set.
Definition: Particle.h:879
double getPDGLifetime() const
Returns particle nominal lifetime.
Definition: Particle.cc:617
void setMomentumVertexErrorMatrix(const TMatrixFSym &errMatrix)
Sets 7x7 error matrix.
Definition: Particle.cc:397
ROOT::Math::XYZVector getMomentum() const
Returns momentum vector.
Definition: Particle.h:526
void setPValue(double pValue)
Sets chi^2 probability of fit.
Definition: Particle.h:338
EParticleSourceObject getParticleSource() const
Returns particle source as defined with enum EParticleSourceObject.
Definition: Particle.h:450
void updateMass(const int pdgCode)
Updates particle mass with the mass of the particle corresponding to the given PDG.
Definition: Particle.cc:601
double m_z
position component z
Definition: Particle.h:1023
void print() const
Prints the contents of a Particle object to standard output.
Definition: Particle.cc:1174
void updateMomentum(const ROOT::Math::PxPyPzEVector &p4, const ROOT::Math::XYZVector &vertex, const TMatrixFSym &errMatrix, double pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
Definition: Particle.h:358
EFlavorType
describes flavor type, see getFlavorType().
Definition: Particle.h:94
@ c_Unflavored
Is its own antiparticle or we don't know whether it is a particle/antiparticle.
Definition: Particle.h:95
@ c_Flavored
Is either particle or antiparticle.
Definition: Particle.h:96
int m_pdgCode
PDG code.
Definition: Particle.h:1012
TMatrixFSym getMomentumVertexErrorMatrix() const
Returns 7x7 error matrix.
Definition: Particle.cc:424
PropertyFlags
Flags that describe the particle property, which are used in the MC matching.
Definition: Particle.h:116
@ c_IsIgnoreNeutrino
Is the particle MC matched with the ignore missing neutrino flag set?
Definition: Particle.h:122
@ c_IsIgnoreRadiatedPhotons
Is the particle MC matched with the ignore radiated photon flag set?
Definition: Particle.h:119
@ c_IsIgnoreGamma
Is the particle MC matched with the ignore missing gamma flag set?
Definition: Particle.h:123
@ c_IsUnspecified
Ordinary particles.
Definition: Particle.h:118
@ c_IsIgnoreBrems
Is the particle MC matched with the ignore added Brems gamma flag set?
Definition: Particle.h:124
@ c_IsIgnoreDecayInFlight
Is the particle MC matched with the ignore DecayInFlight flag set?
Definition: Particle.h:126
@ c_IsIgnoreMisID
Is the particle MC matched with the ignore MisID flag set?
Definition: Particle.h:125
@ c_IsIgnoreIntermediate
Is the particle MC matched with the ignore intermediate resonances flag set?
Definition: Particle.h:120
@ c_IsIgnoreMassive
Is the particle MC matched with the ignore missing massive particle flag set?
Definition: Particle.h:121
double getECLClusterEnergy() const
Returns the energy of the ECLCluster for the particle.
Definition: Particle.cc:923
const TrackFitResult * getTrackFitResult() const
Returns the pointer to the TrackFitResult that was used to create this Particle (ParticleType == c_Tr...
Definition: Particle.cc:858
double getMomentumMagnitude() const
Returns momentum magnitude.
Definition: Particle.h:535
Particle(const ROOT::Math::PxPyPzEVector &momentum, const int pdgCode, EFlavorType flavorType, const EParticleSourceObject particleType, const unsigned mdstIndex)
Constructor for final state particles.
void removeDaughter(const Particle *daughter, const bool updateType=true)
Removes index of daughter from daughters index array.
Definition: Particle.cc:692
double getP() const
Returns momentum magnitude (same as getMomentumMagnitude but with shorter name)
Definition: Particle.h:544
double m_y
position component y
Definition: Particle.h:1022
void setPDGCode(const int pdg)
Sets PDG code.
Definition: Particle.h:262
unsigned int getExtraInfoSize() const
Return the size of the extra info array.
Definition: Particle.h:885
int m_identifier
Identifier that can be used to identify whether the particle is unique or is a copy or representation...
Definition: Particle.h:1041
const std::vector< int > & getDaughterProperties() const
Returns a vector of properties of daughter particles.
Definition: Particle.h:703
int m_properties
particle property
Definition: Particle.h:1031
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:635
ClassDefOverride(Particle, 16)
Class to store reconstructed particles.
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1290
ECLCluster::EHypothesisBit getECLClusterEHypothesisBit() const
Returns the ECLCluster EHypothesisBit for this Particle.
Definition: Particle.h:959
Particle(const ROOT::Math::PxPyPzEVector &momentum, const int pdgCode, EFlavorType flavorType, const std::vector< int > &daughterIndices, TClonesArray *arrayPointer=nullptr)
Constructor for composite particles.
bool replaceDaughterRecursively(const Particle *oldDaughter, Particle *newDaughter)
Apply replaceDaughter to all Particles in the decay tree by looping recursively through it,...
Definition: Particle.cc:727
double getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:479
Defines interface for accessing relations of objects in StoreArray.
TClonesArray * getArrayPointer() const
Returns the pointer to the raw DataStore array holding this object (protected since these arrays are ...
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
Object holding information for V0s.
Definition: V0.h:34
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.