CMS 3D CMS Logo

SiStripHitEfficiencyHarvester.cc
Go to the documentation of this file.
1 // user includes
24 
25 //system includes
26 #include <sstream>
27 #include <numeric> // for std::accumulate
28 
29 // ROOT includes
30 #include "TGraphAsymmErrors.h"
31 #include "TCanvas.h"
32 #include "TLegend.h"
33 #include "TStyle.h"
34 #include "TEfficiency.h"
35 #include "TTree.h"
36 
37 // custom made printout
38 #define LOGPRINT edm::LogPrint("SiStripHitEfficiencyHarvester")
39 
41 public:
43  ~SiStripHitEfficiencyHarvester() override = default;
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46  void endRun(edm::Run const&, edm::EventSetup const&) override;
48 
49 private:
52  const bool isAtPCL_;
54  const bool doStoreOnTree_;
56  const unsigned int nTEClayers_;
57  const double threshold_;
58  const int nModsMin_;
59  const float effPlotMin_;
60  const double tkMapMin_;
62 
67 
68  std::unique_ptr<TrackerTopology> tTopo_;
69  std::unique_ptr<TkDetMap> tkDetMap_;
70  std::unique_ptr<SiStripQuality> stripQuality_;
71  std::vector<DetId> stripDetIds_;
72 
73  int goodlayertotal[bounds::k_END_OF_LAYS_AND_RINGS];
74  int goodlayerfound[bounds::k_END_OF_LAYS_AND_RINGS];
75  int alllayertotal[bounds::k_END_OF_LAYS_AND_RINGS];
76  int alllayerfound[bounds::k_END_OF_LAYS_AND_RINGS];
77 
78  // information for the TTree
79  TTree* tree;
80  unsigned int t_DetId, t_found, t_total;
81  unsigned char t_layer;
83  float t_threshold;
84 
85  void writeBadStripPayload(const SiStripQuality& quality) const;
86  void printTotalStatistics(const std::array<long, bounds::k_END_OF_LAYERS>& layerFound,
87  const std::array<long, bounds::k_END_OF_LAYERS>& layerTotal) const;
88  void printAndWriteBadModules(const SiStripQuality& quality, const SiStripDetInfo& detInfo) const;
89  bool checkMapsValidity(const std::vector<MonitorElement*>& maps, const std::string& type) const;
90  unsigned int countTotalHits(const std::vector<MonitorElement*>& maps); /* to check if TK was ON */
91  void makeSummary(DQMStore::IGetter& getter, TFileService& fs) const;
92  void makeSummaryVsBX(DQMStore::IGetter& getter, TFileService& fs) const;
93  void makeSummaryVsLumi(DQMStore::IGetter& getter) const;
94  void makeSummaryVsCM(DQMStore::IGetter& getter, TFileService& fs) const;
95 };
96 
98  : inputFolder_(conf.getParameter<std::string>("inputFolder")),
99  isAtPCL_(conf.getParameter<bool>("isAtPCL")),
100  autoIneffModTagging_(conf.getUntrackedParameter<bool>("AutoIneffModTagging", false)),
101  doStoreOnDB_(conf.getParameter<bool>("doStoreOnDB")),
102  doStoreOnTree_(conf.getUntrackedParameter<bool>("doStoreOnTree")),
103  showRings_(conf.getUntrackedParameter<bool>("ShowRings", false)),
104  showEndcapSides_(conf.getUntrackedParameter<bool>("ShowEndcapSides", true)),
105  showTOB6TEC9_(conf.getUntrackedParameter<bool>("ShowTOB6TEC9", false)),
106  showOnlyGoodModules_(conf.getUntrackedParameter<bool>("ShowOnlyGoodModules", false)),
107  nTEClayers_(showRings_ ? 7 : 9), // number of rings or wheels
108  threshold_(conf.getParameter<double>("Threshold")),
109  nModsMin_(conf.getParameter<int>("nModsMin")),
110  effPlotMin_(conf.getUntrackedParameter<double>("EffPlotMin", 0.9)),
111  tkMapMin_(conf.getUntrackedParameter<double>("TkMapMin", 0.9)),
112  title_(conf.getParameter<std::string>("Title")),
113  record_(conf.getParameter<std::string>("Record")),
114  tTopoToken_(esConsumes<edm::Transition::EndRun>()),
115  tkDetMapToken_(esConsumes<edm::Transition::EndRun>()),
116  stripQualityToken_(esConsumes<edm::Transition::EndRun>()),
117  tkGeomToken_(esConsumes<edm::Transition::EndRun>()) {
118  // zero in all counts
119  for (int l = 0; l < bounds::k_END_OF_LAYS_AND_RINGS; l++) {
120  goodlayertotal[l] = 0;
121  goodlayerfound[l] = 0;
122  alllayertotal[l] = 0;
123  alllayerfound[l] = 0;
124  }
125 }
126 
128  if (!tTopo_) {
129  tTopo_ = std::make_unique<TrackerTopology>(iSetup.getData(tTopoToken_));
130  }
131  if (!tkDetMap_) {
132  tkDetMap_ = std::make_unique<TkDetMap>(iSetup.getData(tkDetMapToken_));
133  }
134  if (!stripQuality_) {
135  stripQuality_ = std::make_unique<SiStripQuality>(iSetup.getData(stripQualityToken_));
136  }
137  if (stripDetIds_.empty()) {
138  const auto& tkGeom = iSetup.getData(tkGeomToken_);
139  for (const auto& det : tkGeom.detUnits()) {
140  if (dynamic_cast<const StripGeomDetUnit*>(det)) {
141  stripDetIds_.push_back(det->geographicalId());
142  }
143  }
144  }
145 }
146 
147 bool SiStripHitEfficiencyHarvester::checkMapsValidity(const std::vector<MonitorElement*>& maps,
148  const std::string& type) const {
149  std::vector<bool> isThere;
150  isThere.reserve(maps.size());
151  std::transform(maps.begin() + 1, maps.end(), std::back_inserter(isThere), [](auto& x) { return !(x == nullptr); });
152 
153  int count{0};
154  for (const auto& it : isThere) {
155  count++;
156  LogDebug("SiStripHitEfficiencyHarvester") << " layer: " << count << " " << it << std::endl;
157  if (it)
158  LogDebug("SiStripHitEfficiencyHarvester") << "resolving to " << maps[count]->getName() << std::endl;
159  }
160 
161  // check on the input TkHistoMap
162  bool areMapsAvailable{true};
163  int layerCount{0};
164  for (const auto& it : isThere) {
165  layerCount++;
166  if (!it) {
167  edm::LogError("SiStripHitEfficiencyHarvester")
168  << type << " TkHistoMap for layer " << layerCount << " was not found.\n -> Aborting!";
169  areMapsAvailable = false;
170  break;
171  }
172  }
173  return areMapsAvailable;
174 }
175 
176 unsigned int SiStripHitEfficiencyHarvester::countTotalHits(const std::vector<MonitorElement*>& maps) {
177  return std::accumulate(maps.begin() + 1, maps.end(), 0, [](unsigned int total, MonitorElement* item) {
178  return total + item->getEntries();
179  });
180 }
181 
183  if (!isAtPCL_) {
185 
186  if (doStoreOnTree_) {
187  // store information per DetId in the output tree
188  tree = fs->make<TTree>("ModEff", "ModEff");
189  tree->Branch("DetId", &t_DetId, "DetId/i");
190  tree->Branch("Layer", &t_layer, "Layer/b");
191  tree->Branch("FoundHits", &t_found, "FoundHits/i");
192  tree->Branch("AllHits", &t_total, "AllHits/i");
193  tree->Branch("IsTaggedIneff", &t_isTaggedIneff, "IsTaggedIneff/O");
194  tree->Branch("TagThreshold", &t_threshold, "TagThreshold/F");
195  }
196  }
197 
199  LOGPRINT << "A module is bad if efficiency < " << threshold_ << " and has at least " << nModsMin_ << " nModsMin.";
200  else
201  LOGPRINT << "A module is bad if the upper limit on the efficiency is < to the avg in the layer - " << threshold_
202  << " and has at least " << nModsMin_ << " nModsMin.";
203 
204  auto h_module_total = std::make_unique<TkHistoMap>(tkDetMap_.get());
205  h_module_total->loadTkHistoMap(fmt::format("{}/TkDetMaps", inputFolder_), "perModule_total");
206  auto h_module_found = std::make_unique<TkHistoMap>(tkDetMap_.get());
207  h_module_found->loadTkHistoMap(fmt::format("{}/TkDetMaps", inputFolder_), "perModule_found");
208 
209  // collect how many layers are missing
210  const auto& totalMaps = h_module_total->getAllMaps();
211  const auto& foundMaps = h_module_found->getAllMaps();
212 
213  LogDebug("SiStripHitEfficiencyHarvester")
214  << "totalMaps.size(): " << totalMaps.size() << " foundMaps.size() " << foundMaps.size() << std::endl;
215 
216  // check on the input TkHistoMaps
217  bool isTotalMapAvailable = this->checkMapsValidity(totalMaps, std::string("Total"));
218  bool isFoundMapAvailable = this->checkMapsValidity(foundMaps, std::string("Found"));
219 
220  LogDebug("SiStripHitEfficiencyHarvester")
221  << "isTotalMapAvailable: " << isTotalMapAvailable << " isFoundMapAvailable " << isFoundMapAvailable << std::endl;
222 
223  // no input TkHistoMaps -> early return
224  if (!isTotalMapAvailable or !isFoundMapAvailable)
225  return;
226 
227  LogDebug("SiStripHitEfficiencyHarvester")
228  << "Entries in total TkHistoMap for layer 3: " << h_module_total->getMap(3)->getEntries() << ", found "
229  << h_module_found->getMap(3)->getEntries();
230 
231  // count how many hits in the denominator we have
232  const unsigned int totalHits = this->countTotalHits(totalMaps);
233 
234  // set colz
235  for (size_t i = 1; i < totalMaps.size(); i++) {
236  h_module_total->getMap(i)->setOption("colz");
237  h_module_found->getMap(i)->setOption("colz");
238  }
239 
240  // come back to the main folder
242 
243  std::vector<MonitorElement*> hEffInLayer(std::size_t(1), nullptr);
244  hEffInLayer.reserve(bounds::k_END_OF_LAYERS);
245  for (std::size_t i = 1; i != bounds::k_END_OF_LAYERS; ++i) {
246  const auto lyrName = ::layerName(i, showRings_, nTEClayers_);
247  hEffInLayer.push_back(booker.book1D(
248  Form("eff_layer%i", int(i)), Form("Module efficiency in layer %s", lyrName.c_str()), 201, 0, 1.005));
249  }
250  std::array<long, bounds::k_END_OF_LAYERS> layerTotal{};
251  std::array<long, bounds::k_END_OF_LAYERS> layerFound{};
252  layerTotal.fill(0);
253  layerFound.fill(0);
254 
256  // Effiency calculation, bad module tagging, and tracker maps //
258 
259  TrackerMap tkMap{" Detector Inefficiency "};
260  TrackerMap tkMapBad{" Inefficient Modules "};
261  TrackerMap tkMapEff{title_};
262  TrackerMap tkMapNum{" Detector numerator "};
263  TrackerMap tkMapDen{" Detector denominator "};
264  std::map<unsigned int, double> badModules;
265 
266  // load the FEDError map
267  const auto& EventStats = getter.get(fmt::format("{}/EventInfo/EventStats", inputFolder_));
268  const int totalEvents = EventStats->getBinContent(1., 1.); // first bin contains info on number of events run
269  calibData_.FEDErrorOccupancy = std::make_unique<TkHistoMap>(tkDetMap_.get());
270  calibData_.FEDErrorOccupancy->loadTkHistoMap(fmt::format("{}/FEDErrorTkDetMaps", inputFolder_),
271  "perModule_FEDErrors");
272 
273  // tag as bad from FEDErrors the modules that have an error on 75% of the events
274  calibData_.fillMapFromTkMap(totalEvents, 0.75, stripDetIds_);
275 
276  for (const auto& [badId, fraction] : calibData_.fedErrorCounts) {
277  LogDebug("SiStripHitEfficiencyHarvester")
278  << __PRETTY_FUNCTION__ << " bad module from FEDError " << badId << "," << fraction << std::endl;
279  }
280 
281  for (auto det : stripDetIds_) {
282  auto layer = ::checkLayer(det, tTopo_.get());
283  const auto num = h_module_found->getValue(det);
284  const auto denom = h_module_total->getValue(det);
285  if (denom) {
286  // use only the "good" modules
287  if (stripQuality_->getBadApvs(det) == 0 && calibData_.checkFedError(det)) {
288  const auto eff = num / denom;
289  hEffInLayer[layer]->Fill(eff);
290  if (!autoIneffModTagging_) {
291  if ((denom >= nModsMin_) && (eff < threshold_)) {
292  // We have a bad module, put it in the list!
293  badModules[det] = eff;
294  tkMapBad.fillc(det, 255, 0, 0);
295  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
296  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom;
297  } else {
298  //Fill the bad list with empty results for every module
299  tkMapBad.fillc(det, 255, 255, 255);
300  }
301  if (eff < threshold_)
302  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
303  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom;
304 
305  if (denom < nModsMin_) {
306  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
307  << det.rawId() << " is under occupancy at " << denom;
308  }
309 
310  if (doStoreOnTree_ && !isAtPCL_) {
311  t_DetId = det.rawId();
312  t_layer = layer;
313  t_found = num;
314  t_total = denom;
315  t_isTaggedIneff = false;
316  t_threshold = 0;
317  tree->Fill();
318  }
319  }
320 
321  //Put any module into the TKMap
322  tkMap.fill(det, 1. - eff);
323  tkMapEff.fill(det, eff);
324  tkMapNum.fill(det, num);
325  tkMapDen.fill(det, denom);
326 
327  layerTotal[layer] += denom;
328  layerFound[layer] += num;
329 
330  // for the summary
331  // Have to do the decoding for which side to go on (ugh)
332  if (layer <= bounds::k_LayersAtTOBEnd) {
335  } else if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
336  if (tTopo_->tidSide(det) == 1) {
339  } else if (tTopo_->tidSide(det) == 2) {
340  goodlayerfound[layer + 3] += num;
341  goodlayertotal[layer + 3] += denom;
342  }
343  } else if (layer > bounds::k_LayersAtTIDEnd && layer <= bounds::k_LayersAtTECEnd) {
344  if (tTopo_->tecSide(det) == 1) {
345  goodlayerfound[layer + 3] += num;
346  goodlayertotal[layer + 3] += denom;
347  } else if (tTopo_->tecSide(det) == 2) {
350  }
351  }
352  } // if the module is good!
353 
354  //Do the one where we don't exclude bad modules!
355  if (layer <= bounds::k_LayersAtTOBEnd) {
356  alllayerfound[layer] += num;
358  } else if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
359  if (tTopo_->tidSide(det) == 1) {
360  alllayerfound[layer] += num;
362  } else if (tTopo_->tidSide(det) == 2) {
363  alllayerfound[layer + 3] += num;
364  alllayertotal[layer + 3] += denom;
365  }
366  } else if (layer > bounds::k_LayersAtTIDEnd && layer <= bounds::k_LayersAtTECEnd) {
367  if (tTopo_->tecSide(det) == 1) {
368  alllayerfound[layer + 3] += num;
369  alllayertotal[layer + 3] += denom;
370  } else if (tTopo_->tecSide(det) == 2) {
373  }
374  }
375 
376  } // if denom
377  } // loop on DetIds
378 
379  if (autoIneffModTagging_) {
380  for (Long_t i = 1; i <= k_LayersAtTECEnd; i++) {
381  //Compute threshold to use for each layer
382  hEffInLayer[i]->getTH1()->GetXaxis()->SetRange(
383  3, hEffInLayer[i]->getNbinsX() + 1); // Remove from the avg modules below 1%
384  const double layer_min_eff = hEffInLayer[i]->getMean() - std::max(2.5 * hEffInLayer[i]->getRMS(), threshold_);
385  LOGPRINT << "Layer " << i << " threshold for bad modules: <" << layer_min_eff
386  << " (layer mean: " << hEffInLayer[i]->getMean() << " rms: " << hEffInLayer[i]->getRMS() << ")";
387 
388  hEffInLayer[i]->getTH1()->GetXaxis()->SetRange(1, hEffInLayer[i]->getNbinsX() + 1);
389 
390  for (auto det : stripDetIds_) {
391  // use only the "good" modules
392  if (stripQuality_->getBadApvs(det) == 0 && calibData_.checkFedError(det)) {
393  const auto layer = ::checkLayer(det, tTopo_.get());
394  if (layer == i) {
395  const auto num = h_module_found->getValue(det);
396  const auto denom = h_module_total->getValue(det);
397  if (denom) {
398  const auto eff = num / denom;
399  const auto eff_up = TEfficiency::Bayesian(denom, num, .99, 1, 1, true);
400 
401  if ((denom >= nModsMin_) && (eff_up < layer_min_eff)) {
402  //We have a bad module, put it in the list!
403  badModules[det] = eff;
404  tkMapBad.fillc(det, 255, 0, 0);
405  if (!isAtPCL_ && doStoreOnTree_) {
406  t_isTaggedIneff = true;
407  }
408  } else {
409  //Fill the bad list with empty results for every module
410  tkMapBad.fillc(det, 255, 255, 255);
411  if (!isAtPCL_ && doStoreOnTree_) {
412  t_isTaggedIneff = false;
413  }
414  }
415  if (eff_up < layer_min_eff + 0.08) {
416  // printing message also for modules sligthly above (8%) the limit
417  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
418  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom
419  << " , upper limit: " << eff_up;
420  }
421  if (denom < nModsMin_) {
422  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
423  << det.rawId() << " layer " << layer << " is under occupancy at " << denom;
424  }
425 
426  if (!isAtPCL_ && doStoreOnTree_) {
427  t_DetId = det.rawId();
428  t_layer = layer;
429  t_found = num;
430  t_total = denom;
431  t_threshold = layer_min_eff;
432  tree->Fill();
433  } // if storing tree
434  } // if denom
435  } // layer = i
436  } // if there are no bad APVs
437  } // loop on detids
438  } // loop on layers
439  } // if auto tagging
440 
441  tkMap.save(true, 0, 0, "SiStripHitEffTKMap_NEW.png");
442  tkMapBad.save(true, 0, 0, "SiStripHitEffTKMapBad_NEW.png");
443  tkMapEff.save(true, tkMapMin_, 1., "SiStripHitEffTKMapEff_NEW.png");
444  tkMapNum.save(true, 0, 0, "SiStripHitEffTKMapNum_NEW.png");
445  tkMapDen.save(true, 0, 0, "SiStripHitEffTKMapDen_NEW.png");
446 
447  const auto detInfo =
449  SiStripQuality pQuality{detInfo};
450  //This is the list of the bad strips, use to mask out entire APVs
451  //Now simply go through the bad hit list and mask out things that
452  //are bad!
453  for (const auto it : badModules) {
454  const auto det = it.first;
455  std::vector<unsigned int> badStripList;
456  //We need to figure out how many strips are in this particular module
457  //To Mask correctly!
458  const auto nStrips = detInfo.getNumberOfApvsAndStripLength(det).first * sistrip::STRIPS_PER_APV;
459  LOGPRINT << "Number of strips module " << det << " is " << nStrips;
460  badStripList.push_back(pQuality.encode(0, nStrips, 0));
461  //Now compact into a single bad module
462  LOGPRINT << "ID1 shoudl match list of modules above " << det;
463  pQuality.compact(det, badStripList);
464  pQuality.put(det, SiStripQuality::Range(badStripList.begin(), badStripList.end()));
465  }
466  pQuality.fillBadComponents();
467  if (doStoreOnDB_) {
468  if (totalHits > 0u) {
469  writeBadStripPayload(pQuality);
470  } else {
471  edm::LogPrint("SiStripHitEfficiencyHarvester")
472  << __PRETTY_FUNCTION__ << " There are no SiStrip hits for a valid measurement, skipping!";
473  }
474  } else {
475  edm::LogInfo("SiStripHitEfficiencyHarvester") << "Will not produce payload!";
476  }
477 
478  printTotalStatistics(layerFound, layerTotal); // statistics by layer and subdetector
479  //LOGPRINT << "\n-----------------\nNew IOV starting from run " << e.id().run() << " event " << e.id().event()
480  // << " lumiBlock " << e.luminosityBlock() << " time " << e.time().value() << "\n-----------------\n";
481  printAndWriteBadModules(pQuality, detInfo); // TODO
482 
483  if (!isAtPCL_) {
485  makeSummary(getter, *fs); // TODO
486  makeSummaryVsBX(getter, *fs); // TODO
487  makeSummaryVsCM(getter, *fs); // TODO
488  }
489 
490  makeSummaryVsLumi(getter); // TODO
491 }
492 
494  const std::array<long, bounds::k_END_OF_LAYERS>& layerFound,
495  const std::array<long, bounds::k_END_OF_LAYERS>& layerTotal) const {
496  //Calculate the statistics by layer
497  int totalfound = 0;
498  int totaltotal = 0;
499  double layereff;
500  int subdetfound[5];
501  int subdettotal[5];
502 
503  for (Long_t i = 1; i < 5; i++) {
504  subdetfound[i] = 0;
505  subdettotal[i] = 0;
506  }
507 
508  for (Long_t i = 1; i <= bounds::k_LayersAtTECEnd; i++) {
509  layereff = double(layerFound[i]) / double(layerTotal[i]);
510  LOGPRINT << "Layer " << i << " (" << ::layerName(i, showRings_, nTEClayers_) << ") has total efficiency "
511  << layereff << " " << layerFound[i] << "/" << layerTotal[i];
512  totalfound += layerFound[i];
513  totaltotal += layerTotal[i];
514  if (i <= bounds::k_LayersAtTIBEnd) {
515  subdetfound[1] += layerFound[i];
516  subdettotal[1] += layerTotal[i];
517  }
518  if (i > bounds::k_LayersAtTIBEnd && i <= bounds::k_LayersAtTOBEnd) {
519  subdetfound[2] += layerFound[i];
520  subdettotal[2] += layerTotal[i];
521  }
522  if (i > bounds::k_LayersAtTOBEnd && i <= bounds::k_LayersAtTIDEnd) {
523  subdetfound[3] += layerFound[i];
524  subdettotal[3] += layerTotal[i];
525  }
526  if (i > bounds::k_LayersAtTIDEnd) {
527  subdetfound[4] += layerFound[i];
528  subdettotal[4] += layerTotal[i];
529  }
530  }
531 
532  LOGPRINT << "The total efficiency is " << double(totalfound) / double(totaltotal);
533  LOGPRINT << " TIB: " << double(subdetfound[1]) / subdettotal[1] << " " << subdetfound[1] << "/"
534  << subdettotal[1];
535  LOGPRINT << " TOB: " << double(subdetfound[2]) / subdettotal[2] << " " << subdetfound[2] << "/"
536  << subdettotal[2];
537  LOGPRINT << " TID: " << double(subdetfound[3]) / subdettotal[3] << " " << subdetfound[3] << "/"
538  << subdettotal[3];
539  LOGPRINT << " TEC: " << double(subdetfound[4]) / subdettotal[4] << " " << subdetfound[4] << "/"
540  << subdettotal[4];
541 }
542 
544  SiStripBadStrip pBadStrip{};
545  const auto pQdvBegin = quality.getDataVectorBegin();
546  for (auto rIt = quality.getRegistryVectorBegin(); rIt != quality.getRegistryVectorEnd(); ++rIt) {
547  const auto range = SiStripBadStrip::Range(pQdvBegin + rIt->ibegin, pQdvBegin + rIt->iend);
548  if (!pBadStrip.put(rIt->detid, range))
549  edm::LogError("SiStripHitEfficiencyHarvester") << "detid already exists in SiStripBadStrip";
550  }
552  if (poolDbService.isAvailable()) {
553  poolDbService->writeOneIOV(pBadStrip, poolDbService->currentTime(), record_);
554  } else {
555  throw cms::Exception("PoolDBService required");
556  }
557 }
558 
560  // use goodlayer_total/found and alllayer_total/found, collapse side and/or ring if needed
561 
562  unsigned int nLayers = 34;
563  if (showRings_)
564  nLayers = 30;
565  if (!showEndcapSides_) {
566  if (!showRings_)
567  nLayers = 22;
568  else
569  nLayers = 20;
570  }
571 
572  TH1F* found = fs.make<TH1F>("found", "found", nLayers + 1, 0, nLayers + 1);
573  TH1F* all = fs.make<TH1F>("all", "all", nLayers + 1, 0, nLayers + 1);
574  TH1F* found2 = fs.make<TH1F>("found2", "found2", nLayers + 1, 0, nLayers + 1);
575  TH1F* all2 = fs.make<TH1F>("all2", "all2", nLayers + 1, 0, nLayers + 1);
576  // first bin only to keep real data off the y axis so set to -1
577  found->SetBinContent(0, -1);
578  all->SetBinContent(0, 1);
579 
580  // new ROOT version: TGraph::Divide don't handle null or negative values
581  for (unsigned int i = 1; i < nLayers + 2; ++i) {
582  found->SetBinContent(i, 1e-6);
583  all->SetBinContent(i, 1);
584  found2->SetBinContent(i, 1e-6);
585  all2->SetBinContent(i, 1);
586  }
587 
588  TCanvas* c7 = new TCanvas("c7", " test ", 10, 10, 800, 600);
589  c7->SetFillColor(0);
590  c7->SetGrid();
591 
592  unsigned int nLayers_max = nLayers + 1; // barrel+endcap
593  if (!showEndcapSides_)
594  nLayers_max = 11; // barrel
595  for (unsigned int i = 1; i < nLayers_max; ++i) {
596  LOGPRINT << "Fill only good modules layer " << i << ": S = " << goodlayerfound[i]
597  << " B = " << goodlayertotal[i];
598  if (goodlayertotal[i] > 5) {
599  found->SetBinContent(i, goodlayerfound[i]);
600  all->SetBinContent(i, goodlayertotal[i]);
601  }
602 
603  LOGPRINT << "Filling all modules layer " << i << ": S = " << alllayerfound[i] << " B = " << alllayertotal[i];
604  if (alllayertotal[i] > 5) {
605  found2->SetBinContent(i, alllayerfound[i]);
606  all2->SetBinContent(i, alllayertotal[i]);
607  }
608  }
609 
610  // endcap - merging sides
611  if (!showEndcapSides_) {
612  for (unsigned int i = 11; i < 14; ++i) { // TID disks
613  LOGPRINT << "Fill only good modules layer " << i << ": S = " << goodlayerfound[i] + goodlayerfound[i + 3]
614  << " B = " << goodlayertotal[i] + goodlayertotal[i + 3];
615  if (goodlayertotal[i] + goodlayertotal[i + 3] > 5) {
616  found->SetBinContent(i, goodlayerfound[i] + goodlayerfound[i + 3]);
617  all->SetBinContent(i, goodlayertotal[i] + goodlayertotal[i + 3]);
618  }
619  LOGPRINT << "Filling all modules layer " << i << ": S = " << alllayerfound[i] + alllayerfound[i + 3]
620  << " B = " << alllayertotal[i] + alllayertotal[i + 3];
621  if (alllayertotal[i] + alllayertotal[i + 3] > 5) {
622  found2->SetBinContent(i, alllayerfound[i] + alllayerfound[i + 3]);
623  all2->SetBinContent(i, alllayertotal[i] + alllayertotal[i + 3]);
624  }
625  }
626  for (unsigned int i = 17; i < 17 + nTEClayers_; ++i) { // TEC disks
627  LOGPRINT << "Fill only good modules layer " << i - 3
628  << ": S = " << goodlayerfound[i] + goodlayerfound[i + nTEClayers_]
629  << " B = " << goodlayertotal[i] + goodlayertotal[i + nTEClayers_];
630  if (goodlayertotal[i] + goodlayertotal[i + nTEClayers_] > 5) {
631  found->SetBinContent(i - 3, goodlayerfound[i] + goodlayerfound[i + nTEClayers_]);
632  all->SetBinContent(i - 3, goodlayertotal[i] + goodlayertotal[i + nTEClayers_]);
633  }
634  LOGPRINT << "Filling all modules layer " << i - 3
635  << ": S = " << alllayerfound[i] + alllayerfound[i + nTEClayers_]
636  << " B = " << alllayertotal[i] + alllayertotal[i + nTEClayers_];
637  if (alllayertotal[i] + alllayertotal[i + nTEClayers_] > 5) {
638  found2->SetBinContent(i - 3, alllayerfound[i] + alllayerfound[i + nTEClayers_]);
639  all2->SetBinContent(i - 3, alllayertotal[i] + alllayertotal[i + nTEClayers_]);
640  }
641  }
642  }
643 
644  found->Sumw2();
645  all->Sumw2();
646 
647  found2->Sumw2();
648  all2->Sumw2();
649 
650  TGraphAsymmErrors* gr = fs.make<TGraphAsymmErrors>(nLayers + 1);
651  gr->SetName("eff_good");
652  gr->BayesDivide(found, all);
653 
654  TGraphAsymmErrors* gr2 = fs.make<TGraphAsymmErrors>(nLayers + 1);
655  gr2->SetName("eff_all");
656  gr2->BayesDivide(found2, all2);
657 
658  for (unsigned int j = 0; j < nLayers + 1; j++) {
659  gr->SetPointError(j, 0., 0., gr->GetErrorYlow(j), gr->GetErrorYhigh(j));
660  gr2->SetPointError(j, 0., 0., gr2->GetErrorYlow(j), gr2->GetErrorYhigh(j));
661  }
662 
663  gr->GetXaxis()->SetLimits(0, nLayers);
664  gr->SetMarkerColor(2);
665  gr->SetMarkerSize(1.2);
666  gr->SetLineColor(2);
667  gr->SetLineWidth(4);
668  gr->SetMarkerStyle(20);
669  gr->SetMinimum(effPlotMin_);
670  gr->SetMaximum(1.001);
671  gr->GetYaxis()->SetTitle("Efficiency");
672  gStyle->SetTitleFillColor(0);
673  gStyle->SetTitleBorderSize(0);
674  gr->SetTitle(title_.c_str());
675 
676  gr2->GetXaxis()->SetLimits(0, nLayers);
677  gr2->SetMarkerColor(1);
678  gr2->SetMarkerSize(1.2);
679  gr2->SetLineColor(1);
680  gr2->SetLineWidth(4);
681  gr2->SetMarkerStyle(21);
682  gr2->SetMinimum(effPlotMin_);
683  gr2->SetMaximum(1.001);
684  gr2->GetYaxis()->SetTitle("Efficiency");
685  gr2->SetTitle(title_.c_str());
686 
687  for (Long_t k = 1; k < nLayers + 1; k++) {
688  TString label;
689  if (showEndcapSides_)
690  label = ::layerSideName(k, showRings_, nTEClayers_);
691  else
693  if (!showTOB6TEC9_) {
694  if (k == 10)
695  label = "";
696  if (!showRings_ && k == nLayers)
697  label = "";
698  if (!showRings_ && showEndcapSides_ && k == 25)
699  label = "";
700  }
701  if (!showRings_) {
702  if (showEndcapSides_) {
703  gr->GetXaxis()->SetBinLabel(((k + 1) * 100 + 2) / (nLayers)-4, label);
704  gr2->GetXaxis()->SetBinLabel(((k + 1) * 100 + 2) / (nLayers)-4, label);
705  } else {
706  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-6, label);
707  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-6, label);
708  }
709  } else {
710  if (showEndcapSides_) {
711  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-4, label);
712  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-4, label);
713  } else {
714  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-7, label);
715  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-7, label);
716  }
717  }
718  }
719 
720  gr->Draw("AP");
721  gr->GetXaxis()->SetNdivisions(36);
722 
723  c7->cd();
724  TPad* overlay = new TPad("overlay", "", 0, 0, 1, 1);
725  overlay->SetFillStyle(4000);
726  overlay->SetFillColor(0);
727  overlay->SetFrameFillStyle(4000);
728  overlay->Draw("same");
729  overlay->cd();
731  gr2->Draw("AP");
732 
733  TLegend* leg = new TLegend(0.70, 0.27, 0.88, 0.40);
734  leg->AddEntry(gr, "Good Modules", "p");
736  leg->AddEntry(gr2, "All Modules", "p");
737  leg->SetTextSize(0.020);
738  leg->SetFillColor(0);
739  leg->Draw("same");
740 
741  c7->SaveAs("Summary.png");
742  c7->SaveAs("Summary.root");
743 }
744 
746  // use found/totalVsBx_layer%i [0,23)
747 }
748 
750  for (unsigned int iLayer = 1; iLayer != (showRings_ ? 20 : 22); ++iLayer) {
751  auto hfound = getter.get(fmt::format("{}/VsLumi/layerfound_vsLumi_layer_{}", inputFolder_, iLayer))->getTH1F();
752  auto htotal = getter.get(fmt::format("{}/VsLumi/layertotal_vsLumi_layer_{}", inputFolder_, iLayer))->getTH1F();
753 
754  if (hfound == nullptr or htotal == nullptr) {
755  if (hfound == nullptr)
756  edm::LogError("SiStripHitEfficiencyHarvester")
757  << fmt::format("{}/VsLumi/layerfound_vsLumi_layer_{}", inputFolder_, iLayer) << " was not found!";
758  if (htotal == nullptr)
759  edm::LogError("SiStripHitEfficiencyHarvester")
760  << fmt::format("{}/VsLumi/layertotal_vsLumi_layer_{}", inputFolder_, iLayer) << " was not found!";
761  // no input histograms -> continue in the loop
762  continue;
763  }
764 
765  if (!hfound->GetSumw2())
766  hfound->Sumw2();
767  if (!htotal->GetSumw2())
768  htotal->Sumw2();
769  for (Long_t i = 0; i != hfound->GetNbinsX() + 1; ++i) {
770  if (hfound->GetBinContent(i) == 0)
771  hfound->SetBinContent(i, 1e-6);
772  if (htotal->GetBinContent(i) == 0)
773  htotal->SetBinContent(i, 1);
774  }
775  LogDebug("SiStripHitEfficiencyHarvester")
776  << "Total hits for layer " << iLayer << " (vs lumi): " << htotal->GetEntries() << ", found "
777  << hfound->GetEntries();
778  }
779  // continue
780 }
781 
783 
784 namespace {
785  void setBadComponents(int i,
786  int comp,
788  std::stringstream ssV[4][19],
789  int nBad[4][19][4],
790  int nAPV) {
791  ssV[i][comp] << "\n\t\t " << bc.detid << " \t " << bc.BadModule << " \t " << ((bc.BadFibers) & 0x1) << " ";
792  if (nAPV == 4)
793  ssV[i][comp] << "x " << ((bc.BadFibers >> 1) & 0x1);
794 
795  if (nAPV == 6)
796  ssV[i][comp] << ((bc.BadFibers >> 1) & 0x1) << " " << ((bc.BadFibers >> 2) & 0x1);
797  ssV[i][comp] << " \t " << ((bc.BadApvs) & 0x1) << " " << ((bc.BadApvs >> 1) & 0x1) << " ";
798  if (nAPV == 4)
799  ssV[i][comp] << "x x " << ((bc.BadApvs >> 2) & 0x1) << " " << ((bc.BadApvs >> 3) & 0x1);
800  if (nAPV == 6)
801  ssV[i][comp] << ((bc.BadApvs >> 2) & 0x1) << " " << ((bc.BadApvs >> 3) & 0x1) << " " << ((bc.BadApvs >> 4) & 0x1)
802  << " " << ((bc.BadApvs >> 5) & 0x1) << " ";
803 
804  if (bc.BadApvs) {
805  nBad[i][0][2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
806  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
807  nBad[i][comp][2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
808  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
809  }
810  if (bc.BadFibers) {
811  nBad[i][0][1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
812  nBad[i][comp][1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
813  }
814  if (bc.BadModule) {
815  nBad[i][0][0]++;
816  nBad[i][comp][0]++;
817  }
818  }
819 } // namespace
820 
822  const SiStripDetInfo& detInfo) const {
824  //try to write out what's in the quality record
826  int nTkBadComp[4]; //k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
827  int nBadComp[4][19][4];
828  //legend: nBadComp[i][j][k]= SubSystem i, layer/disk/wheel j, BadModule/Fiber/Apv k
829  // i: 0=TIB, 1=TID, 2=TOB, 3=TEC
830  // k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
831  std::stringstream ssV[4][19];
832 
833  for (int i = 0; i < 4; ++i) {
834  nTkBadComp[i] = 0;
835  for (int j = 0; j < 19; ++j) {
836  ssV[i][j].str("");
837  for (int k = 0; k < 4; ++k)
838  nBadComp[i][j][k] = 0;
839  }
840  }
841 
842  for (const auto& bc : quality.getBadComponentList()) {
843  // Full Tk
844  if (bc.BadModule)
845  nTkBadComp[0]++;
846  if (bc.BadFibers)
847  nTkBadComp[1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
848  if (bc.BadApvs)
849  nTkBadComp[2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
850  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
851  // single subsystem
852  DetId det(bc.detid);
853  if ((det.subdetId() >= SiStripSubdetector::TIB) && (det.subdetId() <= SiStripSubdetector::TEC)) {
854  const auto nAPV = detInfo.getNumberOfApvsAndStripLength(det).first;
855  switch (det.subdetId()) {
857  setBadComponents(0, tTopo_->tibLayer(det), bc, ssV, nBadComp, nAPV);
858  break;
861  (tTopo_->tidSide(det) == 2 ? tTopo_->tidWheel(det) : tTopo_->tidWheel(det) + 3),
862  bc,
863  ssV,
864  nBadComp,
865  nAPV);
866  break;
868  setBadComponents(2, tTopo_->tobLayer(det), bc, ssV, nBadComp, nAPV);
869  break;
872  (tTopo_->tecSide(det) == 2 ? tTopo_->tecWheel(det) : tTopo_->tecWheel(det) + 9),
873  bc,
874  ssV,
875  nBadComp,
876  nAPV);
877  break;
878  default:
879  break;
880  }
881  }
882  }
883  // single strip info
884  for (auto rp = quality.getRegistryVectorBegin(); rp != quality.getRegistryVectorEnd(); ++rp) {
885  DetId det{rp->detid};
886  int subdet = -999;
887  int component = -999;
888  switch (det.subdetId()) {
890  subdet = 0;
891  component = tTopo_->tibLayer(det);
892  break;
894  subdet = 1;
895  component = tTopo_->tidSide(det) == 2 ? tTopo_->tidWheel(det) : tTopo_->tidWheel(det) + 3;
896  break;
898  subdet = 2;
899  component = tTopo_->tobLayer(det);
900  break;
902  subdet = 3;
903  component = tTopo_->tecSide(det) == 2 ? tTopo_->tecWheel(det) : tTopo_->tecWheel(det) + 9;
904  break;
905  default:
906  break;
907  }
908 
909  const auto pQdvBegin = quality.getDataVectorBegin();
910  const auto sqrange = SiStripQuality::Range(pQdvBegin + rp->ibegin, pQdvBegin + rp->iend);
911  float percentage = 0;
912  for (int it = 0; it < sqrange.second - sqrange.first; it++) {
913  unsigned int range = quality.decode(*(sqrange.first + it)).range;
914  nTkBadComp[3] += range;
915  nBadComp[subdet][0][3] += range;
916  nBadComp[subdet][component][3] += range;
917  percentage += range;
918  }
919  if (percentage != 0)
920  percentage /= (sistrip::STRIPS_PER_APV * detInfo.getNumberOfApvsAndStripLength(det).first);
921  if (percentage > 1)
922  edm::LogError("SiStripHitEfficiencyHarvester") << "PROBLEM detid " << det.rawId() << " value " << percentage;
923  }
924 
925  // printout
926  std::ostringstream ss;
927  ss << "\n-----------------\nGlobal Info\n-----------------";
928  ss << "\nBadComp \t Modules \tFibers "
929  "\tApvs\tStrips\n----------------------------------------------------------------";
930  ss << "\nTracker:\t\t" << nTkBadComp[0] << "\t" << nTkBadComp[1] << "\t" << nTkBadComp[2] << "\t" << nTkBadComp[3];
931  ss << "\nTIB:\t\t\t" << nBadComp[0][0][0] << "\t" << nBadComp[0][0][1] << "\t" << nBadComp[0][0][2] << "\t"
932  << nBadComp[0][0][3];
933  ss << "\nTID:\t\t\t" << nBadComp[1][0][0] << "\t" << nBadComp[1][0][1] << "\t" << nBadComp[1][0][2] << "\t"
934  << nBadComp[1][0][3];
935  ss << "\nTOB:\t\t\t" << nBadComp[2][0][0] << "\t" << nBadComp[2][0][1] << "\t" << nBadComp[2][0][2] << "\t"
936  << nBadComp[2][0][3];
937  ss << "\nTEC:\t\t\t" << nBadComp[3][0][0] << "\t" << nBadComp[3][0][1] << "\t" << nBadComp[3][0][2] << "\t"
938  << nBadComp[3][0][3];
939  ss << "\n";
940 
941  for (int i = 1; i < 5; ++i)
942  ss << "\nTIB Layer " << i << " :\t\t" << nBadComp[0][i][0] << "\t" << nBadComp[0][i][1] << "\t" << nBadComp[0][i][2]
943  << "\t" << nBadComp[0][i][3];
944  ss << "\n";
945  for (int i = 1; i < 4; ++i)
946  ss << "\nTID+ Disk " << i << " :\t\t" << nBadComp[1][i][0] << "\t" << nBadComp[1][i][1] << "\t" << nBadComp[1][i][2]
947  << "\t" << nBadComp[1][i][3];
948  for (int i = 4; i < 7; ++i)
949  ss << "\nTID- Disk " << i - 3 << " :\t\t" << nBadComp[1][i][0] << "\t" << nBadComp[1][i][1] << "\t"
950  << nBadComp[1][i][2] << "\t" << nBadComp[1][i][3];
951  ss << "\n";
952  for (int i = 1; i < 7; ++i)
953  ss << "\nTOB Layer " << i << " :\t\t" << nBadComp[2][i][0] << "\t" << nBadComp[2][i][1] << "\t" << nBadComp[2][i][2]
954  << "\t" << nBadComp[2][i][3];
955  ss << "\n";
956  for (int i = 1; i < 10; ++i)
957  ss << "\nTEC+ Disk " << i << " :\t\t" << nBadComp[3][i][0] << "\t" << nBadComp[3][i][1] << "\t" << nBadComp[3][i][2]
958  << "\t" << nBadComp[3][i][3];
959  for (int i = 10; i < 19; ++i)
960  ss << "\nTEC- Disk " << i - 9 << " :\t\t" << nBadComp[3][i][0] << "\t" << nBadComp[3][i][1] << "\t"
961  << nBadComp[3][i][2] << "\t" << nBadComp[3][i][3];
962  ss << "\n";
963 
964  ss << "\n----------------------------------------------------------------\n\t\t Detid \tModules Fibers "
965  "Apvs\n----------------------------------------------------------------";
966  for (int i = 1; i < 5; ++i)
967  ss << "\nTIB Layer " << i << " :" << ssV[0][i].str();
968  ss << "\n";
969  for (int i = 1; i < 4; ++i)
970  ss << "\nTID+ Disk " << i << " :" << ssV[1][i].str();
971  for (int i = 4; i < 7; ++i)
972  ss << "\nTID- Disk " << i - 3 << " :" << ssV[1][i].str();
973  ss << "\n";
974  for (int i = 1; i < 7; ++i)
975  ss << "\nTOB Layer " << i << " :" << ssV[2][i].str();
976  ss << "\n";
977  for (int i = 1; i < 10; ++i)
978  ss << "\nTEC+ Disk " << i << " :" << ssV[3][i].str();
979  for (int i = 10; i < 19; ++i)
980  ss << "\nTEC- Disk " << i - 9 << " :" << ssV[3][i].str();
981 
982  LOGPRINT << ss.str();
983 
984  // store also bad modules in log file
985  std::ofstream badModules;
986  badModules.open("BadModules_NEW.log");
987  badModules << "\n----------------------------------------------------------------\n\t\t Detid \tModules Fibers "
988  "Apvs\n----------------------------------------------------------------";
989  for (int i = 1; i < 5; ++i)
990  badModules << "\nTIB Layer " << i << " :" << ssV[0][i].str();
991  badModules << "\n";
992  for (int i = 1; i < 4; ++i)
993  badModules << "\nTID+ Disk " << i << " :" << ssV[1][i].str();
994  for (int i = 4; i < 7; ++i)
995  badModules << "\nTID- Disk " << i - 3 << " :" << ssV[1][i].str();
996  badModules << "\n";
997  for (int i = 1; i < 7; ++i)
998  badModules << "\nTOB Layer " << i << " :" << ssV[2][i].str();
999  badModules << "\n";
1000  for (int i = 1; i < 10; ++i)
1001  badModules << "\nTEC+ Disk " << i << " :" << ssV[3][i].str();
1002  for (int i = 10; i < 19; ++i)
1003  badModules << "\nTEC- Disk " << i - 9 << " :" << ssV[3][i].str();
1004  badModules.close();
1005 }
1006 
1009  desc.add<std::string>("inputFolder", "AlCaReco/SiStripHitEfficiency");
1010  desc.add<bool>("isAtPCL", false);
1011  desc.add<bool>("doStoreOnDB", false);
1012  desc.add<std::string>("Record", "SiStripBadStrip");
1013  desc.add<double>("Threshold", 0.1);
1014  desc.add<std::string>("Title", "Hit Efficiency");
1015  desc.add<int>("nModsMin", 5);
1016  desc.addUntracked<bool>("doStoreOnTree", false);
1017  desc.addUntracked<bool>("AutoIneffModTagging", false);
1018  desc.addUntracked<double>("TkMapMin", 0.9);
1019  desc.addUntracked<double>("EffPlotMin", 0.9);
1020  desc.addUntracked<bool>("ShowRings", false);
1021  desc.addUntracked<bool>("ShowEndcapSides", true);
1022  desc.addUntracked<bool>("ShowTOB6TEC9", false);
1023  desc.addUntracked<bool>("ShowOnlyGoodModules", false);
1024  descriptions.addWithDefaultLabel(desc);
1025 }
1026 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
constexpr char const * layerName[numberOfLayers]
int goodlayerfound[bounds::k_END_OF_LAYS_AND_RINGS]
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
int goodlayertotal[bounds::k_END_OF_LAYS_AND_RINGS]
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
std::unique_ptr< TrackerTopology > tTopo_
void endRun(edm::Run const &, edm::EventSetup const &) override
int alllayertotal[bounds::k_END_OF_LAYS_AND_RINGS]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void printAndWriteBadModules(const SiStripQuality &quality, const SiStripDetInfo &detInfo) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void makeSummaryVsLumi(DQMStore::IGetter &getter) const
Log< level::Error, false > LogError
std::unordered_map< uint32_t, int > fedErrorCounts
const edm::ESGetToken< SiStripQuality, SiStripQualityRcd > stripQualityToken_
nStrips
1.2 is to make the matching window safely the two nearest strips 0.35 is the size of an ME0 chamber i...
std::unique_ptr< TkHistoMap > FEDErrorOccupancy
constexpr std::array< uint8_t, layerIndexSize > layer
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
char const * label
def overlay(hists, ytitle, header, addon)
Definition: compare.py:122
const bool checkFedError(const DetId det)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
void printTotalStatistics(const std::array< long, bounds::k_END_OF_LAYERS > &layerFound, const std::array< long, bounds::k_END_OF_LAYERS > &layerTotal) const
unsigned int countTotalHits(const std::vector< MonitorElement *> &maps)
void writeBadStripPayload(const SiStripQuality &quality) const
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
Transition
Definition: Transition.h:12
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
void makeSummaryVsBX(DQMStore::IGetter &getter, TFileService &fs) const
SiStripDetInfo read(std::string filePath)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void makeSummaryVsCM(DQMStore::IGetter &getter, TFileService &fs) const
Log< level::Warning, true > LogPrint
SiStripHitEfficiencyHarvester(const edm::ParameterSet &)
Log< level::Info, false > LogInfo
Definition: DetId.h:17
void setBadComponents(int i, int component, const SiStripQuality::BadComponent &BC, int NBadComponent[4][19][4])
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
void fillMapFromTkMap(const int nevents, const float threshold, const std::vector< DetId > &stripDetIds)
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
Constants and enumerated types for FED/FEC systems.
virtual TH1F * getTH1F() const
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
~SiStripHitEfficiencyHarvester() override=default
HLT enums.
static const uint16_t STRIPS_PER_APV
std::pair< ContainerIterator, ContainerIterator > Range
std::unique_ptr< SiStripQuality > stripQuality_
int alllayerfound[bounds::k_END_OF_LAYS_AND_RINGS]
static constexpr char const *const kDefaultFile
bool isAvailable() const
Definition: Service.h:40
Definition: tree.py:1
const edm::ESGetToken< TkDetMap, TrackerTopologyRcd > tkDetMapToken_
string quality
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
bool checkMapsValidity(const std::vector< MonitorElement *> &maps, const std::string &type) const
#define str(s)
void makeSummary(DQMStore::IGetter &getter, TFileService &fs) const
Definition: Run.h:45
#define LogDebug(id)
virtual double getBinContent(int binx) const
get content of bin (1-D)
unsigned transform(const HcalDetId &id, unsigned transformCode)