CMS 3D CMS Logo

SiPixelFEDChannelContainer_PayloadInspector.cc
Go to the documentation of this file.
1 
11 
19 
28 
29 // the data format of the condition to be inspected
35 
36 #include <memory>
37 #include <sstream>
38 #include <iostream>
39 #include <fmt/printf.h>
40 
41 // include ROOT
42 #include "TH2F.h"
43 #include "TLegend.h"
44 #include "TCanvas.h"
45 #include "TLine.h"
46 #include "TGraph.h"
47 #include "TGaxis.h"
48 #include "TStyle.h"
49 #include "TLatex.h"
50 #include "TPave.h"
51 #include "TPaveStats.h"
52 
53 namespace {
54 
55  using namespace cond::payloadInspector;
56 
57  /************************************************
58  1d histogram of SiPixelFEDChannelContainer of 1 IOV
59  *************************************************/
60 
61  template <SiPixelPI::DetType myType>
62  class SiPixelFEDChannelContainerMap : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
63  public:
64  SiPixelFEDChannelContainerMap()
65  : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>("SiPixelFEDChannelContainer scenarios count"),
67  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
68  // for inputs
69  PlotBase::addInputParam("Scenarios");
70 
71  // hardcoded connection to the MC cabling tag, though luck
72  m_condDbCabling = "frontier://FrontierProd/CMS_CONDITIONS";
73  m_CablingTagName = "SiPixelFedCablingMap_phase1_v7";
74 
75  m_connectionPool.setParameters(m_connectionPset);
76  m_connectionPool.configure();
77  }
78 
79  bool fill() override {
80  std::vector<std::string> the_scenarios = {};
81 
82  auto paramValues = PlotBase::inputParamValues();
83  auto ip = paramValues.find("Scenarios");
84  if (ip != paramValues.end()) {
85  auto input = ip->second;
86  typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
87  boost::char_separator<char> sep{","};
88  tokenizer tok{input, sep};
89  for (const auto& t : tok) {
90  the_scenarios.push_back(t);
91  }
92  } else {
93  edm::LogWarning("SiPixelFEDChannelContainerTest")
94  << "\n WARNING!!!! \n The needed parameter Scenarios has not been passed. Will use all the scenarios in "
95  "the file!"
96  << "\n Buckle your seatbelts... this might take a while... \n\n";
97  the_scenarios.push_back("all");
98  }
99 
100  Phase1PixelROCMaps theROCMap("");
101 
102  auto tag = PlotBase::getTag<0>();
103  auto tagname = tag.name;
104  auto iov = tag.iovs.front();
105 
106  // open db session for the cabling map
107  edm::LogPrint("SiPixelFEDChannelContainerTest") << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
108  << "Query the condition database " << m_condDbCabling;
109 
110  cond::persistency::Session condDbSession = m_connectionPool.createSession(m_condDbCabling);
111  condDbSession.transaction().start(true);
112 
113  // query the database
114  edm::LogPrint("SiPixelFEDChannelContainerTest") << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
115  << "Reading IOVs from tag " << m_CablingTagName;
116 
119 
120  // get the list of payloads for the Cabling Map
121  std::vector<std::tuple<cond::Time_t, cond::Hash>> m_cabling_iovs;
122  condDbSession.readIov(m_CablingTagName).selectRange(MIN_VAL, MAX_VAL, m_cabling_iovs);
123 
124  std::vector<unsigned int> listOfCablingIOVs;
125  std::transform(m_cabling_iovs.begin(),
126  m_cabling_iovs.end(),
127  std::back_inserter(listOfCablingIOVs),
128  [](std::tuple<cond::Time_t, cond::Hash> myIOV2) -> unsigned int { return std::get<0>(myIOV2); });
129 
130  edm::LogPrint("SiPixelFEDChannelContainerTest")
131  << " Number of SiPixelFedCablngMap payloads: " << listOfCablingIOVs.size() << std::endl;
132 
133  auto it = std::find(
134  listOfCablingIOVs.begin(), listOfCablingIOVs.end(), closest_from_below(listOfCablingIOVs, std::get<0>(iov)));
135  int index = std::distance(listOfCablingIOVs.begin(), it);
136 
137  edm::LogPrint("SiPixelFEDChannelContainerTest")
138  << " using the SiPixelFedCablingMap with hash: " << std::get<1>(m_cabling_iovs.at(index)) << std::endl;
139 
140  auto theCablingMap = condDbSession.fetchPayload<SiPixelFedCablingMap>(std::get<1>(m_cabling_iovs.at(index)));
141  theCablingMap->initializeRocs();
142  // auto theCablingTree = (*theCablingMap).cablingTree();
143 
144  //auto map = theCablingMap->det2fedMap();
145  //for (const auto &element : map){
146  // std::cout << element.first << " " << element.second << std::endl;
147  //}
148 
149  std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
150  const auto& scenarioMap = payload->getScenarioMap();
151 
152  auto pIndexConverter = PixelIndices(numColumns, numRows);
153 
154  for (const auto& scenario : scenarioMap) {
155  std::string scenName = scenario.first;
156 
157  if (std::find_if(the_scenarios.begin(), the_scenarios.end(), compareKeys(scenName)) != the_scenarios.end() ||
158  the_scenarios[0] == "all") {
159  edm::LogPrint("SiPixelFEDChannelContainerTest") << "\t Found Scenario: " << scenName << " ==> dumping it";
160  } else {
161  continue;
162  }
163 
164  //if (strcmp(scenName.c_str(),"320824_103") != 0) continue;
165 
166  const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
167  for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
168  const auto t_detid = disabledChannels.detId();
169  int subid = DetId(t_detid).subdetId();
170  LogDebug("SiPixelFEDChannelContainerTest") << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
171 
172  std::bitset<16> badRocsFromFEDChannels;
173 
174  for (const auto& ch : disabledChannels) {
175  std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
176  ch.fed,
177  ch.link,
178  ch.roc_first,
179  ch.roc_last);
180 
181  LogDebug("SiPixelFEDChannelContainerTest") << toOut_ << std::endl;
182  const std::vector<sipixelobjects::CablingPathToDetUnit>& path =
183  theCablingMap->pathToDetUnit(disabledChannels.detId());
184  for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
185  for (const auto p : path) {
186  const sipixelobjects::PixelROC* myroc = theCablingMap->findItem(p);
187  if (myroc->idInDetUnit() == static_cast<unsigned int>(i_roc)) {
188  sipixelobjects::LocalPixel::RocRowCol local = {39, 25}; //corresponding to center of ROC row,col
190  int chipIndex(0), colROC(0), rowROC(0);
191 
192  pIndexConverter.transformToROC(global.col, global.row, chipIndex, colROC, rowROC);
193 
194  LogDebug("SiPixelFEDChannelContainerTest")
195  << " => i_roc:" << i_roc << " " << global.col << "-" << global.row << " | => " << chipIndex
196  << " : (" << colROC << "," << rowROC << ")" << std::endl;
197 
198  badRocsFromFEDChannels[chipIndex] = true;
199  }
200  }
201  }
202  }
203 
204  LogDebug("SiPixelFEDChannelContainerTest") << badRocsFromFEDChannels << std::endl;
205 
206  auto myDetId = DetId(t_detid);
207 
208  if (subid == PixelSubdetector::PixelBarrel) {
209  theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
210  } // if it's barrel
211  else if (subid == PixelSubdetector::PixelEndcap) {
212  theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
213  } // if it's endcap
214  else {
215  throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
216  } // else nonsense
217  } // loop on the channels
218  } // loop on the scenarios
219 
220  gStyle->SetOptStat(0);
221  //=========================
222  TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
223  canvas.cd();
224 
225  auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
226 
227  std::string IOVstring = (unpacked.first == 0)
228  ? std::to_string(unpacked.second)
229  : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
230 
231  const auto headerText = fmt::sprintf("#color[4]{%s}, IOV: #color[4]{%s}", tagname, IOVstring);
232 
233  switch (myType) {
234  case SiPixelPI::t_barrel:
235  theROCMap.drawBarrelMaps(canvas, headerText);
236  break;
238  theROCMap.drawForwardMaps(canvas, headerText);
239  break;
240  case SiPixelPI::t_all:
241  theROCMap.drawMaps(canvas, headerText);
242  break;
243  default:
244  throw cms::Exception("SiPixelQualityMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
245  }
246 
247  std::string fileName(m_imageFileName);
248  canvas.SaveAs(fileName.c_str());
249 
250  // close the DB session
251  condDbSession.transaction().commit();
252 
253  return true;
254  }
255 
256  public:
257  inline unsigned int closest_from_above(std::vector<unsigned int> const& vec, unsigned int value) {
258  auto const it = std::lower_bound(vec.begin(), vec.end(), value);
259  return vec.at(it - vec.begin() - 1);
260  }
261 
262  inline unsigned int closest_from_below(std::vector<unsigned int> const& vec, unsigned int value) {
263  auto const it = std::upper_bound(vec.begin(), vec.end(), value);
264  return vec.at(it - vec.begin() - 1);
265  }
266 
267  // auxilliary check
268  struct compareKeys {
270  compareKeys(std::string const& i) : key(i) {}
271 
272  bool operator()(std::string const& i) { return (key == i); }
273  };
274 
275  private:
276  // tough luck, we can only do phase-1...
277  static constexpr int numColumns = 416;
278  static constexpr int numRows = 160;
279  static constexpr int n_rings = 2;
280  static constexpr int n_layers = 4;
281 
282  // graphics
283  static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
284 
285  TrackerTopology m_trackerTopo;
286  edm::ParameterSet m_connectionPset;
287  cond::persistency::ConnectionPool m_connectionPool;
288  std::string m_CablingTagName;
289  std::string m_condDbCabling;
290  };
291 
292  using SiPixelBPixFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_barrel>;
293  using SiPixelFPixFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_forward>;
294  using SiPixelFullFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_all>;
295 
296  /************************************************
297  1d histogram of SiPixelFEDChannelContainer of 1 IOV
298  *************************************************/
299 
300  class SiPixelFEDChannelContainerScenarios : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
301  public:
302  SiPixelFEDChannelContainerScenarios()
303  : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>("SiPixelFEDChannelContainer scenarios count") {}
304  bool fill() override {
305  auto tag = PlotBase::getTag<0>();
306  auto tagname = tag.name;
307  auto iov = tag.iovs.front();
308 
309  TGaxis::SetMaxDigits(3);
310 
311  std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
312  std::vector<std::string> scenarios = payload->getScenarioList();
313  sort(scenarios.begin(), scenarios.end());
314 
315  TCanvas canvas("Canv", "Canv", 1200, 1000);
316  canvas.cd();
317  canvas.SetGrid();
318  auto h1 = std::make_unique<TH1F>("Count",
319  "SiPixelFEDChannelContainer Bad Roc count;Scenario index;n. of bad ROCs",
320  scenarios.size(),
321  1,
322  scenarios.size());
323  h1->SetStats(false);
324 
325  canvas.SetTopMargin(0.06);
326  canvas.SetBottomMargin(0.12);
327  canvas.SetLeftMargin(0.12);
328  canvas.SetRightMargin(0.05);
329  canvas.Modified();
330 
331  int scenarioIndex = 0;
332  for (const auto& scenario : scenarios) {
333  scenarioIndex++;
334  int badRocCount = 0;
335  LogDebug("SiPixelFEDChannelContainerScenarios") << scenario << std::endl;
336  auto badChannelCollection = payload->getDetSetBadPixelFedChannels(scenario);
337  for (const auto& disabledChannels : *badChannelCollection) {
338  for (const auto& ch : disabledChannels) {
339  int local_bad_rocs = ch.roc_last - ch.roc_first;
340  badRocCount += local_bad_rocs;
341  } // loop on the channels
342  } // loop on the DetSetVector
343 
344  h1->SetBinContent(scenarioIndex, badRocCount);
345  } // loop on scenarios
346 
347  TGaxis::SetExponentOffset(-0.1, 0.01, "y"); // Y offset
348  TGaxis::SetExponentOffset(-0.03, -0.10, "x"); // Y and Y offset for X axis
349 
350  h1->SetTitle("");
351  h1->GetYaxis()->SetRangeUser(0., h1->GetMaximum() * 1.30);
352  h1->SetFillColor(kRed);
353  h1->SetMarkerStyle(20);
354  h1->SetMarkerSize(1);
355  h1->Draw("bar2");
356 
358 
359  canvas.Update();
360 
361  TLegend legend = TLegend(0.30, 0.88, 0.95, 0.94);
362  //legend.SetHeader(("#splitline{Payload hash: #bf{" + (std::get<1>(iov)) + "}}{Total Scenarios:"+std::to_string(scenarioIndex)+"}").c_str(),"C"); // option "C" allows to center the header
363 
364  legend.SetHeader(fmt::sprintf("Payload hash: #bf{%s}", std::get<1>(iov)).c_str(), "C");
365  legend.AddEntry(h1.get(), fmt::sprintf("total scenarios: #bf{%s}", std::to_string(scenarioIndex)).c_str(), "F");
366  legend.SetTextSize(0.025);
367  legend.Draw("same");
368 
369  auto ltx = TLatex();
370  ltx.SetTextFont(62);
371  //ltx.SetTextColor(kBlue);
372  //ltx.SetTextAlign(11);
373  ltx.SetTextSize(0.040);
374  ltx.DrawLatexNDC(
375  gPad->GetLeftMargin(),
376  1 - gPad->GetTopMargin() + 0.01,
377  fmt::sprintf("#color[4]{%s} IOV: #color[4]{%s}", tagname, std::to_string(std::get<0>(iov))).c_str());
378 
379  std::string fileName(m_imageFileName);
380  canvas.SaveAs(fileName.c_str());
381 
382  return true;
383 
384  } // fill
385  };
386 
387 } // namespace
388 
389 // Register the classes as boost python plugin
391  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerMap);
392  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerMap);
393  PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerMap);
394  PAYLOAD_INSPECTOR_CLASS(SiPixelFEDChannelContainerScenarios);
395 }
const TimeTypeSpecs timeTypeSpecs[]
Definition: Time.cc:16
Time_t beginValue
Definition: Time.h:41
void start(bool readOnly=true)
Definition: Session.cc:18
const Time_t MIN_VAL(0)
std::string to_string(const V &value)
Definition: OMSAccess.h:71
std::unique_ptr< T > fetchPayload(const cond::Hash &payloadHash)
Definition: Session.h:213
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
Transaction & transaction()
Definition: Session.cc:52
scenario
Definition: constants.h:173
identify pixel inside single ROC
Definition: LocalPixel.h:7
static std::string const input
Definition: EdmProvDump.cc:47
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
std::pair< unsigned int, unsigned int > unpack(cond::Time_t since)
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:55
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
IOVProxy readIov(const std::string &tag)
Definition: Session.cc:63
const std::map< std::string, std::string > & inputParamValues() const
void addInputParam(const std::string &paramName)
Definition: DetId.h:17
void makeNicePlotStyle(TH1 *hist)
row and collumn in ROC representation
Definition: LocalPixel.h:13
IOVArray selectRange(const cond::Time_t &begin, const cond::Time_t &end)
Definition: IOVProxy.cc:197
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
Log< level::Warning, false > LogWarning
Time_t endValue
Definition: Time.h:42
const Time_t MAX_VAL(std::numeric_limits< Time_t >::max())
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:37
#define LogDebug(id)
unsigned transform(const HcalDetId &id, unsigned transformCode)