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_SelectorUtils_interface_PFJetIDSelectionFunctor_h
2 #define PhysicsTools_SelectorUtils_interface_PFJetIDSelectionFunctor_h
3 
15 #ifndef __GCCXML__
17 #endif
20 
22 
23 #include <TMath.h>
24 class PFJetIDSelectionFunctor : public Selector<pat::Jet> {
25 
26  public: // interface
27 
30 
32 
33 #ifndef __GCCXML__
36  {}
37 #endif
38 
40  {
41  std::string versionStr = params.getParameter<std::string>("version");
42  std::string qualityStr = params.getParameter<std::string>("quality");
43 
44  if ( versionStr == "FIRSTDATA" )
46  else if( versionStr == "RUNIISTARTUP" )
48  else version_ = RUNIISTARTUP;//set RUNII STARTUP 50 ns points as default
49 
50  if ( qualityStr == "LOOSE") quality_ = LOOSE;
51  else if ( qualityStr == "TIGHT") quality_ = TIGHT;
52  else quality_ = LOOSE;
53 
54  push_back("CHF" );
55  push_back("NHF" );
56  push_back("CEF" );
57  push_back("NEF" );
58  push_back("NCH" );
59  push_back("nConstituents");
60  if(version_ == RUNIISTARTUP){
61  push_back("NEF_FW");
62  push_back("nNeutrals_FW");
63  }
64 
65 
66  // Set some default cuts for LOOSE, TIGHT
67  if ( quality_ == LOOSE ) {
68  set("CHF", 0.0);
69  set("NHF", 0.99);
70  set("CEF", 0.99);
71  set("NEF", 0.99);
72  set("NCH", 0);
73  set("nConstituents", 1);
74  if(version_ == RUNIISTARTUP){
75  set("NEF_FW",0.90);
76  set("nNeutrals_FW",10);
77  }
78  } else if ( quality_ == TIGHT ) {
79  set("CHF", 0.0);
80  set("NHF", 0.9);
81  set("CEF", 0.99);
82  set("NEF", 0.9);
83  set("NCH", 0);
84  set("nConstituents", 1);
85  if(version_ == RUNIISTARTUP){
86  set("NEF_FW",0.90);
87  set("nNeutrals_FW",10);
88  }
89  }
90 
91 
92  // Now check the configuration to see if the user changed anything
93  if ( params.exists("CHF") ) set("CHF", params.getParameter<double>("CHF") );
94  if ( params.exists("NHF") ) set("NHF", params.getParameter<double>("NHF") );
95  if ( params.exists("CEF") ) set("CEF", params.getParameter<double>("CEF") );
96  if ( params.exists("NEF") ) set("NEF", params.getParameter<double>("NEF") );
97  if ( params.exists("NCH") ) set("NCH", params.getParameter<int> ("NCH") );
98  if ( params.exists("nConstituents") ) set("nConstituents", params.getParameter<int> ("nConstituents") );
99  if(version_ == RUNIISTARTUP){
100  if ( params.exists("NEF_FW") ) set("NEF_FW", params.getParameter<double> ("NEF_FW") );
101  if ( params.exists("nNeutrals_FW") ) set("nNeutrals_FW", params.getParameter<int> ("nNeutrals_FW") );
102  }
103 
104  if ( params.exists("cutsToIgnore") )
105  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
106 
107 
108  indexNConstituents_ = index_type (&bits_, "nConstituents");
109  indexNEF_ = index_type (&bits_, "NEF");
110  indexNHF_ = index_type (&bits_, "NHF");
111  indexCEF_ = index_type (&bits_, "CEF");
112  indexCHF_ = index_type (&bits_, "CHF");
113  indexNCH_ = index_type (&bits_, "NCH");
114  if(version_ == RUNIISTARTUP){
115  indexNEF_FW_ = index_type (&bits_, "NEF_FW");
116  indexNNeutrals_FW_ = index_type (&bits_, "nNeutrals_FW");
117  }
118 
120 
121  }
122 
123 
125  Quality_t quality ) :
126  version_(version), quality_(quality)
127  {
128 
129  push_back("CHF" );
130  push_back("NHF" );
131  push_back("CEF" );
132  push_back("NEF" );
133  push_back("NCH" );
134  push_back("nConstituents");
135  if(version_ == RUNIISTARTUP){
136  push_back("NEF_FW");
137  push_back("nNeutrals_FW");
138  }
139 
140  // Set some default cuts for LOOSE, TIGHT
141  if ( quality_ == LOOSE ) {
142  set("CHF", 0.0);
143  set("NHF", 0.99);
144  set("CEF", 0.99);
145  set("NEF", 0.99);
146  set("NCH", 0);
147  set("nConstituents", 1);
148  if(version_ == RUNIISTARTUP){
149  set("NEF_FW",0.90);
150  set("nNeutrals_FW",10);
151  }
152  } else if ( quality_ == TIGHT ) {
153  set("CHF", 0.0);
154  set("NHF", 0.9);
155  set("CEF", 0.99);
156  set("NEF", 0.9);
157  set("NCH", 0);
158  set("nConstituents", 1);
159  if(version_ == RUNIISTARTUP){
160  set("NEF_FW",0.90);
161  set("nNeutrals_FW",10);
162  }
163  }
164 
165 
166  indexNConstituents_ = index_type (&bits_, "nConstituents");
167  indexNEF_ = index_type (&bits_, "NEF");
168  indexNHF_ = index_type (&bits_, "NHF");
169  indexCEF_ = index_type (&bits_, "CEF");
170  indexCHF_ = index_type (&bits_, "CHF");
171  indexNCH_ = index_type (&bits_, "NCH");
172  if(version_ == RUNIISTARTUP){
173  indexNEF_FW_ = index_type (&bits_, "NEF_FW");
174  indexNNeutrals_FW_ = index_type (&bits_, "nNeutrals_FW");
175  }
176 
178  }
179 
180 
181  //
182  // Accessor from PAT jets
183  //
185  {
186  if ( version_ == FIRSTDATA || version_ == RUNIISTARTUP ) {
187  if ( jet.currentJECLevel() == "Uncorrected" || !jet.jecSetsAvailable() )
188  return firstDataCuts( jet, ret, version_);
189  else
190  return firstDataCuts( jet.correctedJet("Uncorrected"), ret, version_ );
191  }
192  else {
193  return false;
194  }
195  }
197 
198  //
199  // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID.
200  // This can be used with reco quantities.
201  //
203  {
204  if ( version_ == FIRSTDATA || version_ == RUNIISTARTUP ){ return firstDataCuts( jet, ret, version_);
205  }
206  else {
207  return false;
208  }
209  }
210 
211  bool operator()( const reco::PFJet & jet )
212  {
213  retInternal_.set(false);
214  operator()(jet, retInternal_);
216  return (bool)retInternal_;
217  }
218 
219  //
220  // cuts based on craft 08 analysis.
221  //
222  bool firstDataCuts( reco::Jet const & jet,
224  {
225  ret.set(false);
226 
227  // cache some variables
228  double chf = 0.0;
229  double nhf = 0.0;
230  double cef = 0.0;
231  double nef = 0.0;
232  int nch = 0;
233  int nconstituents = 0;
234  int nneutrals = 0;
235 
236  // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
237  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
238  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
239  reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
240 
241  if ( patJet != 0 ) {
242  if ( patJet->isPFJet() ) {
243  chf = patJet->chargedHadronEnergyFraction();
244  nhf = patJet->neutralHadronEnergyFraction();
245  cef = patJet->chargedEmEnergyFraction();
246  nef = patJet->neutralEmEnergyFraction();
247  nch = patJet->chargedMultiplicity();
248  nconstituents = patJet->numberOfDaughters();
249  nneutrals = patJet->neutralMultiplicity();
250  }
251  // Handle the special case where this is a composed jet for
252  // subjet analyses
253  else if ( patJet->isBasicJet() ) {
254  double e_chf = 0.0;
255  double e_nhf = 0.0;
256  double e_cef = 0.0;
257  double e_nef = 0.0;
258  nch = 0;
259  nconstituents = 0;
260  nneutrals = 0;
261 
262  for ( reco::Jet::const_iterator ibegin = patJet->begin(),
263  iend = patJet->end(), isub = ibegin;
264  isub != iend; ++isub ) {
265  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
266  e_chf += pfsub->chargedHadronEnergy();
267  e_nhf += pfsub->neutralHadronEnergy();
268  e_cef += pfsub->chargedEmEnergy();
269  e_nef += pfsub->neutralEmEnergy();
270  nch += pfsub->chargedMultiplicity();
271  nconstituents += pfsub->numberOfDaughters();
272  nneutrals += pfsub->neutralMultiplicity();
273  }
274  double e = patJet->energy();
275  if ( e > 0.000001 ) {
276  chf = e_chf / e;
277  nhf = e_nhf / e;
278  cef = e_cef / e;
279  nef = e_nef / e;
280  } else {
281  chf = nhf = cef = nef = 0.0;
282  }
283  }
284  } // end if pat jet
285  else if ( pfJet != 0 ) {
286  // CV: need to compute energy fractions in a way that works for corrected as well as for uncorrected PFJets
287  double jetEnergyUncorrected =
288  pfJet->chargedHadronEnergy()
289  + pfJet->neutralHadronEnergy()
290  + pfJet->photonEnergy()
291  + pfJet->electronEnergy()
292  + pfJet->muonEnergy()
293  + pfJet->HFEMEnergy();
294  if ( jetEnergyUncorrected > 0. ) {
295  chf = pfJet->chargedHadronEnergy() / jetEnergyUncorrected;
296  nhf = pfJet->neutralHadronEnergy() / jetEnergyUncorrected;
297  cef = pfJet->chargedEmEnergy() / jetEnergyUncorrected;
298  nef = pfJet->neutralEmEnergy() / jetEnergyUncorrected;
299  }
300  nch = pfJet->chargedMultiplicity();
301  nconstituents = pfJet->numberOfDaughters();
302  nneutrals = pfJet->neutralMultiplicity();
303  } // end if PF jet
304  // Handle the special case where this is a composed jet for
305  // subjet analyses
306  else if ( basicJet != 0 ) {
307  double e_chf = 0.0;
308  double e_nhf = 0.0;
309  double e_cef = 0.0;
310  double e_nef = 0.0;
311  nch = 0;
312  nconstituents = 0;
313  for ( reco::Jet::const_iterator ibegin = basicJet->begin(),
314  iend = patJet->end(), isub = ibegin;
315  isub != iend; ++isub ) {
316  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
317  e_chf += pfsub->chargedHadronEnergy();
318  e_nhf += pfsub->neutralHadronEnergy();
319  e_cef += pfsub->chargedEmEnergy();
320  e_nef += pfsub->neutralEmEnergy();
321  nch += pfsub->chargedMultiplicity();
322  nconstituents += pfsub->numberOfDaughters();
323  nneutrals += pfsub->neutralMultiplicity();
324  }
325  double e = basicJet->energy();
326  if ( e > 0.000001 ) {
327  chf = e_chf / e;
328  nhf = e_nhf / e;
329  cef = e_cef / e;
330  nef = e_nef / e;
331  }
332  } // end if basic jet
333 
334 
335 
336  // Cuts for |eta| < 2.4 for FIRSTDATA and RUNIISTARTUP
337  if ( ignoreCut(indexCEF_) || ( cef < cut(indexCEF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCEF_);
338  if ( ignoreCut(indexCHF_) || ( chf > cut(indexCHF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCHF_);
339  if ( ignoreCut(indexNCH_) || ( nch > cut(indexNCH_, int()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexNCH_);
340 
341  if(version_ == FIRSTDATA){// Cuts for all eta for FIRSTDATA
342  if ( ignoreCut(indexNConstituents_) || ( nconstituents > cut(indexNConstituents_, int()) ) ) passCut( ret, indexNConstituents_);
343  if ( ignoreCut(indexNEF_) || ( nef < cut(indexNEF_, double()) ) ) passCut( ret, indexNEF_);
344  if ( ignoreCut(indexNHF_) || ( nhf < cut(indexNHF_, double()) ) ) passCut( ret, indexNHF_);
345  }else if(version_ == RUNIISTARTUP){
346  // Cuts for |eta| <= 3.0 for RUNIISTARTUP scenario
347  if ( ignoreCut(indexNConstituents_) || ( nconstituents > cut(indexNConstituents_, int()) || std::abs(jet.eta()) > 3.0 ) ) passCut( ret, indexNConstituents_);
348  if ( ignoreCut(indexNEF_) || ( nef < cut(indexNEF_, double()) || std::abs(jet.eta()) > 3.0 ) ) passCut( ret, indexNEF_);
349  if ( ignoreCut(indexNHF_) || ( nhf < cut(indexNHF_, double()) || std::abs(jet.eta()) > 3.0 ) ) passCut( ret, indexNHF_);
350  // Cuts for |eta| > 3.0 for RUNIISTARTUP scenario
351  if ( ignoreCut(indexNEF_FW_) || ( nef < cut(indexNEF_FW_, double()) || std::abs(jet.eta()) <= 3.0 ) ) passCut( ret, indexNEF_FW_);
352  if ( ignoreCut(indexNNeutrals_FW_) || ( nneutrals > cut(indexNNeutrals_FW_, int()) || std::abs(jet.eta()) <= 3.0 ) ) passCut( ret, indexNNeutrals_FW_);
353  }
354  //std::cout << "<PFJetIDSelectionFunctor::firstDataCuts>:" << std::endl;
355  //std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi() << std::endl;
356  //ret.print(std::cout);
357 
358  setIgnored( ret );
359  return (bool)ret;
360  }
361 
362  private: // member variables
363 
366 
373 
376 
377 };
378 
379 #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:105
float photonEnergy() const
photonEnergy
Definition: PFJet.h:106
tuple ret
prodAgent to be discontinued
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:368
bool jecSetsAvailable() const
Definition: Jet.h:128
float muonEnergy() const
muonEnergy
Definition: PFJet.h:114
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:142
virtual double energy() const final
energy
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:366
PFJetIDSelectionFunctor(edm::ParameterSet const &params, edm::ConsumesCollector &iC)
Base class for all types of Jets.
Definition: Jet.h:20
float chargedEmEnergyFraction() const
chargedEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:370
bool exists(std::string const &parameterName) const
checks if a parameter exists
pat::strbitset::index_type index_type
Definition: Selector.h:29
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:224
int neutralMultiplicity() const
neutralMultiplicity
Definition: Jet.h:416
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:287
Jets made from CaloTowers.
Definition: BasicJet.h:20
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
Jets made from PFObjects.
Definition: PFJet.h:21
PFJetIDSelectionFunctor(Version_t version, Quality_t quality)
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:150
bool operator()(const reco::PFJet &jet, pat::strbitset &ret)
float electronEnergy() const
electronEnergy
Definition: PFJet.h:110
PFJetIDSelectionFunctor(edm::ParameterSet const &params)
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:286
virtual size_t numberOfDaughters() const
number of daughters
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:176
float HFEMEnergy() const
HFEMEnergy.
Definition: PFJet.h:122
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:159
const_iterator end() const
last daughter const_iterator
Definition: Candidate.h:146
bool operator()(const pat::Jet &jet, pat::strbitset &ret)
This provides the interface for base classes to select objects.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:46
bool operator()(const reco::PFJet &jet)
Functor that operates on &lt;T&gt;
Definition: Selector.h:24
PF Jet selector for pat::Jets.
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:251
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:157
virtual size_t numberOfDaughters() const
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:138
Analysis-level calorimeter jet class.
Definition: Jet.h:77
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:212
bool isBasicJet() const
check to see if the jet is no more than a reco::BasicJet
Definition: Jet.h:253
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:144
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:102
virtual double eta() const final
momentum pseudorapidity
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:167
float neutralEmEnergyFraction() const
neutralEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:372
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
int chargedMultiplicity() const
chargedMultiplicity
Definition: Jet.h:657
bool firstDataCuts(reco::Jet const &jet, pat::strbitset &ret, Version_t version_)
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:98
int cut(index_type const &i, int val) const
Access the int cut values at index &quot;s&quot;.
Definition: Selector.h:194