CMS 3D CMS Logo

PFJet.h
Go to the documentation of this file.
1 #ifndef JetReco_PFJet_h
2 #define JetReco_PFJet_h
3 
18 
19 namespace reco {
20  class PFJet : public Jet {
21  public:
24 
25  struct Specific {
29  mPhotonEnergy(0),
30  mElectronEnergy(0),
31  mMuonEnergy(0),
32  mHFHadronEnergy(0),
33  mHFEMEnergy(0),
34 
42 
46 
49 
50  mHOEnergy(0) {}
55  float mMuonEnergy;
57  float mHFEMEnergy;
58 
66 
67  //old (deprecated) data members
68  //kept only for backwards compatibility:
74 
75  float mHOEnergy;
76  };
77 
79  PFJet() {}
80 
82  PFJet(const LorentzVector& fP4,
83  const Point& fVertex,
84  const Specific& fSpecific,
85  const Jet::Constituents& fConstituents);
86 
87  PFJet(const LorentzVector& fP4, const Point& fVertex, const Specific& fSpecific);
88 
90  PFJet(const LorentzVector& fP4, const Specific& fSpecific, const Jet::Constituents& fConstituents);
91 
92  ~PFJet() override{};
93 
97  float chargedHadronEnergyFraction() const { return chargedHadronEnergy() / energy(); }
103  float photonEnergy() const { return m_specific.mPhotonEnergy; }
105  float photonEnergyFraction() const { return photonEnergy() / energy(); }
107  float electronEnergy() const { return m_specific.mElectronEnergy; }
109  float electronEnergyFraction() const { return electronEnergy() / energy(); }
111  float muonEnergy() const { return m_specific.mMuonEnergy; }
113  float muonEnergyFraction() const { return muonEnergy() / energy(); }
115  float HFHadronEnergy() const { return m_specific.mHFHadronEnergy; }
117  float HFHadronEnergyFraction() const { return HFHadronEnergy() / energy(); }
119  float HFEMEnergy() const { return m_specific.mHFEMEnergy; }
121  float HFEMEnergyFraction() const { return HFEMEnergy() / energy(); }
122 
137 
139  float chargedEmEnergy() const { return m_specific.mChargedEmEnergy; }
141  float chargedEmEnergyFraction() const { return chargedEmEnergy() / energy(); }
143  float chargedMuEnergy() const { return m_specific.mChargedMuEnergy; }
145  float chargedMuEnergyFraction() const { return chargedMuEnergy() / energy(); }
147  float neutralEmEnergy() const { return m_specific.mNeutralEmEnergy; }
149  float neutralEmEnergyFraction() const { return neutralEmEnergy() / energy(); }
150 
155 
157  float hoEnergy() const { return m_specific.mHOEnergy; }
159  float hoEnergyFraction() const { return hoEnergy() / energy(); }
160 
162  virtual reco::PFCandidatePtr getPFConstituent(unsigned fIndex) const;
163 
165  virtual std::vector<reco::PFCandidatePtr> getPFConstituents() const;
166 
171 
172  // block accessors
173 
174  const Specific& getSpecific() const { return m_specific; }
175 
177  PFJet* clone() const override;
178 
180  std::string print() const override;
181 
182  private:
184  bool overlap(const Candidate&) const override;
185 
186  //Variables specific to to the PFJet class
188  };
189 
190  // streamer
191  std::ostream& operator<<(std::ostream& out, const reco::PFJet& jet);
192 } // namespace reco
193 // temporary fix before include_checcker runs globally
194 #include "DataFormats/JetReco/interface/PFJetCollection.h" //INCLUDECHECKER:SKIP
195 #endif
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:141
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
reco::PFCandidatePtr ConstituentTypePtr
Definition: PFJet.h:22
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:109
int mChargedMultiplicity
Definition: PFJet.h:72
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
int mPhotonMultiplicity
Definition: PFJet.h:61
float mNeutralHadronEnergy
Definition: PFJet.h:52
float mChargedMuEnergy
Definition: PFJet.h:70
float HFHadronEnergyFraction() const
HFHadronEnergyFraction.
Definition: PFJet.h:117
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:154
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:128
Jets made from PFObjects.
Definition: PFJet.h:20
float mChargedHadronEnergy
Definition: PFJet.h:51
reco::TrackRefVector getTrackRefs() const
Definition: PFJet.cc:48
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:41
float mChargedEmEnergy
Definition: PFJet.h:69
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:105
float mElectronEnergy
Definition: PFJet.h:54
Definition: Jet.py:1
float electronEnergy() const
electronEnergy
Definition: PFJet.h:107
int mChargedHadronMultiplicity
Definition: PFJet.h:59
virtual reco::PFCandidatePtr getPFConstituent(unsigned fIndex) const
get specific constituent
Definition: PFJet.cc:27
float mPhotonEnergy
Definition: PFJet.h:53
float mNeutralEmEnergy
Definition: PFJet.h:71
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
math::XYZTLorentzVector LorentzVector
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:97
float HFHadronEnergy() const
HFHadronEnergy.
Definition: PFJet.h:115
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:126
int mElectronMultiplicity
Definition: PFJet.h:62
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:132
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:149
float HFEMEnergy() const
HFEMEnergy.
Definition: PFJet.h:119
int mNeutralHadronMultiplicity
Definition: PFJet.h:60
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:113
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:101
float HFEMEnergyFraction() const
HFEMEnergyFraction.
Definition: PFJet.h:121
int mNeutralMultiplicity
Definition: PFJet.h:73
const Specific & getSpecific() const
Definition: PFJet.h:174
float muonEnergy() const
muonEnergy
Definition: PFJet.h:111
float hoEnergy() const
hoEnergy
Definition: PFJet.h:157
float chargedMuEnergy() const
chargedMuEnergy
Definition: PFJet.h:143
fixed size matrix
float mHFHadronEnergy
Definition: PFJet.h:56
Structure Point Contains parameters of Gaussian fits to DMRs.
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:95
reco::PFCandidateFwdPtr ConstituentTypeFwdPtr
Definition: PFJet.h:23
Specific m_specific
Definition: PFJet.h:187
float photonEnergy() const
photonEnergy
Definition: PFJet.h:103
float hoEnergyFraction() const
hoEnergyFraction
Definition: PFJet.h:159
int mHFHadronMultiplicity
Definition: PFJet.h:64
double energy() const final
energy
float chargedMuEnergyFraction() const
chargedMuEnergyFraction
Definition: PFJet.h:145
~PFJet() override
Definition: PFJet.h:92