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
48 
49  if ( qualityStr == "LOOSE") quality_ = LOOSE;
50  else if ( qualityStr == "TIGHT") quality_ = TIGHT;
51  else quality_ = LOOSE;
52 
53  push_back("CHF" );
54  push_back("NHF" );
55  push_back("CEF" );
56  push_back("NEF" );
57  push_back("NCH" );
58  push_back("nConstituents");
59 
60 
61  // Set some default cuts for LOOSE, TIGHT
62  if ( quality_ == LOOSE ) {
63  set("CHF", 0.0);
64  set("NHF", 0.99);
65  set("CEF", 0.99);
66  set("NEF", 0.99);
67  set("NCH", 0);
68  set("nConstituents", 1);
69  } else if ( quality_ == TIGHT ) {
70  set("CHF", 0.0);
71  set("NHF", 0.9);
72  set("CEF", 0.99);
73  set("NEF", 0.9);
74  set("NCH", 0);
75  set("nConstituents", 1);
76  }
77 
78 
79  // Now check the configuration to see if the user changed anything
80  if ( params.exists("CHF") ) set("CHF", params.getParameter<double>("CHF") );
81  if ( params.exists("NHF") ) set("NHF", params.getParameter<double>("NHF") );
82  if ( params.exists("CEF") ) set("CEF", params.getParameter<double>("CEF") );
83  if ( params.exists("NEF") ) set("NEF", params.getParameter<double>("NEF") );
84  if ( params.exists("NCH") ) set("NCH", params.getParameter<int> ("NCH") );
85  if ( params.exists("nConstuents") ) set("nConstituents", params.getParameter<int> ("nConstituents") );
86 
87  if ( params.exists("cutsToIgnore") )
88  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
89 
90 
91  indexNConstituents_ = index_type (&bits_, "nConstituents");
92  indexNEF_ = index_type (&bits_, "NEF");
93  indexNHF_ = index_type (&bits_, "NHF");
94  indexCEF_ = index_type (&bits_, "CEF");
95  indexCHF_ = index_type (&bits_, "CHF");
96  indexNCH_ = index_type (&bits_, "NCH");
97 
99 
100  }
101 
102 
104  Quality_t quality ) :
105  version_(version), quality_(quality)
106  {
107 
108  push_back("CHF" );
109  push_back("NHF" );
110  push_back("CEF" );
111  push_back("NEF" );
112  push_back("NCH" );
113  push_back("nConstituents");
114 
115 
116  // Set some default cuts for LOOSE, TIGHT
117  if ( quality_ == LOOSE ) {
118  set("CHF", 0.0);
119  set("NHF", 0.99);
120  set("CEF", 0.99);
121  set("NEF", 0.99);
122  set("NCH", 0);
123  set("nConstituents", 1);
124  } else if ( quality_ == TIGHT ) {
125  set("CHF", 0.0);
126  set("NHF", 0.9);
127  set("CEF", 0.99);
128  set("NEF", 0.9);
129  set("NCH", 0);
130  set("nConstituents", 1);
131  }
132 
133 
134  indexNConstituents_ = index_type (&bits_, "nConstituents");
135  indexNEF_ = index_type (&bits_, "NEF");
136  indexNHF_ = index_type (&bits_, "NHF");
137  indexCEF_ = index_type (&bits_, "CEF");
138  indexCHF_ = index_type (&bits_, "CHF");
139  indexNCH_ = index_type (&bits_, "NCH");
140 
142  }
143 
144 
145  //
146  // Accessor from PAT jets
147  //
149  {
150  if ( version_ == FIRSTDATA ) {
151  if ( jet.currentJECLevel() == "Uncorrected" || !jet.jecSetsAvailable() )
152  return firstDataCuts( jet, ret );
153  else
154  return firstDataCuts( jet.correctedJet("Uncorrected"), ret );
155  }
156  else {
157  return false;
158  }
159  }
161 
162  //
163  // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID.
164  // This can be used with reco quantities.
165  //
167  {
168  if ( version_ == FIRSTDATA ) return firstDataCuts( jet, ret );
169  else {
170  return false;
171  }
172  }
173 
174  bool operator()( const reco::PFJet & jet )
175  {
176  retInternal_.set(false);
177  operator()(jet, retInternal_);
179  return (bool)retInternal_;
180  }
181 
182  //
183  // cuts based on craft 08 analysis.
184  //
185  bool firstDataCuts( reco::Jet const & jet,
187  {
188  ret.set(false);
189 
190  // cache some variables
191  double chf = 0.0;
192  double nhf = 0.0;
193  double cef = 0.0;
194  double nef = 0.0;
195  int nch = 0;
196  int nconstituents = 0;
197 
198  // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
199  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
200  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
201  reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
202 
203  if ( patJet != 0 ) {
204  if ( patJet->isPFJet() ) {
205  chf = patJet->chargedHadronEnergyFraction();
206  nhf = patJet->neutralHadronEnergyFraction();
207  cef = patJet->chargedEmEnergyFraction();
208  nef = patJet->neutralEmEnergyFraction();
209  nch = patJet->chargedMultiplicity();
210  nconstituents = patJet->numberOfDaughters();
211  }
212  // Handle the special case where this is a composed jet for
213  // subjet analyses
214  else if ( patJet->isBasicJet() ) {
215  double e_chf = 0.0;
216  double e_nhf = 0.0;
217  double e_cef = 0.0;
218  double e_nef = 0.0;
219  nch = 0;
220  nconstituents = 0;
221 
222  for ( reco::Jet::const_iterator ibegin = patJet->begin(),
223  iend = patJet->end(), isub = ibegin;
224  isub != iend; ++isub ) {
225  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
226  e_chf += pfsub->chargedHadronEnergy();
227  e_nhf += pfsub->neutralHadronEnergy();
228  e_cef += pfsub->chargedEmEnergy();
229  e_nef += pfsub->neutralEmEnergy();
230  nch += pfsub->chargedMultiplicity();
231  nconstituents += pfsub->numberOfDaughters();
232  }
233  double e = patJet->energy();
234  if ( e > 0.000001 ) {
235  chf = e_chf / e;
236  nhf = e_nhf / e;
237  cef = e_cef / e;
238  nef = e_nef / e;
239  } else {
240  chf = nhf = cef = nef = 0.0;
241  }
242  }
243  } // end if pat jet
244  else if ( pfJet != 0 ) {
245  // CV: need to compute energy fractions in a way that works for corrected as well as for uncorrected PFJets
246  double jetEnergyUncorrected =
247  pfJet->chargedHadronEnergy()
248  + pfJet->neutralHadronEnergy()
249  + pfJet->photonEnergy()
250  + pfJet->electronEnergy()
251  + pfJet->muonEnergy()
252  + pfJet->HFEMEnergy();
253  if ( jetEnergyUncorrected > 0. ) {
254  chf = pfJet->chargedHadronEnergy() / jetEnergyUncorrected;
255  nhf = pfJet->neutralHadronEnergy() / jetEnergyUncorrected;
256  cef = pfJet->chargedEmEnergy() / jetEnergyUncorrected;
257  nef = pfJet->neutralEmEnergy() / jetEnergyUncorrected;
258  }
259  nch = pfJet->chargedMultiplicity();
260  nconstituents = pfJet->numberOfDaughters();
261  } // end if PF jet
262  // Handle the special case where this is a composed jet for
263  // subjet analyses
264  else if ( basicJet != 0 ) {
265  double e_chf = 0.0;
266  double e_nhf = 0.0;
267  double e_cef = 0.0;
268  double e_nef = 0.0;
269  nch = 0;
270  nconstituents = 0;
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();
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:105
float photonEnergy() const
photonEnergy
Definition: PFJet.h:106
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
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
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: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
virtual double eta() const
momentum pseudorapidity
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)
virtual double energy() const
energy
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
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
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
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