CMS 3D CMS Logo

B2GHadronicHLTValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTriggerOffline/B2G
4 // Class: B2GHadronicHLTValidation
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 
39 #include "TString.h"
40 //
41 // member functions
42 //
43 
44 // ------------ method called for each event ------------
46  using namespace edm;
47 
48  isAll_ = false;
49  isSel_ = false;
50 
51  // Jets
53  if (!iEvent.getByToken(tokJets_, jets))
54  edm::LogWarning("B2GHadronicHLTValidation") << "Jets collection not found \n";
55  unsigned int nGoodJ = 0;
56  double ht = 0.0;
57  // Check to see if we want asymmetric jet pt cuts
58  if (ptJets0_ > 0.0 || ptJets1_ > 0.0) {
59  if (ptJets0_ > 0.0) {
60  if (!jets->empty() && jets->at(0).pt() > ptJets0_) {
61  nGoodJ++;
62  jet_ = jets->ptrAt(0);
63  }
64  }
65  if (ptJets1_ > 0.0) {
66  if (jets->size() > 1 && jets->at(1).pt() > ptJets1_) {
67  nGoodJ++;
68  jet_ = jets->ptrAt(1);
69  }
70  }
71  } else if (minJets_ > 0 || htMin_ > 0) {
72  for (edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j) {
73  if (j->pt() < ptJets_)
74  continue;
75  if (fabs(j->eta()) > etaJets_)
76  continue;
77  nGoodJ++;
78  ht += j->pt();
79  if (nGoodJ == minJets_)
80  jet_ = jets->ptrAt(j - jets->begin());
81  }
82  }
83 
84  if (nGoodJ >= minJets_ || ht > htMin_)
85  isAll_ = true;
86 
87  // Trigger
88  Handle<edm::TriggerResults> triggerTable;
89  if (!iEvent.getByToken(tokTrigger_, triggerTable))
90  edm::LogWarning("B2GHadronicHLTValidation") << "Trigger collection not found \n";
91  const edm::TriggerNames &triggerNames = iEvent.triggerNames(*triggerTable);
92  bool isInteresting = false;
93  for (unsigned int i = 0; i < triggerNames.triggerNames().size(); ++i) {
94  TString name = triggerNames.triggerNames()[i].c_str();
95  for (unsigned int j = 0; j < vsPaths_.size(); j++) {
96  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
97  isInteresting = true;
98  break;
99  }
100  }
101  if (isInteresting)
102  break;
103  }
104 
105  if (isAll_ && isInteresting)
106  isSel_ = true;
107  else
108  isSel_ = false;
109 
110  // Histos
111  if (isAll_) {
112  if (jet_.isNonnull()) {
113  hDenJetPt->Fill(jet_->pt());
114  hDenJetEta->Fill(jet_->eta());
115  }
116  for (unsigned int idx = 0; idx < vsPaths_.size(); ++idx) {
117  hDenTriggerMon->Fill(idx + 0.5);
118  }
119  }
120  if (isSel_) {
121  hNumJetPt->Fill(jet_->pt());
122  hNumJetEta->Fill(jet_->eta());
123  for (unsigned int i = 0; i < triggerNames.triggerNames().size(); ++i) {
124  TString name = triggerNames.triggerNames()[i].c_str();
125  for (unsigned int j = 0; j < vsPaths_.size(); j++) {
126  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
127  hNumTriggerMon->Fill(j + 0.5);
128  }
129  }
130  }
131  }
132 }
133 
134 // ------------ booking histograms -----------
136  dbe.setCurrentFolder(sDir_);
137  hDenJetPt = dbe.book1D("PtLastJetAll", "PtLastJetAll", 60, 0., 3000.);
138  hDenJetEta = dbe.book1D("EtaLastJetAll", "EtaLastJetAll", 30, -3., 3.);
139  hNumJetPt = dbe.book1D("PtLastJetSel", "PtLastJetSel", 60, 0., 3000.);
140  hNumJetEta = dbe.book1D("EtaLastJetSel", "EtaLastJetSel", 30, -3., 3.);
141  // determine number of bins for trigger monitoring
142  unsigned int nPaths = vsPaths_.size();
143  // monitored trigger occupancy for single lepton triggers
144  hNumTriggerMon = dbe.book1D("TriggerMonSel", "TriggerMonSel", nPaths, 0., nPaths);
145  hDenTriggerMon = dbe.book1D("TriggerMonAll", "TriggerMonAll", nPaths, 0., nPaths);
146  // set bin labels for trigger monitoring
148 }
149 
150 // ------------ method fills 'descriptions' with the allowed parameters for the
151 // module ------------
153  // The following says we do not know what parameters are allowed so do no
154  // validation
155  // Please change this to state exactly what you do use, even if it is no
156  // parameters
158  desc.setUnknown();
159  descriptions.addDefault(desc);
160 }
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
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
void triggerBinLabels(const std::vector< std::string > &labels)
set configurable labels for trigger monitoring histograms
vector< PseudoJet > jets
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
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
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
std::vector< std::string > vsPaths_
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:256
Definition: Run.h:45