Skip to content

Commit b03913c

Browse files
committed
add SiPixelFEDChannelContainerMapWeigthed to SiPixelFEDChannelContainer_PayloadInspector
- displays aggregate map of the masked components for all scenarios, weighted on the probability per PU unit from SiPixelQualityProbabilities assuming a flat PU profile in the range encoded in SiPixelQualityProbabilities The SiPixelQualityProbabilities tag comes from user input.
1 parent c88bc2e commit b03913c

File tree

1 file changed

+216
-0
lines changed

1 file changed

+216
-0
lines changed

CondCore/SiPixelPlugins/plugins/SiPixelFEDChannelContainer_PayloadInspector.cc

Lines changed: 216 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -28,6 +28,7 @@
2828

2929
// the data format of the condition to be inspected
3030
#include "CondFormats/SiPixelObjects/interface/SiPixelFEDChannelContainer.h"
31+
#include "CondFormats/SiPixelObjects/interface/SiPixelQualityProbabilities.h" // to display aggregate probability
3132
#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
3233
#include "DataFormats/DetId/interface/DetId.h"
3334
#include "FWCore/MessageLogger/interface/MessageLogger.h"
@@ -449,6 +450,218 @@ namespace {
449450
using SiPixelFPixFEDChannelContainerMap = SiPixelFEDChannelContainerMapSimple<SiPixelPI::t_forward>;
450451
using SiPixelFullFEDChannelContainerMap = SiPixelFEDChannelContainerMapSimple<SiPixelPI::t_all>;
451452

453+
/*
454+
Produces an aggregate map of the masked components for all scenarios,
455+
weighted on the probability per PU unit from SiPixelQualityProbabilities
456+
assuming a flat PU profile in the range encoded in SiPixelQualityProbabilities
457+
The SiPixelQualityProbabilities tag comes from user input
458+
*/
459+
template <SiPixelPI::DetType myType>
460+
class SiPixelFEDChannelContainerMapWeigthed : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
461+
public:
462+
SiPixelFEDChannelContainerMapWeigthed()
463+
: PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>(
464+
"SiPixelFEDChannelContainer Pixel Track Map of one (or more scenarios)") {
465+
// for inputs
466+
PlotBase::addInputParam("SiPixelQualityProbabilitiesTag");
467+
468+
// hardcoded connection to the SiPixelQualityProbability tag, though luck
469+
m_condSiPixelProb = "frontier://FrontierProd/CMS_CONDITIONS";
470+
m_connectionPool.setParameters(m_connectionPset);
471+
m_connectionPool.configure();
472+
}
473+
474+
bool fill() override {
475+
auto paramValues = PlotBase::inputParamValues();
476+
auto ip = paramValues.find("SiPixelQualityProbabilitiesTag");
477+
if (ip != paramValues.end()) {
478+
m_SiPixelProbTagName = ip->second;
479+
} else {
480+
edm::LogWarning(k_ClassName) << "\n WARNING!!!! \n The needed parameter SiPixelQualityProbabilitiesTag was not "
481+
"inputed from the user \n Display will be aborted \n\n";
482+
return false;
483+
}
484+
485+
Phase1PixelROCMaps theROCMap("", "Masking Probability [%]");
486+
487+
auto tag = PlotBase::getTag<0>();
488+
auto tagname = tag.name;
489+
auto iov = tag.iovs.front();
490+
491+
// open db session for the cabling map
492+
edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
493+
<< "Query the condition database " << m_condSiPixelProb;
494+
495+
cond::persistency::Session condDbSession = m_connectionPool.createSession(m_condSiPixelProb);
496+
condDbSession.transaction().start(true);
497+
498+
// query the database
499+
edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
500+
<< "Reading IOVs from tag " << m_SiPixelProbTagName;
501+
502+
const auto MIN_VAL = cond::timeTypeSpecs[cond::runnumber].beginValue;
503+
const auto MAX_VAL = cond::timeTypeSpecs[cond::runnumber].endValue;
504+
505+
// get the list of payloads for the Cabling Map
506+
std::vector<std::tuple<cond::Time_t, cond::Hash>> m_pixelProb_iovs;
507+
condDbSession.readIov(m_SiPixelProbTagName).selectRange(MIN_VAL, MAX_VAL, m_pixelProb_iovs);
508+
509+
// in MC there should be only 1 IOV, oh well...
510+
edm::LogPrint(k_ClassName) << " using the SiPixelQualityProbabilities with hash: "
511+
<< std::get<1>(m_pixelProb_iovs.front()) << std::endl;
512+
513+
auto probabilitiesPayload =
514+
condDbSession.fetchPayload<SiPixelQualityProbabilities>(std::get<1>(m_pixelProb_iovs.front()));
515+
516+
const auto& PUbins = probabilitiesPayload->getPileUpBins();
517+
518+
SiPixelQualityProbabilities::probabilityMap m_probabilities = probabilitiesPayload->getProbability_Map();
519+
520+
// find the PU-averaged (assuming flat PU in the range) probabilities for each scenario
521+
std::map<std::string, float> puAvgedProbabilities;
522+
for (const auto& [PUbin, ProbMap] : m_probabilities) {
523+
float totProbInPUbin = 0.f;
524+
for (const auto& [scenName, prob] : ProbMap) {
525+
totProbInPUbin += prob;
526+
if (puAvgedProbabilities.find(scenName) == puAvgedProbabilities.end()) {
527+
puAvgedProbabilities[scenName] += prob;
528+
} else {
529+
puAvgedProbabilities.insert({scenName, prob});
530+
}
531+
}
532+
LogDebug(k_ClassName) << "PU bin: " << PUbin << " tot probability " << totProbInPUbin << std::endl;
533+
}
534+
535+
std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
536+
const auto& scenarioMap = payload->getScenarioMap();
537+
538+
float totProb{0.f};
539+
for (const auto& [scenName, prob] : puAvgedProbabilities) {
540+
// only sum up the scenarios that are in the SiPixelFEDChannelContainer payload!
541+
if (scenarioMap.find(scenName) != scenarioMap.end()) {
542+
LogDebug(k_ClassName) << scenName << " : " << prob << std::endl;
543+
totProb += prob;
544+
}
545+
}
546+
547+
LogDebug(k_ClassName) << "Total probability to normalize to: " << totProb << std::endl;
548+
549+
//normalize the probabilities per scenario to the toal probability
550+
for (auto& pair : puAvgedProbabilities) {
551+
pair.second /= totProb;
552+
}
553+
554+
for (const auto& scenario : scenarioMap) {
555+
std::string scenName = scenario.first;
556+
LogDebug(k_ClassName) << "\t Found Scenario: " << scenName << " ==> dumping it";
557+
558+
// calculate the weight
559+
float w_frac = 0.f;
560+
if (puAvgedProbabilities.find(scenName) != puAvgedProbabilities.end()) {
561+
w_frac = puAvgedProbabilities[scenName];
562+
}
563+
564+
// if scenario is not in the probability payload, continue
565+
if (w_frac == 0.f)
566+
continue;
567+
568+
LogDebug(k_ClassName) << "scen: " << scenName << " weight: " << w_frac << " log(weight):" << log10(w_frac)
569+
<< std::endl;
570+
571+
const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
572+
for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
573+
const auto t_detid = disabledChannels.detId();
574+
int subid = DetId(t_detid).subdetId();
575+
LogDebug(k_ClassName) << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
576+
577+
std::bitset<16> badRocsFromFEDChannels;
578+
579+
for (const auto& ch : disabledChannels) {
580+
std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
581+
ch.fed,
582+
ch.link,
583+
ch.roc_first,
584+
ch.roc_last);
585+
586+
LogDebug(k_ClassName) << toOut_ << std::endl;
587+
for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
588+
badRocsFromFEDChannels.set(i_roc);
589+
}
590+
}
591+
592+
LogDebug(k_ClassName) << badRocsFromFEDChannels << std::endl;
593+
594+
const auto& myDetId = DetId(t_detid);
595+
596+
if (subid == PixelSubdetector::PixelBarrel) {
597+
theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, w_frac * 100);
598+
} // if it's barrel
599+
else if (subid == PixelSubdetector::PixelEndcap) {
600+
theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, w_frac * 100);
601+
} // if it's endcap
602+
else {
603+
throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
604+
} // else nonsense
605+
} // loop on the channels
606+
} // loop on the scenarios
607+
608+
gStyle->SetOptStat(0);
609+
//=========================
610+
TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
611+
canvas.cd();
612+
613+
auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
614+
615+
std::string IOVstring = (unpacked.first == 0)
616+
? std::to_string(unpacked.second)
617+
: (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
618+
619+
const auto headerText =
620+
fmt::sprintf("#bf{#scale[0.6]{#color[2]{%s}, #color[4]{%s}}}", tagname, m_SiPixelProbTagName);
621+
622+
switch (myType) {
623+
case SiPixelPI::t_barrel:
624+
theROCMap.drawBarrelMaps(canvas, headerText);
625+
break;
626+
case SiPixelPI::t_forward:
627+
theROCMap.drawForwardMaps(canvas, headerText);
628+
break;
629+
case SiPixelPI::t_all:
630+
theROCMap.drawMaps(canvas, headerText);
631+
break;
632+
default:
633+
throw cms::Exception("LogicError") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
634+
}
635+
636+
// add list of scenarios watermark
637+
canvas.cd();
638+
auto ltx = TLatex();
639+
ltx.SetTextFont(62);
640+
//ltx.SetTextColor(kMagenta);
641+
ltx.SetTextSize(0.023);
642+
ltx.DrawLatexNDC(gPad->GetLeftMargin() - 0.09, gPad->GetBottomMargin() - 0.09, "");
643+
644+
std::string fileName(m_imageFileName);
645+
canvas.SaveAs(fileName.c_str());
646+
return true;
647+
}
648+
649+
private:
650+
// graphics
651+
static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
652+
static constexpr const char* k_ClassName = "SiPixelFEDChannelContainerMapWeigthed";
653+
654+
// parameters for auxilliary DB connection
655+
edm::ParameterSet m_connectionPset;
656+
cond::persistency::ConnectionPool m_connectionPool;
657+
std::string m_SiPixelProbTagName;
658+
std::string m_condSiPixelProb;
659+
};
660+
661+
using SiPixelBPixFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_barrel>;
662+
using SiPixelFPixFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_forward>;
663+
using SiPixelFullFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_all>;
664+
452665
/************************************************
453666
1d histogram of number of SiPixelFEDChannelContainer scenarios
454667
*************************************************/
@@ -547,5 +760,8 @@ PAYLOAD_INSPECTOR_MODULE(SiPixelFEDChannelContainer) {
547760
PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerMap);
548761
PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerMap);
549762
PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerMap);
763+
PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerWeightedMap);
764+
PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerWeightedMap);
765+
PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerWeightedMap);
550766
PAYLOAD_INSPECTOR_CLASS(SiPixelFEDChannelContainerScenarios);
551767
}

0 commit comments

Comments
 (0)