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