CMS 3D CMS Logo

SiPixelDynamicInefficiency_PayloadInspector.cc
Go to the documentation of this file.
1 
10 
18 
19 // the data format of the condition to be inspected
26 
27 #include <memory>
28 #include <sstream>
29 #include <iostream>
30 
31 // include ROOT
32 #include "TCanvas.h"
33 #include "TGraph.h"
34 #include "TH2F.h"
35 #include "TLatex.h"
36 #include "TLegend.h"
37 #include "TLine.h"
38 #include "TPave.h"
39 #include "TPaveStats.h"
40 #include "TStyle.h"
41 
42 namespace {
43 
44  using namespace cond::payloadInspector;
45  namespace SiPixDynIneff {
46 
47  // different types of geometrical inefficiency factors
48  enum factor { geom = 0, colgeom = 1, chipgeom = 2, pu = 3, INVALID = 4 };
49  const std::array<std::string, 5> factorString = {
50  {"pixel geometry", "column geometry", "chip geometry", "PU", "invalid"}};
51 
52  using FactorMap = std::map<unsigned int, double>;
53  using PUFactorMap = std::map<unsigned int, std::vector<double> >;
54 
55  // constants for ROC level simulation for Phase1
56  enum shiftEnumerator { FPixRocIdShift = 3, BPixRocIdShift = 6 };
57  const int rocIdMaskBits = 0x1F;
58 
59  struct packedBadRocFraction {
60  std::vector<int> badRocNumber;
61  std::vector<float> badRocFrac;
62  };
63 
64  using BRFractions = std::unordered_map<uint32_t, packedBadRocFraction>;
65 
66  //_________________________________________________
67  BRFractions pbrf(std::shared_ptr<SiPixelDynamicInefficiency> payload) {
68  BRFractions f;
69  const std::map<uint32_t, double>& PixelGeomFactorsDBIn = payload->getPixelGeomFactors();
70 
71  // first fill
72  for (const auto db_factor : PixelGeomFactorsDBIn) {
73  int subid = DetId(db_factor.first).subdetId();
74  int shift = (subid == static_cast<int>(PixelSubdetector::PixelBarrel)) ? BPixRocIdShift : FPixRocIdShift;
75  unsigned int rocMask = rocIdMaskBits << shift;
76  unsigned int rocId = (((db_factor.first) & rocMask) >> shift);
77  uint32_t rawid = db_factor.first & (~rocMask);
78 
79  if (f.find(rawid) == f.end()) {
80  packedBadRocFraction p;
81  f.insert(std::make_pair(rawid, p));
82  }
83 
84  if (rocId != 0) {
85  rocId--;
86  double factor = db_factor.second;
87  double badFraction = 1 - factor;
88 
89  f.at(rawid).badRocNumber.emplace_back(rocId);
90  f.at(rawid).badRocFrac.emplace_back(badFraction);
91  }
92  }
93  return f;
94  }
95 
96  //_________________________________________________
97  bool isPhase0(const BRFractions& fractions) {
100  const auto& p0detIds = reader.getAllDetIds();
101  std::vector<uint32_t> ownDetIds;
102 
103  std::transform(fractions.begin(),
104  fractions.end(),
105  std::back_inserter(ownDetIds),
106  [](std::pair<uint32_t, packedBadRocFraction> d) -> uint32_t { return d.first; });
107 
108  for (const auto& det : ownDetIds) {
109  // if found at least one phase-0 detId early return
110  if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
111  return true;
112  }
113  }
114  return false;
115  }
116 
117  //_________________________________________________
118  double getMatchingGeomFactor(const DetId& detid,
119  const std::map<unsigned int, double>& map_geomfactor,
120  const std::vector<uint32_t>& detIdmasks) {
121  double geomfactor_db = 1;
122  for (auto map_element : map_geomfactor) {
123  const DetId mapid = DetId(map_element.first);
124  if (mapid.subdetId() != detid.subdetId())
125  continue;
126  size_t __i = 0;
127  for (; __i < detIdmasks.size(); __i++) {
128  DetId maskid = DetId(detIdmasks.at(__i));
129  if (maskid.subdetId() != mapid.subdetId())
130  continue;
131  if ((detid.rawId() & maskid.rawId()) != (mapid.rawId() & maskid.rawId()) &&
132  (mapid.rawId() & maskid.rawId()) != DetId(mapid.det(), mapid.subdetId()).rawId())
133  break;
134  }
135  if (__i != detIdmasks.size())
136  continue;
137  geomfactor_db *= map_element.second;
138  }
139  return geomfactor_db;
140  }
141 
142  //_________________________________________________
143  std::vector<double> getMatchingPUFactors(const DetId& detid,
144  const std::map<unsigned int, std::vector<double> >& map_pufactory,
145  const std::vector<uint32_t>& detIdmasks) {
146  std::vector<double> pufactors_db;
147  for (const auto& map_element : map_pufactory) {
148  const DetId mapid = DetId(map_element.first);
149  if (mapid.subdetId() != detid.subdetId())
150  continue;
151  size_t __i = 0;
152  for (; __i < detIdmasks.size(); __i++) {
153  DetId maskid = DetId(detIdmasks.at(__i));
154  if (maskid.subdetId() != mapid.subdetId())
155  continue;
156  if ((detid.rawId() & maskid.rawId()) != (mapid.rawId() & maskid.rawId()) &&
157  (mapid.rawId() & maskid.rawId()) != DetId(mapid.det(), mapid.subdetId()).rawId())
158  break;
159  }
160  if (__i != detIdmasks.size())
161  continue;
162  pufactors_db = map_element.second;
163  }
164  return pufactors_db;
165  }
166 
167  /*
168  //(Not used for the moment)
169  //_________________________________________________
170  bool matches(const DetId& detid, const DetId& db_id, const std::vector<uint32_t>& DetIdmasks) {
171  if (detid.subdetId() != db_id.subdetId())
172  return false;
173  for (size_t i = 0; i < DetIdmasks.size(); ++i) {
174  DetId maskid = DetId(DetIdmasks.at(i));
175  if (maskid.subdetId() != db_id.subdetId())
176  continue;
177  if ((detid.rawId() & maskid.rawId()) != (db_id.rawId() & maskid.rawId()) &&
178  (db_id.rawId() & maskid.rawId()) != DetId(db_id.det(), db_id.subdetId()).rawId())
179  return false;
180  }
181  return true;
182  }
183  */
184 
185  //_________________________________________________
186  bool checkPhase(const SiPixelPI::phase phase, const std::vector<uint32_t>& masks_db) {
187  const char* inputFile;
188  switch (phase) {
190  inputFile = "Geometry/TrackerCommonData/data/trackerParameters.xml";
191  break;
193  inputFile = "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
194  break;
196  inputFile = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
197  break;
198  }
199 
200  // create the standalone tracker topology
201  const auto& tkTopo =
203 
204  // Check what masks we would get using the current geometry
205  // It has to match what is in the db content!!
206 
207  std::vector<uint32_t> masks_geom;
209 
210  masks_geom.push_back(tkTopo.pxbDetId(max, 0, 0).rawId());
211  masks_geom.push_back(tkTopo.pxbDetId(0, max, 0).rawId());
212  masks_geom.push_back(tkTopo.pxbDetId(0, 0, max).rawId());
213  masks_geom.push_back(tkTopo.pxfDetId(max, 0, 0, 0, 0).rawId());
214  masks_geom.push_back(tkTopo.pxfDetId(0, max, 0, 0, 0).rawId());
215  masks_geom.push_back(tkTopo.pxfDetId(0, 0, max, 0, 0).rawId());
216  masks_geom.push_back(tkTopo.pxfDetId(0, 0, 0, max, 0).rawId());
217  masks_geom.push_back(tkTopo.pxfDetId(0, 0, 0, 0, max).rawId());
218 
219  return (masks_geom.size() == masks_db.size() &&
220  std::equal(masks_geom.begin(), masks_geom.end(), masks_db.begin()));
221  }
222  } // namespace SiPixDynIneff
223 
224  /************************************************
225  test class
226  *************************************************/
227 
228  class SiPixelDynamicInefficiencyTest : public Histogram1D<SiPixelDynamicInefficiency, SINGLE_IOV> {
229  public:
230  SiPixelDynamicInefficiencyTest()
232  "SiPixelDynamicInefficiency test", "SiPixelDynamicInefficiency test", 1, 0.0, 1.0) {}
233 
234  bool fill() override {
235  auto tag = PlotBase::getTag<0>();
236  for (auto const& iov : tag.iovs) {
237  std::shared_ptr<SiPixelDynamicInefficiency> payload = Base::fetchPayload(std::get<1>(iov));
238  if (payload.get()) {
239  fillWithValue(1.);
240 
241  std::map<unsigned int, double> map_pixelgeomfactor = payload->getPixelGeomFactors();
242  std::map<unsigned int, double> map_colgeomfactor = payload->getColGeomFactors();
243  std::map<unsigned int, double> map_chipgeomfactor = payload->getChipGeomFactors();
244  std::map<unsigned int, std::vector<double> > map_pufactor = payload->getPUFactors();
245  std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
246  double theInstLumiScaleFactor_db = payload->gettheInstLumiScaleFactor_();
247 
248  edm::LogPrint("SiPixelDynamicInefficiencyTest") << "-------------------------------------------------------";
249  edm::LogPrint("SiPixelDynamicInefficiencyTest") << "Printing out DB content:\n";
250 
251  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " PixelGeomFactors:";
252  for (auto pixel : map_pixelgeomfactor)
253  edm::LogPrint("SiPixelDynamicInefficiencyTest")
254  << " MapID = " << pixel.first << "\tFactor = " << pixel.second;
255  edm::LogPrint("SiPixelDynamicInefficiencyTest");
256 
257  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " ColGeomFactors:";
258  for (auto col : map_colgeomfactor)
259  edm::LogPrint("SiPixelDynamicInefficiencyTest")
260  << " MapID = " << col.first << "\tFactor = " << col.second;
261  edm::LogPrint("SiPixelDynamicInefficiencyTest");
262 
263  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " ChipGeomFactors:";
264  for (auto chip : map_chipgeomfactor)
265  edm::LogPrint("SiPixelDynamicInefficiencyTest")
266  << " MapID = " << chip.first << "\tFactor = " << chip.second;
267  edm::LogPrint("SiPixelDynamicInefficiencyTest");
268 
269  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " PUFactors:";
270  for (auto pu : map_pufactor) {
271  edm::LogPrint("SiPixelDynamicInefficiencyTest")
272  << " MapID = " << pu.first << "\t Factor" << (pu.second.size() > 1 ? "s" : "") << " = ";
273  for (size_t i = 0, n = pu.second.size(); i < n; ++i)
274  edm::LogPrint("SiPixelDynamicInefficiencyTest") << pu.second[i] << ((i == n - 1) ? "\n" : ", ");
275  }
276  edm::LogPrint("SiPixelDynamicInefficiencyTest");
277 
278  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " DetIdmasks:";
279  for (auto mask : detIdmasks_db)
280  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " MaskID = " << mask;
281  edm::LogPrint("SiPixelDynamicInefficiencyTest");
282 
283  edm::LogPrint("SiPixelDynamicInefficiencyTest") << " theInstLumiScaleFactor = " << theInstLumiScaleFactor_db;
284 
285  } // payload
286  } // iovs
287  return true;
288  } // fill
289  };
290 
291  /************************************************
292  occupancy style map whole Pixel of inefficient ROCs
293  *************************************************/
294  template <SiPixelPI::DetType myType>
295  class SiPixelIneffROCfromDynIneffMap : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
296  public:
297  SiPixelIneffROCfromDynIneffMap()
298  : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixel Inefficient ROC from Dyn Ineff Pixel Map"),
300  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
301 
302  bool fill() override {
303  auto tag = PlotBase::getTag<0>();
304  auto iov = tag.iovs.front();
305  auto tagname = tag.name;
306  std::shared_ptr<SiPixelDynamicInefficiency> payload = fetchPayload(std::get<1>(iov));
307 
308  const auto fr = SiPixDynIneff::pbrf(payload);
309 
310  if (SiPixDynIneff::isPhase0(fr)) {
311  edm::LogError("SiPixelDynamicInefficiency_PayloadInspector")
312  << "SiPixelIneffROCfromDynIneff maps are not supported for non-Phase1 Pixel geometries !";
313  TCanvas canvas("Canv", "Canv", 1200, 1000);
315  std::string fileName(m_imageFileName);
316  canvas.SaveAs(fileName.c_str());
317  return false;
318  }
319 
320  Phase1PixelROCMaps theMap("", "bad pixel fraction in ROC [%]");
321 
322  for (const auto& element : fr) {
323  auto rawid = element.first;
324  int subid = DetId(rawid).subdetId();
325  auto packedinfo = element.second;
326  auto badRocs = packedinfo.badRocNumber;
327  auto badRocsF = packedinfo.badRocFrac;
328 
329  for (size_t i = 0; i < badRocs.size(); i++) {
330  std::bitset<16> rocToMark;
331  rocToMark.set(badRocs[i]);
332  if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
333  (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
334  (myType == SiPixelPI::t_all)) {
335  theMap.fillSelectedRocs(rawid, rocToMark, badRocsF[i] * 100.f);
336  }
337  }
338  }
339 
340  gStyle->SetOptStat(0);
341  //=========================
342  TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
343  canvas.cd();
344 
345  auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
346 
347  std::string IOVstring = (unpacked.first == 0)
348  ? std::to_string(unpacked.second)
349  : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
350 
351  const auto headerText = fmt::sprintf("#color[4]{%s}, IOV: #color[4]{%s}", tagname, IOVstring);
352 
353  switch (myType) {
354  case SiPixelPI::t_barrel:
355  theMap.drawBarrelMaps(canvas, headerText);
356  break;
358  theMap.drawForwardMaps(canvas, headerText);
359  break;
360  case SiPixelPI::t_all:
361  theMap.drawMaps(canvas, headerText);
362  break;
363  default:
364  throw cms::Exception("SiPixelIneffROCfromDynIneffMap") << "\nERROR: unrecognized Pixel Detector part ";
365  }
366 
367  std::string fileName(m_imageFileName);
368  canvas.SaveAs(fileName.c_str());
369 #ifdef MMDEBUG
370  canvas.SaveAs("outAll.root");
371 #endif
372 
373  return true;
374  }
375 
376  private:
377  static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
378  TrackerTopology m_trackerTopo;
379  };
380 
381  using SiPixelBPixIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_barrel>;
382  using SiPixelFPixIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_forward>;
383  using SiPixelFullIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_all>;
384 
385  /************************************************
386  occupancy style map whole Pixel, difference of payloads
387  *************************************************/
388  template <SiPixelPI::DetType myType, IOVMultiplicity nIOVs, int ntags>
389  class SiPixelIneffROCComparisonBase : public PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags> {
390  public:
391  SiPixelIneffROCComparisonBase()
392  : PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags>(
393  Form("SiPixelDynamicInefficiency %s Pixel Map", SiPixelPI::DetNames[myType].c_str())),
395  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
396 
397  bool fill() override {
398  // trick to deal with the multi-ioved tag and two tag case at the same time
399  auto theIOVs = PlotBase::getTag<0>().iovs;
400  auto f_tagname = PlotBase::getTag<0>().name;
401  std::string l_tagname = "";
402  auto firstiov = theIOVs.front();
403  std::tuple<cond::Time_t, cond::Hash> lastiov;
404 
405  // we don't support (yet) comparison with more than 2 tags
406  assert(this->m_plotAnnotations.ntags < 3);
407 
408  if (this->m_plotAnnotations.ntags == 2) {
409  auto tag2iovs = PlotBase::getTag<1>().iovs;
410  l_tagname = PlotBase::getTag<1>().name;
411  lastiov = tag2iovs.front();
412  } else {
413  lastiov = theIOVs.back();
414  }
415 
416  std::shared_ptr<SiPixelDynamicInefficiency> last_payload = this->fetchPayload(std::get<1>(lastiov));
417  std::shared_ptr<SiPixelDynamicInefficiency> first_payload = this->fetchPayload(std::get<1>(firstiov));
418 
419  const auto fp = SiPixDynIneff::pbrf(last_payload);
420  const auto lp = SiPixDynIneff::pbrf(first_payload);
421 
422  if (SiPixDynIneff::isPhase0(fp) || SiPixDynIneff::isPhase0(lp)) {
423  edm::LogError("SiPixelDynamicInefficiency_PayloadInspector")
424  << "SiPixelDynamicInefficiency comparison maps are not supported for non-Phase1 Pixel geometries !";
425  TCanvas canvas("Canv", "Canv", 1200, 1000);
427  std::string fileName(this->m_imageFileName);
428  canvas.SaveAs(fileName.c_str());
429  return false;
430  }
431 
432  Phase1PixelROCMaps theMap("", "#Delta payload A - payload B");
433 
434  gStyle->SetOptStat(0);
435  //=========================
436  TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
437  canvas.cd();
438 
439  auto f_unpacked = SiPixelPI::unpack(std::get<0>(firstiov));
440  auto l_unpacked = SiPixelPI::unpack(std::get<0>(lastiov));
441 
442  std::string f_IOVstring = (f_unpacked.first == 0)
443  ? std::to_string(f_unpacked.second)
444  : (std::to_string(f_unpacked.first) + "," + std::to_string(f_unpacked.second));
445 
446  std::string l_IOVstring = (l_unpacked.first == 0)
447  ? std::to_string(l_unpacked.second)
448  : (std::to_string(l_unpacked.first) + "," + std::to_string(l_unpacked.second));
449 
450  std::string headerText;
451 
452  if (this->m_plotAnnotations.ntags == 2) {
453  headerText =
454  fmt::sprintf("#color[2]{A: %s, %s} - #color[4]{B: %s, %s}", f_tagname, f_IOVstring, l_tagname, l_IOVstring);
455  } else {
456  headerText = fmt::sprintf("%s,IOV #color[2]{A: %s} - #color[4]{B: %s} ", f_tagname, f_IOVstring, l_IOVstring);
457  }
458 
459  switch (myType) {
460  case SiPixelPI::t_barrel:
461  theMap.drawBarrelMaps(canvas, headerText);
462  break;
464  theMap.drawForwardMaps(canvas, headerText);
465  break;
466  case SiPixelPI::t_all:
467  theMap.drawMaps(canvas, headerText);
468  break;
469  default:
470  throw cms::Exception("SiPixelDynamicInefficiencyMapComparison")
471  << "\nERROR: unrecognized Pixel Detector part ";
472  }
473 
474  // first loop on the first payload (newest)
475  fillTheMapFromPayload(theMap, fp, false);
476 
477  // then loop on the second payload (oldest)
478  fillTheMapFromPayload(theMap, lp, true); // true will subtract
479 
480  std::string fileName(this->m_imageFileName);
481  canvas.SaveAs(fileName.c_str());
482 #ifdef MMDEBUG
483  canvas.SaveAs("outAll.root");
484 #endif
485 
486  return true;
487  }
488 
489  private:
490  static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
491  TrackerTopology m_trackerTopo;
492 
493  //____________________________________________________________________________________________
494  void fillTheMapFromPayload(Phase1PixelROCMaps& theMap, const SiPixDynIneff::BRFractions& fr, bool subtract) {
495  for (const auto& element : fr) {
496  auto rawid = element.first;
497  int subid = DetId(rawid).subdetId();
498  auto packedinfo = element.second;
499  auto badRocs = packedinfo.badRocNumber;
500  auto badRocsF = packedinfo.badRocFrac;
501 
502  for (size_t i = 0; i < badRocs.size(); i++) {
503  std::bitset<16> rocToMark;
504  rocToMark.set(badRocs[i]);
505  if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
506  (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
507  (myType == SiPixelPI::t_all)) {
508  theMap.fillSelectedRocs(rawid, rocToMark, badRocsF[i] * (subtract ? -1. : 1.));
509  }
510  }
511  }
512  }
513  };
514 
515  /*
516  These are not implemented for the time being, since the SiPixelDynamicInefficiency is a condition
517  used only in simulation, hence there is no such thing as a multi-IoV Dynamic Inefficiency tag
518  */
519 
520  using SiPixelBPixIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_barrel, MULTI_IOV, 1>;
521  using SiPixelFPixIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_forward, MULTI_IOV, 1>;
522  using SiPixelFullIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_all, MULTI_IOV, 1>;
523 
524  using SiPixelBPixIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_barrel, SINGLE_IOV, 2>;
525  using SiPixelFPixIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_forward, SINGLE_IOV, 2>;
526  using SiPixelFullIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_all, SINGLE_IOV, 2>;
527 
528  /************************************************
529  Full Pixel Tracker Map class (for geometrical factors)
530  *************************************************/
531  template <SiPixDynIneff::factor theFactor>
532  class SiPixelDynamicInefficiencyFullPixelMap : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
533  public:
534  SiPixelDynamicInefficiencyFullPixelMap()
535  : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixelDynamicInefficiency Map") {
536  label_ = "SiPixelDynamicInefficiencyFullPixelMap";
537  payloadString = fmt::sprintf("%s Dynamic Inefficiency", SiPixDynIneff::factorString[theFactor]);
538  }
539 
540  bool fill() override {
541  gStyle->SetPalette(1);
542  auto tag = PlotBase::getTag<0>();
543  auto iov = tag.iovs.front();
544  std::shared_ptr<SiPixelDynamicInefficiency> payload = this->fetchPayload(std::get<1>(iov));
545 
546  if (payload.get()) {
547  Phase1PixelSummaryMap fullMap("", fmt::sprintf("%s", payloadString), fmt::sprintf("%s", payloadString));
548  fullMap.createTrackerBaseMap();
549 
550  SiPixDynIneff::FactorMap theMap{};
551  switch (theFactor) {
552  case SiPixDynIneff::geom:
553  theMap = payload->getPixelGeomFactors();
554  break;
555  case SiPixDynIneff::colgeom:
556  theMap = payload->getColGeomFactors();
557  break;
558  case SiPixDynIneff::chipgeom:
559  theMap = payload->getChipGeomFactors();
560  break;
561  default:
562  throw cms::Exception(label_) << "\nERROR: unrecognized type of geometry factor ";
563  }
564 
565  std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
566 
567  if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, detIdmasks_db)) {
568  edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
569  TCanvas canvas("Canv", "Canv", 1200, 1000);
571  std::string fileName(m_imageFileName);
572  canvas.SaveAs(fileName.c_str());
573  return false;
574  }
575 
578  const auto& p1detIds = reader.getAllDetIds();
579  for (const auto& det : p1detIds) {
580  const auto& value = SiPixDynIneff::getMatchingGeomFactor(det, theMap, detIdmasks_db);
581  fullMap.fillTrackerMap(det, value);
582  }
583 
584  const auto& range = fullMap.getZAxisRange();
585  if (range.first == range.second) {
586  // in case the map is completely filled with one value;
587  // set the z-axis to be meaningful
588  fullMap.setZAxisRange(range.first - 0.01, range.second + 0.01);
589  }
590 
591  TCanvas canvas("Canv", "Canv", 3000, 2000);
592  fullMap.printTrackerMap(canvas);
593 
594  auto ltx = TLatex();
595  ltx.SetTextFont(62);
596  ltx.SetTextSize(0.025);
597  ltx.SetTextAlign(11);
598  ltx.DrawLatexNDC(
599  gPad->GetLeftMargin() + 0.01,
600  gPad->GetBottomMargin() + 0.01,
601  ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
602 
603  std::string fileName(this->m_imageFileName);
604  canvas.SaveAs(fileName.c_str());
605  }
606  return true;
607  }
608 
609  protected:
610  std::string payloadString;
611  std::string label_;
612  };
613 
614  using SiPixelDynamicInefficiencyGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::geom>;
615  using SiPixelDynamicInefficiencyColGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::colgeom>;
616  using SiPixelDynamicInefficiencyChipGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::chipgeom>;
617 
618  /************************************************
619  Full Pixel Tracker Map class (for PU factors)
620  *************************************************/
621  class SiPixelDynamicInefficiencyPUPixelMaps : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
622  public:
623  SiPixelDynamicInefficiencyPUPixelMaps()
624  : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixelDynamicInefficiency Map") {
625  label_ = "SiPixelDynamicInefficiencyFullPixelMap";
626  payloadString = fmt::sprintf("%s Dynamic Inefficiency", SiPixDynIneff::factorString[SiPixDynIneff::pu]);
627  }
628 
629  bool fill() override {
630  gStyle->SetPalette(1);
631  auto tag = PlotBase::getTag<0>();
632  auto iov = tag.iovs.front();
633  std::shared_ptr<SiPixelDynamicInefficiency> payload = this->fetchPayload(std::get<1>(iov));
634 
635  if (payload.get()) {
636  std::vector<Phase1PixelSummaryMap> maps;
637 
638  SiPixDynIneff::PUFactorMap theMap = payload->getPUFactors();
639  std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
640 
641  if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, detIdmasks_db)) {
642  edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
643  TCanvas canvas("Canv", "Canv", 1200, 1000);
645  std::string fileName(m_imageFileName);
646  canvas.SaveAs(fileName.c_str());
647  return false;
648  }
649 
650  unsigned int depth = maxDepthOfPUArray(theMap);
651 
652  // create the maps
653  for (unsigned int i = 0; i < depth; i++) {
654  maps.emplace_back(
655  "", fmt::sprintf("%s, factor %i", payloadString, i), fmt::sprintf("%s, factor %i", payloadString, i));
656  maps[i].createTrackerBaseMap();
657  }
658 
659  // retrieve the list of phase1 detids
660  const auto& reader =
662  const auto& p1detIds = reader.getAllDetIds();
663 
664  // fill the maps
665  for (const auto& det : p1detIds) {
666  const auto& values = SiPixDynIneff::getMatchingPUFactors(det, theMap, detIdmasks_db);
667  int index = 0;
668  for (const auto& value : values) {
669  maps[index].fillTrackerMap(det, value);
670  index++;
671  }
672  }
673 
674  // in case the map is completely filled with one value;
675  // set the z-axis to be meaningful
676  for (unsigned int i = 0; i < depth; i++) {
677  const auto& range = maps[i].getZAxisRange();
678  if (range.first == range.second) {
679  maps[i].setZAxisRange(range.first - 0.01, range.second + 0.01);
680  }
681  }
682 
683  // determine how the plot will be paginated
684  auto sides = getClosestFactors(depth);
685  TCanvas canvas("Canv", "Canv", sides.second * 900, sides.first * 600);
686  canvas.Divide(sides.second, sides.first);
687 
688  // print the sub-canvases
689  for (unsigned int i = 0; i < depth; i++) {
690  maps[i].printTrackerMap(canvas, 0.035, i + 1);
691  auto ltx = TLatex();
692  ltx.SetTextFont(62);
693  ltx.SetTextSize(0.025);
694  ltx.SetTextAlign(11);
695  ltx.DrawLatexNDC(
696  gPad->GetLeftMargin() + 0.01,
697  gPad->GetBottomMargin() + 0.01,
698  ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
699  }
700 
701  std::string fileName(this->m_imageFileName);
702  canvas.SaveAs(fileName.c_str());
703  }
704  return true;
705  }
706 
707  protected:
708  std::string payloadString;
709  std::string label_;
710 
711  private:
712  unsigned int maxDepthOfPUArray(const std::map<unsigned int, std::vector<double> >& map_pufactor) {
713  unsigned int size{0};
714  for (const auto& [id, vec] : map_pufactor) {
715  if (vec.size() > size)
716  size = vec.size();
717  }
718  return size;
719  }
720 
721  std::pair<int, int> getClosestFactors(int input) {
722  if ((input % 2 != 0) && input > 1) {
723  input += 1;
724  }
725 
726  int testNum = (int)sqrt(input);
727  while (input % testNum != 0) {
728  testNum--;
729  }
730  return std::make_pair(testNum, input / testNum);
731  }
732  };
733 
734 } // namespace
735 
736 // Register the classes as boost python plugin
738  PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyTest);
739  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixIneffROCfromDynIneffMap);
740  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixIneffROCfromDynIneffMap);
741  PAYLOAD_INSPECTOR_CLASS(SiPixelFullIneffROCfromDynIneffMap);
742  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixIneffROCsMapCompareTwoTags);
743  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixIneffROCsMapCompareTwoTags);
744  PAYLOAD_INSPECTOR_CLASS(SiPixelFullIneffROCsMapCompareTwoTags);
745  PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyGeomFactorMap);
746  PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyColGeomFactorMap);
747  PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyChipGeomFactorMap);
748  PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyPUPixelMaps);
749 }
size
Write out results.
std::string to_string(const V &value)
Definition: OMSAccess.h:71
static constexpr char const *const kPh1DefaultFile
reader
Definition: DQM.py:105
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr uint32_t mask
Definition: gpuClustering.h:26
bool equal(const T &first, const T &second)
Definition: Equal.h:32
static std::string const input
Definition: EdmProvDump.cc:50
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
std::pair< unsigned int, unsigned int > unpack(cond::Time_t since)
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Log< level::Warning, true > LogPrint
d
Definition: ztail.py:151
Definition: DetId.h:17
const std::array< std::string, 3 > DetNames
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void displayNotSupported(TCanvas &canv, const unsigned int size)
void fillSelectedRocs(const uint32_t &detid, const std::bitset< 16 > &theROCs, double value)
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
col
Definition: cuy.py:1009
def canvas(sub, attr)
Definition: svgfig.py:482
static unsigned int const shift
static constexpr char const *const kPh0DefaultFile
unsigned transform(const HcalDetId &id, unsigned transformCode)