CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTScalers.cc
Go to the documentation of this file.
1 //
2 // Revision 1.30 2011/03/30 21:44:03 fwyzard
3 // make sure HLTConfigProvider is used only if succesfully initialized
4 //
5 // Revision 1.29 2011/03/30 21:35:40 fwyzard
6 // make sure all members are initialized
7 //
8 // Revision 1.28 2011/03/29 09:46:03 rekovic
9 // clean vector pairPDPaths in beginRun and tidy up
10 //
11 // Revision 1.27 2011/03/24 18:35:38 rekovic
12 // Change name for pd histo
13 //
14 // Revision 1.26 2011/03/24 18:25:45 rekovic
15 // Add single 1D plot of streamA content
16 //
17 // Revision 1.25 2010/07/20 02:58:27 wmtan
18 // Add missing #include files
19 //
20 // Revision 1.24 2010/03/17 20:54:51 wittich
21 // add scalers that I manually reset on beginLumi
22 //
23 // Revision 1.23 2010/02/25 17:34:01 wdd
24 // Central migration of TriggerNames class interface
25 //
26 // Revision 1.22 2010/02/24 17:43:47 wittich
27 // - keep trying to get path names if it doesn't work first time
28 // - move the Bx histograms out of raw to the toplevel directory.
29 //
30 // Revision 1.21 2010/02/11 23:54:28 wittich
31 // modify how the monitoring histo is filled
32 //
33 // Revision 1.20 2010/02/11 00:11:08 wmtan
34 // Adapt to moved framework header
35 //
36 // Revision 1.19 2010/02/02 13:53:05 wittich
37 // fix duplicate histogram name
38 //
39 // Revision 1.18 2010/02/02 11:42:53 wittich
40 // new diagnostic histograms
41 //
42 // Revision 1.17 2009/11/20 00:39:12 lorenzo
43 // fixes
44 //
45 // Revision 1.16 2008/09/03 13:59:06 wittich
46 // make HLT DQM path configurable via python parameter,
47 // which defaults to HLT/HLTScalers_EvF
48 //
49 // Revision 1.15 2008/09/03 02:13:47 wittich
50 // - bug fix in L1Scalers
51 // - configurable dqm directory in L1SCalers
52 // - other minor tweaks in HLTScalers
53 //
54 
55 #include <iostream>
56 
57 // FW
60 
63 
66 
67 // HLT
71 
76 
77 using namespace edm;
78 
80  : folderName_(ps.getUntrackedParameter<std::string>("dqmFolder",
81  "HLT/HLTScalers_EvF")),
82  processname_(ps.getParameter<std::string>("processname")),
83  pairPDPaths_(),
84  trigResultsSource_(consumes<TriggerResults>(
85  ps.getParameter<edm::InputTag>("triggerResults"))),
86  scalersN_(0),
87  scalersException_(0),
88  hltCorrelations_(0),
89  detailedScalers_(0),
90  nProc_(0),
91  nLumiBlock_(0),
92  hltBx_(0),
93  hltBxVsPath_(0),
94  hltOverallScaler_(0),
95  hltOverallScalerN_(0),
96  diagnostic_(0),
97  sentPaths_(false),
98  monitorDaemon_(ps.getUntrackedParameter<bool>("MonitorDaemon", false)),
99  nev_(0),
100  nLumi_(0) {
101  LogDebug("HLTScalers") << "HLTScalers: constructor....";
102 }
103 
105  LogDebug("HLTScalers") << "HLTScalers::beginRun, run " << run.id();
106 
107  // HLT config does not change within runs!
108  bool changed = false;
109 
110  // clear vector pairPDPaths_
111  pairPDPaths_.clear();
112 
113  if (not hltConfig_.init(run, c, processname_, changed)) {
114  edm::LogError("TrigXMonitor") << "HLTConfigProvider failed to initialize.";
115  } else {
116  // check if trigger name in (new) config
117  // cout << "Available TriggerNames are: " << endl;
118  // hltConfig_.dump("Triggers");
119 
120  if (hltConfig_.streamIndex("A") < hltConfig_.streamNames().size()) {
121  // get hold of PD names and constituent path names
122  const std::vector<std::string>& PD = hltConfig_.streamContent("A");
123 
124  for (unsigned int i = 0; i < PD.size(); i++) {
125  const std::vector<std::string>& datasetPaths =
127  pairPDPaths_.push_back(make_pair(PD[i], datasetPaths));
128  }
129 
130  // push stream A and its PDs
131  pairPDPaths_.push_back(make_pair("A", PD));
132 
133  } else {
134  LogDebug("HLTScalers")
135  << "HLTScalers::beginRun, steamm A not in the HLT menu ";
136  }
137  }
138 }
139 
141  edm::EventSetup const&) {
142  std::string rawdir(folderName_ + "/raw");
143  iBooker.setCurrentFolder(rawdir);
144 
145  nProc_ = iBooker.bookInt("nProcessed");
146  nLumiBlock_ = iBooker.bookInt("nLumiBlock");
147  diagnostic_ = iBooker.book1D("hltMerge", "HLT merging diagnostic", 1, 0.5, 1.5);
148  // fill for ever accepted event
150  iBooker.book1D("hltOverallScaler", "HLT Overall Scaler", 1, 0.5, 1.5);
151  hltOverallScalerN_ = iBooker.book1D("hltOverallScalerN",
152  "Reset HLT Overall Scaler", 1, 0.5, 1.5);
153 
154  // DQM: Previously the number of trigger paths was determined on the first
155  // event, by taking the size of the htlResults to book the histogram.
156  // Now we use the size of the hltConfig instead.
157  int npath = hltConfig_.size();
158  unsigned int nPD = pairPDPaths_.size();
159 
160  // need to get maxModules dynamically
161  int maxModules = 200;
162 
163  scalersPD_ = iBooker.book1D("pdScalers", "PD scalers (stream A)", nPD, -0.5,
164  nPD - 0.5);
166  iBooker.book2D("detailedHltScalers", "HLT Scalers", npath, -0.5,
167  npath - 0.5, maxModules, 0, maxModules - 1);
168  scalers_ =
169  iBooker.book1D("hltScalers", "HLT scalers", npath, -0.5, npath - 0.5);
170  scalersN_ = iBooker.book1D("hltScalersN", "Reset HLT scalers", npath, -0.5,
171  npath - 0.5);
172  scalersException_ = iBooker.book1D("hltExceptions", "HLT Exception scalers",
173  npath, -0.5, npath - 0.5);
175  iBooker.book2D("hltCorrelations", "HLT Scalers", npath, -0.5, npath - 0.5,
176  npath, -0.5, npath - 0.5);
177 
178  // these two belong in top-level
179  iBooker.setCurrentFolder(folderName_);
180  hltBxVsPath_ = iBooker.book2D("hltBxVsPath", "HLT Accept vs Bunch Number",
181  3600, -0.5, 3599.5, npath, -0.5, npath - 0.5);
182  hltBx_ =
183  iBooker.book1D("hltBx", "Bx of HLT Accepted Events ", 3600, -0.5, 3599.5);
184 }
185 
187  const edm::EventSetup& c) {
188  LogDebug("HLTScalers") << "Start of luminosity block.";
189  // reset the N guys
190  if (scalersN_) scalersN_->Reset();
192 }
193 
195  nProc_->Fill(++nev_);
196  diagnostic_->setBinContent(1, 1); // this ME is never touched -
197  // it just tells you how the merging is doing.
198 
200  bool b = e.getByToken(trigResultsSource_, hltResults);
201  if ( !b ) {
202  Labels l;
204 
205  edm::LogInfo("HLTScalers") << "getByLabel for TriggerResults failed"
206  << " with label " << l.module;
207  return;
208  }
209  int npath = hltResults->size();
210  unsigned int nPD = pairPDPaths_.size();
212 
213  const edm::TriggerNames& trigNames = e.triggerNames(*hltResults);
214  // for some reason this doesn't appear to work on the first event sometimes
215  if (!sentPaths_) {
216  const edm::TriggerNames& names = e.triggerNames(*hltResults);
217 
218  // save path names in DQM-accessible format
219  int q = 0;
220  for (TriggerNames::Strings::const_iterator j = names.triggerNames().begin();
221  j != names.triggerNames().end(); ++j) {
222  LogDebug("HLTScalers") << q << ": " << *j;
223  ++q;
224  scalers_->getTH1()->GetXaxis()->SetBinLabel(q, j->c_str());
225  }
226 
227  for (unsigned int i = 0; i < nPD; i++) {
228  LogDebug("HLTScalers") << i << ": " << pairPDPaths_[i].first << std::endl;
229  scalersPD_->getTH1()->GetXaxis()->SetBinLabel(
230  i + 1, pairPDPaths_[i].first.c_str());
231  }
232 
233  sentPaths_ = true;
234  }
235 
236  bool accept = false;
237  int bx = e.bunchCrossing();
238  for (int i = 0; i < npath; ++i) {
239  // state returns 0 on ready, 1 on accept, 2 on fail, 3 on exception.
240  // these are defined in HLTEnums.h
241  for (unsigned int j = 0; j < hltResults->index(i); ++j) {
243  }
244  if (hltResults->state(i) == hlt::Pass) {
245  scalers_->Fill(i);
246  scalersN_->Fill(i);
247  hltBxVsPath_->Fill(bx, i);
248  accept = true;
249  for (int j = i + 1; j < npath; ++j) {
250  if (hltResults->state(j) == hlt::Pass) {
251  hltCorrelations_->Fill(i, j); // fill
253  }
254  }
255  } else if (hltResults->state(i) == hlt::Exception) {
257  }
258  }
259  if (accept) {
260  hltOverallScaler_->Fill(1.0);
261  hltOverallScalerN_->Fill(1.0);
262  hltBx_->Fill(int(bx));
263  }
264 
265  bool anyGroupPassed = false;
266  for (unsigned int mi = 0; mi < pairPDPaths_.size(); mi++) {
267  bool groupPassed = false;
268 
269  for (unsigned int i = 0; i < pairPDPaths_[mi].second.size(); i++) {
270  // string hltPathName = hist_2d->GetXaxis()->GetBinLabel(i);
271  std::string hltPathName = pairPDPaths_[mi].second[i];
272 
273  // check if this is hlt path name
274  // unsigned int pathByIndex = triggerNames.triggerIndex(hltPathName);
275  unsigned int pathByIndex =
276  trigNames.triggerIndex(pairPDPaths_[mi].second[i]);
277  if (pathByIndex >= hltResults->size()) continue;
278 
279  // check if its L1 passed
280  // comment out below but set groupL1Passed to true always
281  // if(hasL1Passed(hltPathName,triggerNames)) groupL1Passed = true;
282  // groupL1Passed = true;
283 
284  // Fill HLTPassed Matrix and HLTPassFail Matrix
285  // --------------------------------------------------------
286 
287  if (hltResults->accept(pathByIndex)) {
288  groupPassed = true;
289  break;
290  }
291  }
292 
293  if (groupPassed) {
294  scalersPD_->Fill(mi);
295  anyGroupPassed = true;
296  }
297  }
298 
299  if (anyGroupPassed) scalersPD_->Fill(pairPDPaths_.size() - 1);
300 }
301 
303  const edm::EventSetup& c) {
304  // put this in as a first-pass for figuring out the rate
305  // each lumi block is 23 seconds in length
306  nLumiBlock_->Fill(lumiSeg.id().luminosityBlock());
307 
308  LogDebug("HLTScalers") << "End of luminosity block.";
309 }
310 
312  LogDebug("HLTScalers") << "HLTScalers::endRun , run " << run.id();
313 }
#define LogDebug(id)
unsigned int size() const
number of trigger paths in trigger table
LuminosityBlockID id() const
std::vector< std::pair< std::string, std::vector< std::string > > > pairPDPaths_
Definition: HLTScalers.h:86
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:220
void setBinContent(int binx, double content)
set content of bin (1-D)
RunID const & id() const
Definition: RunBase.h:41
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: HLTScalers.cc:140
static const HistoName names[]
Definition: hltDiff.cc:316
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
HLTScalers(const edm::ParameterSet &ps)
Definition: HLTScalers.cc:79
MonitorElement * scalersN_
Definition: HLTScalers.h:91
MonitorElement * detailedScalers_
Definition: HLTScalers.h:94
std::string processname_
Definition: HLTScalers.h:85
MonitorElement * bookInt(Args &&...args)
Definition: DQMStore.h:103
int bunchCrossing() const
Definition: EventBase.h:66
MonitorElement * hltOverallScaler_
Definition: HLTScalers.h:98
void endLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
DQM Client Diagnostic should be performed here:
Definition: HLTScalers.cc:302
MonitorElement * scalers_
Definition: HLTScalers.h:90
unsigned int streamIndex(const std::string &stream) const
index of stream with name
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
MonitorElement * diagnostic_
Definition: HLTScalers.h:100
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
void Fill(long long x)
U second(std::pair< T, U > const &p)
MonitorElement * hltBxVsPath_
Definition: HLTScalers.h:97
void endRun(const edm::Run &run, const edm::EventSetup &c)
Definition: HLTScalers.cc:311
const std::vector< std::string > & streamNames() const
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
const std::vector< std::string > & streamContent(unsigned int stream) const
names of datasets in stream with index i
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
int j
Definition: DBlmapReader.cc:9
TH1 * getTH1(void) const
MonitorElement * hltCorrelations_
Definition: HLTScalers.h:93
HLTConfigProvider hltConfig_
Definition: HLTScalers.h:83
char const * module
Definition: ProductLabels.h:5
void beginLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
Definition: HLTScalers.cc:186
const std::vector< std::string > & datasetContent(unsigned int dataset) const
names of trigger paths in dataset with index i
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
MonitorElement * scalersException_
Definition: HLTScalers.h:92
MonitorElement * nLumiBlock_
Definition: HLTScalers.h:96
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * hltOverallScalerN_
Definition: HLTScalers.h:99
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
double b
Definition: hdecay.h:120
LuminosityBlockNumber_t luminosityBlock() const
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
bool sentPaths_
Definition: HLTScalers.h:102
std::string folderName_
Definition: HLTScalers.h:84
MonitorElement * hltBx_
Definition: HLTScalers.h:97
MonitorElement * nProc_
Definition: HLTScalers.h:95
volatile std::atomic< bool > shutdown_flag false
MonitorElement * scalersPD_
Definition: HLTScalers.h:89
void analyze(const edm::Event &e, const edm::EventSetup &c)
Definition: HLTScalers.cc:194
void dqmBeginRun(const edm::Run &run, const edm::EventSetup &c)
Definition: HLTScalers.cc:104
void Reset(void)
reset ME (ie. contents, errors, etc)
Definition: Run.h:43
edm::EDGetTokenT< edm::TriggerResults > trigResultsSource_
Definition: HLTScalers.h:87