CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopSingleLeptonHLTValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTriggerOffline/Top
4 // Class: TopSingleLeptonHLTValidation
5 //
15 //
16 // Original Author: Elvire Bouvier
17 // Created: Thu, 16 Jan 2014 16:27:35 GMT
18 //
19 //
21 
22 // system include files
23 #include <memory>
24 
25 // user include files
27 
30 
33 
37 #include "TString.h"
39 //
40 // member functions
41 //
42 
43 // ------------ method called for each event ------------
44  void
46 {
47  using namespace edm;
48 
49  isAll_ = false; isSel_ = false;
50 
51  // Electrons
53  if (!iEvent.getByToken(tokElectrons_,electrons))
54  edm::LogWarning("TopSingleLeptonHLTValidation") << "Electrons collection not found \n";
55  unsigned int nGoodE = 0;
57  if (e->pt() < ptElectrons_) continue;
58  if (fabs(e->eta()) > etaElectrons_) continue;
59  if ((e->dr03TkSumPt()+e->dr03EcalRecHitSumEt()+e->dr03HcalTowerSumEt())/e->pt() > isoElectrons_ ) continue;
60  nGoodE++;
61  if (nGoodE == 1) elec_ = &(*e);
62  }
63  // Muons
65  if (!iEvent.getByToken(tokMuons_,muons))
66  edm::LogWarning("TopSingleLeptonHLTValidation") << "Muons collection not found \n";
67  unsigned int nGoodM = 0;
68  for(edm::View<reco::Muon>::const_iterator m = muons->begin(); m != muons->end(); ++m){
69  if (!m->isPFMuon() || !m->isGlobalMuon()) continue;
70  if (m->pt() < ptMuons_) continue;
71  if (fabs(m->eta()) > etaMuons_) continue;
72  if (((m->pfIsolationR04()).sumChargedHadronPt+(m->pfIsolationR04()).sumPhotonEt+(m->pfIsolationR04()).sumNeutralHadronEt)/m->pt() > isoMuons_ ) continue;
73  nGoodM++;
74  if (nGoodM == 1) mu_ = &(*m);
75  }
76  // Jets
78  if (!iEvent.getByToken(tokJets_,jets))
79  edm::LogWarning("TopSingleLeptonHLTValidation") << "Jets collection not found \n";
80  unsigned int nGoodJ = 0;
81  if (minJets_ == 4) {
82  for(edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j){
83  if (nGoodJ < 1 && j->pt() < 55) continue;
84  if (nGoodJ < 2 && j->pt() < 45) continue;
85  if (nGoodJ < 3 && j->pt() < 35) continue;
86  if (nGoodJ >= 3 && j->pt() < 20) continue;
87  if (fabs(j->eta()) > etaJets_) continue;
88  nGoodJ++;
89  if (nGoodJ == minJets_) jet_ = &(*j);
90  }
91  }
92  else {
93  for(edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j){
94  if (j->pt() < ptJets_) continue;
95  if (fabs(j->eta()) > etaJets_) continue;
96  nGoodJ++;
97  if (nGoodJ == minJets_) jet_ = &(*j);
98  }
99  }
100 
101  if (nGoodE >= minElectrons_ && nGoodM >= minMuons_ && nGoodJ >= minJets_) isAll_ = true;
102 
103 
104  //Trigger
105  Handle<edm::TriggerResults> triggerTable;
106  if (!iEvent.getByToken(tokTrigger_,triggerTable))
107  edm::LogWarning("TopSingleLeptonHLTValidation") << "Trigger collection not found \n";
108  const edm::TriggerNames& triggerNames = iEvent.triggerNames(*triggerTable);
109  for (unsigned int i=0; i<triggerNames.triggerNames().size(); ++i) {
110 
111  TString name = triggerNames.triggerNames()[i].c_str();
112  bool isInteresting = false;
113  for (unsigned int j=0; j<vsPaths_.size(); j++) {
114  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) isInteresting = true;
115  }
116 
117  if (isAll_ && isInteresting) isSel_ = true;
118  else isSel_ = false;
119 
120  //Histos
121  if (isAll_) {
122  if (minElectrons_ > 0) {
123  hDenLeptonPt->Fill(elec_->pt());
125  }
126  if (minMuons_ > 0) {
127  hDenLeptonPt->Fill(mu_->pt());
128  hDenLeptonEta->Fill(mu_->eta());
129  }
130  hDenJetPt->Fill(jet_->pt());
131  hDenJetEta->Fill(jet_->eta());
132  }
133  if (isSel_) {
134  if (minElectrons_ > 0) {
135  hNumLeptonPt->Fill(elec_->pt());
137  }
138  if (minMuons_ > 0) {
139  hNumLeptonPt->Fill(mu_->pt());
140  hNumLeptonEta->Fill(mu_->eta());
141  }
142  hNumJetPt->Fill(jet_->pt());
143  hNumJetEta->Fill(jet_->eta());
144  }
145  }
146 }
147 
148 
149 // ------------ booking histograms -----------
150  void
152 {
153  dbe.setCurrentFolder(sDir_);
154  hDenLeptonPt = dbe.book1D("PtLeptonAll", "PtLeptonAll", 50, 0., 250.);
155  hDenLeptonEta = dbe.book1D("EtaLeptonAll", "EtaLeptonAll", 30, -3. , 3.);
156  hDenJetPt = dbe.book1D("PtLastJetAll", "PtLastJetAll", 60, 0., 300.);
157  hDenJetEta = dbe.book1D("EtaLastJetAll", "EtaLastJetAll", 30, -3., 3.);
158  hNumLeptonPt = dbe.book1D("PtLeptonSel", "PtLeptonSel", 50, 0., 250.);
159  hNumLeptonEta = dbe.book1D("EtaLeptonSel", "EtaLeptonSel", 30, -3. , 3.);
160  hNumJetPt = dbe.book1D("PtLastJetSel", "PtLastJetSel", 60, 0., 300.);
161  hNumJetEta = dbe.book1D("EtaLastJetSel", "EtaLastJetSel", 30, -3., 3.);
162 }
163 
164 
165 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
166 void
168  //The following says we do not know what parameters are allowed so do no validation
169  // Please change this to state exactly what you do use, even if it is no parameters
171  desc.setUnknown();
172  descriptions.addDefault(desc);
173 }
174 
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:204
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
void Fill(long long x)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< edm::TriggerResults > tokTrigger_
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
vector< PseudoJet > jets
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
int j
Definition: DBlmapReader.cc:9
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
edm::EDGetTokenT< edm::View< reco::Muon > > tokMuons_
tuple muons
Definition: patZpeak.py:38
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual float pt() const GCC11_FINAL
transverse momentum
edm::EDGetTokenT< edm::View< reco::GsfElectron > > tokElectrons_
Definition: Run.h:41