CMS 3D CMS Logo

HiggsToWW2LeptonsSkim.cc
Go to the documentation of this file.
1 
12 
14 
15 // Muons:
17 
18 // Electrons
20 
22 
23 #include <iostream>
24 
26 
27 using namespace edm;
28 using namespace std;
29 using namespace reco;
30 
31 HiggsToWW2LeptonsSkim::HiggsToWW2LeptonsSkim(const edm::ParameterSet& iConfig) : nEvents_(0), nAccepted_(0) {
32  // Reconstructed objects
33  theGLBMuonToken = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("GlobalMuonCollectionLabel"));
34  theGsfEToken = consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("ElectronCollectionLabel"));
35 
36  singleTrackPtMin_ = iConfig.getParameter<double>("SingleTrackPtMin");
37  diTrackPtMin_ = iConfig.getParameter<double>("DiTrackPtMin");
38  etaMin_ = iConfig.getParameter<double>("etaMin");
39  etaMax_ = iConfig.getParameter<double>("etaMax");
40 }
41 
43 
45  edm::LogVerbatim("HiggsToWW2LeptonsSkim")
46  << "Events read " << nEvents_ << " Events accepted " << nAccepted_ << "\nEfficiency "
47  << ((double)nAccepted_) / ((double)nEvents_) << std::endl;
48 }
49 
51  nEvents_++;
52  bool accepted = false;
53  bool accepted1 = false;
54  int nTrackOver2ndCut = 0;
55 
56  // Handle<CandidateCollection> tracks;
57 
59 
60  // Get the muon track collection from the event
62  event.getByToken(theGLBMuonToken, muTracks);
63 
64  if (muTracks.isValid()) {
65  reco::TrackCollection::const_iterator muons;
66 
67  // Loop over muon collections and count how many muons there are,
68  // and how many are above threshold
69  for (muons = muTracks->begin(); muons != muTracks->end(); ++muons) {
70  if (muons->eta() > etaMin_ && muons->eta() < etaMax_) {
71  if (muons->pt() > singleTrackPtMin_)
72  accepted1 = true;
73  if (muons->pt() > diTrackPtMin_)
74  nTrackOver2ndCut++;
75  }
76  }
77  }
78 
79  // Now look at electrons:
80 
81  // Get the electron track collection from the event
83 
84  event.getByToken(theGsfEToken, pTracks);
85 
86  if (pTracks.isValid()) {
87  const reco::GsfElectronCollection* eTracks = pTracks.product();
88 
89  reco::GsfElectronCollection::const_iterator electrons;
90 
91  // Loop over electron collections and count how many muons there are,
92  // and how many are above threshold
93  for (electrons = eTracks->begin(); electrons != eTracks->end(); ++electrons) {
94  if (electrons->eta() > etaMin_ && electrons->eta() < etaMax_) {
95  if (electrons->pt() > singleTrackPtMin_)
96  accepted1 = true;
97  if (electrons->pt() > diTrackPtMin_)
98  nTrackOver2ndCut++;
99  }
100  }
101  }
102 
103  if (accepted1 && nTrackOver2ndCut >= 2)
104  accepted = true;
105 
106  if (accepted)
107  nAccepted_++;
108 
109  return accepted;
110 }
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
MessageLogger.h
HiggsToWW2LeptonsSkim::endJob
void endJob() override
Definition: HiggsToWW2LeptonsSkim.cc:44
edm::Handle::product
T const * product() const
Definition: Handle.h:70
edm
HLT enums.
Definition: AlignableModifier.h:19
HiggsToWW2LeptonsSkim::singleTrackPtMin_
double singleTrackPtMin_
Definition: HiggsToWW2LeptonsSkim.h:38
HiggsToWW2LeptonsSkim::HiggsToWW2LeptonsSkim
HiggsToWW2LeptonsSkim(const edm::ParameterSet &)
Definition: HiggsToWW2LeptonsSkim.cc:31
reco::GsfElectronCollection
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
Definition: GsfElectronFwd.h:14
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle< reco::TrackCollection >
HiggsToWW2LeptonsSkim::theGLBMuonToken
edm::EDGetTokenT< reco::TrackCollection > theGLBMuonToken
Definition: HiggsToWW2LeptonsSkim.h:46
HiggsToWW2LeptonsSkim.h
Track.h
cms::dd::accepted
bool accepted(std::vector< std::regex > const &, std::string_view)
GsfElectron.h
HiggsToWW2LeptonsSkim::~HiggsToWW2LeptonsSkim
~HiggsToWW2LeptonsSkim() override
Definition: HiggsToWW2LeptonsSkim.cc:42
HiggsToWW2LeptonsSkim::diTrackPtMin_
double diTrackPtMin_
Definition: HiggsToWW2LeptonsSkim.h:39
edm::ParameterSet
Definition: ParameterSet.h:36
HiggsToWW2LeptonsSkim::theGsfEToken
edm::EDGetTokenT< reco::GsfElectronCollection > theGsfEToken
Definition: HiggsToWW2LeptonsSkim.h:47
edm::LogVerbatim
Definition: MessageLogger.h:297
edm::EventSetup
Definition: EventSetup.h:57
HiggsToWW2LeptonsSkim::etaMin_
double etaMin_
Definition: HiggsToWW2LeptonsSkim.h:40
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
std
Definition: JetResolutionObject.h:76
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
HiggsToWW2LeptonsSkim::etaMax_
double etaMax_
Definition: HiggsToWW2LeptonsSkim.h:41
HiggsToWW2LeptonsSkim::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Definition: HiggsToWW2LeptonsSkim.cc:50
Candidate.h
HepMCProduct.h
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
HiggsToWW2LeptonsSkim::nAccepted_
unsigned int nAccepted_
Definition: HiggsToWW2LeptonsSkim.h:43
HiggsToWW2LeptonsSkim::nEvents_
unsigned int nEvents_
Definition: HiggsToWW2LeptonsSkim.h:42