Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:36

0001 // Implementation based on FastTimerService.
0002 // Target usage is to generate fake DAQ histograms for performance tests on RUBU machines (effect of running fastHadd) and to include in unit tests.
0003 
0004 // CMSSW headers
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/LuminosityBlock.h"
0007 #include "FWCore/Framework/interface/Run.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
0013 
0014 namespace {
0015 
0016   struct RunBasedHistograms {
0017     // overall event count and event types
0018     dqm::reco::MonitorElement *events_processed;
0019     std::vector<dqm::reco::MonitorElement *> element_array;
0020 
0021     RunBasedHistograms()
0022         :  // overall event count and event types
0023           events_processed(nullptr),
0024           element_array() {}
0025   };
0026 }  // namespace
0027 
0028 class DaqTestHistograms : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
0029 public:
0030   explicit DaqTestHistograms(edm::ParameterSet const &);
0031   ~DaqTestHistograms() override = default;
0032 
0033   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0034 
0035 private:
0036   void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override;
0037   void bookHistograms(DQMStore::IBooker &,
0038                       edm::Run const &,
0039                       edm::EventSetup const &,
0040                       RunBasedHistograms &) const override;
0041   void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override;
0042 
0043   // module configuration
0044   const std::string m_dqm_path;
0045   const uint32_t m_lumisections_range;
0046   const uint32_t m_num_histograms;
0047 };
0048 
0049 void DaqTestHistograms::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0050   edm::ParameterSetDescription desc;
0051   desc.addUntracked<std::string>("dqmPath", "DAQTEST/Test");
0052   desc.addUntracked<uint32_t>("lumisectionRange", 25);
0053   desc.addUntracked<uint32_t>("numberOfHistograms", 10);
0054   descriptions.add("dqmHLTTestMonitor", desc);
0055 }
0056 
0057 DaqTestHistograms::DaqTestHistograms(edm::ParameterSet const &config)
0058     :  // module configuration
0059       m_dqm_path(config.getUntrackedParameter<std::string>("dqmPath")),
0060       m_lumisections_range(config.getUntrackedParameter<uint32_t>("lumisectionRange")),
0061       m_num_histograms(config.getUntrackedParameter<uint32_t>("numberOfHistograms")) {}
0062 
0063 void DaqTestHistograms::dqmBeginRun(edm::Run const &run,
0064                                     edm::EventSetup const &setup,
0065                                     RunBasedHistograms &histograms) const {
0066   histograms.element_array.resize(m_num_histograms);
0067 }
0068 
0069 void DaqTestHistograms::bookHistograms(DQMStore::IBooker &booker,
0070                                        edm::Run const &run,
0071                                        edm::EventSetup const &setup,
0072                                        RunBasedHistograms &histograms) const {
0073   // book the overall event count and event types histograms
0074   booker.setCurrentFolder(m_dqm_path);
0075   histograms.events_processed = booker.book1D(
0076       "events", "Processed events vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0077   for (size_t i = 0; i < histograms.element_array.size(); i++) {
0078     std::stringstream strs;
0079     strs << "element " << i;
0080     std::stringstream strs2;
0081     strs2 << "e vs ls " << i;
0082     histograms.element_array[i] =
0083         booker.book1D(strs.str(), strs2.str(), m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0084   }
0085 }
0086 
0087 void DaqTestHistograms::dqmAnalyze(edm::Event const &event,
0088                                    edm::EventSetup const &setup,
0089                                    RunBasedHistograms const &histograms) const {
0090   unsigned int lumisection = event.luminosityBlock();
0091 
0092   histograms.events_processed->Fill(lumisection);
0093   for (size_t i = 0; i < histograms.element_array.size(); i++) {
0094     histograms.element_array[i]->Fill(lumisection);
0095   }
0096 }
0097 
0098 //define this as a plug-in
0099 #include "FWCore/Framework/interface/MakerMacros.h"
0100 DEFINE_FWK_MODULE(DaqTestHistograms);