CMS 3D CMS Logo

Public Member Functions | Protected Attributes

WSelector Class Reference

Example class of an EventSelector to apply a simple W Boson selection. More...

#include <PhysicsTools/FWLite/interface/WSelector.h>

Inheritance diagram for WSelector:
EventSelector EventSelector

List of all members.

Public Member Functions

pat::MET const & met () const
 return MET of W boson
pat::MET const & met () const
 return MET of W boson
virtual bool operator() (edm::EventBase const &event, pat::strbitset &ret)
 here is where the selection occurs
virtual bool operator() (edm::EventBase const &event, pat::strbitset &ret)
 here is where the selection occurs
pat::Muon const & wMuon () const
 return muon candidate of W boson
pat::Muon const & wMuon () const
 return muon candidate of W boson
 WSelector (edm::ParameterSet const &params)
 constructor
 WSelector (edm::ParameterSet const &params)
 constructor
virtual ~WSelector ()
 destructor
virtual ~WSelector ()
 destructor

Protected Attributes

pat::MET const * met_
 MET from W boson.
index_type metIndex_
 index for MET cut
edm::InputTag metSrc_
 met input
index_type muonPtIndex_
 index for muon Pt cut
edm::InputTag muonSrc_
 muon input
pat::Muon const * wMuon_
 muon candidate from W boson

Detailed Description

Example class of an EventSelector to apply a simple W Boson selection.

Example class for of an EventSelector as defined in the SelectorUtils package. EventSelectors may be used to facilitate cutflows and to implement selections independent from the event loop in FWLite or full framework.

Definition at line 15 of file WSelector.h.


Constructor & Destructor Documentation

WSelector::WSelector ( edm::ParameterSet const &  params) [inline]

constructor

Definition at line 19 of file WSelector.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), met_, and wMuon_.

                                           :
    muonSrc_(params.getParameter<edm::InputTag>("muonSrc")),
    metSrc_ (params.getParameter<edm::InputTag>("metSrc")) 
  {
    double muonPtMin = params.getParameter<double>("muonPtMin");
    double metMin    = params.getParameter<double>("metMin");
    push_back("Muon Pt", muonPtMin );
    push_back("MET"    , metMin    );
    set("Muon Pt"); set("MET");
    wMuon_ = 0; met_ = 0;
    if ( params.exists("cutsToIgnore") ){
      setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
    }
    retInternal_ = getBitTemplate();
  }
virtual WSelector::~WSelector ( ) [inline, virtual]

destructor

Definition at line 35 of file WSelector.h.

{}
WSelector::WSelector ( edm::ParameterSet const &  params) [inline]

constructor

Definition at line 19 of file WSelectorFast.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), met_, and wMuon_.

                                           :
    muonSrc_(params.getParameter<edm::InputTag>("muonSrc")),
    metSrc_ (params.getParameter<edm::InputTag>("metSrc")) 
  {
    double muonPtMin = params.getParameter<double>("muonPtMin");
    double metMin    = params.getParameter<double>("metMin");
    push_back("Muon Pt", muonPtMin );
    push_back("MET"    , metMin    );
    set("Muon Pt"); set("MET");
    wMuon_ = 0; met_ = 0;
    if ( params.exists("cutsToIgnore") ){
      setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
    }
    retInternal_ = getBitTemplate();
  }
virtual WSelector::~WSelector ( ) [inline, virtual]

destructor

Definition at line 35 of file WSelectorFast.h.

{}

Member Function Documentation

pat::MET const& WSelector::met ( ) const [inline]

return MET of W boson

Definition at line 39 of file WSelector.h.

References met_.

Referenced by operator()().

{ return *met_;  }
pat::MET const& WSelector::met ( ) const [inline]

return MET of W boson

Definition at line 39 of file WSelectorFast.h.

References met_.

{ return *met_;  }
virtual bool WSelector::operator() ( edm::EventBase const &  event,
pat::strbitset ret 
) [inline, virtual]

here is where the selection occurs

Definition at line 42 of file WSelectorFast.h.

References GOODCOLL_filter_cfg::cut, met(), met_, metIndex_, metSrc_, muonPtIndex_, patZpeak::muons, muonSrc_, reco::LeafCandidate::pt(), run_regression::ret, pat::strbitset::set(), and wMuon_.

                                                                          {
    ret.set(false);
    // Handle to the muon collection
    edm::Handle<std::vector<pat::Muon> > muons;    
    // Handle to the MET collection
    edm::Handle<std::vector<pat::MET> > met;
    // get the objects from the event
    bool gotMuons = event.getByLabel(muonSrc_, muons);
    bool gotMET   = event.getByLabel(metSrc_, met   );
    // get the MET, require to be > minimum
    if( gotMET ){
      met_ = &met->at(0);
      if( met_->pt() > cut(metIndex_, double()) || ignoreCut(metIndex_) ) 
        passCut(ret, metIndex_);
    }
    // get the highest pt muon, require to have pt > minimum
    if( gotMuons ){
      if( !ignoreCut(muonPtIndex_) ){
        if( muons->size() > 0 ){
          wMuon_ = &muons->at(0);
          if ( wMuon_->pt() > cut(muonPtIndex_, double()) || ignoreCut(muonPtIndex_) ) 
            passCut(ret, muonPtIndex_);
        }
      } 
      else{
        passCut( ret, muonPtIndex_);
      }
    }
    setIgnored(ret);
    return (bool)ret;
  }
virtual bool WSelector::operator() ( edm::EventBase const &  event,
pat::strbitset ret 
) [inline, virtual]

here is where the selection occurs

Definition at line 42 of file WSelector.h.

References GOODCOLL_filter_cfg::cut, met(), met_, metSrc_, patZpeak::muons, muonSrc_, reco::LeafCandidate::pt(), run_regression::ret, pat::strbitset::set(), and wMuon_.

                                                                          {
    ret.set(false);
    // Handle to the muon collection
    edm::Handle<std::vector<pat::Muon> > muons;    
    // Handle to the MET collection
    edm::Handle<std::vector<pat::MET> > met;
    // get the objects from the event
    bool gotMuons = event.getByLabel(muonSrc_, muons);
    bool gotMET   = event.getByLabel(metSrc_, met   );
    // get the MET, require to be > minimum
    if( gotMET ){
      met_ = &met->at(0);
      if( met_->pt() > cut("MET", double()) || ignoreCut("MET") ) 
        passCut(ret, "MET");
    }
    // get the highest pt muon, require to have pt > minimum
    if( gotMuons ){
      if( !ignoreCut("Muon Pt") ){
        if( muons->size() > 0 ){
          wMuon_ = &muons->at(0);
          if( wMuon_->pt() > cut("Muon Pt", double()) || ignoreCut("Muon Pt") ) 
            passCut(ret, "Muon Pt");
        }
      } 
      else{
        passCut( ret, "Muon Pt");
      }
    }
    setIgnored(ret);
    return (bool)ret;
  }
pat::Muon const& WSelector::wMuon ( ) const [inline]

return muon candidate of W boson

Definition at line 37 of file WSelectorFast.h.

References wMuon_.

{ return *wMuon_;}
pat::Muon const& WSelector::wMuon ( ) const [inline]

return muon candidate of W boson

Definition at line 37 of file WSelector.h.

References wMuon_.

{ return *wMuon_;}

Member Data Documentation

pat::MET const * WSelector::met_ [protected]

MET from W boson.

Definition at line 82 of file WSelector.h.

Referenced by met(), operator()(), and WSelector().

index_type WSelector::metIndex_ [protected]

index for MET cut

Definition at line 86 of file WSelectorFast.h.

Referenced by operator()().

met input

Definition at line 78 of file WSelector.h.

Referenced by operator()().

index_type WSelector::muonPtIndex_ [protected]

index for muon Pt cut

Definition at line 84 of file WSelectorFast.h.

Referenced by operator()().

muon input

Definition at line 76 of file WSelector.h.

Referenced by operator()().

pat::Muon const * WSelector::wMuon_ [protected]

muon candidate from W boson

Definition at line 80 of file WSelector.h.

Referenced by operator()(), wMuon(), and WSelector().