CMS 3D CMS Logo

LeptonRecoSkim.cc
Go to the documentation of this file.
1 
2 // -*- C++ -*-
3 //
4 // Package: LeptonRecoSkim
5 // Class: LeptonRecoSkim
6 //
14 //
15 // Original Author: Massimiliano Chiorboli,40 4-A01,+41227671535,
16 // Created: Wed Mar 31 21:49:08 CEST 2010
17 // $Id: LeptonRecoSkim.cc,v 1.1 2010/11/05 18:37:51 torimoto Exp $
18 //
19 //
20 
22 
23 using namespace edm;
24 using namespace reco;
25 using namespace std;
26 
27 //
28 // constructors and destructor
29 //
31  : m_CaloGeoToken(esConsumes()),
32  m_CaloTopoToken(esConsumes()),
33  hltLabel(iConfig.getParameter<edm::InputTag>("HltLabel")),
34  filterName(iConfig.getParameter<std::string>("@module_label")),
35  gsfElectronCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("electronCollection"))),
36  pfCandidateCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("pfElectronCollection"))),
37  muonCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("muonCollection"))),
38  caloJetCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("caloJetCollection"))),
39  pfJetCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("PFJetCollection"))),
40  ebRecHitCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("ecalBarrelRecHitsCollection"))),
41  eeRecHitCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("ecalEndcapRecHitsCollection"))),
42  useElectronSelection(iConfig.getParameter<bool>("UseElectronSelection")),
43  usePfElectronSelection(iConfig.getParameter<bool>("UsePfElectronSelection")),
44  useMuonSelection(iConfig.getParameter<bool>("UseMuonSelection")),
45  useHtSelection(iConfig.getParameter<bool>("UseHtSelection")),
46  usePFHtSelection(iConfig.getParameter<bool>("UsePFHtSelection")),
47  ptElecMin(iConfig.getParameter<double>("electronPtMin")),
48  ptPfElecMin(iConfig.getParameter<double>("pfElectronPtMin")),
49  nSelectedElectrons(iConfig.getParameter<int>("electronN")),
50  nSelectedPfElectrons(iConfig.getParameter<int>("pfElectronN")),
51  ptGlobalMuonMin(iConfig.getParameter<double>("globalMuonPtMin")),
52  ptTrackerMuonMin(iConfig.getParameter<double>("trackerMuonPtMin")),
53  nSelectedMuons(iConfig.getParameter<int>("muonN")),
54  htMin(iConfig.getParameter<double>("HtMin")),
55  pfHtMin(iConfig.getParameter<double>("PFHtMin")),
56  htJetThreshold(iConfig.getParameter<double>("HtJetThreshold")),
57  pfHtJetThreshold(iConfig.getParameter<double>("PFHtJetThreshold")) {
58  NeventsTotal = 0;
59  NeventsFiltered = 0;
60  NHltMu9 = 0;
61  NHltDiMu3 = 0;
62 
63  NtotalElectrons = 0;
64  NmvaElectrons = 0;
65 }
66 
68  // do anything here that needs to be done at desctruction time
69  // (e.g. close files, deallocate resources etc.)
70 }
71 
72 //
73 // member functions
74 //
75 
76 // ------------ method called on each new Event ------------
78  bool accept = false;
79 
80  NeventsTotal++;
81 
82  ElectronCutPassed = false;
83  PfElectronCutPassed = false;
84  MuonCutPassed = false;
85  HtCutPassed = false;
86  PFHtCutPassed = false;
87 
88  // edm::Handle<TriggerResults> trhv;
89  // iEvent.getByLabel(hltLabel,trhv);
90 
91  // const edm::TriggerNames& triggerNames_ = iEvent.triggerNames(*trhv);
92 
93  // if(trhv->at(triggerNames_.triggerIndex("HLT_Mu9")).accept()) NHltMu9++;
94  // if(trhv->at(triggerNames_.triggerIndex("HLT_DoubleMu3")).accept()) NHltDiMu3++;
95 
96  this->handleObjects(iEvent, iSetup);
97 
99  int nElecPassingCut = 0;
100  for (unsigned int i = 0; i < theElectronCollection->size(); i++) {
101  GsfElectron electron = (*theElectronCollection)[i];
102  // if (electron.ecalDrivenSeed()) {
103  float elpt = electron.pt();
104  // if (electron.sigmaIetaIeta() < 0.002) continue;
105  if (elpt > ptElecMin) {
106  NtotalElectrons++;
107  nElecPassingCut++;
108  if (electron.mva_e_pi() > -0.1)
109  NmvaElectrons++;
110  }
111  LogDebug("LeptonRecoSkim") << "elpt = " << elpt << endl;
112  // } // closes if (electron.ecalDrivenSeed()) {
113  }
114  if (nElecPassingCut >= nSelectedElectrons)
115  ElectronCutPassed = true;
116  } else
117  ElectronCutPassed = true;
118 
120  int nPfElecPassingCut = 0;
121  for (unsigned int i = 0; i < thePfCandidateCollection->size(); i++) {
122  const reco::PFCandidate& thePfCandidate = (*thePfCandidateCollection)[i];
123  if (thePfCandidate.particleId() != reco::PFCandidate::e)
124  continue;
125  if (thePfCandidate.gsfTrackRef().isNull())
126  continue;
127  float pfelpt = thePfCandidate.pt();
128  // if (electron.sigmaIetaIeta() < 0.002) continue;
129  if (pfelpt > ptPfElecMin)
130  nPfElecPassingCut++;
131  LogDebug("LeptonRecoSkim") << "pfelpt = " << pfelpt << endl;
132  }
133  if (nPfElecPassingCut >= nSelectedPfElectrons)
134  PfElectronCutPassed = true;
135  } else
136  PfElectronCutPassed = true;
137 
138  if (useMuonSelection) {
139  int nMuonPassingCut = 0;
140  for (unsigned int i = 0; i < theMuonCollection->size(); i++) {
141  Muon muon = (*theMuonCollection)[i];
142  if (!(muon.isGlobalMuon() || muon.isTrackerMuon()))
143  continue;
144  const TrackRef siTrack = muon.innerTrack();
145  const TrackRef globalTrack = muon.globalTrack();
146  float muonpt;
147  float ptMuonMin;
148  // if (siTrack.isNonnull() && globalTrack.isNonnull()) {
149  if (muon.isGlobalMuon() && muon.isTrackerMuon()) {
150  muonpt = max(siTrack->pt(), globalTrack->pt());
151  ptMuonMin = ptGlobalMuonMin;
152  } else if (muon.isGlobalMuon() && !(muon.isTrackerMuon())) {
153  muonpt = globalTrack->pt();
154  ptMuonMin = ptGlobalMuonMin;
155  } else if (muon.isTrackerMuon() && !(muon.isGlobalMuon())) {
156  muonpt = siTrack->pt();
157  ptMuonMin = ptTrackerMuonMin;
158  } else {
159  muonpt = 0; // if we get here this is a STA only muon
160  ptMuonMin = 999999.;
161  }
162  if (muonpt > ptMuonMin)
163  nMuonPassingCut++;
164  LogDebug("RecoSelectorCuts") << "muonpt = " << muonpt << endl;
165  }
166  if (nMuonPassingCut >= nSelectedMuons)
167  MuonCutPassed = true;
168  } else
169  MuonCutPassed = true;
170 
171  if (useHtSelection) {
172  double Ht = 0;
173  for (unsigned int i = 0; i < theCaloJetCollection->size(); i++) {
174  // if((*theCaloJetCollection)[i].eta()<2.6 && (*theCaloJetCollection)[i].emEnergyFraction() <= 0.01) continue;
175  if ((*theCaloJetCollection)[i].pt() > htJetThreshold)
176  Ht += (*theCaloJetCollection)[i].pt();
177  }
178  if (Ht > htMin)
179  HtCutPassed = true;
180  } else
181  HtCutPassed = true;
182 
183  if (usePFHtSelection) {
184  double PFHt = 0;
185  for (unsigned int i = 0; i < thePFJetCollection->size(); i++) {
186  if ((*thePFJetCollection)[i].pt() > pfHtJetThreshold)
187  PFHt += (*thePFJetCollection)[i].pt();
188  }
189  if (PFHt > pfHtMin)
190  PFHtCutPassed = true;
191  } else
192  PFHtCutPassed = true;
193 
195  accept = true;
196 
197  if (accept)
198  NeventsFiltered++;
199 
200  firstEvent = false;
201  return accept;
202 }
203 
204 // ------------ method called once each job just before starting event loop ------------
206 
207 // ------------ method called once each job just after ending the event loop ------------
209  edm::LogPrint("LeptonRecoSkim") << "Filter Name = " << filterName << endl;
210  edm::LogPrint("LeptonRecoSkim") << "Total number of events = " << NeventsTotal << endl;
211  edm::LogPrint("LeptonRecoSkim") << "Total HLT_Mu9 = " << NHltMu9 << endl;
212  edm::LogPrint("LeptonRecoSkim") << "Total HLT_DoubleMu3 = " << NHltDiMu3 << endl;
213  edm::LogPrint("LeptonRecoSkim") << "Filtered events = " << NeventsFiltered << endl;
214  edm::LogPrint("LeptonRecoSkim") << "Filter Efficiency = " << (float)NeventsFiltered / (float)NeventsTotal
215  << endl;
216  edm::LogPrint("LeptonRecoSkim") << endl;
217  edm::LogPrint("LeptonRecoSkim") << "N total electrons = " << NtotalElectrons << endl;
218  edm::LogPrint("LeptonRecoSkim") << "N mva>-0.1 electrons = " << NmvaElectrons << endl;
219  edm::LogPrint("LeptonRecoSkim") << endl;
220  edm::LogPrint("LeptonRecoSkim") << endl;
221 }
222 
224  //Get the electrons
226 
227  //Get the pf electrons
229 
230  //Get the Muons
232 
233  //Get the CaloJets
235 
236  //Get the PfJets
238 
239  // Get the ECAL rechhits to clean the spikes
240  // Get EB RecHits
242  // Get EE RecHits
244 
245  // Get topology for spike cleaning
248 }
249 
250 //define this as a plug-in
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:462
const reco::CaloJetCollection * theCaloJetCollection
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::string filterName
edm::EDGetTokenT< reco::GsfElectronCollection > gsfElectronCollectionToken_
double pt() const final
transverse momentum
edm::EDGetTokenT< reco::CaloJetCollection > caloJetCollectionToken_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > m_CaloGeoToken
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollectionToken_
void handleObjects(const edm::Event &, const edm::EventSetup &iSetup)
void endJob() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool filter(edm::Event &, const edm::EventSetup &) override
bool usePfElectronSelection
const reco::PFCandidateCollection * thePfCandidateCollection
edm::EDGetTokenT< EcalRecHitCollection > eeRecHitCollectionToken_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
double ptTrackerMuonMin
int iEvent
Definition: GenABIO.cc:224
Definition: Muon.py:1
const reco::GsfElectronCollection * theElectronCollection
void beginJob() override
bool useElectronSelection
edm::EDGetTokenT< reco::MuonCollection > muonCollectionToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
bool isNull() const
Checks for null.
Definition: Ref.h:235
Log< level::Warning, true > LogPrint
edm::EDGetTokenT< EcalRecHitCollection > ebRecHitCollectionToken_
const reco::MuonCollection * theMuonCollection
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
const EcalRecHitCollection * theEcalBarrelCollection
HLT enums.
const edm::ESGetToken< CaloTopology, CaloTopologyRecord > m_CaloTopoToken
~LeptonRecoSkim() override
LeptonRecoSkim(const edm::ParameterSet &)
const reco::PFJetCollection * thePFJetCollection
double ptGlobalMuonMin
edm::EDGetTokenT< reco::PFJetCollection > pfJetCollectionToken_
const CaloGeometry * theCaloGeometry
double pfHtJetThreshold
const CaloTopology * theCaloTopology
const EcalRecHitCollection * theEcalEndcapCollection
#define LogDebug(id)
virtual ParticleType particleId() const
Definition: PFCandidate.h:392