CMS 3D CMS Logo

PFJet.cc
Go to the documentation of this file.
1 // PFJet.cc
2 // Fedor Ratnikov UMd
3 #include <sstream>
4 #include <typeinfo>
5 
9 
10 //Own header file
12 
13 using namespace reco;
14 
16  const Point& fVertex,
17  const Specific& fSpecific,
18  const Jet::Constituents& fConstituents)
19  : Jet(fP4, fVertex, fConstituents), m_specific(fSpecific) {}
20 
21 PFJet::PFJet(const LorentzVector& fP4, const Point& fVertex, const Specific& fSpecific)
22  : Jet(fP4, fVertex), m_specific(fSpecific) {}
23 
24 PFJet::PFJet(const LorentzVector& fP4, const Specific& fSpecific, const Jet::Constituents& fConstituents)
25  : Jet(fP4, Point(0, 0, 0), fConstituents), m_specific(fSpecific) {}
26 
28  Constituent dau = daughterPtr(fIndex);
29  if (dau.isNonnull() && dau.isAvailable()) {
30  const PFCandidate* pfCandidate = dynamic_cast<const PFCandidate*>(dau.get());
31  if (pfCandidate) {
32  return edm::Ptr<PFCandidate>(dau.id(), pfCandidate, dau.key());
33  } else {
34  throw cms::Exception("Invalid Constituent") << "PFJet constituent is not of PFCandidate type";
35  }
36  } else {
37  return PFCandidatePtr();
38  }
39 }
40 
41 std::vector<reco::PFCandidatePtr> PFJet::getPFConstituents() const {
42  std::vector<PFCandidatePtr> result;
43  for (unsigned i = 0; i < numberOfDaughters(); i++)
44  result.push_back(getPFConstituent(i));
45  return result;
46 }
47 
49  // result will contain chargedMultiplicity() elements
51  result.reserve(chargedMultiplicity());
52  for (unsigned i = 0; i < numberOfDaughters(); i++) {
54  reco::TrackRef trackref = pfcand->trackRef();
55  if (trackref.isNonnull()) {
56  result.push_back(trackref);
57  }
58  }
59 
60  return result;
61 }
62 
63 PFJet* PFJet::clone() const { return new PFJet(*this); }
64 
65 bool PFJet::overlap(const Candidate&) const { return false; }
66 
68  std::ostringstream out;
69  out << Jet::print() // generic jet info
70  << " PFJet specific:" << std::endl
71  << " charged hadron energy/multiplicity: " << chargedHadronEnergy() << '/' << chargedHadronMultiplicity()
72  << std::endl
73  << " neutral hadron energy/multiplicity: " << neutralHadronEnergy() << '/' << neutralHadronMultiplicity()
74  << std::endl
75  << " photon energy/multiplicity: " << photonEnergy() << '/' << photonMultiplicity() << std::endl
76  << " electron energy/multiplicity: " << electronEnergy() << '/' << electronMultiplicity() << std::endl
77  << " muon energy/multiplicity: " << muonEnergy() << '/' << muonMultiplicity() << std::endl
78  << " HF Hadron energy/multiplicity: " << HFHadronEnergy() << '/' << HFHadronMultiplicity() << std::endl
79  << " HF EM particle energy/multiplicity: " << HFEMEnergy() << '/' << HFEMMultiplicity() << std::endl
80  << " charged/neutral hadrons energy: " << chargedHadronEnergy() << '/' << neutralHadronEnergy() << std::endl
81  << " charged/neutral em energy: " << chargedEmEnergy() << '/' << neutralEmEnergy() << std::endl
82  << " charged muon energy: " << chargedMuEnergy() << '/' << std::endl
83  << " charged/neutral multiplicity: " << chargedMultiplicity() << '/' << neutralMultiplicity() << std::endl
84  << " HO energy: " << hoEnergy() << std::endl;
85  out << " PFCandidate constituents:" << std::endl;
86  std::vector<PFCandidatePtr> constituents = getPFConstituents();
87  for (unsigned i = 0; i < constituents.size(); ++i) {
88  if (constituents[i].get()) {
89  out << " #" << i << " " << *(constituents[i]) << std::endl;
90  } else {
91  out << " #" << i << " PFCandidate is not available in the event" << std::endl;
92  }
93  }
94  return out.str();
95 }
96 
97 std::ostream& reco::operator<<(std::ostream& out, const reco::PFJet& jet) {
98  if (out) {
99  out << "PFJet "
100  << "(pt, eta, phi) = " << jet.pt() << "," << jet.eta() << "," << jet.phi()
101  << " (Rch,Rnh,Rgamma,Re,Rmu,RHFHad,RHFEM) = " << jet.chargedHadronEnergyFraction() << ","
102  << jet.neutralHadronEnergyFraction() << "," << jet.photonEnergyFraction() << "," << jet.electronEnergyFraction()
103  << "," << jet.muonEnergyFraction() << "," << jet.HFHadronEnergyFraction() << "," << jet.HFEMEnergyFraction();
104  }
105  return out;
106 }
virtual CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:99
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:152
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:139
std::string print() const override
Print object in details.
Definition: PFJet.cc:67
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:124
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:147
PFJet * clone() const override
Polymorphic clone.
Definition: PFJet.cc:63
PFJet()
Definition: PFJet.h:79
int HFHadronMultiplicity() const
HFHadronMultiplicity.
Definition: PFJet.h:134
std::vector< Constituent > Constituents
Definition: Jet.h:23
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:154
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:128
Jets made from PFObjects.
Definition: PFJet.h:20
reco::TrackRefVector getTrackRefs() const
Definition: PFJet.cc:48
daughters dau
collection of references to daughters
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:41
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
Definition: Jet.py:1
float electronEnergy() const
electronEnergy
Definition: PFJet.h:107
virtual reco::PFCandidatePtr getPFConstituent(unsigned fIndex) const
get specific constituent
Definition: PFJet.cc:27
virtual std::string print() const
Print object.
int HFEMMultiplicity() const
HFEMMultiplicity.
Definition: PFJet.h:136
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:130
bool overlap(const Candidate &) const override
Polymorphic overlap.
Definition: PFJet.cc:65
float HFHadronEnergy() const
HFHadronEnergy.
Definition: PFJet.h:115
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:126
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:132
size_t numberOfDaughters() const override
number of daughters
float HFEMEnergy() const
HFEMEnergy.
Definition: PFJet.h:119
float muonEnergy() const
muonEnergy
Definition: PFJet.h:111
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
float hoEnergy() const
hoEnergy
Definition: PFJet.h:157
float chargedMuEnergy() const
chargedMuEnergy
Definition: PFJet.h:143
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:95
float photonEnergy() const
photonEnergy
Definition: PFJet.h:103
edm::Ptr< PFCandidate > PFCandidatePtr
persistent Ptr to a PFCandidate