1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
|
#include <string>
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/Histograms/interface/MonitorElementCollection.h"
template <typename BOOKERLIKE, typename ME, bool DOLUMI = false>
class BookerFiller {
private:
struct PolygonDef {
int nPoints;
std::vector<double> x;
std::vector<double> y;
};
std::vector<PolygonDef> m_polygons;
public:
BookerFiller(std::string folder, int howmany) {
this->howmany = howmany;
this->folder = folder;
m_polygons = {
{6, {37.5, 25.0, -0.5, -0.5, 25.0, 37.5}, {5.0, -0.5, -0.5, 10.0, 10.0, 5.0}}, // polygon-1: n, x, y
{5, {37.5, 25.0, 75.0, 62.5, 37.5}, {5.0, 10.0, 10.0, 5.0, 5.0}}, // polygon-2: n, x, y
{5, {37.5, 25.0, 75.0, 62.5, 37.5}, {5.0, -0.5, -0.5, 5.0, 5.0}}, // polygon-3: n, x, y
{6, {62.5, 75.0, 100.0, 100.0, 75.0, 62.5}, {5.0, 10.0, 10.0, -0.5, -0.5, 5.0}} // polygon-3: n, x, y
};
}
BookerFiller() {}
void bookall(BOOKERLIKE& ibooker) {
mes_1D.clear();
mes_2D.clear();
mes_3D.clear();
for (int i = 0; i < howmany; i++) {
ibooker.setCurrentFolder(folder);
auto num = std::to_string(i);
// this does *not* test most of the booking API, just one ME of each kind.
mes_1D.push_back(ibooker.bookFloat("float" + num));
mes_1D.push_back(ibooker.bookInt("int" + num));
mes_1D.push_back(ibooker.book1D("th1f" + num, "1D Float Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1S("th1s" + num, "1D Short Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1I("th1i" + num, "1D Integer Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1DD("th1d" + num, "1D Double Histogram " + num, 101, -0.5, 100.5));
mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
auto th2poly_main =
ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5);
int nCells = th2poly_main->getNcells();
// Only add bins if they don't exist yet (nCells<=9)
if (nCells <= 9) {
for (const auto& poly : m_polygons) {
th2poly_main->addBin(poly.nPoints, poly.x.data(), poly.y.data());
}
}
mes_2D.push_back(th2poly_main);
mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(
ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_3D.push_back(
ibooker.book3D("th3f" + num, "3D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
mes_3D.push_back(ibooker.bookProfile2D(
"thprofile2d" + num, "2D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
if (DOLUMI) {
auto scope = typename BOOKERLIKE::UseLumiScope(ibooker);
ibooker.setCurrentFolder(folder + "/lumi");
mes_1D.push_back(ibooker.bookFloat("float" + num));
mes_1D.push_back(ibooker.bookInt("int" + num));
mes_1D.push_back(ibooker.book1D("th1f" + num, "1D Float Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1S("th1s" + num, "1D Short Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1I("th1i" + num, "1D Integer Histogram " + num, 101, -0.5, 100.5));
mes_1D.push_back(ibooker.book1DD("th1d" + num, "1D Double Histogram " + num, 101, -0.5, 100.5));
mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
auto th2poly_lumi =
ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5);
int nCells_lumi = th2poly_lumi->getNcells();
// Only add bins if they don't exist yet (nCells<=9)
if (nCells_lumi <= 9) {
for (const auto& poly : m_polygons) {
th2poly_lumi->addBin(poly.nPoints, poly.x.data(), poly.y.data());
}
}
mes_2D.push_back(th2poly_lumi);
mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_2D.push_back(
ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
mes_3D.push_back(
ibooker.book3D("th3f" + num, "3D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
mes_3D.push_back(ibooker.bookProfile2D(
"thprofile2d" + num, "2D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
}
}
}
void fillall(double x, double y, double z) const {
for (auto me : mes_1D) {
me->Fill(x);
}
for (auto me : mes_2D) {
me->Fill(x, y);
}
for (auto me : mes_3D) {
me->Fill(x, y, z);
}
}
std::vector<ME*> mes_1D;
std::vector<ME*> mes_2D;
std::vector<ME*> mes_3D;
std::string folder;
int howmany;
};
#include "DQMServices/Core/interface/DQMEDAnalyzer.h"
class TestDQMEDAnalyzer : public DQMEDAnalyzer {
public:
explicit TestDQMEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Normal/test")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("test", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
mymes_.bookall(ibooker);
}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
BookerFiller<DQMStore::IBooker, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMEDAnalyzer);
#include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
class TestDQMOneEDAnalyzer : public DQMOneEDAnalyzer<> {
public:
explicit TestDQMOneEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMOneEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "One/testone")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testone", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
mymes_.bookall(ibooker);
}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
BookerFiller<DQMStore::IBooker, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMOneEDAnalyzer);
class TestDQMOneFillRunEDAnalyzer : public DQMOneEDAnalyzer<> {
public:
explicit TestDQMOneFillRunEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMOneFillRunEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "One/testonefillrun")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testonefillrun", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
mymes_.bookall(ibooker);
}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
void dqmEndRun(edm::Run const& run, edm::EventSetup const&) override { mymes_.fillall(0, run.run(), myvalue_); }
BookerFiller<DQMStore::IBooker, MonitorElement> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMOneFillRunEDAnalyzer);
class TestDQMOneLumiEDAnalyzer : public DQMOneLumiEDAnalyzer<> {
public:
explicit TestDQMOneLumiEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMOneLumiEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "One/testonelumi")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testonelumi", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
mymes_.bookall(ibooker);
}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
BookerFiller<DQMStore::IBooker, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMOneLumiEDAnalyzer);
class TestDQMOneLumiFillLumiEDAnalyzer : public DQMOneLumiEDAnalyzer<> {
public:
explicit TestDQMOneLumiFillLumiEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMOneLumiFillLumiEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "One/testonelumifilllumi")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testonelumifilllumi", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
mymes_.bookall(ibooker);
}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
void dqmBeginLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {
mymes_.fillall(lumi.luminosityBlock(), lumi.run(), myvalue_);
}
BookerFiller<DQMStore::IBooker, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMOneLumiFillLumiEDAnalyzer);
#include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
typedef BookerFiller<dqm::reco::DQMStore::IBooker, dqm::reco::MonitorElement> TestHistograms;
class TestDQMGlobalEDAnalyzer : public DQMGlobalEDAnalyzer<TestHistograms> {
public:
explicit TestDQMGlobalEDAnalyzer(const edm::ParameterSet& iConfig)
: folder_(iConfig.getParameter<std::string>("folder")),
howmany_(iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMGlobalEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Global/testglobal")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testglobal", desc);
}
private:
void bookHistograms(DQMStore::IBooker& ibooker,
edm::Run const&,
edm::EventSetup const&,
TestHistograms& h) const override {
h = TestHistograms(this->folder_, this->howmany_);
h.bookall(ibooker);
}
void dqmAnalyze(edm::Event const& iEvent, edm::EventSetup const&, TestHistograms const& h) const override {
h.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
private:
std::string folder_;
int howmany_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMGlobalEDAnalyzer);
class TestDQMGlobalRunSummaryEDAnalyzer : public DQMGlobalRunSummaryEDAnalyzer<TestHistograms, int> {
public:
explicit TestDQMGlobalRunSummaryEDAnalyzer(const edm::ParameterSet& iConfig)
: folder_(iConfig.getParameter<std::string>("folder")),
howmany_(iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {}
~TestDQMGlobalRunSummaryEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Global/testglobalrunsummary")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testglobalrunsummary", desc);
}
private:
std::shared_ptr<int> globalBeginRunSummary(edm::Run const&, edm::EventSetup const&) const override {
return std::make_shared<int>(0);
}
void bookHistograms(DQMStore::IBooker& ibooker,
edm::Run const&,
edm::EventSetup const&,
TestHistograms& h) const override {
h = TestHistograms(this->folder_, this->howmany_);
h.bookall(ibooker);
}
void dqmAnalyze(edm::Event const& iEvent, edm::EventSetup const&, TestHistograms const& h) const override {
h.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
void streamEndRunSummary(edm::StreamID, edm::Run const&, edm::EventSetup const&, int* runSummaryCache) const override {
(*runSummaryCache) += 1;
}
void globalEndRunSummary(edm::Run const&, edm::EventSetup const&, int* runSummaryCache) const override {}
void dqmEndRun(edm::Run const& run,
edm::EventSetup const& setup,
TestHistograms const& h,
int const& runSummaryCache) const override {}
private:
std::string folder_;
int howmany_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestDQMGlobalRunSummaryEDAnalyzer);
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
class TestLegacyEDAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
public:
typedef dqm::legacy::DQMStore DQMStore;
typedef dqm::legacy::MonitorElement MonitorElement;
explicit TestLegacyEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {
usesResource("DQMStore");
}
~TestLegacyEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Legacy/testlegacy")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testlegacy", desc);
}
private:
void beginRun(edm::Run const&, edm::EventSetup const&) override {
edm::Service<DQMStore> store;
mymes_.bookall(*store);
}
void endRun(edm::Run const&, edm::EventSetup const&) override {}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
}
BookerFiller<DQMStore, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestLegacyEDAnalyzer);
class TestLegacyFillRunEDAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
public:
typedef dqm::legacy::DQMStore DQMStore;
typedef dqm::legacy::MonitorElement MonitorElement;
explicit TestLegacyFillRunEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {
usesResource("DQMStore");
}
~TestLegacyFillRunEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Legacy/testlegacyfillrun")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testlegacyfillrun", desc);
}
private:
void beginRun(edm::Run const& run, edm::EventSetup const&) override {
edm::Service<DQMStore> store;
mymes_.bookall(*store);
mymes_.fillall(0, run.run(), myvalue_);
}
void endRun(edm::Run const&, edm::EventSetup const&) override {}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
BookerFiller<DQMStore, MonitorElement> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestLegacyFillRunEDAnalyzer);
class TestLegacyFillLumiEDAnalyzer
: public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns, edm::one::WatchLuminosityBlocks> {
public:
typedef dqm::legacy::DQMStore DQMStore;
typedef dqm::legacy::MonitorElement MonitorElement;
explicit TestLegacyFillLumiEDAnalyzer(const edm::ParameterSet& iConfig)
: mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
myvalue_(iConfig.getParameter<double>("value")) {
usesResource("DQMStore");
}
~TestLegacyFillLumiEDAnalyzer() override {}
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folder", "Legacy/testlegacyfilllumi")->setComment("Where to put all the histograms");
desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
descriptions.add("testlegacyfilllumi", desc);
}
private:
void beginRun(edm::Run const&, edm::EventSetup const&) override {
edm::Service<DQMStore> store;
mymes_.bookall(*store);
}
void endRun(edm::Run const&, edm::EventSetup const&) override {}
void beginLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {
mymes_.fillall(lumi.luminosityBlock(), lumi.run(), myvalue_);
}
void endLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {}
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
BookerFiller<DQMStore, MonitorElement, /* DOLUMI */ true> mymes_;
double myvalue_;
};
DEFINE_FWK_MODULE(TestLegacyFillLumiEDAnalyzer);
|