File indexing completed on 2024-04-06 12:10:07
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0015 #include "DQMServices/Core/interface/DQMStore.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
0018 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
0019 #include "DataFormats/L1GlobalTrigger/interface/L1GtFdlWord.h"
0020 #include "DataFormats/Scalers/interface/DcsStatus.h"
0021 #include "DataFormats/Scalers/interface/DcsStatus.h"
0022 #include "DataFormats/Scalers/interface/Level1TriggerRates.h"
0023 #include "DataFormats/Scalers/interface/Level1TriggerRates.h"
0024 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
0025 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
0026 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
0027 #include "DataFormats/Scalers/interface/LumiScalers.h"
0028 #include "DataFormats/Scalers/interface/LumiScalers.h"
0029 #include "DataFormats/Scalers/interface/ScalersRaw.h"
0030 #include "DataFormats/Scalers/interface/ScalersRaw.h"
0031 #include "DataFormats/Scalers/interface/TimeSpec.h"
0032 #include "DataFormats/Scalers/interface/TimeSpec.h"
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 #include "FWCore/Framework/interface/Run.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 #include "FWCore/ParameterSet/interface/Registry.h"
0038 #include "FWCore/ServiceRegistry/interface/Service.h"
0039
0040 class DQMScalInfo : public DQMEDAnalyzer {
0041 public:
0042
0043 DQMScalInfo(const edm::ParameterSet& ps);
0044
0045
0046 ~DQMScalInfo() override = default;
0047
0048 protected:
0049
0050 void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0051 void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0052
0053 private:
0054 void makeL1Scalars(const edm::Event& e);
0055 void makeLumiScalars(const edm::Event& e);
0056
0057 edm::ParameterSet parameters_;
0058 std::string scalfolder_;
0059 edm::EDGetTokenT<L1GlobalTriggerReadoutRecord> gtCollection_;
0060 edm::EDGetTokenT<DcsStatusCollection> dcsStatusCollection_;
0061 edm::EDGetTokenT<Level1TriggerScalersCollection> l1tscollectionToken_;
0062 edm::EDGetTokenT<LumiScalersCollection> lumicollectionToken_;
0063
0064
0065 MonitorElement* hlresync_;
0066 MonitorElement* hlOC0_;
0067 MonitorElement* hlTE_;
0068 MonitorElement* hlstart_;
0069 MonitorElement* hlEC0_;
0070 MonitorElement* hlHR_;
0071 MonitorElement* hphysTrig_;
0072
0073 MonitorElement* hinstLumi_;
0074 };
0075
0076 using namespace std;
0077 using namespace edm;
0078
0079
0080
0081 DQMScalInfo::DQMScalInfo(const edm::ParameterSet& ps) {
0082 parameters_ = ps;
0083
0084 scalfolder_ = parameters_.getUntrackedParameter<std::string>("dqmScalFolder", "Scal");
0085 gtCollection_ = consumes<L1GlobalTriggerReadoutRecord>(
0086 parameters_.getUntrackedParameter<edm::InputTag>("gtCollection", edm::InputTag("gtDigis")));
0087 dcsStatusCollection_ = consumes<DcsStatusCollection>(
0088 parameters_.getUntrackedParameter<edm::InputTag>("dcsStatusCollection", edm::InputTag("scalersRawToDigi")));
0089 l1tscollectionToken_ = consumes<Level1TriggerScalersCollection>(
0090 parameters_.getUntrackedParameter<edm::InputTag>("l1TSCollection", edm::InputTag("scalersRawToDigi")));
0091 lumicollectionToken_ = consumes<LumiScalersCollection>(
0092 parameters_.getUntrackedParameter<edm::InputTag>("lumiCollection", edm::InputTag("scalersRawToDigi")));
0093 }
0094
0095 void DQMScalInfo::bookHistograms(DQMStore::IBooker& ibooker,
0096 edm::Run const& ,
0097 edm::EventSetup const& ) {
0098 const int maxNbins = 2001;
0099
0100
0101 ibooker.cd();
0102 ibooker.setCurrentFolder(scalfolder_ + "/L1TriggerScalers/");
0103 const int fracLS = 16;
0104 const int maxLS = 250;
0105 hlresync_ = ibooker.book1D("lresync", "Orbit of last resync", fracLS * maxLS, 0, maxLS * 262144);
0106 hlOC0_ = ibooker.book1D("lOC0", "Orbit of last OC0", fracLS * maxLS, 0, maxLS * 262144);
0107 hlTE_ = ibooker.book1D("lTE", "Orbit of last TestEnable", fracLS * maxLS, 0, maxLS * 262144);
0108 hlstart_ = ibooker.book1D("lstart", "Orbit of last Start", fracLS * maxLS, 0, maxLS * 262144);
0109 hlEC0_ = ibooker.book1D("lEC0", "Orbit of last EC0", fracLS * maxLS, 0, maxLS * 262144);
0110 hlHR_ = ibooker.book1D("lHR", "Orbit of last HardReset", fracLS * maxLS, 0, maxLS * 262144);
0111
0112 hphysTrig_ = ibooker.book1D("Physics_Triggers", "Physics Triggers", maxNbins, -0.5, double(maxNbins) - 0.5);
0113 hphysTrig_->setAxisTitle("Lumi Section", 1);
0114
0115 ibooker.cd();
0116 ibooker.setCurrentFolder(scalfolder_ + "/LumiScalers/");
0117 hinstLumi_ = ibooker.book1D("Instant_Lumi", "Instant Lumi", maxNbins, -0.5, double(maxNbins) - 0.5);
0118 }
0119
0120 void DQMScalInfo::analyze(const edm::Event& e, const edm::EventSetup& c) {
0121 makeL1Scalars(e);
0122 makeLumiScalars(e);
0123 return;
0124 }
0125
0126 void DQMScalInfo::makeL1Scalars(const edm::Event& e) {
0127 edm::Handle<Level1TriggerScalersCollection> l1ts;
0128 e.getByToken(l1tscollectionToken_, l1ts);
0129 edm::Handle<LumiScalersCollection> lumiScalers;
0130 e.getByToken(lumicollectionToken_, lumiScalers);
0131
0132 auto it = l1ts->begin();
0133
0134 if (l1ts->empty())
0135 return;
0136 hlresync_->Fill((*l1ts)[0].lastResync());
0137 hlOC0_->Fill((*l1ts)[0].lastOrbitCounter0());
0138 hlTE_->Fill((*l1ts)[0].lastTestEnable());
0139 hlstart_->Fill((*l1ts)[0].lastStart());
0140 hlEC0_->Fill((*l1ts)[0].lastEventCounter0());
0141 hlHR_->Fill((*l1ts)[0].lastHardReset());
0142
0143 unsigned int lumisection = it->lumiSegmentNr();
0144 if (lumisection) {
0145 hphysTrig_->setBinContent(lumisection + 1, it->l1AsPhysics());
0146 }
0147
0148 return;
0149 }
0150
0151 void DQMScalInfo::makeLumiScalars(const edm::Event& e) {
0152 edm::Handle<LumiScalersCollection> lumiScalers;
0153 e.getByToken(lumicollectionToken_, lumiScalers);
0154
0155 auto it = lumiScalers->begin();
0156
0157 if (!lumiScalers->empty()) {
0158 unsigned int lumisection = it->sectionNumber();
0159 if (lumisection) {
0160 hinstLumi_->setBinContent(lumisection + 1, it->instantLumi());
0161 }
0162 }
0163
0164 return;
0165 }
0166
0167 #include "FWCore/PluginManager/interface/ModuleDef.h"
0168 #include "FWCore/Framework/interface/MakerMacros.h"
0169 DEFINE_FWK_MODULE(DQMScalInfo);