CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFJetIDSelectionFunctor.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
2 #define PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
3 
4 
24 
26 
27 #include <TMath.h>
28 class PFJetIDSelectionFunctor : public Selector<pat::Jet> {
29 
30  public: // interface
31 
34 
36 
38  {
39  std::string versionStr = params.getParameter<std::string>("version");
40  std::string qualityStr = params.getParameter<std::string>("quality");
41 
42  if ( versionStr == "FIRSTDATA" )
44  else
46 
47  if ( qualityStr == "LOOSE") quality_ = LOOSE;
48  else if ( qualityStr == "TIGHT") quality_ = TIGHT;
49  else quality_ = LOOSE;
50 
51  push_back("CHF" );
52  push_back("NHF" );
53  push_back("CEF" );
54  push_back("NEF" );
55  push_back("NCH" );
56  push_back("nConstituents");
57 
58 
59  // Set some default cuts for LOOSE, TIGHT
60  if ( quality_ == LOOSE ) {
61  set("CHF", 0.0);
62  set("NHF", 0.99);
63  set("CEF", 0.99);
64  set("NEF", 0.99);
65  set("NCH", 0);
66  set("nConstituents", 1);
67  } else if ( quality_ == TIGHT ) {
68  set("CHF", 0.0);
69  set("NHF", 0.9);
70  set("CEF", 0.99);
71  set("NEF", 0.9);
72  set("NCH", 0);
73  set("nConstituents", 1);
74  }
75 
76 
77  // Now check the configuration to see if the user changed anything
78  if ( params.exists("CHF") ) set("CHF", params.getParameter<double>("CHF") );
79  if ( params.exists("NHF") ) set("NHF", params.getParameter<double>("NHF") );
80  if ( params.exists("CEF") ) set("CEF", params.getParameter<double>("CEF") );
81  if ( params.exists("NEF") ) set("NEF", params.getParameter<double>("NEF") );
82  if ( params.exists("NCH") ) set("NCH", params.getParameter<int> ("NCH") );
83  if ( params.exists("nConstuents") ) set("nConstituents", params.getParameter<int> ("nConstituents") );
84 
85  if ( params.exists("cutsToIgnore") )
86  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
87 
88 
89  indexNConstituents_ = index_type (&bits_, "nConstituents");
90  indexNEF_ = index_type (&bits_, "NEF");
91  indexNHF_ = index_type (&bits_, "NHF");
92  indexCEF_ = index_type (&bits_, "CEF");
93  indexCHF_ = index_type (&bits_, "CHF");
94  indexNCH_ = index_type (&bits_, "NCH");
95 
97 
98  }
99 
100 
102  Quality_t quality ) :
103  version_(version), quality_(quality)
104  {
105 
106  push_back("CHF" );
107  push_back("NHF" );
108  push_back("CEF" );
109  push_back("NEF" );
110  push_back("NCH" );
111  push_back("nConstituents");
112 
113 
114  // Set some default cuts for LOOSE, TIGHT
115  if ( quality_ == LOOSE ) {
116  set("CHF", 0.0);
117  set("NHF", 0.99);
118  set("CEF", 0.99);
119  set("NEF", 0.99);
120  set("NCH", 0);
121  set("nConstituents", 1);
122  } else if ( quality_ == TIGHT ) {
123  set("CHF", 0.0);
124  set("NHF", 0.9);
125  set("CEF", 0.99);
126  set("NEF", 0.9);
127  set("NCH", 0);
128  set("nConstituents", 1);
129  }
130 
131 
132  indexNConstituents_ = index_type (&bits_, "nConstituents");
133  indexNEF_ = index_type (&bits_, "NEF");
134  indexNHF_ = index_type (&bits_, "NHF");
135  indexCEF_ = index_type (&bits_, "CEF");
136  indexCHF_ = index_type (&bits_, "CHF");
137  indexNCH_ = index_type (&bits_, "NCH");
138 
140  }
141 
142 
143  //
144  // Accessor from PAT jets
145  //
147  {
148  if ( version_ == FIRSTDATA ) {
149  if ( jet.currentJECLevel() == "Uncorrected" || !jet.jecSetsAvailable() )
150  return firstDataCuts( jet, ret );
151  else
152  return firstDataCuts( jet.correctedJet("Uncorrected"), ret );
153  }
154  else {
155  return false;
156  }
157  }
159 
160  //
161  // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID.
162  // This can be used with reco quantities.
163  //
165  {
166  if ( version_ == FIRSTDATA ) return firstDataCuts( jet, ret );
167  else {
168  return false;
169  }
170  }
171 
172  bool operator()( const reco::PFJet & jet )
173  {
174  retInternal_.set(false);
175  operator()(jet, retInternal_);
177  return (bool)retInternal_;
178  }
179 
180  //
181  // cuts based on craft 08 analysis.
182  //
183  bool firstDataCuts( reco::Jet const & jet,
184  pat::strbitset & ret)
185  {
186  ret.set(false);
187 
188  // cache some variables
189  double chf = 0.0;
190  double nhf = 0.0;
191  double cef = 0.0;
192  double nef = 0.0;
193  int nch = 0;
194  int nconstituents = 0;
195 
196  // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
197  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
198  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
199  reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
200 
201  if ( patJet != 0 ) {
202  if ( patJet->isPFJet() ) {
203  chf = patJet->chargedHadronEnergyFraction();
204  nhf = ( patJet->neutralHadronEnergy() + patJet->HFHadronEnergy() ) / patJet->energy();
205  cef = patJet->chargedEmEnergyFraction();
206  nef = patJet->neutralEmEnergyFraction();
207  nch = patJet->chargedMultiplicity();
208  nconstituents = patJet->numberOfDaughters();
209  }
210  // Handle the special case where this is a composed jet for
211  // subjet analyses
212  else if ( patJet->isBasicJet() ) {
213  double e_chf = 0.0;
214  double e_nhf = 0.0;
215  double e_cef = 0.0;
216  double e_nef = 0.0;
217  nch = 0;
218  nconstituents = 0;
219 
220  for ( reco::Jet::const_iterator ibegin = patJet->begin(),
221  iend = patJet->end(), isub = ibegin;
222  isub != iend; ++isub ) {
223  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
224  e_chf += pfsub->chargedHadronEnergy();
225  e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
226  e_cef += pfsub->chargedEmEnergy();
227  e_nef += pfsub->neutralEmEnergy();
228  nch += pfsub->chargedMultiplicity();
229  nconstituents += pfsub->numberOfDaughters();
230  }
231  double e = patJet->energy();
232  if ( e > 0.000001 ) {
233  chf = e_chf / e;
234  nhf = e_nhf / e;
235  cef = e_cef / e;
236  nef = e_nef / e;
237  } else {
238  chf = nhf = cef = nef = 0.0;
239  }
240  }
241  } // end if pat jet
242  else if ( pfJet != 0 ) {
243  // CV: need to compute energy fractions in a way that works for corrected as well as for uncorrected PFJets
244  double jetEnergyUncorrected =
245  pfJet->chargedHadronEnergy()
246  + pfJet->neutralHadronEnergy()
247  + pfJet->photonEnergy()
248  + pfJet->electronEnergy()
249  + pfJet->muonEnergy()
250  + pfJet->HFHadronEnergy()
251  + pfJet->HFEMEnergy();
252  if ( jetEnergyUncorrected > 0. ) {
253  chf = pfJet->chargedHadronEnergy() / jetEnergyUncorrected;
254  nhf = ( pfJet->neutralHadronEnergy() + pfJet->HFHadronEnergy() ) / jetEnergyUncorrected;
255  cef = pfJet->chargedEmEnergy() / jetEnergyUncorrected;
256  nef = pfJet->neutralEmEnergy() / jetEnergyUncorrected;
257  }
258  nch = pfJet->chargedMultiplicity();
259  nconstituents = pfJet->numberOfDaughters();
260  } // end if PF jet
261  // Handle the special case where this is a composed jet for
262  // subjet analyses
263  else if ( basicJet != 0 ) {
264  double e_chf = 0.0;
265  double e_nhf = 0.0;
266  double e_cef = 0.0;
267  double e_nef = 0.0;
268  nch = 0;
269  nconstituents = 0;
270 
271  for ( reco::Jet::const_iterator ibegin = basicJet->begin(),
272  iend = patJet->end(), isub = ibegin;
273  isub != iend; ++isub ) {
274  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
275  e_chf += pfsub->chargedHadronEnergy();
276  e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
277  e_cef += pfsub->chargedEmEnergy();
278  e_nef += pfsub->neutralEmEnergy();
279  nch += pfsub->chargedMultiplicity();
280  nconstituents += pfsub->numberOfDaughters();
281  }
282  double e = basicJet->energy();
283  if ( e > 0.000001 ) {
284  chf = e_chf / e;
285  nhf = e_nhf / e;
286  cef = e_cef / e;
287  nef = e_nef / e;
288  }
289  } // end if basic jet
290 
291 
292  // Cuts for all |eta|:
293  if ( ignoreCut(indexNConstituents_) || nconstituents > cut(indexNConstituents_, int() ) ) passCut( ret, indexNConstituents_);
294  if ( ignoreCut(indexNEF_) || ( nef < cut(indexNEF_, double()) ) ) passCut( ret, indexNEF_);
295  if ( ignoreCut(indexNHF_) || ( nhf < cut(indexNHF_, double()) ) ) passCut( ret, indexNHF_);
296  // Cuts for |eta| < 2.4:
297  if ( ignoreCut(indexCEF_) || ( cef < cut(indexCEF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCEF_);
298  if ( ignoreCut(indexCHF_) || ( chf > cut(indexCHF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCHF_);
299  if ( ignoreCut(indexNCH_) || ( nch > cut(indexNCH_, int()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexNCH_);
300 
301  //std::cout << "<PFJetIDSelectionFunctor::firstDataCuts>:" << std::endl;
302  //std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi() << std::endl;
303  //ret.print(std::cout);
304 
305  setIgnored( ret );
306  return (bool)ret;
307  }
308 
309  private: // member variables
310 
313 
320 
321 };
322 
323 #endif
T getParameter(std::string const &) const
void set(std::string const &s, bool val=true)
Set a given selection cut, on or off.
Definition: Selector.h:106
float photonEnergy() const
photonEnergy
Definition: PFJet.h:99
bool jecSetsAvailable() const
Definition: Jet.h:117
float muonEnergy() const
muonEnergy
Definition: PFJet.h:107
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:135
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:350
Base class for all types of Jets.
Definition: Jet.h:21
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.cc:277
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: Jet.h:590
float chargedEmEnergyFraction() const
chargedEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:354
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define abs(x)
Definition: mlp_lapack.h:159
pat::strbitset::index_type index_type
Definition: Selector.h:30
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:225
bool firstDataCuts(reco::Jet const &jet, pat::strbitset &ret)
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:288
float HFHadronEnergy() const
HFHadronEnergy.
Definition: PFJet.h:111
Jets made from CaloTowers.
Definition: BasicJet.h:21
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:148
Jets made from PFObjects.
Definition: PFJet.h:22
virtual double eta() const
momentum pseudorapidity
PFJetIDSelectionFunctor(Version_t version, Quality_t quality)
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:143
bool operator()(const reco::PFJet &jet, pat::strbitset &ret)
float electronEnergy() const
electronEnergy
Definition: PFJet.h:103
PFJetIDSelectionFunctor(edm::ParameterSet const &params)
virtual double energy() const
energy
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:287
virtual size_t numberOfDaughters() const
number of daughters
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:177
float HFHadronEnergy() const
HFHadronEnergy.
Definition: Jet.h:372
float HFEMEnergy() const
HFEMEnergy.
Definition: PFJet.h:115
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:160
bool operator()(const pat::Jet &jet, pat::strbitset &ret)
This provides the interface for base classes to select objects.
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:47
bool operator()(const reco::PFJet &jet)
Functor that operates on &lt;T&gt;
Definition: Selector.h:25
PF Jet selector for pat::Jets.
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:235
virtual size_t numberOfDaughters() const
Definition: Jet.h:439
virtual const_iterator begin() const
first daughter const_iterator
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:144
std::string currentJECLevel() const
return the name of the current step of jet energy corrections
Definition: Jet.h:127
Analysis-level calorimeter jet class.
Definition: Jet.h:72
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:213
virtual const_iterator end() const
last daughter const_iterator
bool isBasicJet() const
check to see if the jet is no more than a reco::BasicJet
Definition: Jet.h:237
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:95
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:168
float neutralEmEnergyFraction() const
neutralEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:356
int chargedMultiplicity() const
chargedMultiplicity
Definition: Jet.h:618
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:91
int cut(index_type const &i, int val) const
Access the int cut values at index &quot;s&quot;.
Definition: Selector.h:195