CMS 3D CMS Logo

B2GSingleLeptonHLTValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTriggerOffline/B2G
4 // Class: B2GSingleLeptonHLTValidation
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("B2GSingleLeptonHLTValidation") << "Electrons collection not found \n";
55  unsigned int nGoodE = 0;
56  for(edm::View<reco::GsfElectron>::const_iterator e = electrons->begin(); e != electrons->end(); ++e){
57  if (e->pt() < ptElectrons_) continue;
58  if (fabs(e->eta()) > etaElectrons_) continue;
59  nGoodE++;
60  if (nGoodE == 1) elec_ = electrons->ptrAt( e - electrons->begin());
61  }
62  // Muons
64  if (!iEvent.getByToken(tokMuons_,muons))
65  edm::LogWarning("B2GSingleLeptonHLTValidation") << "Muons collection not found \n";
66  unsigned int nGoodM = 0;
67  for(edm::View<reco::Muon>::const_iterator m = muons->begin(); m != muons->end(); ++m){
68  if (!m->isPFMuon() || !m->isGlobalMuon()) continue;
69  if (m->pt() < ptMuons_) continue;
70  if (fabs(m->eta()) > etaMuons_) continue;
71  nGoodM++;
72  if (nGoodM == 1) mu_ = muons->ptrAt( m - muons->begin() );
73  }
74  // Jets
76  if (!iEvent.getByToken(tokJets_,jets))
77  edm::LogWarning("B2GSingleLeptonHLTValidation") << "Jets collection not found \n";
78  unsigned int nGoodJ = 0;
79 
80  // Check to see if we want asymmetric jet pt cuts
81  if ( ptJets0_ > 0.0 || ptJets1_ > 0.0 ) {
82  if ( ptJets0_ > 0.0 ) {
83  if ( !jets->empty() && jets->at(0).pt() > ptJets0_ ) {
84  nGoodJ++;
85  jet_ = jets->ptrAt(0) ;
86  }
87  }
88  if ( ptJets1_ > 0.0 ) {
89  if ( jets->size() > 1 && jets->at(1).pt() > ptJets1_ ) {
90  nGoodJ++;
91  jet_ = jets->ptrAt(1) ;
92  }
93  }
94  } else {
95  for(edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j){
96  if (j->pt() < ptJets_) continue;
97  if (fabs(j->eta()) > etaJets_) continue;
98  nGoodJ++;
99  if (nGoodJ == minJets_) jet_ = jets->ptrAt( j -jets->begin() );
100  }
101  }
102 
103  if (nGoodE >= minElectrons_ && nGoodM >= minMuons_ && nGoodJ >= minJets_) isAll_ = true;
104 
105 
106  //Trigger
107  Handle<edm::TriggerResults> triggerTable;
108  if (!iEvent.getByToken(tokTrigger_,triggerTable))
109  edm::LogWarning("B2GSingleLeptonHLTValidation") << "Trigger collection not found \n";
110  const edm::TriggerNames& triggerNames = iEvent.triggerNames(*triggerTable);
111  bool isInteresting = false;
112  for (unsigned int i=0; i<triggerNames.triggerNames().size(); ++i) {
113 
114  TString name = triggerNames.triggerNames()[i].c_str();
115  for (unsigned int j=0; j<vsPaths_.size(); j++) {
116  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
117  isInteresting = true;
118  break;
119  }
120  }
121  if (isInteresting) break;
122  }
123 
124  if (isAll_ && isInteresting) isSel_ = true;
125  else isSel_ = false;
126 
127  //Histos
128  if (isAll_) {
129  if (minElectrons_ > 0 && elec_.isNonnull() ) {
130  hDenLeptonPt->Fill(elec_->pt());
132  }
133  if (minMuons_ > 0 && mu_.isNonnull() ) {
134  hDenLeptonPt->Fill(mu_->pt());
135  hDenLeptonEta->Fill(mu_->eta());
136  }
137  if ( jet_.isNonnull() ) {
138  hDenJetPt->Fill(jet_->pt());
139  hDenJetEta->Fill(jet_->eta());
140  }
141  for(unsigned int idx=0; idx<vsPaths_.size(); ++idx){
142  hDenTriggerMon->Fill(idx+0.5);
143  }
144 
145  }
146  if (isSel_) {
147  if (minElectrons_ > 0 && elec_.isNonnull() ) {
148  hNumLeptonPt->Fill(elec_->pt());
150  }
151  if (minMuons_ > 0 && mu_.isNonnull() ) {
152  hNumLeptonPt->Fill(mu_->pt());
153  hNumLeptonEta->Fill(mu_->eta());
154  }
155  if ( jet_.isNonnull() ) {
156  hNumJetPt->Fill(jet_->pt());
157  hNumJetEta->Fill(jet_->eta());
158  }
159  for (unsigned int i=0; i<triggerNames.triggerNames().size(); ++i) {
160  TString name = triggerNames.triggerNames()[i].c_str();
161  for (unsigned int j=0; j<vsPaths_.size(); j++) {
162  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
163  hNumTriggerMon->Fill(j+0.5 );
164  }
165  }
166  }
167  }
168 }
169 
170 
171 // ------------ booking histograms -----------
172  void
174 {
175  dbe.setCurrentFolder(sDir_);
176  hDenLeptonPt = dbe.book1D("PtLeptonAll", "PtLeptonAll", 50, 0., 2500.);
177  hDenLeptonEta = dbe.book1D("EtaLeptonAll", "EtaLeptonAll", 30, -3. , 3.);
178  hDenJetPt = dbe.book1D("PtLastJetAll", "PtLastJetAll", 60, 0., 3000.);
179  hDenJetEta = dbe.book1D("EtaLastJetAll", "EtaLastJetAll", 30, -3., 3.);
180  hNumLeptonPt = dbe.book1D("PtLeptonSel", "PtLeptonSel", 50, 0., 2500.);
181  hNumLeptonEta = dbe.book1D("EtaLeptonSel", "EtaLeptonSel", 30, -3. , 3.);
182  hNumJetPt = dbe.book1D("PtLastJetSel", "PtLastJetSel", 60, 0., 3000.);
183  hNumJetEta = dbe.book1D("EtaLastJetSel", "EtaLastJetSel", 30, -3., 3.);
184  // determine number of bins for trigger monitoring
185  unsigned int nPaths = vsPaths_.size();
186  // monitored trigger occupancy for single lepton triggers
187  hNumTriggerMon = dbe.book1D("TriggerMonSel" , "TriggerMonSel", nPaths, 0., nPaths);
188  hDenTriggerMon = dbe.book1D("TriggerMonAll" , "TriggerMonAll", nPaths, 0., nPaths);
189  // set bin labels for trigger monitoring
191 }
192 
193 
194 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
195 void
197  //The following says we do not know what parameters are allowed so do no validation
198  // Please change this to state exactly what you do use, even if it is no parameters
200  desc.setUnknown();
201  descriptions.addDefault(desc);
202 }
203 
edm::EDGetTokenT< edm::View< reco::GsfElectron > > tokElectrons_
double eta() const final
momentum pseudorapidity
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
double pt() const final
transverse momentum
edm::Ptr< reco::GsfElectron > elec_
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::TriggerResults > tokTrigger_
void addDefault(ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
vector< PseudoJet > jets
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void triggerBinLabels(const std::vector< std::string > &labels)
set configurable labels for trigger monitoring histograms
edm::EDGetTokenT< edm::View< reco::Muon > > tokMuons_
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
Definition: Run.h:44