CMS 3D CMS Logo

SiPixelQuality_PayloadInspector.cc
Go to the documentation of this file.
1 
10 
18 
19 // the data format of the condition to be inspected
25 
26 #include <memory>
27 #include <sstream>
28 #include <iostream>
29 
30 // include ROOT
31 #include "TH2F.h"
32 #include "TLegend.h"
33 #include "TCanvas.h"
34 #include "TLine.h"
35 #include "TGraph.h"
36 #include "TStyle.h"
37 #include "TLatex.h"
38 #include "TPave.h"
39 #include "TPaveStats.h"
40 
41 namespace {
42 
43  using namespace cond::payloadInspector;
44 
45  /************************************************
46  test class
47  *************************************************/
48 
49  class SiPixelQualityTest : public Histogram1D<SiPixelQuality, SINGLE_IOV> {
50  public:
51  SiPixelQualityTest()
52  : Histogram1D<SiPixelQuality, SINGLE_IOV>("SiPixelQuality test", "SiPixelQuality test", 10, 0.0, 10.0) {}
53 
54  bool fill() override {
55  auto tag = PlotBase::getTag<0>();
56  for (auto const& iov : tag.iovs) {
57  std::shared_ptr<SiPixelQuality> payload = Base::fetchPayload(std::get<1>(iov));
58  if (payload.get()) {
59  fillWithValue(1.);
60 
61  auto theDisabledModules = payload->getBadComponentList();
62  for (const auto& mod : theDisabledModules) {
63  int BadRocCount(0);
64  for (unsigned short n = 0; n < 16; n++) {
65  unsigned short mask = 1 << n; // 1 << n = 2^{n} using bitwise shift
66  if (mod.BadRocs & mask)
67  BadRocCount++;
68  }
69  COUT << "detId:" << mod.DetID << " error type:" << mod.errorType << " BadRocs:" << BadRocCount << std::endl;
70  }
71  } // payload
72  } // iovs
73  return true;
74  } // fill
75  };
76 
77  /************************************************
78  summary class
79  *************************************************/
80 
81  class SiPixelQualityBadRocsSummary : public PlotImage<SiPixelQuality, MULTI_IOV, 1> {
82  public:
83  SiPixelQualityBadRocsSummary() : PlotImage<SiPixelQuality, MULTI_IOV, 1>("SiPixel Quality Summary") {}
84 
85  bool fill() override {
86  auto tag = PlotBase::getTag<0>();
87  for (const auto& iov : tag.iovs) {
88  std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));
89  auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
90 
91  COUT << "======================= " << unpacked.first << " : " << unpacked.second << std::endl;
92  auto theDisabledModules = payload->getBadComponentList();
93  for (const auto& mod : theDisabledModules) {
94  COUT << "detId: " << mod.DetID << " |error type: " << mod.errorType << " |BadRocs: " << mod.BadRocs
95  << std::endl;
96  }
97  }
98 
99  //=========================
100  TCanvas canv("Partition summary", "partition summary", 1200, 1000);
101  canv.SetBottomMargin(0.11);
102  canv.SetLeftMargin(0.13);
103  canv.SetRightMargin(0.05);
104  canv.cd();
106 
107  std::string fileName(m_imageFileName);
108  canv.SaveAs(fileName.c_str());
109 
110  return true;
111  }
112  };
113 
114  /************************************************
115  time history class
116  *************************************************/
117 
118  class SiPixelQualityBadRocsTimeHistory : public TimeHistoryPlot<SiPixelQuality, std::pair<double, double> > {
119  public:
120  SiPixelQualityBadRocsTimeHistory()
121  : TimeHistoryPlot<SiPixelQuality, std::pair<double, double> >("bad ROCs count vs time", "bad ROCs count") {}
122 
123  std::pair<double, double> getFromPayload(SiPixelQuality& payload) override {
124  return std::make_pair(extractBadRocCount(payload), 0.);
125  }
126 
127  unsigned int extractBadRocCount(SiPixelQuality& payload) {
128  unsigned int BadRocCount(0);
129  auto theDisabledModules = payload.getBadComponentList();
130  for (const auto& mod : theDisabledModules) {
131  for (unsigned short n = 0; n < 16; n++) {
132  unsigned short mask = 1 << n; // 1 << n = 2^{n} using bitwise shift
133  if (mod.BadRocs & mask)
134  BadRocCount++;
135  }
136  }
137  return BadRocCount;
138  }
139  };
140 
141  /************************************************
142  occupancy style map whole Pixel
143  *************************************************/
144  template <SiPixelPI::DetType myType>
145  class SiPixelQualityMap : public PlotImage<SiPixelQuality, SINGLE_IOV> {
146  public:
147  SiPixelQualityMap()
148  : PlotImage<SiPixelQuality, SINGLE_IOV>("SiPixelQuality Pixel Map"),
150  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
151 
152  bool fill() override {
153  auto tag = PlotBase::getTag<0>();
154  auto iov = tag.iovs.front();
155  auto tagname = tag.name;
156  std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));
157 
158  Phase1PixelROCMaps theMap("");
159 
160  auto theDisabledModules = payload->getBadComponentList();
161  if (this->isPhase0(theDisabledModules)) {
162  edm::LogError("SiPixelQuality_PayloadInspector")
163  << "SiPixelQuality maps are not supported for non-Phase1 Pixel geometries !";
164  TCanvas canvas("Canv", "Canv", 1200, 1000);
166  std::string fileName(m_imageFileName);
167  canvas.SaveAs(fileName.c_str());
168  return false;
169  }
170 
171  for (const auto& mod : theDisabledModules) {
172  int subid = DetId(mod.DetID).subdetId();
173 
174  if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
175  (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
176  (myType == SiPixelPI::t_all)) {
177  std::bitset<16> bad_rocs(mod.BadRocs);
178  if (payload->IsModuleBad(mod.DetID)) {
179  theMap.fillWholeModule(mod.DetID, 1.);
180  } else {
181  theMap.fillSelectedRocs(mod.DetID, bad_rocs, 1.);
182  }
183  }
184  }
185 
186  gStyle->SetOptStat(0);
187  //=========================
188  TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
189  canvas.cd();
190 
191  auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
192 
193  std::string IOVstring = (unpacked.first == 0)
194  ? std::to_string(unpacked.second)
195  : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
196 
197  const auto headerText = fmt::sprintf("#color[4]{%s}, IOV: #color[4]{%s}", tagname, IOVstring);
198 
199  switch (myType) {
200  case SiPixelPI::t_barrel:
201  theMap.drawBarrelMaps(canvas, headerText);
202  break;
204  theMap.drawForwardMaps(canvas, headerText);
205  break;
206  case SiPixelPI::t_all:
207  theMap.drawMaps(canvas, headerText);
208  break;
209  default:
210  throw cms::Exception("SiPixelQualityMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
211  }
212 
213  std::string fileName(m_imageFileName);
214  canvas.SaveAs(fileName.c_str());
215 #ifdef MMDEBUG
216  canvas.SaveAs("outAll.root");
217 #endif
218 
219  return true;
220  }
221 
222  private:
223  static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
224  TrackerTopology m_trackerTopo;
225 
226  //_________________________________________________
227  bool isPhase0(std::vector<SiPixelQuality::disabledModuleType> mods) {
230  const auto& p0detIds = reader.getAllDetIds();
231 
232  std::vector<uint32_t> ownDetIds;
233  std::transform(mods.begin(),
234  mods.end(),
235  std::back_inserter(ownDetIds),
236  [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });
237 
238  for (const auto& det : ownDetIds) {
239  // if found at least one phase-0 detId early return
240  if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
241  return true;
242  }
243  }
244  return false;
245  }
246  };
247 
248  using SiPixelBPixQualityMap = SiPixelQualityMap<SiPixelPI::t_barrel>;
249  using SiPixelFPixQualityMap = SiPixelQualityMap<SiPixelPI::t_forward>;
250  using SiPixelFullQualityMap = SiPixelQualityMap<SiPixelPI::t_all>;
251 
252  /************************************************
253  occupancy style map whole Pixel, difference of payloads
254  *************************************************/
255  template <SiPixelPI::DetType myType, IOVMultiplicity nIOVs, int ntags>
256  class SiPixelQualityMapComparisonBase : public PlotImage<SiPixelQuality, nIOVs, ntags> {
257  public:
258  SiPixelQualityMapComparisonBase()
259  : PlotImage<SiPixelQuality, nIOVs, ntags>(
260  Form("SiPixelQuality %s Pixel Map", SiPixelPI::DetNames[myType].c_str())),
262  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
263 
264  bool fill() override {
265  // trick to deal with the multi-ioved tag and two tag case at the same time
266  auto theIOVs = PlotBase::getTag<0>().iovs;
267  auto f_tagname = PlotBase::getTag<0>().name;
268  std::string l_tagname = "";
269  auto firstiov = theIOVs.front();
270  std::tuple<cond::Time_t, cond::Hash> lastiov;
271 
272  // we don't support (yet) comparison with more than 2 tags
273  assert(this->m_plotAnnotations.ntags < 3);
274 
275  if (this->m_plotAnnotations.ntags == 2) {
276  auto tag2iovs = PlotBase::getTag<1>().iovs;
277  l_tagname = PlotBase::getTag<1>().name;
278  lastiov = tag2iovs.front();
279  } else {
280  lastiov = theIOVs.back();
281  }
282 
283  std::shared_ptr<SiPixelQuality> last_payload = this->fetchPayload(std::get<1>(lastiov));
284  std::shared_ptr<SiPixelQuality> first_payload = this->fetchPayload(std::get<1>(firstiov));
285 
286  if (this->isPhase0(first_payload) || this->isPhase0(last_payload)) {
287  edm::LogError("SiPixelQuality_PayloadInspector")
288  << "SiPixelQuality comparison maps are not supported for non-Phase1 Pixel geometries !";
289  TCanvas canvas("Canv", "Canv", 1200, 1000);
291  std::string fileName(this->m_imageFileName);
292  canvas.SaveAs(fileName.c_str());
293  return false;
294  }
295 
296  Phase1PixelROCMaps theMap("", "#Delta payload A - payload B");
297 
298  gStyle->SetOptStat(0);
299  //=========================
300  TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
301  canvas.cd();
302 
303  auto f_unpacked = SiPixelPI::unpack(std::get<0>(firstiov));
304  auto l_unpacked = SiPixelPI::unpack(std::get<0>(lastiov));
305 
306  std::string f_IOVstring = (f_unpacked.first == 0)
307  ? std::to_string(f_unpacked.second)
308  : (std::to_string(f_unpacked.first) + "," + std::to_string(f_unpacked.second));
309 
310  std::string l_IOVstring = (l_unpacked.first == 0)
311  ? std::to_string(l_unpacked.second)
312  : (std::to_string(l_unpacked.first) + "," + std::to_string(l_unpacked.second));
313 
314  std::string headerText;
315 
316  if (this->m_plotAnnotations.ntags == 2) {
317  headerText = fmt::sprintf(
318  "#Delta #color[2]{A: %s, %s} - #color[4]{B: %s, %s}", f_tagname, f_IOVstring, l_tagname, l_IOVstring);
319  } else {
320  headerText =
321  fmt::sprintf("%s, #Delta IOV #color[2]{A: %s} - #color[4]{B: %s} ", f_tagname, f_IOVstring, l_IOVstring);
322  }
323 
324  switch (myType) {
325  case SiPixelPI::t_barrel:
326  theMap.drawBarrelMaps(canvas, headerText);
327  break;
329  theMap.drawForwardMaps(canvas, headerText);
330  break;
331  case SiPixelPI::t_all:
332  theMap.drawMaps(canvas, headerText);
333  break;
334  default:
335  throw cms::Exception("SiPixelQualityMapComparison")
336  << "\nERROR: unrecognized Pixel Detector part " << std::endl;
337  }
338 
339  // first loop on the first payload (newest)
340  fillTheMapFromPayload(theMap, first_payload, false);
341 
342  // then loop on the second payload (oldest)
343  fillTheMapFromPayload(theMap, last_payload, true); // true will subtract
344 
345  std::string fileName(this->m_imageFileName);
346  canvas.SaveAs(fileName.c_str());
347 #ifdef MMDEBUG
348  canvas.SaveAs("outAll.root");
349 #endif
350 
351  return true;
352  }
353 
354  private:
355  static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
356  TrackerTopology m_trackerTopo;
357 
358  //_________________________________________________
359  bool isPhase0(const std::shared_ptr<SiPixelQuality>& payload) {
360  const auto mods = payload->getBadComponentList();
363  const auto& p0detIds = reader.getAllDetIds();
364 
365  std::vector<uint32_t> ownDetIds;
366  std::transform(mods.begin(),
367  mods.end(),
368  std::back_inserter(ownDetIds),
369  [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });
370 
371  for (const auto& det : ownDetIds) {
372  // if found at least one phase-0 detId early return
373  if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
374  return true;
375  }
376  }
377  return false;
378  }
379 
380  //____________________________________________________________________________________________
381  void fillTheMapFromPayload(Phase1PixelROCMaps& theMap,
382  const std::shared_ptr<SiPixelQuality>& payload,
383  bool subtract) {
384  const auto theDisabledModules = payload->getBadComponentList();
385  for (const auto& mod : theDisabledModules) {
386  int subid = DetId(mod.DetID).subdetId();
387  if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
388  (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
389  (myType == SiPixelPI::t_all)) {
390  std::bitset<16> bad_rocs(mod.BadRocs);
391  if (payload->IsModuleBad(mod.DetID)) {
392  theMap.fillWholeModule(mod.DetID, (subtract ? -1. : 1.));
393  } else {
394  theMap.fillSelectedRocs(mod.DetID, bad_rocs, (subtract ? -1. : 1.));
395  }
396  }
397  }
398  }
399  };
400 
401  using SiPixelBPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, MULTI_IOV, 1>;
402  using SiPixelFPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, MULTI_IOV, 1>;
403  using SiPixelFullQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, MULTI_IOV, 1>;
404  using SiPixelBPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, SINGLE_IOV, 2>;
405  using SiPixelFPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, SINGLE_IOV, 2>;
406  using SiPixelFullQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, SINGLE_IOV, 2>;
407 
408 } // namespace
409 
410 // Register the classes as boost python plugin
412  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityTest);
413  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsSummary);
414  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsTimeHistory);
415  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMap);
416  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMap);
417  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMap);
418  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareSingleTag);
419  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareSingleTag);
420  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareSingleTag);
421  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareTwoTags);
422  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareTwoTags);
423  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareTwoTags);
424 }
SiPixelPI::unpack
std::pair< unsigned int, unsigned int > unpack(cond::Time_t since)
Definition: SiPixelPayloadInspectorHelper.h:114
Phase1PixelROCMaps
Definition: Phase1PixelROCMaps.h:81
svgfig.canvas
def canvas(*sub, **attr)
Definition: svgfig.py:482
SiPixelPI::t_barrel
Definition: SiPixelPayloadInspectorHelper.h:574
PixelSubdetector.h
MessageLogger.h
PixelSubdetector::PixelEndcap
Definition: PixelSubdetector.h:11
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
cond::payloadInspector::Histogram1
Definition: PayloadInspector.h:711
Phase1PixelROCMaps::fillWholeModule
void fillWholeModule(const uint32_t &detid, double value)
Definition: Phase1PixelROCMaps.cc:232
SiPixelPI
Definition: SiPixelPayloadInspectorHelper.h:37
contentValuesFiles.fullPath
fullPath
Definition: contentValuesFiles.py:64
TrackerTopology
Definition: TrackerTopology.h:16
PayloadInspector.h
mod
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
SiPixelQuality::disabledModuleType
Definition: SiPixelQuality.h:29
SiPixelPI::t_all
Definition: SiPixelPayloadInspectorHelper.h:574
cms::cuda::assert
assert(be >=bs)
PAYLOAD_INSPECTOR_CLASS
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
Definition: PayloadInspectorModule.h:10
MillePedeFileConverter_cfg.fileName
fileName
Definition: MillePedeFileConverter_cfg.py:32
SiPixelPI::displayNotSupported
void displayNotSupported(TCanvas &canv, const unsigned int size)
Definition: SiPixelPayloadInspectorHelper.h:782
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
SiPixelDetInfoFileReader::kPh0DefaultFile
constexpr static char const *const kPh0DefaultFile
Definition: SiPixelDetInfoFileReader.h:36
PayloadInspectorModule.h
mps_splice.mods
mods
Definition: mps_splice.py:43
DetId
Definition: DetId.h:17
SiPixelDetInfoFileReader
Definition: SiPixelDetInfoFileReader.h:28
edm::FileInPath
Definition: FileInPath.h:61
StandaloneTrackerTopology.h
DQM.reader
reader
Definition: DQM.py:105
HcalDetIdTransform::transform
unsigned transform(const HcalDetId &id, unsigned transformCode)
Definition: HcalDetIdTransform.cc:7
jets_cff.payload
payload
Definition: jets_cff.py:32
cond::payloadInspector
Definition: PayloadInspector.h:42
electronEcalRecHitIsolationLcone_cfi.subtract
subtract
Definition: electronEcalRecHitIsolationLcone_cfi.py:24
PAYLOAD_INSPECTOR_MODULE
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: PayloadInspectorModule.h:8
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
FileInPath.h
ntuplemaker.fill
fill
Definition: ntuplemaker.py:304
Time.h
SiPixelQuality
Definition: SiPixelQuality.h:27
makeGlobalPositionRcd_cfg.tag
tag
Definition: makeGlobalPositionRcd_cfg.py:6
cond::payloadInspector::MULTI_IOV
Definition: PayloadInspector.h:295
Phase1PixelROCMaps::fillSelectedRocs
void fillSelectedRocs(const uint32_t &detid, const std::bitset< 16 > &theROCs, double value)
Definition: Phase1PixelROCMaps.cc:255
cond::payloadInspector::TimeHistoryPlot
Definition: PayloadInspector.h:627
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
Phase1PixelROCMaps.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
createPayload.tagname
tagname
Definition: createPayload.py:183
SiPixelDetInfoFileReader.h
COUT
#define COUT
Definition: PVValidationHelpers.h:13
SiPixelPI::t_forward
Definition: SiPixelPayloadInspectorHelper.h:574
std
Definition: JetResolutionObject.h:76
DetId.h
Exception
Definition: hltDiff.cc:245
cond::payloadInspector::PlotImage
Definition: PayloadInspector.h:894
SiPixelQuality.h
ztail.d
d
Definition: ztail.py:151
SiPixelPayloadInspectorHelper.h
SiPixelPI::DetNames
const std::array< std::string, 3 > DetNames
Definition: SiPixelPayloadInspectorHelper.h:575
cond::payloadInspector::SINGLE_IOV
Definition: PayloadInspector.h:295
StandaloneTrackerTopology::fromTrackerParametersXMLFile
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
Definition: StandaloneTrackerTopology.cc:168