Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <string>
0002 #include <vector>
0003 #include <map>
0004 
0005 #include "DQM/TrackingMonitor/interface/GetLumi.h"
0006 #include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
0007 #include "DataFormats/L1TGlobal/interface/GlobalAlgBlk.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/ParameterSet/interface/Registry.h"
0016 #include "FWCore/PrescaleService/interface/PrescaleService.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 #include "FWCore/Utilities/interface/EDGetToken.h"
0019 
0020 namespace {
0021   typedef dqm::reco::DQMStore DQMStore;
0022   typedef dqm::reco::MonitorElement MonitorElement;
0023 
0024   struct MEbinning {
0025     int nbins;
0026     double xmin;
0027     double xmax;
0028   };
0029 
0030   struct Histograms {
0031     dqm::reco::MonitorElement* psColumnIndexVsLS;
0032   };
0033 }  // namespace
0034 
0035 //
0036 // class declaration
0037 //
0038 
0039 class PSMonitor : public DQMGlobalEDAnalyzer<Histograms> {
0040 public:
0041   PSMonitor(const edm::ParameterSet&);
0042   ~PSMonitor() override = default;
0043   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0044   static void fillHistoPSetDescription(edm::ParameterSetDescription& pset, int value);
0045 
0046 protected:
0047   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, Histograms&) const override;
0048   void dqmAnalyze(edm::Event const& event, edm::EventSetup const& setup, Histograms const&) const override;
0049 
0050 private:
0051   void getHistoPSet(edm::ParameterSet& pset, MEbinning& mebinning);
0052 
0053   std::string folderName_;
0054 
0055   edm::EDGetTokenT<GlobalAlgBlkBxCollection> ugtBXToken_;
0056 
0057   MEbinning ps_binning_;
0058   MEbinning ls_binning_;
0059 };
0060 
0061 // -----------------------------
0062 //  constructors and destructor
0063 // -----------------------------
0064 
0065 PSMonitor::PSMonitor(const edm::ParameterSet& config)
0066     : folderName_(config.getParameter<std::string>("folderName")),
0067       ugtBXToken_(consumes<GlobalAlgBlkBxCollection>(config.getParameter<edm::InputTag>("ugtBXInputTag"))) {
0068   edm::ParameterSet histoPSet = config.getParameter<edm::ParameterSet>("histoPSet");
0069   edm::ParameterSet psColumnPSet = histoPSet.getParameter<edm::ParameterSet>("psColumnPSet");
0070   edm::ParameterSet lsPSet = histoPSet.getParameter<edm::ParameterSet>("lsPSet");
0071 
0072   getHistoPSet(psColumnPSet, ps_binning_);
0073   getHistoPSet(lsPSet, ls_binning_);
0074 }
0075 
0076 void PSMonitor::getHistoPSet(edm::ParameterSet& pset, MEbinning& mebinning) {
0077   mebinning.nbins = pset.getParameter<int32_t>("nbins");
0078   mebinning.xmin = 0.;
0079   mebinning.xmax = double(pset.getParameter<int32_t>("nbins"));
0080 }
0081 
0082 void PSMonitor::bookHistograms(DQMStore::IBooker& booker,
0083                                edm::Run const& run,
0084                                edm::EventSetup const& setup,
0085                                Histograms& histograms) const {
0086   std::string histname, histtitle;
0087 
0088   std::string currentFolder = folderName_;
0089   booker.setCurrentFolder(currentFolder);
0090 
0091   int nbins;
0092   double xmin, xmax;
0093   std::vector<std::string> labels;
0094   edm::Service<edm::service::PrescaleService> prescaleService;
0095   if (prescaleService.isAvailable() and not prescaleService->getLvl1Labels().empty()) {
0096     labels = prescaleService->getLvl1Labels();
0097     nbins = labels.size();
0098     xmin = 0.;
0099     xmax = double(labels.size());
0100   } else {
0101     nbins = ps_binning_.nbins;
0102     xmin = ps_binning_.xmin;
0103     xmax = ps_binning_.xmax;
0104     labels.resize(nbins, "");
0105   }
0106 
0107   histname = "psColumnIndexVsLS";
0108   histtitle = "PS column index vs LS";
0109   auto me =
0110       booker.book2D(histname, histtitle, ls_binning_.nbins, ls_binning_.xmin, ls_binning_.xmax, nbins, xmin, xmax);
0111   me->setAxisTitle("LS", 1);
0112   me->setAxisTitle("PS column index", 2);
0113   histograms.psColumnIndexVsLS = me;
0114 
0115   int bin = 1;
0116   for (auto const& l : labels) {
0117     histograms.psColumnIndexVsLS->setBinLabel(bin, l, 2);
0118     bin++;
0119   }
0120 }
0121 
0122 void PSMonitor::dqmAnalyze(edm::Event const& event, edm::EventSetup const& setup, Histograms const& histograms) const {
0123   int ls = event.id().luminosityBlock();
0124   int psColumn = -1;
0125 
0126   edm::Handle<GlobalAlgBlkBxCollection> ugtBXhandle;
0127   event.getByToken(ugtBXToken_, ugtBXhandle);
0128   if (ugtBXhandle.isValid() and not ugtBXhandle->isEmpty(0)) {
0129     psColumn = ugtBXhandle->at(0, 0).getPreScColumn();
0130   }
0131   histograms.psColumnIndexVsLS->Fill(ls, psColumn);
0132 }
0133 
0134 void PSMonitor::fillHistoPSetDescription(edm::ParameterSetDescription& pset, int value) {
0135   pset.add<int>("nbins", value);
0136 }
0137 
0138 void PSMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0139   edm::ParameterSetDescription desc;
0140   desc.add<edm::InputTag>("ugtBXInputTag", edm::InputTag("hltGtStage2Digis"));
0141   desc.add<std::string>("folderName", "HLT/PSMonitoring");
0142 
0143   edm::ParameterSetDescription histoPSet;
0144 
0145   edm::ParameterSetDescription psColumnPSet;
0146   fillHistoPSetDescription(psColumnPSet, 8);
0147   histoPSet.add("psColumnPSet", psColumnPSet);
0148 
0149   edm::ParameterSetDescription lsPSet;
0150   fillHistoPSetDescription(lsPSet, 2500);
0151   histoPSet.add("lsPSet", lsPSet);
0152 
0153   desc.add("histoPSet", histoPSet);
0154 
0155   descriptions.add("psMonitoring", desc);
0156 }
0157 
0158 // Define this as a plug-in
0159 #include "FWCore/Framework/interface/MakerMacros.h"
0160 DEFINE_FWK_MODULE(PSMonitor);