Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:14

0001 #include <cassert>
0002 #include <sstream>
0003 
0004 #include "DQM/TrigXMonitorClient/interface/L1ScalersClient.h"
0005 
0006 #include "FWCore/Framework/interface/LuminosityBlock.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ServiceRegistry/interface/Service.h"
0010 
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 
0013 using edm::LogInfo;
0014 using edm::LogWarning;
0015 
0016 #define SECS_PER_LUMI_SECTION 23.31040958083832;
0017 const int kPerHisto = 20;
0018 const int kNumAlgoHistos = MAX_ALGOS / kPerHisto;  // this hasta be w/o remainders
0019 const int kNumTTHistos = MAX_TT / kPerHisto;       // this hasta be w/o remainders
0020 
0021 /// Constructors
0022 L1ScalersClient::L1ScalersClient(const edm::ParameterSet &ps)
0023     : dbe_(nullptr),
0024       nLumi_(0),
0025       l1AlgoCurrentRate_(nullptr),
0026       l1TechTrigCurrentRate_(nullptr),
0027       selected_(nullptr),
0028       bxSelected_(nullptr),
0029       algoSelected_(ps.getUntrackedParameter<std::vector<int>>("algoMonitorBits", std::vector<int>())),
0030       techSelected_(ps.getUntrackedParameter<std::vector<int>>("techMonitorBits", std::vector<int>())),
0031       folderName_(ps.getUntrackedParameter<std::string>("dqmFolder", "L1T/L1Scalers_EvF")),
0032       currentLumiBlockNumber_(0),
0033       first_algo(true),
0034       first_tt(true) {
0035   LogDebug("Status") << "constructor";
0036   // get back-end interface
0037   usesResource("DQMStore");
0038   dbe_ = edm::Service<DQMStore>().operator->();
0039   assert(dbe_ != nullptr);  // blammo!
0040   dbe_->setCurrentFolder(folderName_);
0041 
0042   l1AlgoCurrentRate_ =
0043       dbe_->book1D("algo_cur_rate", "current lumi section rate per Algo Bits", MAX_ALGOS, -0.5, MAX_ALGOS - 0.5);
0044 
0045   l1TechTrigCurrentRate_ =
0046       dbe_->book1D("tt_cur_rate", "current lumi section rate per Tech. Trig.s", MAX_TT, -0.5, MAX_TT - 0.5);
0047   // ----------------------
0048   numSelected_ = algoSelected_.size() + techSelected_.size();
0049   selected_ = dbe_->book1D("l1BitsSel",
0050                            "Selected L1 Algorithm"
0051                            " and tech Bits",
0052                            numSelected_,
0053                            -0.5,
0054                            numSelected_ - 0.5);
0055   bxSelected_ = dbe_->book2D(
0056       "l1BitsBxSel", "Selected L1 Algorithm Bits vs Bx", 3600, -0.5, 3599.5, numSelected_, -0.5, numSelected_ - 0.5);
0057   int j = 1;
0058   for (unsigned int i = 0; i < algoSelected_.size(); ++i) {
0059     char title[256];
0060     snprintf(title, 256, "Algo %d", algoSelected_[i]);
0061     selected_->setBinLabel(j, title);
0062     bxSelected_->setBinLabel(j, title, 2);
0063     ++j;
0064   }
0065   for (unsigned int i = 0; i < techSelected_.size(); ++i) {
0066     char title[256];
0067     snprintf(title, 256, "Tech %d", techSelected_[i]);
0068     selected_->setBinLabel(j, title);
0069     bxSelected_->setBinLabel(j, title, 2);
0070     ++j;
0071   }
0072 
0073   // book individual bit rates vs lumi for algo bits.
0074   totalAlgoRate_ = dbe_->book1D("totAlgoRate", "Total Algo Rate", MAX_LUMI_SEG, -0.5, MAX_LUMI_SEG - 0.5);
0075   totalTtRate_ = dbe_->book1D("totTtRate", "Total Tech Rate", MAX_LUMI_SEG, -0.5, MAX_LUMI_SEG - 0.5);
0076 
0077   totAlgoPrevCount = 0UL;
0078   totTtPrevCount = 0UL;
0079 
0080   std::string algodir = "/AlgoRates";
0081   dbe_->setCurrentFolder(folderName_ + algodir);
0082 
0083   for (int i = 0; i < MAX_ALGOS; ++i) {
0084     l1AlgoScalerCounters_[i] = 0UL;
0085     l1AlgoRateHistories_[i] = nullptr;  // not really needed but ...
0086     char name[256];
0087     snprintf(name, 256, "rate_algobit%03d", i);
0088     LogDebug("Parameter") << "name " << i << " is " << name;
0089     l1AlgoRateHistories_[i] = dbe_->book1D(name, name, MAX_LUMI_SEG, -0.5, MAX_LUMI_SEG - 0.5);
0090   }
0091 
0092   // book individual bit rates vs lumi for technical trigger bits.
0093 
0094   std::string techdir = "/TechRates";
0095   dbe_->setCurrentFolder(folderName_ + techdir);
0096 
0097   for (int i = 0; i < MAX_TT; ++i) {
0098     l1TechTrigScalerCounters_[i] = 0UL;
0099     l1TechTrigRateHistories_[i] = nullptr;  // not really needed but ...
0100     char name[256];
0101     snprintf(name, 256, "rate_ttbit%03d", i);
0102     LogDebug("Parameter") << "name " << i << " is " << name;
0103     l1TechTrigRateHistories_[i] = dbe_->book1D(name, name, MAX_LUMI_SEG, -0.5, MAX_LUMI_SEG - 0.5);
0104   }
0105 
0106   // split l1 scalers up into groups of 20, assuming total of 140 bits
0107   std::string algodir2 = "/AlgoBits";
0108   dbe_->setCurrentFolder(folderName_ + algodir2);
0109 
0110   char metitle1[64];  // histo name
0111   char mename1[64];   // ME name
0112   for (int k = 0; k < kNumAlgoHistos; k++) {
0113     int npath_low = kPerHisto * k;
0114     int npath_high = kPerHisto * (k + 1) - 1;
0115     snprintf(mename1, 64, "L1AlgoBits_%0d", k);
0116     snprintf(metitle1, 64, "L1 rates - Algo Bits %d to %d", npath_low, npath_high);
0117     l1AlgoCurrentRatePerAlgo_[k] = dbe_->book1D(mename1, metitle1, kPerHisto, -0.5 + npath_low, npath_high + 0.5);
0118   }
0119 
0120   // split l1 scalers up into groups of 20, assuming total of 80 technical bits
0121 
0122   std::string techdir2 = "/TechBits";
0123   dbe_->setCurrentFolder(folderName_ + techdir2);
0124 
0125   char metitle2[64];  // histo name
0126   char mename2[64];   // ME name
0127   for (int k = 0; k < kNumTTHistos; k++) {
0128     int npath_low = kPerHisto * k;
0129     int npath_high = kPerHisto * (k + 1) - 1;
0130     snprintf(mename2, 64, "L1TechBits_%0d", k);
0131     snprintf(metitle2, 64, "L1 rates - Tech. Trig. Bits %d to %d", npath_low, npath_high);
0132     l1TechTrigCurrentRatePerAlgo_[k] = dbe_->book1D(mename2, metitle2, kPerHisto, -0.5 + npath_low, npath_high + 0.5);
0133   }
0134 
0135   std::ostringstream params;
0136   params << "Algo: ";
0137   for (unsigned int i = 0; i < algoSelected_.size(); ++i) {
0138     params << algoSelected_[i] << " ";
0139   }
0140   params << ", Tech: ";
0141   for (unsigned int i = 0; i < techSelected_.size(); ++i) {
0142     params << techSelected_[i] << " ";
0143   }
0144   LogDebug("Parameter") << "L1 bits to monitor are " << params.str();
0145 }
0146 
0147 /// BeginJob
0148 void L1ScalersClient::beginJob(void) {
0149   LogDebug("Status") << "beingJob";
0150   if (dbe_) {
0151     dbe_->setCurrentFolder(folderName_);
0152   }
0153 }
0154 
0155 /// BeginRun
0156 void L1ScalersClient::beginRun(const edm::Run &run, const edm::EventSetup &c) {}
0157 
0158 /// EndRun
0159 void L1ScalersClient::endRun(const edm::Run &run, const edm::EventSetup &c) {}
0160 
0161 /// End LumiBlock
0162 /// DQM Client Diagnostic should be performed here
0163 void L1ScalersClient::endLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c) {
0164   nLumi_ = lumiSeg.id().luminosityBlock();
0165 
0166   // get EvF data
0167 
0168   MonitorElement *algoScalers = dbe_->get(folderName_ + std::string("/l1AlgoBits"));
0169   MonitorElement *ttScalers = dbe_->get(folderName_ + std::string("/l1TechBits"));
0170 
0171   if (algoScalers == nullptr || ttScalers == nullptr) {
0172     LogInfo("Status") << "cannot get l1 scalers histogram, bailing out.";
0173     return;
0174   }
0175 
0176   int nalgobits = algoScalers->getNbinsX();
0177   int nttbits = ttScalers->getNbinsX();
0178 
0179   if (nalgobits > MAX_ALGOS)
0180     nalgobits = MAX_ALGOS;  // HARD CODE FOR NOW
0181   if (nttbits > MAX_TT)
0182     nttbits = MAX_TT;  // HARD CODE FOR NOW
0183 
0184   LogDebug("Status") << "I see " << nalgobits << " algo paths. ";
0185   LogDebug("Status") << "I see " << nttbits << " tt paths. ";
0186 
0187   // set the bin labels on the first go-through
0188   if (first_algo) {
0189     for (int i = 0; i < nalgobits; ++i) {
0190       int whichHisto = i / kPerHisto;
0191       int whichBin = i % kPerHisto + 1;
0192       char pname[256];
0193       snprintf(pname, 256, "AlgoBit%03d", i);
0194       l1AlgoCurrentRatePerAlgo_[whichHisto]->setBinLabel(whichBin, pname);
0195       snprintf(pname, 256, "Rate - Algorithm Bit %03d", i);
0196       l1AlgoRateHistories_[i]->setTitle(pname);
0197     }
0198     first_algo = false;
0199   }
0200 
0201   // set the bin labels on the first go-through
0202   if (first_tt) {
0203     for (int i = 0; i < nttbits; ++i) {
0204       int whichHisto = i / kPerHisto;
0205       int whichBin = i % kPerHisto + 1;
0206       char pname[256];
0207       snprintf(pname, 256, "TechBit%03d", i);
0208       l1TechTrigCurrentRatePerAlgo_[whichHisto]->setBinLabel(whichBin, pname);
0209       snprintf(pname, 256, "Rate - Technical Bit %03d", i);
0210       l1TechTrigRateHistories_[i]->setTitle(pname);
0211     }
0212     first_tt = false;
0213   }
0214 
0215   MonitorElement *nLumi = dbe_->get(folderName_ + std::string("nLumiBlock"));
0216 
0217   int testval = (nLumi != nullptr ? nLumi->getIntValue() : -1);
0218   LogDebug("Parameter") << "Lumi Block from DQM: " << testval << ", local is " << nLumi_;
0219 
0220   int nL = (nLumi != nullptr ? nLumi->getIntValue() : nLumi_);
0221   if (nL > MAX_LUMI_SEG) {
0222     LogDebug("Status") << "Too many Lumi segments, " << nL << " is greater than MAX_LUMI_SEG,"
0223                        << " wrapping to " << (nL % MAX_LUMI_SEG);
0224     nL = nL % MAX_LUMI_SEG;
0225   }
0226   float delta_t = (nL - currentLumiBlockNumber_) * SECS_PER_LUMI_SECTION;
0227   if (delta_t < 0) {
0228     LogDebug("Status") << " time is negative ... " << delta_t;
0229     delta_t = -delta_t;
0230   } else if (nL == currentLumiBlockNumber_) {  // divide-by-zero
0231     LogInfo("Status") << "divide by zero: same lumi section 2x " << nL;
0232     return;
0233   }
0234   // selected ---------------------  fill in the rates for th
0235   int currSlot = 1;  // for selected bits histogram
0236   MonitorElement *algoBx = dbe_->get(folderName_ + std::string("/l1AlgoBits_Vs_Bx"));
0237   // selected ---------------------  end
0238   for (int i = 1; i <= nalgobits; ++i) {  // bins start at 1
0239     float current_count = algoScalers->getBinContent(i);
0240     // selected -------------------- start
0241     int bit = i - 1;  //
0242     if (std::find(algoSelected_.begin(), algoSelected_.end(), bit) != algoSelected_.end()) {
0243       selected_->setBinContent(currSlot, current_count);
0244       if (algoBx) {
0245         for (int j = 1; j <= 3600; ++j) {
0246           bxSelected_->setBinContent(j, currSlot, algoBx->getBinContent(j, i));
0247         }
0248       }
0249       ++currSlot;
0250     }
0251     // selected -------------------- end
0252     float rate = (current_count - l1AlgoScalerCounters_[i - 1]) / delta_t;
0253     if (rate > 1E-3) {
0254       LogDebug("Parameter") << "rate path " << i << " is " << rate;
0255     }
0256     l1AlgoCurrentRate_->setBinContent(i, rate);
0257     l1AlgoCurrentRatePerAlgo_[i / kPerHisto]->setBinContent(i % kPerHisto, rate);
0258     // currentRate_->setBinError(i, error);
0259     l1AlgoScalerCounters_[i - 1] = (unsigned long)(current_count);
0260     l1AlgoRateHistories_[i - 1]->setBinContent(nL, rate);
0261   }
0262   // selected ----------------- start
0263   MonitorElement *techBx = dbe_->get(folderName_ + std::string("/l1TechBits_Vs_Bx"));
0264   // selected ----------------- end
0265 
0266   for (int i = 1; i <= nttbits; ++i) {  // bins start at 1
0267     float current_count = ttScalers->getBinContent(i);
0268     // selected -------------------- start
0269     int bit = i - 1;  //
0270     if (std::find(techSelected_.begin(), techSelected_.end(), bit) != techSelected_.end()) {
0271       selected_->setBinContent(currSlot, current_count);
0272       if (techBx) {
0273         for (int j = 1; j <= 3600; ++j) {
0274           bxSelected_->setBinContent(j, currSlot, techBx->getBinContent(j, i));
0275         }
0276       }
0277       ++currSlot;
0278     }
0279     // selected -------------------- end
0280     float rate = (current_count - l1TechTrigScalerCounters_[i - 1]) / delta_t;
0281     if (rate > 1E-3) {
0282       LogDebug("Parameter") << "rate path " << i << " is " << rate;
0283     }
0284     l1TechTrigCurrentRate_->setBinContent(i, rate);
0285     l1TechTrigCurrentRatePerAlgo_[i / kPerHisto]->setBinContent(i % kPerHisto, rate);
0286     // currentRate_->setBinError(i, error);
0287     l1TechTrigScalerCounters_[i - 1] = (unsigned long)(current_count);
0288     l1TechTrigRateHistories_[i - 1]->setBinContent(nL, rate);
0289   }
0290 
0291   //  compute total rate
0292   MonitorElement *l1AlgoCounter = dbe_->get(folderName_ + std::string("/l1AlgoCounter"));
0293   MonitorElement *l1TtCounter = dbe_->get(folderName_ + std::string("/l1TtCounter"));
0294   if (l1AlgoCounter != nullptr && l1TtCounter != nullptr) {
0295     float totAlgoCount = l1AlgoCounter->getIntValue();
0296     float totTtCount = l1TtCounter->getIntValue();
0297     float totAlgRate = (totAlgoCount - totAlgoPrevCount) / delta_t;
0298     float totTtRate = (totTtCount - totTtPrevCount) / delta_t;
0299     totalAlgoRate_->setBinContent(nL, totAlgRate);
0300     totAlgoPrevCount = totAlgoCount;
0301     totalTtRate_->setBinContent(nL, totTtRate);
0302     totTtPrevCount = totTtCount;
0303   }
0304 
0305   currentLumiBlockNumber_ = nL;
0306 }
0307 
0308 // unused
0309 void L1ScalersClient::analyze(const edm::Event &e, const edm::EventSetup &c) {}