CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetIDSelectionFunctor.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_interface_JetIDSelectionFunctor_h
2 #define PhysicsTools_PatUtils_interface_JetIDSelectionFunctor_h
3 
4 
24 
27 
28 #include <TMath.h>
29 class JetIDSelectionFunctor : public Selector<pat::Jet> {
30 
31  public: // interface
32 
35 
37 
38 
40  std::string versionStr = parameters.getParameter<std::string>("version");
41  std::string qualityStr = parameters.getParameter<std::string>("quality");
42  Quality_t quality = N_QUALITY;
43 
44  if ( qualityStr == "MINIMAL" ) quality = MINIMAL;
45  else if ( qualityStr == "LOOSE_AOD" ) quality = LOOSE_AOD;
46  else if ( qualityStr == "LOOSE" ) quality = LOOSE;
47  else if ( qualityStr == "TIGHT" ) quality = TIGHT;
48  else
49  throw cms::Exception("InvalidInput") << "Expect quality to be one of MINIMAL, LOOSE_AOD, LOOSE,TIGHT" << std::endl;
50 
51 
52  if ( versionStr == "CRAFT08" ) {
53  initialize( CRAFT08, quality );
54  }
55  else if ( versionStr == "PURE09" ) {
56  initialize( PURE09, quality );
57  }
58  else if ( versionStr == "DQM09" ) {
59  initialize( DQM09, quality );
60  }
61  else {
62  throw cms::Exception("InvalidInput") << "Expect version to be one of CRAFT08, PURE09, DQM09" << std::endl;
63  }
64  }
65 
67  initialize(version, quality);
68  }
69 
71  {
72 
73  version_ = version;
74  quality_ = quality;
75 
76  push_back("MINIMAL_EMF");
77 
78  push_back("LOOSE_AOD_fHPD");
79  push_back("LOOSE_AOD_N90Hits");
80  push_back("LOOSE_AOD_EMF");
81 
82  push_back("LOOSE_fHPD");
83  push_back("LOOSE_N90Hits");
84  push_back("LOOSE_EMF");
85 
86  push_back("TIGHT_fHPD");
87  push_back("TIGHT_EMF");
88 
89  push_back("LOOSE_nHit");
90  push_back("LOOSE_als");
91  push_back("LOOSE_fls");
92  push_back("LOOSE_foot");
93 
94  push_back("TIGHT_nHit");
95  push_back("TIGHT_als");
96  push_back("TIGHT_fls");
97  push_back("TIGHT_foot");
98  push_back("widths");
99  push_back("EF_N90Hits");
100  push_back("EF_EMF");
101 
102 
103  bool use_09_fwd_id = version_ != CRAFT08; // CRAFT08 predates the 09 forward ID cuts
104  bool use_dqm_09 = version_ == DQM09 && quality_ != LOOSE_AOD;
105 
106  // all appropriate for version and format (AOD complications) are on by default
107  set( "MINIMAL_EMF" );
108  set( "LOOSE_AOD_fHPD" );
109  set( "LOOSE_AOD_N90Hits" );
110  set( "LOOSE_AOD_EMF", ! use_09_fwd_id ); // except in CRAFT08, this devolves into MINIMAL_EMF
111  set( "LOOSE_fHPD" );
112  set( "LOOSE_N90Hits" );
113  set( "LOOSE_EMF", ! use_09_fwd_id ); // except in CRAFT08, this devolves into MINIMAL_EMF
114  set( "TIGHT_fHPD" );
115  set( "TIGHT_EMF" );
116 
117  set( "LOOSE_nHit", use_09_fwd_id );
118  set( "LOOSE_als", use_09_fwd_id );
119  set( "TIGHT_nHit", use_09_fwd_id );
120  set( "TIGHT_als", use_09_fwd_id );
121  set( "widths", use_09_fwd_id );
122  set( "EF_N90Hits", use_09_fwd_id );
123  set( "EF_EMF", use_09_fwd_id );
124 
125  set( "LOOSE_fls", use_dqm_09 );
126  set( "LOOSE_foot", use_dqm_09 );
127  set( "TIGHT_fls", use_dqm_09 );
128  set( "TIGHT_foot", use_dqm_09 );
129 
130 
131 
132 
133  index_MINIMAL_EMF_ = index_type(&bits_, "MINIMAL_EMF");
134 
135  index_LOOSE_AOD_fHPD_ = index_type(&bits_, "LOOSE_AOD_fHPD");
136  index_LOOSE_AOD_N90Hits_ = index_type(&bits_, "LOOSE_AOD_N90Hits");
137  index_LOOSE_AOD_EMF_ = index_type(&bits_, "LOOSE_AOD_EMF");
138 
139  index_LOOSE_fHPD_ = index_type(&bits_, "LOOSE_fHPD");
140  index_LOOSE_N90Hits_ = index_type(&bits_, "LOOSE_N90Hits");
141  index_LOOSE_EMF_ = index_type(&bits_, "LOOSE_EMF");
142 
143  index_TIGHT_fHPD_ = index_type(&bits_, "TIGHT_fHPD");
144  index_TIGHT_EMF_ = index_type(&bits_, "TIGHT_EMF");
145 
146  index_LOOSE_nHit_ = index_type(&bits_, "LOOSE_nHit");
147  index_LOOSE_als_ = index_type(&bits_, "LOOSE_als");
148  index_LOOSE_fls_ = index_type(&bits_, "LOOSE_fls");
149  index_LOOSE_foot_ = index_type(&bits_, "LOOSE_foot");
150 
151  index_TIGHT_nHit_ = index_type(&bits_, "TIGHT_nHit");
152  index_TIGHT_als_ = index_type(&bits_, "TIGHT_als");
153  index_TIGHT_fls_ = index_type(&bits_, "TIGHT_fls");
154  index_TIGHT_foot_ = index_type(&bits_, "TIGHT_foot");
155  index_widths_ = index_type(&bits_, "widths");
156  index_EF_N90Hits_ = index_type(&bits_, "EF_N90Hits");
157  index_EF_EMF_ = index_type(&bits_, "EF_EMF");
158 
159  // now set the return values for the ignored parts
160  bool use_loose_aod = false;
161  bool use_loose = false;
162  bool use_tight = false;
163  bool use_tight_09_fwd_id = false;
164  bool use_loose_09_fwd_id = false;
165  // if ( quality_ == MINIMAL ) nothing to do...
166  if ( quality_ == LOOSE ) {
167  use_loose = true;
168  if( use_09_fwd_id ) use_loose_09_fwd_id = true;
169  }
170  if ( quality_ == LOOSE_AOD ) {
171  use_loose_aod = true;
172  if( use_09_fwd_id ) use_loose_09_fwd_id = true;
173  }
174  if ( quality_ == TIGHT ) {
175  use_tight = true;
176  if( use_09_fwd_id ) use_tight_09_fwd_id = true;
177  }
178 
179  if( ! use_loose_aod ) {
180  set("LOOSE_AOD_fHPD", false );
181  set("LOOSE_AOD_N90Hits", false );
182  set("LOOSE_AOD_EMF", false );
183  }
184 
185  if( ! ( use_loose || use_tight ) ) { // the CRAFT08 cuts are cumulative
186  set("LOOSE_N90Hits", false );
187  set("LOOSE_fHPD", false );
188  set("LOOSE_EMF", false );
189  }
190 
191  if( ! use_tight ) {
192  set("TIGHT_fHPD", false );
193  set("TIGHT_EMF", false );
194  }
195 
196  if( ! use_loose_09_fwd_id ) { // the FWD09 cuts are not
197  set( "LOOSE_nHit", false );
198  set( "LOOSE_als", false );
199  if( use_dqm_09 ) {
200  set( "LOOSE_fls", false );
201  set( "LOOSE_foot", false );
202  }
203  } // not using loose 09 fwd ID
204 
205  if( ! use_tight_09_fwd_id ) {
206  set( "TIGHT_nHit", false );
207  set( "TIGHT_als", false );
208  set( "widths", false );
209  set( "EF_N90Hits", false );
210  set( "EF_EMF", false );
211  if( use_dqm_09 ) {
212  set( "TIGHT_fls", false );
213  set( "TIGHT_foot", false );
214  }
215  } // not using tight 09 fwd ID
216 
218  }
219 
220  // this functionality should be migrated into JetIDHelper in future releases
221  unsigned int count_hits( const std::vector<CaloTowerPtr> & towers )
222  {
223  unsigned int nHit = 0;
224  for ( unsigned int iTower = 0; iTower < towers.size() ; ++iTower ) {
225  const std::vector<DetId>& cellIDs = towers[iTower]->constituents(); // cell == recHit
226  nHit += cellIDs.size();
227  }
228  return nHit;
229  }
230 
231  //
232  // Accessor from PAT jets
233  //
235  {
236  if ( ! jet.isCaloJet() && !jet.isJPTJet() ) {
237  edm::LogWarning( "NYI" )<<"Criteria for pat::Jet-s other than CaloJets and JPTJets are not yet implemented";
238  return false;
239  }
240  if ( version_ == CRAFT08 ) return craft08Cuts( jet.p4(), jet.emEnergyFraction(), jet.jetID(), ret );
241  if ( version_ == PURE09 || version_ == DQM09 ) {
242  unsigned int nHit = count_hits( jet.getCaloConstituents() );
243  if ( jet.currentJECLevel() == "Uncorrected" ) {
244  return fwd09Cuts( jet.p4(), jet.emEnergyFraction(), jet.etaetaMoment(), jet.phiphiMoment(), nHit,
245  jet.jetID(), ret );
246  }
247  else {
248  return fwd09Cuts( jet.correctedP4("Uncorrected"), jet.emEnergyFraction(), jet.etaetaMoment(), jet.phiphiMoment(), nHit,
249  jet.jetID(), ret );
250  }
251  }
252  edm::LogWarning( "BadInput | NYI" )<<"Requested version ("<<version_<<") is unknown";
253  return false;
254  }
255 
256  // accessor from PAT jets without the ret
258 
259  //
260  // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID.
261  // This can be used with reco quantities.
262  //
264  double emEnergyFraction,
265  reco::JetID const & jetID,
266  pat::strbitset & ret )
267  {
268  if ( version_ == CRAFT08 ) return craft08Cuts( correctedP4, emEnergyFraction, jetID, ret );
269  edm::LogWarning( "BadInput | NYI" )<<"Requested version ("<<version_
270  <<") is unknown or doesn't match the depreceated interface used";
271  return false;
272  }
275  double emEnergyFraction,
276  reco::JetID const & jetID )
277  {
278  retInternal_.set(false);
279  operator()(correctedP4,emEnergyFraction, jetID, retInternal_);
281  return (bool)retInternal_;
282  }
283 
284  //
285  // Accessor from CaloJet and Jet ID. Jets MUST BE corrected for craft08, but uncorrected for later cuts...
286  // This can be used with reco quantities.
287  //
288  bool operator()( reco::CaloJet const & jet,
289  reco::JetID const & jetID,
290  pat::strbitset & ret )
291  {
292  if ( version_ == CRAFT08 ) return craft08Cuts( jet.p4(), jet.emEnergyFraction(), jetID, ret );
293  if ( version_ == PURE09 || version_ == DQM09 ) {
294  unsigned int nHit = count_hits( jet.getCaloConstituents() );
295  return fwd09Cuts( jet.p4(), jet.emEnergyFraction(), jet.etaetaMoment(), jet.phiphiMoment(), nHit,
296  jetID, ret );
297  }
298  edm::LogWarning( "BadInput | NYI" )<<"Requested version ("<<version_<<") is unknown";
299  return false;
300  }
301 
303  virtual bool operator()( reco::CaloJet const & jet,
304  reco::JetID const & jetID )
305  {
306  retInternal_.set(false);
307  operator()(jet, jetID, retInternal_);
309  return (bool)retInternal_;
310  }
311 
312  //
313  // cuts based on craft 08 analysis.
314  //
316  double emEnergyFraction,
317  reco::JetID const & jetID,
318  pat::strbitset & ret)
319  {
320 
321  ret.set(false);
322 
323  // cache some variables
324  double abs_eta = TMath::Abs( correctedP4.eta() );
325  double corrPt = correctedP4.pt();
326  double emf = emEnergyFraction;
327 
328  if ( ignoreCut(index_MINIMAL_EMF_) || abs_eta > 2.6 || emf > 0.01 ) passCut( ret, index_MINIMAL_EMF_);
329 
330  if ( quality_ == LOOSE_AOD ) {
331  // loose fhpd cut from aod
333  // loose n90 hits cut
335 
336  // loose EMF Cut from aod
337  bool emf_loose = true;
338  if( abs_eta <= 2.6 ) { // HBHE
339  if( emEnergyFraction <= 0.01 ) emf_loose = false;
340  } else { // HF
341  if( emEnergyFraction <= -0.9 ) emf_loose = false;
342  if( corrPt > 80 && emEnergyFraction >= 1 ) emf_loose = false;
343  }
344  if ( ignoreCut(index_LOOSE_AOD_EMF_) || emf_loose ) passCut(ret, index_LOOSE_AOD_EMF_);
345 
346  }
347  else if ( quality_ == LOOSE || quality_ == TIGHT ) {
348  // loose fhpd cut
349  if ( ignoreCut(index_LOOSE_fHPD_) || jetID.fHPD < 0.98 ) passCut( ret, index_LOOSE_fHPD_);
350  // loose n90 hits cut
352 
353  // loose EMF Cut
354  bool emf_loose = true;
355  if( abs_eta <= 2.6 ) { // HBHE
356  if( emEnergyFraction <= 0.01 ) emf_loose = false;
357  } else { // HF
358  if( emEnergyFraction <= -0.9 ) emf_loose = false;
359  if( corrPt > 80 && emEnergyFraction >= 1 ) emf_loose = false;
360  }
361  if ( ignoreCut(index_LOOSE_EMF_) || emf_loose ) passCut(ret, index_LOOSE_EMF_);
362 
363  if ( quality_ == TIGHT ) {
364  // tight fhpd cut
365  bool tight_fhpd = true;
366  if ( corrPt >= 25 && jetID.fHPD >= 0.95 ) tight_fhpd = false; // this was supposed to use raw pT, see AN2009/087 :-(
367  if ( ignoreCut(index_TIGHT_fHPD_) || tight_fhpd ) passCut(ret, index_TIGHT_fHPD_);
368 
369  // tight emf cut
370  bool tight_emf = true;
371  if( abs_eta >= 1 && corrPt >= 80 && emEnergyFraction >= 1 ) tight_emf = false; // outside HB
372  if( abs_eta >= 2.6 ) { // outside HBHE
373  if( emEnergyFraction <= -0.3 ) tight_emf = false;
374  if( abs_eta < 3.25 ) { // HE-HF transition region
375  if( corrPt >= 50 && emEnergyFraction <= -0.2 ) tight_emf = false;
376  if( corrPt >= 80 && emEnergyFraction <= -0.1 ) tight_emf = false;
377  if( corrPt >= 340 && emEnergyFraction >= 0.95 ) tight_emf = false;
378  } else { // HF
379  if( emEnergyFraction >= 0.9 ) tight_emf = false;
380  if( corrPt >= 50 && emEnergyFraction <= -0.2 ) tight_emf = false;
381  if( corrPt >= 50 && emEnergyFraction >= 0.8 ) tight_emf = false;
382  if( corrPt >= 130 && emEnergyFraction <= -0.1 ) tight_emf = false;
383  if( corrPt >= 130 && emEnergyFraction >= 0.7 ) tight_emf = false;
384 
385  } // end if HF
386  }// end if outside HBHE
387  if ( ignoreCut(index_TIGHT_EMF_) || tight_emf ) passCut(ret, index_TIGHT_EMF_);
388  }// end if tight
389  }// end if loose or tight
390 
391  setIgnored( ret );
392 
393  return (bool)ret;
394  }
395 
396  //
397  // cuts based on craft 08 analysis + forward jet ID based on 09 data
398  //
400  double emEnergyFraction, double etaWidth, double phiWidth, unsigned int nHit,
401  reco::JetID const & jetID,
402  pat::strbitset & ret)
403  {
404  ret.set(false);
405 
406  // cache some variables
407  double abs_eta = TMath::Abs( rawP4.eta() );
408  double rawPt = rawP4.pt();
409  double emf = emEnergyFraction;
410  double fhf = jetID.fLong + jetID.fShort;
411  double lnpt = ( rawPt > 0 ) ? TMath::Log( rawPt ) : -10;
412  double lnE = ( rawP4.energy() > 0 ) ? TMath::Log( rawP4.energy() ) : -10;
413 
414  bool HB = abs_eta < 1.0;
415  bool EF = 2.6 <= abs_eta && abs_eta < 3.4 && 0.1 <= fhf && fhf < 0.9;
416  bool HBHE = abs_eta < 2.6 || ( abs_eta < 3.4 && fhf < 0.1 );
417  bool HF = 3.4 <= abs_eta || ( 2.6 <= abs_eta && 0.9 <= fhf );
418  bool HFa = HF && abs_eta < 3.8;
419  bool HFb = HF && ! HFa;
420 
421  // HBHE cuts as in CRAFT08
422  // - but using raw pTs
423  // ========================
424 
425  if ( (!HBHE) || ignoreCut(index_MINIMAL_EMF_) || emf > 0.01 ) passCut( ret, index_MINIMAL_EMF_);
426 
427  // loose fhpd cut from AOD
428  if ( (!HBHE) || ignoreCut(index_LOOSE_AOD_fHPD_) || jetID.approximatefHPD < 0.98 ) passCut( ret, index_LOOSE_AOD_fHPD_);
429  // loose n90 hits cut from AOD
430  if ( (!HBHE) || ignoreCut(index_LOOSE_AOD_N90Hits_) || jetID.hitsInN90 > 1 ) passCut( ret, index_LOOSE_AOD_N90Hits_);
431 
432  // loose fhpd cut
433  if ( (!HBHE) || ignoreCut(index_LOOSE_fHPD_) || jetID.fHPD < 0.98 ) passCut( ret, index_LOOSE_fHPD_);
434  // loose n90 hits cut
435  if ( (!HBHE) || ignoreCut(index_LOOSE_N90Hits_) || jetID.n90Hits > 1 ) passCut( ret, index_LOOSE_N90Hits_);
436 
437  // tight fhpd cut
438  if ( (!HBHE) || ignoreCut(index_TIGHT_fHPD_) || rawPt < 25 || jetID.fHPD < 0.95 ) passCut(ret, index_TIGHT_fHPD_);
439 
440  // tight emf cut
441  if ( (!HBHE) || ignoreCut(index_TIGHT_EMF_) || HB || rawPt < 55 || emf < 1 ) passCut(ret, index_TIGHT_EMF_);
442 
443 
444  // EF - these cuts are only used in "tight", but there's no need for this test here.
445 
446  if( (!EF) || ignoreCut( index_EF_N90Hits_ )
447  || jetID.n90Hits > 1 + 1.5 * TMath::Max( 0., lnpt - 1.5 ) )
448  passCut( ret, index_EF_N90Hits_ );
449 
450  if( (!EF) || ignoreCut( index_EF_EMF_ )
451  || emf > TMath::Max( -0.9, -0.1 - 0.05 * TMath::Power( TMath::Max( 0., 5 - lnpt ), 2. ) ) )
452  passCut( ret, index_EF_EMF_ );
453 
454  // both EF and HF
455 
456  if( ( !( EF || HF ) ) || ignoreCut( index_TIGHT_fls_ )
457  || ( EF && jetID.fLS < TMath::Min( 0.8, 0.1 + 0.016 * TMath::Power( TMath::Max( 0., 6 - lnpt ), 2.5 ) ) )
458  || ( HFa && jetID.fLS < TMath::Min( 0.6, 0.05 + 0.045 * TMath::Power( TMath::Max( 0., 7.5 - lnE ), 2.2 ) ) )
459  || ( HFb && jetID.fLS < TMath::Min( 0.1, 0.05 + 0.07 * TMath::Power( TMath::Max( 0., 7.8 - lnE ), 2. ) ) ) )
460  passCut( ret, index_TIGHT_fls_ );
461 
462  if( ( !( EF || HF ) ) || ignoreCut( index_widths_ )
463  || ( 1E-10 < etaWidth && etaWidth < 0.12 &&
464  1E-10 < phiWidth && phiWidth < 0.12 ) )
465  passCut( ret, index_widths_ );
466 
467  // HF cuts
468 
469  if( (!HF) || ignoreCut( index_LOOSE_nHit_ )
470  || ( HFa && nHit > 1 + 2.4*( lnpt - 1. ) )
471  || ( HFb && nHit > 1 + 3.*( lnpt - 1. ) ) )
472  passCut( ret, index_LOOSE_nHit_ );
473 
474  if( (!HF) || ignoreCut( index_LOOSE_als_ )
475  || ( emf < 0.6 + 0.05 * TMath::Power( TMath::Max( 0., 9 - lnE ), 1.5 ) &&
476  emf > -0.2 - 0.041 * TMath::Power( TMath::Max( 0., 7.5 - lnE ), 2.2 ) ) )
477  passCut( ret, index_LOOSE_als_ );
478 
479  if( (!HF) || ignoreCut( index_LOOSE_fls_ )
480  || ( HFa && jetID.fLS < TMath::Min( 0.9, 0.1 + 0.05 * TMath::Power( TMath::Max( 0., 7.5 - lnE ), 2.2 ) ) )
481  || ( HFb && jetID.fLS < TMath::Min( 0.6, 0.1 + 0.065 * TMath::Power( TMath::Max( 0., 7.5 - lnE ), 2.2 ) ) ) )
482  passCut( ret, index_LOOSE_fls_ );
483 
484  if( (!HF) || ignoreCut( index_LOOSE_foot_ )
485  || jetID.fHFOOT < 0.9 )
486  passCut( ret, index_LOOSE_foot_ );
487 
488  if( (!HF) || ignoreCut( index_TIGHT_nHit_ )
489  || ( HFa && nHit > 1 + 2.7*( lnpt - 0.8 ) )
490  || ( HFb && nHit > 1 + 3.5*( lnpt - 0.8 ) ) )
491  passCut( ret, index_TIGHT_nHit_ );
492 
493  if( (!HF) || ignoreCut( index_TIGHT_als_ )
494  || ( emf < 0.5 + 0.057 * TMath::Power( TMath::Max( 0., 9 - lnE ), 1.5 ) &&
495  emf > TMath::Max( -0.6, -0.1 - 0.026 * TMath::Power( TMath::Max( 0., 8 - lnE ), 2.2 ) ) ) )
496  passCut( ret, index_TIGHT_als_ );
497 
498  if( (!HF) || ignoreCut( index_TIGHT_foot_ )
499  || jetID.fLS < 0.5 )
500  passCut( ret, index_TIGHT_foot_ );
501 
502  setIgnored( ret );
503 
504  return (bool)ret;
505  }
506 
507  private: // member variables
508 
511 
513 
517 
521 
524 
529 
537 
538 };
539 
540 #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
dictionary parameters
Definition: Parameters.py:2
short hitsInN90
Definition: JetID.h:57
void initialize(Version_t version, Quality_t quality)
JetIDSelectionFunctor(Version_t version, Quality_t quality)
Jets made from CaloTowers.
Definition: CaloJet.h:30
float approximatefHPD
Definition: JetID.h:55
unsigned int count_hits(const std::vector< CaloTowerPtr > &towers)
float fHPD
Definition: JetID.h:45
virtual bool operator()(reco::Candidate::LorentzVector const &correctedP4, double emEnergyFraction, reco::JetID const &jetID)
accessor like previous, without the ret
bool operator()(const pat::Jet &jet, pat::strbitset &ret)
This provides the interface for base classes to select objects.
float fLS
Definition: JetID.h:64
pat::strbitset::index_type index_type
Definition: Selector.h:30
float fLong
Definition: JetID.h:63
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:225
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:288
const LorentzVector & correctedP4(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.h:144
JetIDSelectionFunctor(edm::ParameterSet const &parameters)
Jet ID object.
Definition: JetID.h:16
virtual std::vector< CaloTowerPtr > getCaloConstituents() const
get all constituents
Definition: CaloJet.cc:94
bool operator()(reco::Candidate::LorentzVector const &correctedP4, double emEnergyFraction, reco::JetID const &jetID, pat::strbitset &ret)
bool craft08Cuts(reco::Candidate::LorentzVector const &correctedP4, double emEnergyFraction, reco::JetID const &jetID, pat::strbitset &ret)
bool isCaloJet() const
check to see if the jet is a reco::CaloJet
Definition: Jet.h:231
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:287
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:177
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:160
float fShort
Definition: JetID.h:63
bool operator()(reco::CaloJet const &jet, reco::JetID const &jetID, pat::strbitset &ret)
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:47
Functor that operates on &lt;T&gt;
Definition: Selector.h:25
reco::Candidate::LorentzVector correctedP4(const T &rawMEt, const CorrMETData &correction)
bool fwd09Cuts(reco::Candidate::LorentzVector const &rawP4, double emEnergyFraction, double etaWidth, double phiWidth, unsigned int nHit, reco::JetID const &jetID, pat::strbitset &ret)
float emEnergyFraction() const
returns the jet electromagnetic energy fraction
Definition: Jet.h:278
std::vector< CaloTowerPtr > const & getCaloConstituents() const
Definition: Jet.cc:141
bool isJPTJet() const
check to see if the jet is a reco::JPTJet
Definition: Jet.h:233
float etaetaMoment() const
eta-eta second moment, ET weighted
Definition: Jet.cc:224
Jet selector for pat::Jets and for CaloJets.
reco::JetID const & jetID() const
accessing Jet ID information
Definition: Jet.h:458
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
float phiphiMoment() const
phi-phi second moment, ET weighted
Definition: Jet.cc:241
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:213
virtual bool operator()(reco::CaloJet const &jet, reco::JetID const &jetID)
accessor like previous, without the ret
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
float fHFOOT
Definition: JetID.h:64
short n90Hits
Definition: JetID.h:47
float emEnergyFraction() const
Definition: CaloJet.h:98