#include <HiggsToWW2LeptonsSkim.h>
Public Member Functions | |
virtual void | endJob () |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HiggsToWW2LeptonsSkim (const edm::ParameterSet &) | |
~HiggsToWW2LeptonsSkim () | |
Private Attributes | |
double | diTrackPtMin_ |
double | etaMax_ |
double | etaMin_ |
unsigned int | nAccepted_ |
unsigned int | nEvents_ |
edm::InputTag | recTrackLabel |
double | singleTrackPtMin_ |
edm::InputTag | theGLBMuonLabel |
edm::InputTag | theGsfELabel |
This class is an EDFilter for HWW events
Definition at line 28 of file HiggsToWW2LeptonsSkim.h.
HiggsToWW2LeptonsSkim::HiggsToWW2LeptonsSkim | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 37 of file HiggsToWW2LeptonsSkim.cc.
References diTrackPtMin_, etaMax_, etaMin_, edm::ParameterSet::getParameter(), recTrackLabel, singleTrackPtMin_, theGLBMuonLabel, and theGsfELabel.
: nEvents_(0), nAccepted_(0) { // Reconstructed objects recTrackLabel = iConfig.getParameter<edm::InputTag>("RecoTrackLabel"); theGLBMuonLabel = iConfig.getParameter<edm::InputTag>("GlobalMuonCollectionLabel"); theGsfELabel = iConfig.getParameter<edm::InputTag>("ElectronCollectionLabel"); singleTrackPtMin_ = iConfig.getParameter<double>("SingleTrackPtMin"); diTrackPtMin_ = iConfig.getParameter<double>("DiTrackPtMin"); etaMin_ = iConfig.getParameter<double>("etaMin"); etaMax_ = iConfig.getParameter<double>("etaMax"); }
HiggsToWW2LeptonsSkim::~HiggsToWW2LeptonsSkim | ( | ) |
Definition at line 53 of file HiggsToWW2LeptonsSkim.cc.
{ }
void HiggsToWW2LeptonsSkim::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDFilter.
Definition at line 57 of file HiggsToWW2LeptonsSkim.cc.
References nAccepted_, and nEvents_.
{ edm::LogVerbatim("HiggsToWW2LeptonsSkim") << "Events read " << nEvents_ << " Events accepted " << nAccepted_ << "\nEfficiency " << ((double)nAccepted_)/((double)nEvents_) << std::endl; }
bool HiggsToWW2LeptonsSkim::filter | ( | edm::Event & | event, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDFilter.
Definition at line 67 of file HiggsToWW2LeptonsSkim.cc.
References diTrackPtMin_, etaMax_, etaMin_, edm::HandleBase::isValid(), edm::InputTag::label(), ExpressReco_HICollisions_FallBack::muons, nAccepted_, nEvents_, edm::Handle< T >::product(), singleTrackPtMin_, theGLBMuonLabel, and theGsfELabel.
{ nEvents_++; bool accepted = false; bool accepted1 = false; int nTrackOver2ndCut = 0; // Handle<CandidateCollection> tracks; using reco::TrackCollection; // Get the muon track collection from the event edm::Handle<reco::TrackCollection> muTracks; event.getByLabel(theGLBMuonLabel.label(), muTracks); if ( muTracks.isValid() ) { reco::TrackCollection::const_iterator muons; // Loop over muon collections and count how many muons there are, // and how many are above threshold for ( muons = muTracks->begin(); muons != muTracks->end(); ++muons ) { if ( muons->eta() > etaMin_ && muons->eta() < etaMax_ ) { if ( muons->pt() > singleTrackPtMin_ ) accepted1 = true; if ( muons->pt() > diTrackPtMin_ ) nTrackOver2ndCut++; } } } // Now look at electrons: // Get the electron track collection from the event edm::Handle<reco::GsfElectronCollection> pTracks; event.getByLabel(theGsfELabel.label(),pTracks); if ( pTracks.isValid() ) { const reco::GsfElectronCollection* eTracks = pTracks.product(); reco::GsfElectronCollection::const_iterator electrons; // Loop over electron collections and count how many muons there are, // and how many are above threshold for ( electrons = eTracks->begin(); electrons != eTracks->end(); ++electrons ) { if ( electrons->eta() > etaMin_ && electrons->eta() < etaMax_ ) { if ( electrons->pt() > singleTrackPtMin_ ) accepted1 = true; if ( electrons->pt() > diTrackPtMin_ ) nTrackOver2ndCut++; } } } if ( accepted1 && nTrackOver2ndCut >= 2 ) accepted = true; if ( accepted ) nAccepted_++; return accepted; }
double HiggsToWW2LeptonsSkim::diTrackPtMin_ [private] |
Definition at line 38 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().
double HiggsToWW2LeptonsSkim::etaMax_ [private] |
Definition at line 40 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().
double HiggsToWW2LeptonsSkim::etaMin_ [private] |
Definition at line 39 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().
unsigned int HiggsToWW2LeptonsSkim::nAccepted_ [private] |
Definition at line 42 of file HiggsToWW2LeptonsSkim.h.
unsigned int HiggsToWW2LeptonsSkim::nEvents_ [private] |
Definition at line 41 of file HiggsToWW2LeptonsSkim.h.
Definition at line 45 of file HiggsToWW2LeptonsSkim.h.
Referenced by HiggsToWW2LeptonsSkim().
double HiggsToWW2LeptonsSkim::singleTrackPtMin_ [private] |
Definition at line 37 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().
Definition at line 46 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().
Definition at line 47 of file HiggsToWW2LeptonsSkim.h.
Referenced by filter(), and HiggsToWW2LeptonsSkim().