CMS 3D CMS Logo

HiggsToZZ4LeptonsSkim.cc
Go to the documentation of this file.
1 
2 /* \class HiggsTo4LeptonsSkim
3  *
4  * Consult header file for description
5  *
6  * author: Dominique Fortin - UC Riverside
7  *
8  */
9 
10 // system include files
12 
13 // User include files
15 
16 // Muons:
18 
19 // Electrons
21 
22 // C++
23 #include <iostream>
24 #include <vector>
25 
26 using namespace std;
27 using namespace edm;
28 using namespace reco;
29 
30 // Constructor
32  // Local Debug flag
33  debug = pset.getParameter<bool>("DebugHiggsToZZ4LeptonsSkim");
34 
35  // Reconstructed objects
36  theGLBMuonToken = consumes<reco::TrackCollection>(pset.getParameter<edm::InputTag>("GlobalMuonCollectionLabel"));
37  theGsfEToken = consumes<reco::GsfElectronCollection>(pset.getParameter<edm::InputTag>("ElectronCollectionLabel"));
38 
39  // Minimum Pt for leptons for skimming
40  stiffMinPt = pset.getParameter<double>("stiffMinimumPt");
41  softMinPt = pset.getParameter<double>("softMinimumPt");
42  nStiffLeptonMin = pset.getParameter<int>("nStiffLeptonMinimum");
43  nLeptonMin = pset.getParameter<int>("nLeptonMinimum");
44 
45  nEvents = 0;
46  nSelectedEvents = 0;
47 }
48 
49 // Destructor
51  edm::LogVerbatim("HiggsToZZ4LeptonsSkim")
52  << " Number_events_read " << nEvents << " Number_events_kept " << nSelectedEvents << " Efficiency "
53  << ((double)nSelectedEvents) / ((double)nEvents + 0.01) << std::endl;
54 }
55 
56 // Filter event
58  nEvents++;
59 
61 
62  bool keepEvent = false;
63  int nStiffLeptons = 0;
64  int nLeptons = 0;
65 
66  // First look at muons:
67 
68  // Get the muon track collection from the event
70  event.getByToken(theGLBMuonToken, muTracks);
71 
72  if (muTracks.isValid()) {
73  reco::TrackCollection::const_iterator muons;
74 
75  // Loop over muon collections and count how many muons there are,
76  // and how many are above threshold
77  for (muons = muTracks->begin(); muons != muTracks->end(); ++muons) {
78  if (muons->pt() > stiffMinPt)
79  nStiffLeptons++;
80  if (muons->pt() > softMinPt)
81  nLeptons++;
82  }
83  }
84 
85  // Now look at electrons:
86 
87  // Get the electron track collection from the event
89  event.getByToken(theGsfEToken, pTracks);
90 
91  if (pTracks.isValid()) {
92  const reco::GsfElectronCollection* eTracks = pTracks.product();
93 
94  reco::GsfElectronCollection::const_iterator electrons;
95 
96  // Loop over electron collections and count how many muons there are,
97  // and how many are above threshold
98  for (electrons = eTracks->begin(); electrons != eTracks->end(); ++electrons) {
99  float pt_e = electrons->pt();
100  if (pt_e > stiffMinPt)
101  nStiffLeptons++;
102  if (pt_e > softMinPt)
103  nLeptons++;
104  }
105  }
106 
107  // Make decision:
108  if (nStiffLeptons >= nStiffLeptonMin && nLeptons >= nLeptonMin)
109  keepEvent = true;
110 
111  if (keepEvent)
112  nSelectedEvents++;
113 
114  return keepEvent;
115 }
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
HiggsToZZ4LeptonsSkim::HiggsToZZ4LeptonsSkim
HiggsToZZ4LeptonsSkim(const edm::ParameterSet &)
Definition: HiggsToZZ4LeptonsSkim.cc:31
edm::Handle::product
T const * product() const
Definition: Handle.h:70
edm
HLT enums.
Definition: AlignableModifier.h:19
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 >
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
debug
#define debug
Definition: HDRShower.cc:19
Track.h
GsfElectron.h
edm::ParameterSet
Definition: ParameterSet.h:47
HiggsToZZ4LeptonsSkim::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Get event properties to send to builder to fill seed collection.
Definition: HiggsToZZ4LeptonsSkim.cc:57
edm::EventSetup
Definition: EventSetup.h:57
std
Definition: JetResolutionObject.h:76
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HiggsToZZ4LeptonsSkim.h
ParameterSet.h
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
event
Definition: event.py:1
nEvents
UInt_t nEvents
Definition: hcalCalib.cc:41
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
HiggsToZZ4LeptonsSkim::~HiggsToZZ4LeptonsSkim
~HiggsToZZ4LeptonsSkim() override
Definition: HiggsToZZ4LeptonsSkim.cc:50
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27