CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauPiZero.cc
Go to the documentation of this file.
3 
4 namespace reco {
5 
7 {
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) ++nGammas;
12  }
13  return nGammas;
14 }
15 
17 {
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) ++nElectrons;
22  }
23  return nElectrons;
24 }
25 
27 {
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 
39 {
40  double maxDEta = 0;
41  size_t nDaughters = numberOfDaughters();
42  for(size_t i = 0; i < nDaughters; ++i) {
43  double dEta = std::fabs(eta() - daughter(i)->eta());
44  if(dEta > maxDEta)
45  maxDEta = dEta;
46  }
47  return maxDEta;
48 }
49 
51  return algoName_;
52 }
53 
55  return (algoName_ == algo);
56 }
57 
58 namespace {
59 std::ostream& operator<<(std::ostream& out, const reco::Candidate::LorentzVector& p4)
60 {
61  out << "(mass/pt/eta/phi) (" << std::setiosflags(std::ios::fixed) << std::setprecision(2)
62  << p4.mass() << "/" << std::setprecision(1) << p4.pt() << "/" << std::setprecision(2) << p4.eta()
63  << "/" << std::setprecision(2) << p4.phi() << ")";
64  return out;
65 }
66 }
67 
68 void RecoTauPiZero::print(std::ostream& out) const {
69  if (!out) return;
70 
71  out << "RecoTauPiZero: " << this->p4() <<
72  " nDaughters: " << this->numberOfDaughters() <<
73  " (gamma/e) (" << this->numberOfGammas() << "/" << this->numberOfElectrons() << ")" <<
74  " maxDeltaPhi: " << std::setprecision(3) << maxDeltaPhi() <<
75  " maxDeltaEta: " << std::setprecision(3) << maxDeltaEta() <<
76  " algo: " << algo() <<
77  std::endl;
78 
79  for(size_t i = 0; i < this->numberOfDaughters(); ++i)
80  {
81  out << "--- daughter " << i << ": " << daughterPtr(i)->p4() <<
82  " key: " << daughterPtr(i).key() << 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 i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
unsigned int nGammas(const GenJet &jet)
key_type key() const
Definition: Ptr.h:169
CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
void print(std::ostream &out=std::cout) const
double maxDeltaPhi() const
Maximum DeltaPhi between a constituent and the four vector.
double maxDeltaEta() const
Maxmum DeltaEta between a constituent and the four vector.
virtual size_t numberOfDaughters() const
number of daughters
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
virtual float eta() const
momentum pseudorapidity
double p4[4]
Definition: TauolaWrapper.h:92
size_t numberOfElectrons() const
Number of electron constituents.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PiZeroAlgorithm algoName_
Definition: RecoTauPiZero.h:68
tuple out
Definition: dbtoconf.py:99
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
size_t numberOfGammas() const
Number of PFGamma constituents.
Definition: RecoTauPiZero.cc:6
PiZeroAlgorithm algo() const
Algorithm that built this piZero.
bool algoIs(PiZeroAlgorithm algo) const
Check whether a given algo produced this pi zero.
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...