File indexing completed on 2025-03-10 23:53:32
0001 #include <string>
0002
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "DataFormats/Histograms/interface/MonitorElementCollection.h"
0007
0008 template <typename BOOKERLIKE, typename ME, bool DOLUMI = false>
0009 class BookerFiller {
0010 private:
0011 struct PolygonDef {
0012 int nPoints;
0013 std::vector<double> x;
0014 std::vector<double> y;
0015 };
0016 std::vector<PolygonDef> m_polygons;
0017
0018 public:
0019 BookerFiller(std::string folder, int howmany) {
0020 this->howmany = howmany;
0021 this->folder = folder;
0022 m_polygons = {
0023 {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}},
0024 {5, {37.5, 25.0, 75.0, 62.5, 37.5}, {5.0, 10.0, 10.0, 5.0, 5.0}},
0025 {5, {37.5, 25.0, 75.0, 62.5, 37.5}, {5.0, -0.5, -0.5, 5.0, 5.0}},
0026 {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}}
0027 };
0028 }
0029
0030 BookerFiller() {}
0031
0032 void bookall(BOOKERLIKE& ibooker) {
0033 mes_1D.clear();
0034 mes_2D.clear();
0035 mes_3D.clear();
0036
0037 for (int i = 0; i < howmany; i++) {
0038 ibooker.setCurrentFolder(folder);
0039 auto num = std::to_string(i);
0040
0041 mes_1D.push_back(ibooker.bookFloat("float" + num));
0042 mes_1D.push_back(ibooker.bookInt("int" + num));
0043 mes_1D.push_back(ibooker.book1D("th1f" + num, "1D Float Histogram " + num, 101, -0.5, 100.5));
0044 mes_1D.push_back(ibooker.book1S("th1s" + num, "1D Short Histogram " + num, 101, -0.5, 100.5));
0045 mes_1D.push_back(ibooker.book1I("th1i" + num, "1D Integer Histogram " + num, 101, -0.5, 100.5));
0046 mes_1D.push_back(ibooker.book1DD("th1d" + num, "1D Double Histogram " + num, 101, -0.5, 100.5));
0047
0048 mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0049 mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0050 mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0051 auto th2poly_main =
0052 ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5);
0053 int nCells = th2poly_main->getNcells();
0054
0055 if (nCells <= 9) {
0056 for (const auto& poly : m_polygons) {
0057 th2poly_main->addBin(poly.nPoints, poly.x.data(), poly.y.data());
0058 }
0059 }
0060 mes_2D.push_back(th2poly_main);
0061 mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0062 mes_2D.push_back(
0063 ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0064
0065 mes_3D.push_back(
0066 ibooker.book3D("th3f" + num, "3D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
0067 mes_3D.push_back(ibooker.bookProfile2D(
0068 "thprofile2d" + num, "2D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
0069
0070 if (DOLUMI) {
0071 auto scope = typename BOOKERLIKE::UseLumiScope(ibooker);
0072 ibooker.setCurrentFolder(folder + "/lumi");
0073
0074 mes_1D.push_back(ibooker.bookFloat("float" + num));
0075 mes_1D.push_back(ibooker.bookInt("int" + num));
0076 mes_1D.push_back(ibooker.book1D("th1f" + num, "1D Float Histogram " + num, 101, -0.5, 100.5));
0077 mes_1D.push_back(ibooker.book1S("th1s" + num, "1D Short Histogram " + num, 101, -0.5, 100.5));
0078 mes_1D.push_back(ibooker.book1I("th1i" + num, "1D Integer Histogram " + num, 101, -0.5, 100.5));
0079 mes_1D.push_back(ibooker.book1DD("th1d" + num, "1D Double Histogram " + num, 101, -0.5, 100.5));
0080
0081 mes_2D.push_back(ibooker.book2D("th2f" + num, "2D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0082 mes_2D.push_back(ibooker.book2S("th2s" + num, "2D Short Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0083 mes_2D.push_back(ibooker.book2DD("th2d" + num, "2D Double Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0084 auto th2poly_lumi =
0085 ibooker.book2DPoly("th2poly" + num, "2D Polygonal Double Histogram " + num, -0.5, 100.5, -0.5, 10.5);
0086 int nCells_lumi = th2poly_lumi->getNcells();
0087
0088 if (nCells_lumi <= 9) {
0089 for (const auto& poly : m_polygons) {
0090 th2poly_lumi->addBin(poly.nPoints, poly.x.data(), poly.y.data());
0091 }
0092 }
0093 mes_2D.push_back(th2poly_lumi);
0094 mes_2D.push_back(ibooker.book2I("th2i" + num, "2D Integer Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0095 mes_2D.push_back(
0096 ibooker.bookProfile("tprofile" + num, "1D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5));
0097
0098 mes_3D.push_back(
0099 ibooker.book3D("th3f" + num, "3D Float Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
0100 mes_3D.push_back(ibooker.bookProfile2D(
0101 "thprofile2d" + num, "2D Profile Histogram " + num, 101, -0.5, 100.5, 11, -0.5, 10.5, 3, -0.5, 2.5));
0102 }
0103 }
0104 }
0105
0106 void fillall(double x, double y, double z) const {
0107 for (auto me : mes_1D) {
0108 me->Fill(x);
0109 }
0110 for (auto me : mes_2D) {
0111 me->Fill(x, y);
0112 }
0113 for (auto me : mes_3D) {
0114 me->Fill(x, y, z);
0115 }
0116 }
0117
0118 std::vector<ME*> mes_1D;
0119 std::vector<ME*> mes_2D;
0120 std::vector<ME*> mes_3D;
0121 std::string folder;
0122 int howmany;
0123 };
0124
0125 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0126 class TestDQMEDAnalyzer : public DQMEDAnalyzer {
0127 public:
0128 explicit TestDQMEDAnalyzer(const edm::ParameterSet& iConfig)
0129 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0130 myvalue_(iConfig.getParameter<double>("value")) {}
0131
0132 ~TestDQMEDAnalyzer() override {}
0133
0134 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0135 edm::ParameterSetDescription desc;
0136 desc.add<std::string>("folder", "Normal/test")->setComment("Where to put all the histograms");
0137 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0138 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0139 descriptions.add("test", desc);
0140 }
0141
0142 private:
0143 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
0144 mymes_.bookall(ibooker);
0145 }
0146
0147 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0148 mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0149 }
0150
0151 BookerFiller<DQMStore::IBooker, MonitorElement, true> mymes_;
0152 double myvalue_;
0153 };
0154 DEFINE_FWK_MODULE(TestDQMEDAnalyzer);
0155
0156 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0157 class TestDQMOneEDAnalyzer : public DQMOneEDAnalyzer<> {
0158 public:
0159 explicit TestDQMOneEDAnalyzer(const edm::ParameterSet& iConfig)
0160 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0161 myvalue_(iConfig.getParameter<double>("value")) {}
0162
0163 ~TestDQMOneEDAnalyzer() override {}
0164
0165 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0166 edm::ParameterSetDescription desc;
0167 desc.add<std::string>("folder", "One/testone")->setComment("Where to put all the histograms");
0168 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0169 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0170 descriptions.add("testone", desc);
0171 }
0172
0173 private:
0174 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
0175 mymes_.bookall(ibooker);
0176 }
0177
0178 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0179 mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0180 }
0181
0182 BookerFiller<DQMStore::IBooker, MonitorElement, true> mymes_;
0183 double myvalue_;
0184 };
0185 DEFINE_FWK_MODULE(TestDQMOneEDAnalyzer);
0186
0187 class TestDQMOneFillRunEDAnalyzer : public DQMOneEDAnalyzer<> {
0188 public:
0189 explicit TestDQMOneFillRunEDAnalyzer(const edm::ParameterSet& iConfig)
0190 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0191 myvalue_(iConfig.getParameter<double>("value")) {}
0192
0193 ~TestDQMOneFillRunEDAnalyzer() override {}
0194
0195 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0196 edm::ParameterSetDescription desc;
0197 desc.add<std::string>("folder", "One/testonefillrun")->setComment("Where to put all the histograms");
0198 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0199 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0200 descriptions.add("testonefillrun", desc);
0201 }
0202
0203 private:
0204 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
0205 mymes_.bookall(ibooker);
0206 }
0207
0208 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
0209 void dqmEndRun(edm::Run const& run, edm::EventSetup const&) override { mymes_.fillall(0, run.run(), myvalue_); }
0210
0211 BookerFiller<DQMStore::IBooker, MonitorElement> mymes_;
0212 double myvalue_;
0213 };
0214 DEFINE_FWK_MODULE(TestDQMOneFillRunEDAnalyzer);
0215
0216 class TestDQMOneLumiEDAnalyzer : public DQMOneLumiEDAnalyzer<> {
0217 public:
0218 explicit TestDQMOneLumiEDAnalyzer(const edm::ParameterSet& iConfig)
0219 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0220 myvalue_(iConfig.getParameter<double>("value")) {}
0221
0222 ~TestDQMOneLumiEDAnalyzer() override {}
0223
0224 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0225 edm::ParameterSetDescription desc;
0226 desc.add<std::string>("folder", "One/testonelumi")->setComment("Where to put all the histograms");
0227 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0228 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0229 descriptions.add("testonelumi", desc);
0230 }
0231
0232 private:
0233 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
0234 mymes_.bookall(ibooker);
0235 }
0236
0237 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0238 mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0239 }
0240
0241 BookerFiller<DQMStore::IBooker, MonitorElement, true> mymes_;
0242 double myvalue_;
0243 };
0244 DEFINE_FWK_MODULE(TestDQMOneLumiEDAnalyzer);
0245
0246 class TestDQMOneLumiFillLumiEDAnalyzer : public DQMOneLumiEDAnalyzer<> {
0247 public:
0248 explicit TestDQMOneLumiFillLumiEDAnalyzer(const edm::ParameterSet& iConfig)
0249 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0250 myvalue_(iConfig.getParameter<double>("value")) {}
0251
0252 ~TestDQMOneLumiFillLumiEDAnalyzer() override {}
0253
0254 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0255 edm::ParameterSetDescription desc;
0256 desc.add<std::string>("folder", "One/testonelumifilllumi")->setComment("Where to put all the histograms");
0257 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0258 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0259 descriptions.add("testonelumifilllumi", desc);
0260 }
0261
0262 private:
0263 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) override {
0264 mymes_.bookall(ibooker);
0265 }
0266
0267 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
0268 void dqmBeginLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {
0269 mymes_.fillall(lumi.luminosityBlock(), lumi.run(), myvalue_);
0270 }
0271
0272 BookerFiller<DQMStore::IBooker, MonitorElement, true> mymes_;
0273 double myvalue_;
0274 };
0275 DEFINE_FWK_MODULE(TestDQMOneLumiFillLumiEDAnalyzer);
0276
0277 #include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
0278 typedef BookerFiller<dqm::reco::DQMStore::IBooker, dqm::reco::MonitorElement> TestHistograms;
0279
0280 class TestDQMGlobalEDAnalyzer : public DQMGlobalEDAnalyzer<TestHistograms> {
0281 public:
0282 explicit TestDQMGlobalEDAnalyzer(const edm::ParameterSet& iConfig)
0283 : folder_(iConfig.getParameter<std::string>("folder")),
0284 howmany_(iConfig.getParameter<int>("howmany")),
0285 myvalue_(iConfig.getParameter<double>("value")) {}
0286 ~TestDQMGlobalEDAnalyzer() override {}
0287
0288 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0289 edm::ParameterSetDescription desc;
0290 desc.add<std::string>("folder", "Global/testglobal")->setComment("Where to put all the histograms");
0291 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0292 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0293 descriptions.add("testglobal", desc);
0294 }
0295
0296 private:
0297 void bookHistograms(DQMStore::IBooker& ibooker,
0298 edm::Run const&,
0299 edm::EventSetup const&,
0300 TestHistograms& h) const override {
0301 h = TestHistograms(this->folder_, this->howmany_);
0302 h.bookall(ibooker);
0303 }
0304
0305 void dqmAnalyze(edm::Event const& iEvent, edm::EventSetup const&, TestHistograms const& h) const override {
0306 h.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0307 }
0308
0309 private:
0310 std::string folder_;
0311 int howmany_;
0312 double myvalue_;
0313 };
0314 DEFINE_FWK_MODULE(TestDQMGlobalEDAnalyzer);
0315
0316 class TestDQMGlobalRunSummaryEDAnalyzer : public DQMGlobalRunSummaryEDAnalyzer<TestHistograms, int> {
0317 public:
0318 explicit TestDQMGlobalRunSummaryEDAnalyzer(const edm::ParameterSet& iConfig)
0319 : folder_(iConfig.getParameter<std::string>("folder")),
0320 howmany_(iConfig.getParameter<int>("howmany")),
0321 myvalue_(iConfig.getParameter<double>("value")) {}
0322 ~TestDQMGlobalRunSummaryEDAnalyzer() override {}
0323
0324 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0325 edm::ParameterSetDescription desc;
0326 desc.add<std::string>("folder", "Global/testglobalrunsummary")->setComment("Where to put all the histograms");
0327 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0328 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0329 descriptions.add("testglobalrunsummary", desc);
0330 }
0331
0332 private:
0333 std::shared_ptr<int> globalBeginRunSummary(edm::Run const&, edm::EventSetup const&) const override {
0334 return std::make_shared<int>(0);
0335 }
0336
0337 void bookHistograms(DQMStore::IBooker& ibooker,
0338 edm::Run const&,
0339 edm::EventSetup const&,
0340 TestHistograms& h) const override {
0341 h = TestHistograms(this->folder_, this->howmany_);
0342 h.bookall(ibooker);
0343 }
0344
0345 void dqmAnalyze(edm::Event const& iEvent, edm::EventSetup const&, TestHistograms const& h) const override {
0346 h.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0347 }
0348
0349 void streamEndRunSummary(edm::StreamID, edm::Run const&, edm::EventSetup const&, int* runSummaryCache) const override {
0350 (*runSummaryCache) += 1;
0351 }
0352
0353 void globalEndRunSummary(edm::Run const&, edm::EventSetup const&, int* runSummaryCache) const override {}
0354
0355 void dqmEndRun(edm::Run const& run,
0356 edm::EventSetup const& setup,
0357 TestHistograms const& h,
0358 int const& runSummaryCache) const override {}
0359
0360 private:
0361 std::string folder_;
0362 int howmany_;
0363 double myvalue_;
0364 };
0365 DEFINE_FWK_MODULE(TestDQMGlobalRunSummaryEDAnalyzer);
0366
0367 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0368 class TestLegacyEDAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0369 public:
0370 typedef dqm::legacy::DQMStore DQMStore;
0371 typedef dqm::legacy::MonitorElement MonitorElement;
0372
0373 explicit TestLegacyEDAnalyzer(const edm::ParameterSet& iConfig)
0374 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0375 myvalue_(iConfig.getParameter<double>("value")) {
0376 usesResource("DQMStore");
0377 }
0378
0379 ~TestLegacyEDAnalyzer() override {}
0380
0381 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0382 edm::ParameterSetDescription desc;
0383 desc.add<std::string>("folder", "Legacy/testlegacy")->setComment("Where to put all the histograms");
0384 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0385 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0386 descriptions.add("testlegacy", desc);
0387 }
0388
0389 private:
0390 void beginRun(edm::Run const&, edm::EventSetup const&) override {
0391 edm::Service<DQMStore> store;
0392 mymes_.bookall(*store);
0393 }
0394 void endRun(edm::Run const&, edm::EventSetup const&) override {}
0395
0396 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0397 mymes_.fillall(iEvent.luminosityBlock(), iEvent.run(), myvalue_);
0398 }
0399
0400 BookerFiller<DQMStore, MonitorElement, true> mymes_;
0401 double myvalue_;
0402 };
0403 DEFINE_FWK_MODULE(TestLegacyEDAnalyzer);
0404
0405 class TestLegacyFillRunEDAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0406 public:
0407 typedef dqm::legacy::DQMStore DQMStore;
0408 typedef dqm::legacy::MonitorElement MonitorElement;
0409
0410 explicit TestLegacyFillRunEDAnalyzer(const edm::ParameterSet& iConfig)
0411 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0412 myvalue_(iConfig.getParameter<double>("value")) {
0413 usesResource("DQMStore");
0414 }
0415
0416 ~TestLegacyFillRunEDAnalyzer() override {}
0417
0418 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0419 edm::ParameterSetDescription desc;
0420 desc.add<std::string>("folder", "Legacy/testlegacyfillrun")->setComment("Where to put all the histograms");
0421 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0422 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0423 descriptions.add("testlegacyfillrun", desc);
0424 }
0425
0426 private:
0427 void beginRun(edm::Run const& run, edm::EventSetup const&) override {
0428 edm::Service<DQMStore> store;
0429 mymes_.bookall(*store);
0430 mymes_.fillall(0, run.run(), myvalue_);
0431 }
0432 void endRun(edm::Run const&, edm::EventSetup const&) override {}
0433
0434 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
0435
0436 BookerFiller<DQMStore, MonitorElement> mymes_;
0437 double myvalue_;
0438 };
0439 DEFINE_FWK_MODULE(TestLegacyFillRunEDAnalyzer);
0440
0441 class TestLegacyFillLumiEDAnalyzer
0442 : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns, edm::one::WatchLuminosityBlocks> {
0443 public:
0444 typedef dqm::legacy::DQMStore DQMStore;
0445 typedef dqm::legacy::MonitorElement MonitorElement;
0446
0447 explicit TestLegacyFillLumiEDAnalyzer(const edm::ParameterSet& iConfig)
0448 : mymes_(iConfig.getParameter<std::string>("folder"), iConfig.getParameter<int>("howmany")),
0449 myvalue_(iConfig.getParameter<double>("value")) {
0450 usesResource("DQMStore");
0451 }
0452
0453 ~TestLegacyFillLumiEDAnalyzer() override {}
0454
0455 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0456 edm::ParameterSetDescription desc;
0457 desc.add<std::string>("folder", "Legacy/testlegacyfilllumi")->setComment("Where to put all the histograms");
0458 desc.add<int>("howmany", 1)->setComment("How many copies of each ME to put");
0459 desc.add<double>("value", 1)->setComment("Which value to use on the third axis (first two are lumi and run)");
0460 descriptions.add("testlegacyfilllumi", desc);
0461 }
0462
0463 private:
0464 void beginRun(edm::Run const&, edm::EventSetup const&) override {
0465 edm::Service<DQMStore> store;
0466 mymes_.bookall(*store);
0467 }
0468 void endRun(edm::Run const&, edm::EventSetup const&) override {}
0469
0470 void beginLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {
0471 mymes_.fillall(lumi.luminosityBlock(), lumi.run(), myvalue_);
0472 }
0473 void endLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const&) override {}
0474
0475 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {}
0476
0477 BookerFiller<DQMStore, MonitorElement, true> mymes_;
0478 double myvalue_;
0479 };
0480 DEFINE_FWK_MODULE(TestLegacyFillLumiEDAnalyzer);