File indexing completed on 2024-09-12 04:16:24
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018 #include "DQM/SiPixelPhase1Summary/interface/SiPixelPhase1Summary.h"
0019
0020 #include "FWCore/ServiceRegistry/interface/Service.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023
0024 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
0025 #include "DQMServices/Core/interface/DQMStore.h"
0026
0027 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0028 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0029 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0030 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0031 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0032
0033 #include "DataFormats/DetId/interface/DetId.h"
0034 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0035 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0036 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0037 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0038 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0039
0040 #include <string>
0041 #include <cstdlib>
0042 #include <iostream>
0043 #include <fstream>
0044 #include <sstream>
0045
0046 using namespace std;
0047 using namespace edm;
0048
0049 SiPixelPhase1Summary::SiPixelPhase1Summary(const edm::ParameterSet& iConfig)
0050 : DQMEDHarvester(iConfig), conf_(iConfig), firstLumi(true) {
0051 LogInfo("PixelDQM") << "SiPixelPhase1Summary::SiPixelPhase1Summary: Got DQM BackEnd interface" << endl;
0052 topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0053 runOnEndLumi_ = conf_.getParameter<bool>("RunOnEndLumi");
0054 runOnEndJob_ = conf_.getParameter<bool>("RunOnEndJob");
0055
0056 std::vector<edm::ParameterSet> mapPSets = conf_.getParameter<std::vector<edm::ParameterSet> >("SummaryMaps");
0057
0058
0059 for (auto const& mapPSet : mapPSets) {
0060 summaryPlotName_[mapPSet.getParameter<std::string>("MapName")] = mapPSet.getParameter<std::string>("MapHist");
0061 }
0062 deadRocThresholds_ = conf_.getParameter<std::vector<double> >("DeadROCErrorThreshold");
0063 deadRocWarnThresholds_ = conf_.getParameter<std::vector<double> >("DeadROCWarningThreshold");
0064 }
0065
0066 SiPixelPhase1Summary::~SiPixelPhase1Summary() {
0067
0068
0069 LogInfo("PixelDQM") << "SiPixelPhase1Summary::~SiPixelPhase1Summary: Destructor" << endl;
0070 }
0071
0072 void SiPixelPhase1Summary::beginRun(edm::Run const& run, edm::EventSetup const& eSetup) {}
0073
0074 void SiPixelPhase1Summary::dqmEndLuminosityBlock(DQMStore::IBooker& iBooker,
0075 DQMStore::IGetter& iGetter,
0076 const edm::LuminosityBlock& lumiSeg,
0077 edm::EventSetup const& c) {
0078 if (firstLumi) {
0079 bookSummaries(iBooker);
0080 bookTrendPlots(iBooker);
0081 firstLumi = false;
0082 }
0083
0084 if (runOnEndLumi_) {
0085 int lumiSec = lumiSeg.id().luminosityBlock();
0086 fillTrendPlots(iBooker, iGetter, lumiSec);
0087 fillSummaries(iBooker, iGetter, lumiSec);
0088 }
0089
0090
0091 }
0092
0093
0094
0095
0096 void SiPixelPhase1Summary::dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) {
0097 if (firstLumi) {
0098 bookSummaries(iBooker);
0099 bookTrendPlots(iBooker);
0100 firstLumi = false;
0101 }
0102 if (runOnEndJob_) {
0103 if (!runOnEndLumi_)
0104 fillTrendPlots(
0105 iBooker,
0106 iGetter);
0107 fillSummaries(iBooker, iGetter);
0108 }
0109 }
0110
0111
0112
0113
0114 void SiPixelPhase1Summary::bookSummaries(DQMStore::IBooker& iBooker) {
0115 iBooker.cd();
0116
0117 std::vector<std::string> xAxisLabels_ = {"BMO",
0118 "BMI",
0119 "BPO ",
0120 "BPI",
0121 "HCMO_1",
0122 "HCMO_2",
0123 "HCMI_1",
0124 "HCMI_2",
0125 "HCPO_1",
0126 "HCPO_2",
0127 "HCPI_1",
0128 "HCPI_2"};
0129 std::vector<std::string> yAxisLabels_ = {
0130 "1",
0131 "2",
0132 "3",
0133 "4"};
0134
0135 iBooker.setCurrentFolder("PixelPhase1/Summary");
0136
0137 for (const auto& mapInfo : summaryPlotName_) {
0138 auto name = mapInfo.first;
0139 summaryMap_[name] = iBooker.book2D("pixel" + name + "Summary", "Pixel " + name + " Summary", 12, 0, 12, 4, 0, 4);
0140 }
0141
0142 deadROCSummary =
0143 iBooker.book2D("deadROCSummary", "Percentage of filled ROCs (with digis) per layer/ring", 2, 0, 2, 4, 0, 4);
0144 std::vector<std::string> xAxisLabelsReduced_ = {"Barrel", "Forward"};
0145 deadROCSummary->setAxisTitle("Subdetector", 1);
0146 for (unsigned int i = 0; i < xAxisLabelsReduced_.size(); i++) {
0147 deadROCSummary->setBinLabel(i + 1, xAxisLabelsReduced_[i]);
0148 }
0149
0150
0151 iBooker.setCurrentFolder("PixelPhase1/EventInfo");
0152
0153 if (runOnEndLumi_) {
0154
0155 summaryMap_["Grand"] = iBooker.book2D("reportSummaryMap", "Pixel Summary Map", 2, 0, 2, 4, 0, 4);
0156 summaryMap_["Grand"]->setAxisTitle("Subdetector", 1);
0157 for (unsigned int i = 0; i < xAxisLabelsReduced_.size(); i++) {
0158 summaryMap_["Grand"]->setBinLabel(i + 1, xAxisLabelsReduced_[i]);
0159 for (unsigned int j = 0; j < 4; j++) {
0160 summaryMap_["Grand"]->Fill(i, j, -1);
0161 }
0162 }
0163 } else {
0164
0165 summaryMap_["Grand"] = iBooker.book2D("reportSummaryMap", "Pixel Summary Map", 12, 0, 12, 4, 0, 4);
0166 }
0167
0168 reportSummary = iBooker.bookFloat("reportSummary");
0169
0170
0171 for (const auto& summaryMapEntry : summaryMap_) {
0172 if (summaryMapEntry.first == "Grand")
0173 continue;
0174 auto summaryMap = summaryMapEntry.second;
0175 for (unsigned int i = 0; i < xAxisLabels_.size(); i++) {
0176 summaryMap->setBinLabel(i + 1, xAxisLabels_[i], 1);
0177 }
0178 for (unsigned int i = 0; i < yAxisLabels_.size(); i++) {
0179 summaryMap->setBinLabel(i + 1, yAxisLabels_[i], 2);
0180 }
0181 summaryMap->setAxisTitle("Subdetector", 1);
0182 summaryMap->setAxisTitle("Layer/disk", 2);
0183 for (int i = 0; i < summaryMap->getTH1()->GetXaxis()->GetNbins(); i++) {
0184 for (int j = 0; j < summaryMap->getTH1()->GetYaxis()->GetNbins(); j++) {
0185 summaryMap->Fill(i, j, -1.);
0186 }
0187 }
0188 }
0189 reportSummary->Fill(-1.);
0190
0191
0192 iBooker.setCurrentFolder("PixelPhase1/");
0193 }
0194
0195
0196
0197
0198 void SiPixelPhase1Summary::bookTrendPlots(DQMStore::IBooker& iBooker) {
0199
0200 iBooker.setCurrentFolder("PixelPhase1/");
0201 std::vector<string> binAxisLabels = {"Layer 1", "Layer 2", "Layer 3", "Layer 4", "Ring 1", "Ring 2"};
0202 if (runOnEndLumi_) {
0203 std::vector<trendPlots> histoOrder = {layer1, layer2, layer3, layer4, ring1, ring2};
0204 std::vector<string> varName = {"Layer_1", "Layer_2", "Layer_3", "Layer_4", "Ring_1", "Ring_2"};
0205 std::vector<int> yMax = {1536, 3584, 5632, 8192, 4224, 6528};
0206 for (unsigned int i = 0; i < histoOrder.size(); i++) {
0207 string varNameStr = "deadRocTrend" + varName[i];
0208 string varTitle = binAxisLabels[i] + " dead ROC trend";
0209 deadROCTrends_[histoOrder[i]] = iBooker.bookProfile(varNameStr, varTitle, 500, 0., 5000, 0., yMax[i], "");
0210 varNameStr = "ineffRocTrend" + varName[i];
0211 varTitle = binAxisLabels[i] + " inefficient ROC trend";
0212 ineffROCTrends_[histoOrder[i]] = iBooker.bookProfile(varNameStr, varTitle, 500, 0., 5000, 0., yMax[i], "");
0213 deadROCTrends_[histoOrder[i]]->setAxisTitle("Lumisection", 1);
0214 ineffROCTrends_[histoOrder[i]]->setAxisTitle("Lumisection", 1);
0215 }
0216 } else {
0217 deadROCTrends_[offline] = iBooker.bookProfile("deadRocTotal", "N dead ROCs", 6, 0., 6, 0., 8192, "");
0218 ineffROCTrends_[offline] = iBooker.bookProfile("ineffRocTotal", "N inefficient ROCs", 6, 0., 6, 0., 8192, "");
0219 deadROCTrends_[offline]->setAxisTitle("Subdetector", 1);
0220 ineffROCTrends_[offline]->setAxisTitle("Subdetector", 1);
0221 for (unsigned int i = 1; i <= binAxisLabels.size(); i++) {
0222 deadROCTrends_[offline]->setBinLabel(i, binAxisLabels[i - 1]);
0223 ineffROCTrends_[offline]->setBinLabel(i, binAxisLabels[i - 1]);
0224 }
0225 }
0226 }
0227
0228
0229
0230
0231 void SiPixelPhase1Summary::fillSummaries(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter, int lumiSec) {
0232
0233 for (const auto& mapInfo : summaryPlotName_) {
0234 auto name = mapInfo.first;
0235 std::ostringstream histNameStream;
0236 std::string histName;
0237
0238 for (int i = 0; i < 12; i++) {
0239 for (int j = 0; j < 4; j++) {
0240 if (i > 3 && j == 3)
0241 continue;
0242 bool minus = i < 2 || (i > 3 && i < 8);
0243 int iOver2 = floor(i / 2.);
0244 bool outer = (i > 3) ? iOver2 % 2 == 0 : i % 2 == 0;
0245
0246 histNameStream.str("");
0247 histNameStream << topFolderName_.c_str() << "PX" << ((i > 3) ? "Forward" : "Barrel") << "/"
0248 << ((i > 3) ? "HalfCylinder" : "Shell") << "_" << (minus ? "m" : "p") << ((outer) ? "O" : "I")
0249 << "/" << ((i > 3) ? ((i % 2 == 0) ? "PXRing_1/" : "PXRing_2/") : "")
0250 << summaryPlotName_[name].c_str() << "_PX" << ((i > 3) ? "Disk" : "Layer") << "_"
0251 << ((i > 3) ? ((minus) ? "-" : "+") : "") << (j + 1);
0252 histName = histNameStream.str();
0253 MonitorElement* me = iGetter.get(histName);
0254
0255 if (!me) {
0256 edm::LogWarning("SiPixelPhase1Summary") << "ME " << histName << " is not available !!";
0257 continue;
0258 }
0259
0260 if (summaryMap_[name] == nullptr) {
0261 edm::LogWarning("SiPixelPhase1Summary") << "Summary map " << name << " is not available !!";
0262 continue;
0263 }
0264 if (!(me->getQReports()).empty())
0265 summaryMap_[name]->setBinContent(i + 1, j + 1, (me->getQReports())[0]->getQTresult());
0266 else
0267 summaryMap_[name]->setBinContent(i + 1, j + 1, -1);
0268 }
0269 }
0270 }
0271
0272
0273 std::vector<trendPlots> trendOrder = {layer1, layer2, layer3, layer4, ring1, ring2};
0274 std::vector<int> nRocsPerTrend = {1536, 3584, 5632, 8192, 4224, 6528};
0275 std::vector<int> nDisabledRocs = {12, 128, 240, 320, 96, 120};
0276 for (unsigned int i = 0; i < trendOrder.size(); i++) {
0277 int xBin = i < 4 ? 1 : 2;
0278 int yBin = i % 4 + 1;
0279 float numDeadROCs = 0.;
0280 float numFilledROCs = 0.;
0281 float fracFilledROCs = 1.;
0282 unsigned int lastTrendBin = 0;
0283 int numEmptyBins = 0;
0284 float sumOf5Bins = 0.;
0285 std::vector<float> last5TrendBinsValue = {0., 0., 0., 0., 0.};
0286 if (runOnEndLumi_) {
0287 TH1* tempProfile = deadROCTrends_[trendOrder[i]]->getTH1();
0288
0289
0290
0291
0292
0293 lastTrendBin = lumiSec / 10;
0294 if (lastTrendBin >= 5) {
0295
0296 for (unsigned int n = 0; n < 5; n++)
0297 last5TrendBinsValue[n] = tempProfile->GetBinContent(lastTrendBin - 5 + n + 1);
0298 } else {
0299
0300 for (unsigned int n = 0; n < lastTrendBin; n++)
0301 last5TrendBinsValue[n] = tempProfile->GetBinContent(n + 1);
0302 }
0303
0304 for (unsigned int n = 0; n < 5; n++)
0305 if (last5TrendBinsValue[n] == 0)
0306 numEmptyBins++;
0307 if (numEmptyBins == 0) {
0308
0309 std::sort(last5TrendBinsValue.begin(), last5TrendBinsValue.end());
0310 numDeadROCs = last5TrendBinsValue[2];
0311 } else if (numEmptyBins != 5) {
0312
0313 for (unsigned int n = 0; n < 5; n++)
0314 sumOf5Bins += last5TrendBinsValue[n];
0315 numDeadROCs = sumOf5Bins / (5 - numEmptyBins);
0316 }
0317
0318 } else {
0319 TH1* tempProfile = deadROCTrends_[offline]->getTH1();
0320 numDeadROCs = tempProfile->GetBinContent(i + 1);
0321 }
0322
0323 numFilledROCs = nRocsPerTrend[i] - numDeadROCs;
0324
0325 fracFilledROCs = numFilledROCs / (nRocsPerTrend[i] - nDisabledRocs[i]);
0326 if (fracFilledROCs > 1)
0327 fracFilledROCs = 1;
0328 deadROCSummary->setBinContent(xBin, yBin, fracFilledROCs);
0329 deadROCSummary->setBinContent(2, 3, -1);
0330 deadROCSummary->setBinContent(2, 4, -1);
0331 }
0332
0333
0334 float sumOfNonNegBins = 0.;
0335
0336
0337 if (!runOnEndLumi_) {
0338 for (int i = 0; i < 12; i++) {
0339 if (summaryMap_["Grand"] == nullptr) {
0340 edm::LogWarning("SiPixelPhase1Summary") << "Grand summary does not exist!";
0341 break;
0342 }
0343 for (int j = 0; j < 4; j++) {
0344 summaryMap_["Grand"]->setBinContent(
0345 i + 1,
0346 j + 1,
0347 1);
0348 for (auto const& mapInfo : summaryPlotName_) {
0349 auto name = mapInfo.first;
0350 if (summaryMap_[name] == nullptr) {
0351 edm::LogWarning("SiPixelPhase1Summary") << "Summary " << name << " does not exist!";
0352 continue;
0353 }
0354 if (summaryMap_["Grand"]->getBinContent(i + 1, j + 1) > summaryMap_[name]->getBinContent(i + 1, j + 1))
0355 summaryMap_["Grand"]->setBinContent(i + 1, j + 1, summaryMap_[name]->getBinContent(i + 1, j + 1));
0356 }
0357 if (summaryMap_["Grand"]->getBinContent(i + 1, j + 1) > -0.1)
0358 sumOfNonNegBins += summaryMap_["Grand"]->getBinContent(i + 1, j + 1);
0359 }
0360 }
0361 reportSummary->Fill(sumOfNonNegBins / 40.);
0362 }
0363
0364
0365
0366 else {
0367 for (int i = 0; i < 2; i++) {
0368 if (summaryMap_["Grand"] == nullptr) {
0369 edm::LogWarning("SiPixelPhase1Summary") << "Grand summary does not exist!";
0370 break;
0371 }
0372 for (int j = 0; j < 4; j++) {
0373
0374 if (i == 1 && j > 1) {
0375 summaryMap_["Grand"]->setBinContent(i + 1, j + 1, -1);
0376 } else {
0377
0378 summaryMap_["Grand"]->setBinContent(i + 1, j + 1, deadROCSummary->getBinContent(i + 1, j + 1));
0379
0380 sumOfNonNegBins += summaryMap_["Grand"]->getBinContent(i + 1, j + 1);
0381 }
0382 }
0383 }
0384 reportSummary->Fill(sumOfNonNegBins / 6.);
0385 }
0386 }
0387
0388
0389
0390 void SiPixelPhase1Summary::fillTrendPlots(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter, int lumiSec) {
0391
0392 if (lumiSec % 10 != 0)
0393 return;
0394
0395 if (runOnEndLumi_) {
0396 MonitorElement* nClustersAll = iGetter.get("PixelPhase1/Phase1_MechanicalView/num_clusters_per_Lumisection_PXAll");
0397 if (nClustersAll == nullptr) {
0398 edm::LogWarning("SiPixelPhase1Summary") << "All pixel cluster trend plot not available!!";
0399 return;
0400 }
0401 if (nClustersAll->getTH1()->GetBinContent(lumiSec) < 100)
0402 return;
0403 }
0404
0405 std::string histName;
0406
0407
0408 std::vector<trendPlots> trendOrder = {layer1, layer2, layer3, layer4, ring1, ring2};
0409 std::vector<int> nFilledROCs(trendOrder.size(), 0);
0410 std::vector<int> hiEffROCs(trendOrder.size(), 0);
0411 std::vector<int> nRocsPerTrend = {1536, 3584, 5632, 8192, 4224, 6528};
0412 std::vector<string> trendNames = {};
0413
0414 for (auto it : {1, 2, 3, 4}) {
0415 histName = "PXBarrel/digi_occupancy_per_SignedModuleCoord_per_SignedLadderCoord_PXLayer_" + std::to_string(it);
0416 trendNames.push_back(histName);
0417 }
0418 for (auto it : {1, 2}) {
0419 histName = "PXForward/digi_occupancy_per_SignedDiskCoord_per_SignedBladePanelCoord_PXRing_" + std::to_string(it);
0420 trendNames.push_back(histName);
0421 }
0422
0423 for (unsigned int trendIt = 0; trendIt < trendOrder.size(); trendIt++) {
0424 iGetter.cd();
0425 histName = "PixelPhase1/Phase1_MechanicalView/" + trendNames[trendIt];
0426 MonitorElement* tempLayerME = iGetter.get(histName);
0427 if (tempLayerME == nullptr)
0428 continue;
0429 float lowEffValue = 0.25 * (tempLayerME->getTH1()->Integral() / nRocsPerTrend[trendIt]);
0430 for (int i = 1; i <= tempLayerME->getTH1()->GetXaxis()->GetNbins(); i++) {
0431 for (int j = 1; j <= tempLayerME->getTH1()->GetYaxis()->GetNbins(); j++) {
0432 if (tempLayerME->getBinContent(i, j) > 0.)
0433 nFilledROCs[trendIt]++;
0434 if (tempLayerME->getBinContent(i, j) > lowEffValue)
0435 hiEffROCs[trendIt]++;
0436 }
0437 }
0438 if (runOnEndLumi_) {
0439 tempLayerME->Reset();
0440 }
0441 }
0442
0443 if (!runOnEndLumi_) {
0444 for (unsigned int i = 0; i < trendOrder.size(); i++) {
0445 deadROCTrends_[offline]->Fill(i, nRocsPerTrend[i] - nFilledROCs[i]);
0446 ineffROCTrends_[offline]->Fill(i, nFilledROCs[i] - hiEffROCs[i]);
0447 }
0448 } else {
0449 for (unsigned int i = 0; i < trendOrder.size(); i++) {
0450 deadROCTrends_[trendOrder[i]]->Fill(lumiSec - 1, nRocsPerTrend[i] - nFilledROCs[i]);
0451 ineffROCTrends_[trendOrder[i]]->Fill(lumiSec - 1, nFilledROCs[i] - hiEffROCs[i]);
0452 }
0453 }
0454
0455 if (!runOnEndLumi_)
0456 return;
0457
0458 for (auto it : {1, 2, 3, 4}) {
0459 histName = "PixelPhase1/Phase1_MechanicalView/PXBarrel/clusterposition_zphi_PXLayer_" + std::to_string(it);
0460 MonitorElement* toReset = iGetter.get(histName);
0461 if (toReset != nullptr) {
0462 toReset->Reset();
0463 }
0464 histName = "PixelPhase1/FED/Dead Channels per ROC_per_SignedModuleCoord_per_SignedLadderCoord_PXLayer_" +
0465 std::to_string(it);
0466 MonitorElement* twoReset = iGetter.get(histName);
0467 if (twoReset != nullptr) {
0468 twoReset->Reset();
0469 }
0470 }
0471 for (auto it : {"-3", "-2", "-1", "+1", "+2", "+3"}) {
0472 histName = "PixelPhase1/Phase1_MechanicalView/PXForward/clusterposition_xy_PXDisk_" + std::string(it);
0473 MonitorElement* toReset = iGetter.get(histName);
0474 if (toReset != nullptr) {
0475 toReset->Reset();
0476 }
0477 }
0478 for (auto it : {1, 2}) {
0479 histName = "PixelPhase1/FED/Dead Channels per ROC_per_SignedDiskCoord_per_SignedBladePanelCoord_PXRing_" +
0480 std::to_string(it);
0481 MonitorElement* twoReset = iGetter.get(histName);
0482 if (twoReset != nullptr) {
0483 twoReset->Reset();
0484 }
0485 }
0486 }
0487
0488
0489 DEFINE_FWK_MODULE(SiPixelPhase1Summary);