CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
15 PFJet::PFJet (const LorentzVector& fP4, const Point& fVertex,
16  const Specific& fSpecific,
17  const Jet::Constituents& fConstituents)
18  : Jet (fP4, fVertex, fConstituents),
19  m_specific (fSpecific)
20 {}
21 
22 PFJet::PFJet (const LorentzVector& fP4, const Point& fVertex,
23  const Specific& fSpecific)
24  : Jet (fP4, fVertex),
25  m_specific (fSpecific)
26 {}
27 
29  const Specific& fSpecific,
30  const Jet::Constituents& fConstituents)
31  : Jet (fP4, Point(0,0,0), fConstituents),
32  m_specific (fSpecific)
33 {}
34 
36 
37  Constituent dau = daughterPtr (fIndex);
38  if ( dau.isNonnull() && dau.isAvailable() ) {
39  const PFCandidate* pfCandidate = dynamic_cast <const PFCandidate*> (dau.get());
40  if (pfCandidate) {
41  return edm::Ptr<PFCandidate> (dau.id(), pfCandidate, dau.key() );
42  }
43  else {
44  throw cms::Exception("Invalid Constituent") << "PFJet constituent is not of PFCandidate type";
45  }
46  }
47  else {
48  return PFCandidatePtr();
49  }
50 }
51 
52 std::vector <reco::PFCandidatePtr> PFJet::getPFConstituents () const {
53  std::vector <PFCandidatePtr> result;
54  for (unsigned i = 0; i < numberOfDaughters (); i++) result.push_back (getPFConstituent(i));
55  return result;
56 }
57 
58 
60  // result will contain chargedMultiplicity() elements
62  result.reserve( chargedMultiplicity() );
63  for (unsigned i = 0; i < numberOfDaughters (); i++) {
64  const reco::PFCandidatePtr pfcand = getPFConstituent (i);
65  reco::TrackRef trackref = pfcand->trackRef();
66  if( trackref.isNonnull() ) {
67  result.push_back( trackref );
68  }
69  }
70 
71  return result;
72 }
73 
74 
75 PFJet* PFJet::clone () const {
76  return new PFJet (*this);
77 }
78 
79 bool PFJet::overlap( const Candidate & ) const {
80  return false;
81 }
82 
84  std::ostringstream out;
85  out << Jet::print () // generic jet info
86  << " PFJet specific:" << std::endl
87  << " charged hadron energy/multiplicity: " << chargedHadronEnergy () << '/' << chargedHadronMultiplicity () << std::endl
88  << " neutral hadron energy/multiplicity: " << neutralHadronEnergy () << '/' << neutralHadronMultiplicity () << std::endl
89  << " photon energy/multiplicity: " << photonEnergy () << '/' << photonMultiplicity () << std::endl
90  << " electron energy/multiplicity: " << electronEnergy () << '/' << electronMultiplicity () << std::endl
91  << " muon energy/multiplicity: " << muonEnergy () << '/' << muonMultiplicity () << std::endl
92  << " HF Hadron energy/multiplicity: " << HFHadronEnergy () << '/' << HFHadronMultiplicity () << std::endl
93  << " HF EM particle energy/multiplicity: " << HFEMEnergy () << '/' << HFEMMultiplicity () << std::endl
94  << " charged/neutral hadrons energy: " << chargedHadronEnergy () << '/' << neutralHadronEnergy () << std::endl
95  << " charged/neutral em energy: " << chargedEmEnergy () << '/' << neutralEmEnergy () << std::endl
96  << " charged muon energy: " << chargedMuEnergy () << '/' << std::endl
97  << " charged/neutral multiplicity: " << chargedMultiplicity () << '/' << neutralMultiplicity () << std::endl;
98  out << " PFCandidate constituents:" << std::endl;
99  std::vector <PFCandidatePtr> constituents = getPFConstituents ();
100  for (unsigned i = 0; i < constituents.size (); ++i) {
101  if (constituents[i].get()) {
102  out << " #" << i << " " << *(constituents[i]) << std::endl;
103  }
104  else {
105  out << " #" << i << " PFCandidate is not available in the event" << std::endl;
106  }
107  }
108  return out.str ();
109 }
110 
111 std::ostream& reco::operator<<(std::ostream& out, const reco::PFJet& jet) {
112 
113  if(out) {
114  out<<"PFJet "
115  <<"(pt, eta, phi) = "<<jet.pt()<<","<<jet.eta()<<","<<jet.phi()
116  <<" (Rch,Rnh,Rgamma,Re,Rmu,RHFHad,RHFEM) = "
117  <<jet.chargedHadronEnergyFraction()<<","
118  <<jet.neutralHadronEnergyFraction()<<","
119  <<jet.photonEnergyFraction()<<","
120  <<jet.electronEnergyFraction()<<","
121  <<jet.muonEnergyFraction()<<","
122  <<jet.HFHadronEnergyFraction()<<","
123  <<jet.HFEMEnergyFraction();
124  }
125  return out;
126 }
int i
Definition: DBlmapReader.cc:9
float photonEnergy() const
photonEnergy
Definition: PFJet.h:102
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:127
float muonEnergy() const
muonEnergy
Definition: PFJet.h:110
virtual bool overlap(const Candidate &) const
Polymorphic overlap.
Definition: PFJet.cc:79
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:138
virtual float pt() const
transverse momentum
CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
virtual float phi() const
momentum azimuthal angle
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:96
Base class for all types of Jets.
Definition: Jet.h:20
PFJet()
Definition: PFJet.h:78
std::vector< Constituent > Constituents
Definition: Jet.h:23
float HFEMEnergyFraction() const
HFEMEnergyFraction.
Definition: PFJet.h:120
float HFHadronEnergy() const
HFHadronEnergy.
Definition: PFJet.h:114
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:151
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:104
Jets made from PFObjects.
Definition: PFJet.h:21
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:146
virtual reco::PFCandidatePtr getPFConstituent(unsigned fIndex) const
get specific constituent
Definition: PFJet.cc:35
float electronEnergy() const
electronEnergy
Definition: PFJet.h:106
daughters dau
collection of references to daughters
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
virtual size_t numberOfDaughters() const
number of daughters
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:143
float HFEMEnergy() const
HFEMEnergy.
Definition: PFJet.h:118
virtual float eta() const
momentum pseudorapidity
virtual std::string print() const
Print object in details.
Definition: PFJet.cc:83
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:164
tuple result
Definition: query.py:137
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:125
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:100
void reserve(size_type n)
Reserve space for RefVector.
Definition: RefVector.h:95
virtual PFJet * clone() const
Polymorphic clone.
Definition: PFJet.cc:75
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:153
tuple out
Definition: dbtoconf.py:99
key_type key() const
Definition: Ptr.h:169
float HFHadronEnergyFraction() const
HFHadronEnergyFraction.
Definition: PFJet.h:116
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:108
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:123
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:131
int HFEMMultiplicity() const
HFEMMultiplicity.
Definition: PFJet.h:135
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
reco::TrackRefVector getTrackRefs() const
Definition: PFJet.cc:59
bool isAvailable() const
Definition: Ptr.h:158
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:98
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
float chargedMuEnergy() const
chargedMuEnergy
Definition: PFJet.h:142
virtual std::string print() const
Print object.
Definition: Jet.cc:424
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:112
edm::Ptr< PFCandidate > PFCandidatePtr
persistent Ptr to a PFCandidate
int HFHadronMultiplicity() const
HFHadronMultiplicity.
Definition: PFJet.h:133
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:94
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:129