CMS 3D CMS Logo

RecoTauPiZero.h
Go to the documentation of this file.
1 #ifndef DataFormats_TauReco_RecoTauPiZero_h
2 #define DataFormats_TauReco_RecoTauPiZero_h
3 
5 
6 namespace reco {
8  public:
10  // Algorithm where each photon becomes a pi zero
12  kTrivial = 1,
14  kStrips = 3
15  };
16 
18  this->setPdgId(111);
19  }
20 
23  this->setPdgId(111);
24  }
25 
28  const LorentzVector& p4,
29  const Point& vtx = Point(0, 0, 0),
30  int pdgId = 111,
31  int status = 0,
32  bool integerCharge = true,
34  : CompositePtrCandidate(q, p4, vtx, pdgId, status, integerCharge),
36  bendCorrEta_(0.),
37  bendCorrPhi_(0.) {}
38 
41  const PolarLorentzVector& p4,
42  const Point& vtx = Point(0, 0, 0),
43  int pdgId = 111,
44  int status = 0,
45  bool integerCharge = true,
47  : CompositePtrCandidate(q, p4, vtx, pdgId, status, integerCharge),
49  bendCorrEta_(0.),
50  bendCorrPhi_(0.) {}
51 
55  this->setPdgId(111);
56  }
57 
59  ~RecoTauPiZero() override {}
60 
62  size_t numberOfGammas() const;
63 
65  size_t numberOfElectrons() const;
66 
68  double maxDeltaPhi() const;
69 
71  double maxDeltaEta() const;
72 
74  PiZeroAlgorithm algo() const;
75 
77  bool algoIs(PiZeroAlgorithm algo) const;
78 
81  float bendCorrEta() const { return bendCorrEta_; }
82  float bendCorrPhi() const { return bendCorrPhi_; }
85 
86  void print(std::ostream& out = std::cout) const;
87 
88  private:
90 
91  float bendCorrEta_;
92  float bendCorrPhi_;
93  };
94 
95  std::ostream& operator<<(std::ostream& out, const RecoTauPiZero& c);
96 
97 } // namespace reco
98 
99 #endif
PiZeroAlgorithm algo() const
Algorithm that built this piZero.
int Charge
electric charge type
Definition: Candidate.h:34
double maxDeltaEta() const
Maxmum DeltaEta between a constituent and the four vector.
RecoTauPiZero(Charge q, const LorentzVector &p4, const Point &vtx=Point(0, 0, 0), int pdgId=111, int status=0, bool integerCharge=true, PiZeroAlgorithm algoName=kUndefined)
constructor from values
Definition: RecoTauPiZero.h:27
int status() const final
status word
float bendCorrPhi() const
Definition: RecoTauPiZero.h:82
const LorentzVector & p4() const final
four-momentum Lorentz vector
void print(std::ostream &out=std::cout) const
int pdgId() const final
PDG identifier.
RecoTauPiZero(PiZeroAlgorithm algoName)
Definition: RecoTauPiZero.h:21
size_t numberOfGammas() const
Number of PFGamma constituents.
Definition: RecoTauPiZero.cc:7
float bendCorrEta() const
Definition: RecoTauPiZero.h:81
double p() const final
magnitude of momentum vector
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
double maxDeltaPhi() const
Maximum DeltaPhi between a constituent and the four vector.
size_t numberOfElectrons() const
Number of electron constituents.
~RecoTauPiZero() override
destructor
Definition: RecoTauPiZero.h:59
PiZeroAlgorithm algoName_
Definition: RecoTauPiZero.h:89
void setBendCorrEta(float bendCorrEta)
Definition: RecoTauPiZero.h:83
void setBendCorrPhi(float bendCorrPhi)
Definition: RecoTauPiZero.h:84
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
Structure Point Contains parameters of Gaussian fits to DMRs.
bool algoIs(PiZeroAlgorithm algo) const
Check whether a given algo produced this pi zero.
void setPdgId(int pdgId) final
RecoTauPiZero(const Candidate &p, PiZeroAlgorithm algoName=kUndefined)
constructor from a Candidate
Definition: RecoTauPiZero.h:53
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
RecoTauPiZero(Charge q, const PolarLorentzVector &p4, const Point &vtx=Point(0, 0, 0), int pdgId=111, int status=0, bool integerCharge=true, PiZeroAlgorithm algoName=kUndefined)
constructor from values
Definition: RecoTauPiZero.h:40
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38