CMS 3D CMS Logo

WPlusJetsEventSelector.cc
Go to the documentation of this file.
1 
4 
5 #include <iostream>
6 
7 using namespace std;
8 
10  EventSelector(),
11  muonTag_ (params.getParameter<edm::InputTag>("muonSrc") ),
12  electronTag_ (params.getParameter<edm::InputTag>("electronSrc") ),
13  jetTag_ (params.getParameter<edm::InputTag>("jetSrc") ),
14  metTag_ (params.getParameter<edm::InputTag>("metSrc") ),
15  trigTag_ (params.getParameter<edm::InputTag>("trigSrc") ),
16  muTrig_ (params.getParameter<std::string>("muTrig")),
17  eleTrig_ (params.getParameter<std::string>("eleTrig")),
18  pvSelector_ (params.getParameter<edm::ParameterSet>("pvSelector") ),
19  muonIdTight_ (params.getParameter<edm::ParameterSet>("muonIdTight") ),
20  electronIdTight_ (params.getParameter<edm::ParameterSet>("electronIdTight") ),
21  muonIdLoose_ (params.getParameter<edm::ParameterSet>("muonIdLoose") ),
22  electronIdLoose_ (params.getParameter<edm::ParameterSet>("electronIdLoose") ),
23  jetIdLoose_ (params.getParameter<edm::ParameterSet>("jetIdLoose") ),
24  pfjetIdLoose_ (params.getParameter<edm::ParameterSet>("pfjetIdLoose") ),
25  minJets_ (params.getParameter<int> ("minJets") ),
26  muJetDR_ (params.getParameter<double>("muJetDR")),
27  eleJetDR_ (params.getParameter<double>("eleJetDR")),
28  muPlusJets_ (params.getParameter<bool>("muPlusJets") ),
29  ePlusJets_ (params.getParameter<bool>("ePlusJets") ),
30  muPtMin_ (params.getParameter<double>("muPtMin")),
31  muEtaMax_ (params.getParameter<double>("muEtaMax")),
32  eleEtMin_ (params.getParameter<double>("eleEtMin")),
33  eleEtaMax_ (params.getParameter<double>("eleEtaMax")),
34  muPtMinLoose_ (params.getParameter<double>("muPtMinLoose")),
35  muEtaMaxLoose_ (params.getParameter<double>("muEtaMaxLoose")),
36  eleEtMinLoose_ (params.getParameter<double>("eleEtMinLoose")),
37  eleEtaMaxLoose_ (params.getParameter<double>("eleEtaMaxLoose")),
38  jetPtMin_ (params.getParameter<double>("jetPtMin")),
39  jetEtaMax_ (params.getParameter<double>("jetEtaMax")),
40  jetScale_ (params.getParameter<double>("jetScale")),
41  metMin_ (params.getParameter<double>("metMin"))
42 {
43  // make the bitset
44  push_back( "Inclusive" );
45  push_back( "Trigger" );
46  push_back( "PV" );
47  push_back( ">= 1 Lepton" );
48  push_back( "== 1 Tight Lepton" );
49  push_back( "== 1 Tight Lepton, Mu Veto");
50  push_back( "== 1 Lepton" );
51  push_back( "MET Cut" );
52  push_back( "Z Veto" );
53  push_back( "Conversion Veto");
54  push_back( "Cosmic Veto" );
55  push_back( ">=1 Jets" );
56  push_back( ">=2 Jets" );
57  push_back( ">=3 Jets" );
58  push_back( ">=4 Jets" );
59  push_back( ">=5 Jets" );
60 
61 
62  // turn (almost) everything on by default
63  set( "Inclusive" );
64  set( "Trigger" );
65  set( "PV" );
66  set( ">= 1 Lepton" );
67  set( "== 1 Tight Lepton" );
68  set( "== 1 Tight Lepton, Mu Veto");
69  set( "== 1 Lepton" );
70  set( "MET Cut" );
71  set( "Z Veto" );
72  set( "Conversion Veto");
73  set( "Cosmic Veto" );
74  set( ">=1 Jets", minJets_ >= 1);
75  set( ">=2 Jets", minJets_ >= 2);
76  set( ">=3 Jets", minJets_ >= 3);
77  set( ">=4 Jets", minJets_ >= 4);
78  set( ">=5 Jets", minJets_ >= 5);
79 
80 
81  inclusiveIndex_ = index_type(&bits_, std::string("Inclusive" ));
82  triggerIndex_ = index_type(&bits_, std::string("Trigger" ));
83  pvIndex_ = index_type(&bits_, std::string("PV" ));
84  lep1Index_ = index_type(&bits_, std::string(">= 1 Lepton" ));
85  lep2Index_ = index_type(&bits_, std::string("== 1 Tight Lepton" ));
86  lep3Index_ = index_type(&bits_, std::string("== 1 Tight Lepton, Mu Veto"));
87  lep4Index_ = index_type(&bits_, std::string("== 1 Lepton" ));
88  metIndex_ = index_type(&bits_, std::string("MET Cut" ));
89  zvetoIndex_ = index_type(&bits_, std::string("Z Veto" ));
90  conversionIndex_ = index_type(&bits_, std::string("Conversion Veto"));
91  cosmicIndex_ = index_type(&bits_, std::string("Cosmic Veto" ));
92  jet1Index_ = index_type(&bits_, std::string(">=1 Jets"));
93  jet2Index_ = index_type(&bits_, std::string(">=2 Jets"));
94  jet3Index_ = index_type(&bits_, std::string(">=3 Jets"));
95  jet4Index_ = index_type(&bits_, std::string(">=4 Jets"));
96  jet5Index_ = index_type(&bits_, std::string(">=5 Jets"));
97 
98  if ( params.exists("cutsToIgnore") )
99  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
100 
101 
102  retInternal_ = getBitTemplate();
103 }
104 
106 {
107 
108  ret.set(false);
109 
110  selectedJets_.clear();
111  cleanedJets_.clear();
112  selectedMuons_.clear();
113  selectedElectrons_.clear();
114  looseMuons_.clear();
115  looseElectrons_.clear();
116  selectedMETs_.clear();
117 
118 
119  passCut( ret, inclusiveIndex_);
120 
121 
122  bool passTrig = false;
123  if (!ignoreCut(triggerIndex_) ) {
124 
125 
126  edm::Handle<pat::TriggerEvent> triggerEvent;
127  event.getByLabel(trigTag_, triggerEvent);
128 
129  pat::TriggerEvent const * trig = &*triggerEvent;
130 
131  if ( trig->wasRun() && trig->wasAccept() ) {
132 
133  pat::TriggerPath const * muPath = trig->path(muTrig_);
134 
135  pat::TriggerPath const * elePath = trig->path(eleTrig_);
136 
137  if ( muPlusJets_ && muPath != 0 && muPath->wasAccept() ) {
138  passTrig = true;
139  }
140 
141  if ( ePlusJets_ && elePath != 0 && elePath->wasAccept() ) {
142  passTrig = true;
143  }
144  }
145  }
146 
147 
148 
149 
150  if ( ignoreCut(triggerIndex_) ||
151  passTrig ) {
152  passCut(ret, triggerIndex_);
153 
154 
155  bool passPV = false;
156 
157  passPV = pvSelector_( event );
158  if ( ignoreCut(pvIndex_) || passPV ) {
159  passCut(ret, pvIndex_);
160 
161  edm::Handle< vector< pat::Electron > > electronHandle;
162  event.getByLabel (electronTag_, electronHandle);
163 
165  event.getByLabel (muonTag_, muonHandle);
166 
168 
170 
172  event.getByLabel (metTag_, metHandle);
173 
174  int nElectrons = 0;
175  for ( std::vector<pat::Electron>::const_iterator electronBegin = electronHandle->begin(),
176  electronEnd = electronHandle->end(), ielectron = electronBegin;
177  ielectron != electronEnd; ++ielectron ) {
178  ++nElectrons;
179  // Tight cuts
180  if ( ielectron->et() > eleEtMin_ && fabs(ielectron->eta()) < eleEtaMax_ &&
181  electronIdTight_(*ielectron) &&
182  ielectron->electronID( "eidRobustTight" ) > 0 ) {
183  selectedElectrons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Electron>( electronHandle, ielectron - electronBegin ) ) );
184  } else {
185  // Loose cuts
186  if ( ielectron->et() > eleEtMinLoose_ && fabs(ielectron->eta()) < eleEtaMaxLoose_ &&
187  electronIdLoose_(*ielectron) ) {
188  looseElectrons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Electron>( electronHandle, ielectron - electronBegin ) ) );
189  }
190  }
191  }
192 
193 
194  for ( std::vector<pat::Muon>::const_iterator muonBegin = muonHandle->begin(),
195  muonEnd = muonHandle->end(), imuon = muonBegin;
196  imuon != muonEnd; ++imuon ) {
197  if ( !imuon->isGlobalMuon() ) continue;
198 
199  // Tight cuts
200  bool passTight = muonIdTight_(*imuon,event) && imuon->isTrackerMuon() ;
201  if ( imuon->pt() > muPtMin_ && fabs(imuon->eta()) < muEtaMax_ &&
202  passTight ) {
203 
204  selectedMuons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Muon>( muonHandle, imuon - muonBegin ) ) );
205  } else {
206  // Loose cuts
207  if ( imuon->pt() > muPtMinLoose_ && fabs(imuon->eta()) < muEtaMaxLoose_ &&
208  muonIdLoose_(*imuon,event) ) {
209  looseMuons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Muon>( muonHandle, imuon - muonBegin ) ) );
210  }
211  }
212  }
213 
214 
216  metHandle->at(0).charge(),
217  metHandle->at(0).p4() );
218 
219 
220 
221  event.getByLabel (jetTag_, jetHandle);
224  for ( std::vector<pat::Jet>::const_iterator jetBegin = jetHandle->begin(),
225  jetEnd = jetHandle->end(), ijet = jetBegin;
226  ijet != jetEnd; ++ijet ) {
227  reco::ShallowClonePtrCandidate scaledJet ( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Jet>( jetHandle, ijet - jetBegin ),
228  ijet->charge(),
229  ijet->p4() * jetScale_ ) );
230  bool passJetID = false;
231  if ( ijet->isCaloJet() || ijet->isJPTJet() ) passJetID = jetIdLoose_(*ijet, ret1);
232  else passJetID = pfjetIdLoose_(*ijet, ret2);
233  if ( scaledJet.pt() > jetPtMin_ && fabs(scaledJet.eta()) < jetEtaMax_ && passJetID ) {
234  selectedJets_.push_back( scaledJet );
235 
236  if ( muPlusJets_ ) {
237 
238  //Remove some jets
239  bool indeltaR = false;
240  for( std::vector<reco::ShallowClonePtrCandidate>::const_iterator muonBegin = selectedMuons_.begin(),
241  muonEnd = selectedMuons_.end(), imuon = muonBegin;
242  imuon != muonEnd; ++imuon ) {
243  if( reco::deltaR( imuon->eta(), imuon->phi(), scaledJet.eta(), scaledJet.phi() ) < muJetDR_ )
244  { indeltaR = true; }
245  }
246  if( !indeltaR ) {
247  cleanedJets_.push_back( scaledJet );
248  }// end if jet is not within dR of a muon
249  }// end if mu+jets
250  else {
251  //Remove some jets
252  bool indeltaR = false;
253  for( std::vector<reco::ShallowClonePtrCandidate>::const_iterator electronBegin = selectedElectrons_.begin(),
254  electronEnd = selectedElectrons_.end(), ielectron = electronBegin;
255  ielectron != electronEnd; ++ielectron ) {
256  if( reco::deltaR( ielectron->eta(), ielectron->phi(), scaledJet.eta(), scaledJet.phi() ) < eleJetDR_ )
257  { indeltaR = true; }
258  }
259  if( !indeltaR ) {
260  cleanedJets_.push_back( scaledJet );
261  }// end if jet is not within dR of an electron
262  }// end if e+jets
263  }// end if pass id and kin cuts
264  }// end loop over jets
265 
266 
267 
268  int nleptons = 0;
269  if ( muPlusJets_ )
270  nleptons += selectedMuons_.size();
271 
272  if ( ePlusJets_ )
273  nleptons += selectedElectrons_.size();
274 
275  if ( ignoreCut(lep1Index_) ||
276  ( nleptons > 0 ) ){
277  passCut( ret, lep1Index_);
278 
279  if ( ignoreCut(lep2Index_) ||
280  ( nleptons == 1 ) ){
281  passCut( ret, lep2Index_);
282 
283  bool oneMuon =
284  ( selectedMuons_.size() == 1 &&
285  looseMuons_.size() + selectedElectrons_.size() + looseElectrons_.size() == 0
286  );
287  bool oneElectron =
288  ( selectedElectrons_.size() == 1 &&
289  selectedMuons_.size() == 0
290  );
291 
292  bool oneMuonMuVeto =
293  ( selectedMuons_.size() == 1 &&
294  looseMuons_.size() == 0
295  );
296 
297 
298  if ( ignoreCut(lep3Index_) ||
299  ePlusJets_ ||
300  (muPlusJets_ && oneMuonMuVeto)
301  ) {
302  passCut(ret, lep3Index_);
303 
304  if ( ignoreCut(lep4Index_) ||
305  ( (muPlusJets_ && oneMuon) ^ (ePlusJets_ && oneElectron ) )
306  ) {
307  passCut(ret, lep4Index_);
308 
309  bool metCut = met_.pt() > metMin_;
310  if ( ignoreCut(metIndex_) ||
311  metCut ) {
312  passCut( ret, metIndex_ );
313 
314 
315  bool zVeto = true;
316  if ( selectedMuons_.size() == 2 ) {
317  }
318  if ( selectedElectrons_.size() == 2 ) {
319  }
320  if ( ignoreCut(zvetoIndex_) ||
321  zVeto ){
322  passCut(ret, zvetoIndex_);
323 
324 
325  bool conversionVeto = true;
326  if ( ignoreCut(conversionIndex_) ||
327  conversionVeto ) {
328  passCut(ret,conversionIndex_);
329 
330 
331 
332  bool cosmicVeto = true;
333  if ( ignoreCut(cosmicIndex_) ||
334  cosmicVeto ) {
335  passCut(ret,cosmicIndex_);
336 
337  if ( ignoreCut(jet1Index_) ||
338  static_cast<int>(cleanedJets_.size()) >= 1 ){
339  passCut(ret,jet1Index_);
340  } // end if >=1 tight jets
341 
342  if ( ignoreCut(jet2Index_) ||
343  static_cast<int>(cleanedJets_.size()) >= 2 ){
344  passCut(ret,jet2Index_);
345  } // end if >=2 tight jets
346 
347  if ( ignoreCut(jet3Index_) ||
348  static_cast<int>(cleanedJets_.size()) >= 3 ){
349  passCut(ret,jet3Index_);
350  } // end if >=3 tight jets
351 
352  if ( ignoreCut(jet4Index_) ||
353  static_cast<int>(cleanedJets_.size()) >= 4 ){
354  passCut(ret,jet4Index_);
355  } // end if >=4 tight jets
356 
357  if ( ignoreCut(jet5Index_) ||
358  static_cast<int>(cleanedJets_.size()) >= 5 ){
359  passCut(ret,jet5Index_);
360  } // end if >=5 tight jets
361 
362 
363 
364  } // end if cosmic veto
365 
366  } // end if conversion veto
367 
368  } // end if z veto
369 
370  } // end if met cut
371 
372  } // end if == 1 lepton
373 
374  } // end if == 1 tight lepton with a muon veto separately
375 
376  } // end if == 1 tight lepton
377 
378  } // end if >= 1 lepton
379 
380  } // end if PV
381 
382  } // end if trigger
383 
384 
385  setIgnored(ret);
386 
387  return (bool)ret;
388 }
std::vector< reco::ShallowClonePtrCandidate > selectedElectrons_
bool wasAccept() const
Get the success flag.
Definition: TriggerEvent.h:151
T getParameter(std::string const &) const
virtual double pt() const final
transverse momentum
ElectronVPlusJetsIDSelectionFunctor electronIdLoose_
ElectronVPlusJetsIDSelectionFunctor electronIdTight_
Selector< edm::EventBase > EventSelector
Definition: EventSelector.h:20
virtual double eta() const final
momentum pseudorapidity
PFJetIDSelectionFunctor pfjetIdLoose_
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::ShallowClonePtrCandidate met_
JetIDSelectionFunctor jetIdLoose_
virtual double phi() const final
momentum azimuthal angle
virtual int charge() const final
electric charge
Definition: LeafCandidate.h:91
Analysis-level HLTrigger path class.
Definition: TriggerPath.h:39
std::vector< reco::ShallowClonePtrCandidate > selectedJets_
const TriggerPath * path(const std::string &namePath) const
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
std::vector< reco::ShallowClonePtrCandidate > selectedMuons_
Analysis-level trigger event class.
Definition: TriggerEvent.h:42
virtual bool operator()(edm::EventBase const &t, pat::strbitset &ret)
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:144
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:212
std::vector< reco::ShallowClonePtrCandidate > looseElectrons_
HLT enums.
std::vector< reco::ShallowClonePtrCandidate > selectedMETs_
std::vector< reco::ShallowClonePtrCandidate > cleanedJets_
MuonVPlusJetsIDSelectionFunctor muonIdLoose_
MuonVPlusJetsIDSelectionFunctor muonIdTight_
Definition: event.py:1
bool wasRun() const
Get the run flag.
Definition: TriggerEvent.h:149
bool wasAccept() const
Get the success flag.
Definition: TriggerPath.h:121
bool passTrig(const float objEta, float objPhi, const trigger::TriggerEvent &trigEvt, const std::string &filterName, const std::string &processName)
Definition: UtilFuncs.h:14
std::vector< reco::ShallowClonePtrCandidate > looseMuons_