CMS 3D CMS Logo

RecoTauPiZero.cc
Go to the documentation of this file.
4 
5 namespace reco {
6 
8  size_t nGammas = 0;
9  size_t nDaughters = numberOfDaughters();
10  for (size_t i = 0; i < nDaughters; ++i) {
11  if (daughter(i)->pdgId() == 22)
12  ++nGammas;
13  }
14  return nGammas;
15  }
16 
18  size_t nElectrons = 0;
19  size_t nDaughters = numberOfDaughters();
20  for (size_t i = 0; i < nDaughters; ++i) {
21  if (std::abs(daughter(i)->pdgId()) == 11)
22  ++nElectrons;
23  }
24  return nElectrons;
25  }
26 
27  double RecoTauPiZero::maxDeltaPhi() const {
28  double maxDPhi = 0;
29  size_t nDaughters = numberOfDaughters();
30  for (size_t i = 0; i < nDaughters; ++i) {
31  double dPhi = std::fabs(deltaPhi(*this, *daughter(i)));
32  if (dPhi > maxDPhi)
33  maxDPhi = dPhi;
34  }
35  return maxDPhi;
36  }
37 
38  double RecoTauPiZero::maxDeltaEta() const {
39  double maxDEta = 0;
40  size_t nDaughters = numberOfDaughters();
41  for (size_t i = 0; i < nDaughters; ++i) {
42  double dEta = std::fabs(eta() - daughter(i)->eta());
43  if (dEta > maxDEta)
44  maxDEta = dEta;
45  }
46  return maxDEta;
47  }
48 
50 
52 
53  namespace {
54  std::string getPFCandidateType(reco::PFCandidate::ParticleType pfCandidateType) {
55  if (pfCandidateType == reco::PFCandidate::X)
56  return "undefined";
57  else if (pfCandidateType == reco::PFCandidate::h)
58  return "PFChargedHadron";
59  else if (pfCandidateType == reco::PFCandidate::e)
60  return "PFElectron";
61  else if (pfCandidateType == reco::PFCandidate::mu)
62  return "PFMuon";
63  else if (pfCandidateType == reco::PFCandidate::gamma)
64  return "PFGamma";
65  else if (pfCandidateType == reco::PFCandidate::h0)
66  return "PFNeutralHadron";
67  else if (pfCandidateType == reco::PFCandidate::h_HF)
68  return "HF_had";
69  else if (pfCandidateType == reco::PFCandidate::egamma_HF)
70  return "HF_em";
71  else
72  assert(0);
73  }
74  } // namespace
75 
76  void RecoTauPiZero::print(std::ostream& stream) const {
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()
83  << std::endl;
84  }
85  }
86 
87  std::ostream& operator<<(std::ostream& out, const reco::RecoTauPiZero& piZero) {
88  if (!out)
89  return out;
90  piZero.print(out);
91  return out;
92  }
93 } // namespace reco
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
int pdgId() const final
PDG identifier.
ParticleType
particle types
Definition: PFCandidate.h:43
double eta() const final
momentum pseudorapidity
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
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:66
size_t numberOfElectrons() const
Number of electron constituents.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PiZeroAlgorithm algoName_
Definition: RecoTauPiZero.h:89
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
size_t numberOfGammas() const
Number of PFGamma constituents.
Definition: RecoTauPiZero.cc:7
virtual ParticleType particleId() const
Definition: PFCandidate.h:366
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