CMS 3D CMS Logo

WSelectorFast.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_FWLite_WSelectorFast_h
2 #define PhysicsTools_FWLite_WSelectorFast_h
6 
17 class WSelector : public EventSelector {
18 
19 public:
22  muonSrc_(params.getParameter<edm::InputTag>("muonSrc")),
23  metSrc_ (params.getParameter<edm::InputTag>("metSrc"))
24  {
25  double muonPtMin = params.getParameter<double>("muonPtMin");
26  double metMin = params.getParameter<double>("metMin");
27  push_back("Muon Pt", muonPtMin );
28  push_back("MET" , metMin );
29  set("Muon Pt"); set("MET");
30  wMuon_ = 0; met_ = 0;
31  if ( params.exists("cutsToIgnore") ){
32  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
33  }
34  retInternal_ = getBitTemplate();
35  }
37  virtual ~WSelector() {}
39  pat::Muon const& wMuon() const { return *wMuon_;}
41  pat::MET const& met() const { return *met_; }
42 
44  virtual bool operator()( edm::EventBase const & event, pat::strbitset & ret){
45  ret.set(false);
46  // Handle to the muon collection
48  // Handle to the MET collection
50  // get the objects from the event
51  bool gotMuons = event.getByLabel(muonSrc_, muons);
52  bool gotMET = event.getByLabel(metSrc_, met );
53  // get the MET, require to be > minimum
54  if( gotMET ){
55  met_ = &met->at(0);
56  if( met_->pt() > cut(metIndex_, double()) || ignoreCut(metIndex_) )
57  passCut(ret, metIndex_);
58  }
59  // get the highest pt muon, require to have pt > minimum
60  if( gotMuons ){
61  if( !ignoreCut(muonPtIndex_) ){
62  if( muons->size() > 0 ){
63  wMuon_ = &muons->at(0);
64  if ( wMuon_->pt() > cut(muonPtIndex_, double()) || ignoreCut(muonPtIndex_) )
65  passCut(ret, muonPtIndex_);
66  }
67  }
68  else{
69  passCut( ret, muonPtIndex_);
70  }
71  }
72  setIgnored(ret);
73  return (bool)ret;
74  }
75 
76 protected:
82  pat::Muon const* wMuon_;
84  pat::MET const* met_;
86  index_type muonPtIndex_;
88  index_type metIndex_;
89 };
90 #endif // PhysicsTools_FWLite_WSelectorFast_h
Analysis-level MET class.
Definition: MET.h:40
T getParameter(std::string const &) const
index_type metIndex_
index for MET cut
Definition: WSelectorFast.h:88
ret
prodAgent to be discontinued
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::InputTag muonSrc_
muon input
Definition: WSelector.h:75
double pt() const final
transverse momentum
edm::InputTag metSrc_
met input
Definition: WSelector.h:77
pat::Muon const & wMuon() const
return muon candidate of W boson
Definition: WSelectorFast.h:39
WSelector(edm::ParameterSet const &params)
constructor
Definition: WSelectorFast.h:21
pat::MET const & met() const
return MET of W boson
Definition: WSelectorFast.h:41
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:126
A selector of events.
Definition: EventSelector.h:16
virtual ~WSelector()
destructor
Definition: WSelectorFast.h:37
HLT enums.
virtual bool operator()(edm::EventBase const &event, pat::strbitset &ret)
here is where the selection occurs
Definition: WSelectorFast.h:44
pat::MET const * met_
MET from W boson.
Definition: WSelector.h:81
index_type muonPtIndex_
index for muon Pt cut
Definition: WSelectorFast.h:86
pat::Muon const * wMuon_
muon candidate from W boson
Definition: WSelector.h:79
Example class of an EventSelector to apply a simple W Boson selection.
Definition: WSelector.h:16
Analysis-level muon class.
Definition: Muon.h:51
Definition: event.py:1