CMS 3D CMS Logo

TopSingleLeptonHLTValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTriggerOffline/Top
4 // Class: TopSingleLeptonHLTValidation
5 //
17 //
18 // Original Author: Elvire Bouvier
19 // Created: Thu, 16 Jan 2014 16:27:35 GMT
20 //
21 //
23 
24 // system include files
25 #include <memory>
26 
27 // user include files
29 
32 
35 
40 #include "TString.h"
41 //
42 // member functions
43 //
44 
45 // ------------ method called for each event ------------
47  using namespace edm;
48 
49  isAll_ = false;
50  isSel_ = false;
51 
52  // Electrons
54  if (!iEvent.getByToken(tokElectrons_, electrons))
55  edm::LogWarning("TopSingleLeptonHLTValidation") << "Electrons collection not found \n";
56  unsigned int nGoodE = 0;
57  for (edm::View<reco::GsfElectron>::const_iterator e = electrons->begin(); e != electrons->end(); ++e) {
58  if (e->pt() < ptElectrons_)
59  continue;
60  if (fabs(e->eta()) > etaElectrons_)
61  continue;
62  if ((e->dr03TkSumPt() + e->dr03EcalRecHitSumEt() + e->dr03HcalTowerSumEt()) / e->pt() > isoElectrons_)
63  continue;
64  nGoodE++;
65  if (nGoodE == 1)
66  elec_ = &(*e);
67  }
68  // Muons
70  if (!iEvent.getByToken(tokMuons_, muons))
71  edm::LogWarning("TopSingleLeptonHLTValidation") << "Muons collection not found \n";
72  unsigned int nGoodM = 0;
73  for (edm::View<reco::Muon>::const_iterator m = muons->begin(); m != muons->end(); ++m) {
74  if (!m->isPFMuon() || !m->isGlobalMuon())
75  continue;
76  if (m->pt() < ptMuons_)
77  continue;
78  if (fabs(m->eta()) > etaMuons_)
79  continue;
80  if (((m->pfIsolationR04()).sumChargedHadronPt + (m->pfIsolationR04()).sumPhotonEt +
81  (m->pfIsolationR04()).sumNeutralHadronEt) /
82  m->pt() >
83  isoMuons_)
84  continue;
85  nGoodM++;
86  if (nGoodM == 1)
87  mu_ = &(*m);
88  }
89  // Jets
91  if (!iEvent.getByToken(tokJets_, jets))
92  edm::LogWarning("TopSingleLeptonHLTValidation") << "Jets collection not found \n";
93  unsigned int nGoodJ = 0;
94  if (minJets_ == 4) {
95  for (edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j) {
96  if (nGoodJ < 1 && j->pt() < 55)
97  continue;
98  if (nGoodJ < 2 && j->pt() < 45)
99  continue;
100  if (nGoodJ < 3 && j->pt() < 35)
101  continue;
102  if (nGoodJ >= 3 && j->pt() < 20)
103  continue;
104  if (fabs(j->eta()) > etaJets_)
105  continue;
106  nGoodJ++;
107  if (nGoodJ == minJets_)
108  jet_ = &(*j);
109  }
110  } else {
111  for (edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j) {
112  if (j->pt() < ptJets_)
113  continue;
114  if (fabs(j->eta()) > etaJets_)
115  continue;
116  nGoodJ++;
117  if (nGoodJ == minJets_)
118  jet_ = &(*j);
119  }
120  }
121 
122  if (nGoodE >= minElectrons_ && nGoodM >= minMuons_ && nGoodJ >= minJets_)
123  isAll_ = true;
124 
125  // Trigger
126  Handle<edm::TriggerResults> triggerTable;
127  if (!iEvent.getByToken(tokTrigger_, triggerTable))
128  edm::LogWarning("TopSingleLeptonHLTValidation") << "Trigger collection not found \n";
129  const edm::TriggerNames &triggerNames = iEvent.triggerNames(*triggerTable);
130  unsigned int isInteresting = 0;
131  for (unsigned int i = 0; i < triggerNames.triggerNames().size(); ++i) {
132  TString name = triggerNames.triggerNames()[i].c_str();
133  for (unsigned int j = 0; j < vsPaths_.size(); j++) {
134  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
135  if (triggerTable->accept(i)) {
136  isInteresting++;
137  if (isAll_)
138  hNumTriggerMon->Fill(j + 0.5);
139  }
140  }
141  }
142  }
143 
144  if (isAll_ && isInteresting > 0)
145  isSel_ = true;
146  else
147  isSel_ = false;
148 
149  // Histos
150  if (isAll_) {
151  if (minElectrons_ > 0) {
152  hDenLeptonPt->Fill(elec_->pt());
154  }
155  if (minMuons_ > 0) {
156  hDenLeptonPt->Fill(mu_->pt());
157  hDenLeptonEta->Fill(mu_->eta());
158  }
159  hDenJetPt->Fill(jet_->pt());
160  hDenJetEta->Fill(jet_->eta());
161  for (unsigned int idx = 0; idx < vsPaths_.size(); ++idx) {
162  hDenTriggerMon->Fill(idx + 0.5);
163  }
164  }
165  if (isSel_) {
166  if (minElectrons_ > 0) {
167  hNumLeptonPt->Fill(elec_->pt());
169  }
170  if (minMuons_ > 0) {
171  hNumLeptonPt->Fill(mu_->pt());
172  hNumLeptonEta->Fill(mu_->eta());
173  }
174  hNumJetPt->Fill(jet_->pt());
175  hNumJetEta->Fill(jet_->eta());
176  }
177 }
178 
179 // ------------ booking histograms -----------
181  dbe.setCurrentFolder(sDir_);
182  hDenLeptonPt = dbe.book1D("PtLeptonAll", "PtLeptonAll", 50, 0., 250.);
183  hDenLeptonEta = dbe.book1D("EtaLeptonAll", "EtaLeptonAll", 30, -3., 3.);
184  hDenJetPt = dbe.book1D("PtLastJetAll", "PtLastJetAll", 60, 0., 300.);
185  hDenJetEta = dbe.book1D("EtaLastJetAll", "EtaLastJetAll", 30, -3., 3.);
186  hNumLeptonPt = dbe.book1D("PtLeptonSel", "PtLeptonSel", 50, 0., 250.);
187  hNumLeptonEta = dbe.book1D("EtaLeptonSel", "EtaLeptonSel", 30, -3., 3.);
188  hNumJetPt = dbe.book1D("PtLastJetSel", "PtLastJetSel", 60, 0., 300.);
189  hNumJetEta = dbe.book1D("EtaLastJetSel", "EtaLastJetSel", 30, -3., 3.);
190  // determine number of bins for trigger monitoring
191  unsigned int nPaths = vsPaths_.size();
192  // monitored trigger occupancy for single lepton triggers
193  hNumTriggerMon = dbe.book1D("TriggerMonSel", "TriggerMonSel", nPaths, 0., nPaths);
194  hDenTriggerMon = dbe.book1D("TriggerMonAll", "TriggerMonAll", nPaths, 0., nPaths);
195  // set bin labels for trigger monitoring
197 }
198 
199 // ------------ method fills 'descriptions' with the allowed parameters for the
200 // module ------------
202  // The following says we do not know what parameters are allowed so do no
203  // validation
204  // Please change this to state exactly what you do use, even if it is no
205  // parameters
207  desc.setUnknown();
208  descriptions.addDefault(desc);
209 }
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool accept() const
Has at least one path accepted the event?
double pt() const final
transverse momentum
Strings const & triggerNames() const
Definition: TriggerNames.cc:20
void Fill(long long x)
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< edm::TriggerResults > tokTrigger_
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
vector< PseudoJet > jets
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
void triggerBinLabels(const std::vector< std::string > &labels)
set configurable labels for trigger monitoring histograms
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
edm::EDGetTokenT< edm::View< reco::GsfElectron > > tokElectrons_
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< edm::View< reco::Muon > > tokMuons_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:256
Definition: Run.h:45