CMS 3D CMS Logo

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