CMS 3D CMS Logo

RecoTauPiZero.cc
Go to the documentation of this file.
4 
5 namespace reco {
6 
8 {
9  size_t nGammas = 0;
10  size_t nDaughters = numberOfDaughters();
11  for(size_t i = 0; i < nDaughters; ++i) {
12  if(daughter(i)->pdgId() == 22) ++nGammas;
13  }
14  return nGammas;
15 }
16 
18 {
19  size_t nElectrons = 0;
20  size_t nDaughters = numberOfDaughters();
21  for(size_t i = 0; i < nDaughters; ++i) {
22  if(std::abs(daughter(i)->pdgId()) == 11) ++nElectrons;
23  }
24  return nElectrons;
25 }
26 
28 {
29  double maxDPhi = 0;
30  size_t nDaughters = numberOfDaughters();
31  for(size_t i = 0; i < nDaughters; ++i) {
32  double dPhi = std::fabs(deltaPhi(*this, *daughter(i)));
33  if(dPhi > maxDPhi)
34  maxDPhi = dPhi;
35  }
36  return maxDPhi;
37 }
38 
40 {
41  double maxDEta = 0;
42  size_t nDaughters = numberOfDaughters();
43  for(size_t i = 0; i < nDaughters; ++i) {
44  double dEta = std::fabs(eta() - daughter(i)->eta());
45  if(dEta > maxDEta)
46  maxDEta = dEta;
47  }
48  return maxDEta;
49 }
50 
52  return algoName_;
53 }
54 
56  return (algoName_ == algo);
57 }
58 
59 namespace
60 {
61  std::string getPFCandidateType(reco::PFCandidate::ParticleType pfCandidateType)
62  {
63  if ( pfCandidateType == reco::PFCandidate::X ) return "undefined";
64  else if ( pfCandidateType == reco::PFCandidate::h ) return "PFChargedHadron";
65  else if ( pfCandidateType == reco::PFCandidate::e ) return "PFElectron";
66  else if ( pfCandidateType == reco::PFCandidate::mu ) return "PFMuon";
67  else if ( pfCandidateType == reco::PFCandidate::gamma ) return "PFGamma";
68  else if ( pfCandidateType == reco::PFCandidate::h0 ) return "PFNeutralHadron";
69  else if ( pfCandidateType == reco::PFCandidate::h_HF ) return "HF_had";
70  else if ( pfCandidateType == reco::PFCandidate::egamma_HF ) return "HF_em";
71  else assert(0);
72  }
73 }
74 
75 void RecoTauPiZero::print(std::ostream& stream) const
76 {
77  std::cout << "Pt = " << this->pt() << ", eta = " << this->eta() << ", phi = " << this->phi() << std::endl;
78  size_t numDaughters = this->numberOfDaughters();
79  for ( size_t iDaughter = 0; iDaughter < numDaughters; ++iDaughter ) {
80  const reco::PFCandidate* daughter = dynamic_cast<const reco::PFCandidate*>(this->daughterPtr(iDaughter).get());
81  std::cout << " daughter #" << iDaughter << " (" << getPFCandidateType(daughter->particleId()) << "):"
82  << " Pt = " << daughter->pt() << ", eta = " << daughter->eta() << ", phi = " << daughter->phi() << std::endl;
83  }
84 }
85 
86 std::ostream& operator<<(std::ostream& out, const reco::RecoTauPiZero& piZero)
87 {
88  if(!out) return out;
89  piZero.print(out);
90  return out;
91 }
92 }
int pdgId() const final
PDG identifier.
unsigned int nGammas(const GenJet &jet)
ParticleType
particle types
Definition: PFCandidate.h:44
double eta() const final
momentum pseudorapidity
CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
double pt() const final
transverse momentum
void print(std::ostream &out=std::cout) const
size_t numberOfDaughters() const override
number of daughters
double maxDeltaPhi() const
Maximum DeltaPhi between a constituent and the four vector.
double maxDeltaEta() const
Maxmum DeltaEta between a constituent and the four vector.
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
size_t numberOfElectrons() const
Number of electron constituents.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PiZeroAlgorithm algoName_
Definition: RecoTauPiZero.h:90
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
fixed size matrix
size_t numberOfGammas() const
Number of PFGamma constituents.
Definition: RecoTauPiZero.cc:7
virtual ParticleType particleId() const
Definition: PFCandidate.h:373
PiZeroAlgorithm algo() const
Algorithm that built this piZero.
bool algoIs(PiZeroAlgorithm algo) const
Check whether a given algo produced this pi zero.
double phi() const final
momentum azimuthal angle