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> isAvailable;
150  isAvailable.reserve(maps.size());
152  maps.begin() + 1, maps.end(), std::back_inserter(isAvailable), [](auto& x) { return !(x == nullptr); });
153 
154  int count{0};
155  for (const auto& it : isAvailable) {
156  count++;
157  LogDebug("SiStripHitEfficiencyHarvester") << " layer: " << count << " " << it << std::endl;
158  if (it)
159  LogDebug("SiStripHitEfficiencyHarvester") << "resolving to " << maps[count]->getName() << std::endl;
160  }
161 
162  // check on the input TkHistoMap
163  bool areMapsAvailable{true};
164  int layerCount{0};
165  for (const auto& it : isAvailable) {
166  layerCount++;
167  if (!it) {
168  edm::LogError("SiStripHitEfficiencyHarvester")
169  << type << " TkHistoMap for layer " << layerCount << " was not found.\n -> Aborting!";
170  areMapsAvailable = false;
171  break;
172  }
173  }
174  return areMapsAvailable;
175 }
176 
177 unsigned int SiStripHitEfficiencyHarvester::countTotalHits(const std::vector<MonitorElement*>& maps) {
178  return std::accumulate(maps.begin() + 1, maps.end(), 0, [](unsigned int total, MonitorElement* item) {
179  return total + item->getEntries();
180  });
181 }
182 
184  if (!isAtPCL_) {
186 
187  if (doStoreOnTree_) {
188  // store information per DetId in the output tree
189  tree = fs->make<TTree>("ModEff", "ModEff");
190  tree->Branch("DetId", &t_DetId, "DetId/i");
191  tree->Branch("Layer", &t_layer, "Layer/b");
192  tree->Branch("FoundHits", &t_found, "FoundHits/i");
193  tree->Branch("AllHits", &t_total, "AllHits/i");
194  tree->Branch("IsTaggedIneff", &t_isTaggedIneff, "IsTaggedIneff/O");
195  tree->Branch("TagThreshold", &t_threshold, "TagThreshold/F");
196  }
197  }
198 
200  LOGPRINT << "A module is bad if efficiency < " << threshold_ << " and has at least " << nModsMin_ << " nModsMin.";
201  else
202  LOGPRINT << "A module is bad if the upper limit on the efficiency is < to the avg in the layer - " << threshold_
203  << " and has at least " << nModsMin_ << " nModsMin.";
204 
205  auto h_module_total = std::make_unique<TkHistoMap>(tkDetMap_.get());
206  h_module_total->loadTkHistoMap(fmt::format("{}/TkDetMaps", inputFolder_), "perModule_total");
207  auto h_module_found = std::make_unique<TkHistoMap>(tkDetMap_.get());
208  h_module_found->loadTkHistoMap(fmt::format("{}/TkDetMaps", inputFolder_), "perModule_found");
209 
210  // collect how many layers are missing
211  const auto& totalMaps = h_module_total->getAllMaps();
212  const auto& foundMaps = h_module_found->getAllMaps();
213 
214  LogDebug("SiStripHitEfficiencyHarvester")
215  << "totalMaps.size(): " << totalMaps.size() << " foundMaps.size() " << foundMaps.size() << std::endl;
216 
217  // check on the input TkHistoMaps
218  bool isTotalMapAvailable = this->checkMapsValidity(totalMaps, std::string("Total"));
219  bool isFoundMapAvailable = this->checkMapsValidity(foundMaps, std::string("Found"));
220 
221  LogDebug("SiStripHitEfficiencyHarvester")
222  << "isTotalMapAvailable: " << isTotalMapAvailable << " isFoundMapAvailable " << isFoundMapAvailable << std::endl;
223 
224  // no input TkHistoMaps -> early return
225  if (!isTotalMapAvailable or !isFoundMapAvailable)
226  return;
227 
228  LogDebug("SiStripHitEfficiencyHarvester")
229  << "Entries in total TkHistoMap for layer 3: " << h_module_total->getMap(3)->getEntries() << ", found "
230  << h_module_found->getMap(3)->getEntries();
231 
232  // count how many hits in the denominator we have
233  const unsigned int totalHits = this->countTotalHits(totalMaps);
234 
235  // set colz
236  for (size_t i = 1; i < totalMaps.size(); i++) {
237  h_module_total->getMap(i)->setOption("colz");
238  h_module_found->getMap(i)->setOption("colz");
239  }
240 
241  // come back to the main folder
243 
244  std::vector<MonitorElement*> hEffInLayer(std::size_t(1), nullptr);
245  hEffInLayer.reserve(bounds::k_END_OF_LAYERS);
246  for (std::size_t i = 1; i != bounds::k_END_OF_LAYERS; ++i) {
247  const auto lyrName = ::layerName(i, showRings_, nTEClayers_);
248  hEffInLayer.push_back(booker.book1D(
249  Form("eff_layer%i", int(i)), Form("Module efficiency in layer %s", lyrName.c_str()), 201, 0, 1.005));
250  }
251  std::array<long, bounds::k_END_OF_LAYERS> layerTotal{};
252  std::array<long, bounds::k_END_OF_LAYERS> layerFound{};
253  layerTotal.fill(0);
254  layerFound.fill(0);
255 
257  // Effiency calculation, bad module tagging, and tracker maps //
259 
260  TrackerMap tkMap{" Detector Inefficiency "};
261  TrackerMap tkMapBad{" Inefficient Modules "};
262  TrackerMap tkMapEff{title_};
263  TrackerMap tkMapNum{" Detector numerator "};
264  TrackerMap tkMapDen{" Detector denominator "};
265  std::map<unsigned int, double> badModules;
266 
267  for (auto det : stripDetIds_) {
268  auto layer = ::checkLayer(det, tTopo_.get());
269  const auto num = h_module_found->getValue(det);
270  const auto denom = h_module_total->getValue(det);
271  if (denom) {
272  // use only the "good" modules
273  if (stripQuality_->getBadApvs(det) == 0) {
274  const auto eff = num / denom;
275  hEffInLayer[layer]->Fill(eff);
276  if (!autoIneffModTagging_) {
277  if ((denom >= nModsMin_) && (eff < threshold_)) {
278  // We have a bad module, put it in the list!
279  badModules[det] = eff;
280  tkMapBad.fillc(det, 255, 0, 0);
281  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
282  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom;
283  } else {
284  //Fill the bad list with empty results for every module
285  tkMapBad.fillc(det, 255, 255, 255);
286  }
287  if (eff < threshold_)
288  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
289  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom;
290 
291  if (denom < nModsMin_) {
292  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
293  << det.rawId() << " is under occupancy at " << denom;
294  }
295 
296  if (doStoreOnTree_ && !isAtPCL_) {
297  t_DetId = det.rawId();
298  t_layer = layer;
299  t_found = num;
300  t_total = denom;
301  t_isTaggedIneff = false;
302  t_threshold = 0;
303  tree->Fill();
304  }
305  }
306 
307  //Put any module into the TKMap
308  tkMap.fill(det, 1. - eff);
309  tkMapEff.fill(det, eff);
310  tkMapNum.fill(det, num);
311  tkMapDen.fill(det, denom);
312 
313  layerTotal[layer] += denom;
314  layerFound[layer] += num;
315 
316  // for the summary
317  //Have to do the decoding for which side to go on (ugh)
318  if (layer <= bounds::k_LayersAtTOBEnd) {
321  } else if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
322  if (tTopo_->tidSide(det) == 1) {
325  } else if (tTopo_->tidSide(det) == 2) {
326  goodlayerfound[layer + 3] += num;
327  goodlayertotal[layer + 3] += denom;
328  }
329  } else if (layer > bounds::k_LayersAtTIDEnd && layer <= bounds::k_LayersAtTECEnd) {
330  if (tTopo_->tecSide(det) == 1) {
331  goodlayerfound[layer + 3] += num;
332  goodlayertotal[layer + 3] += denom;
333  } else if (tTopo_->tecSide(det) == 2) {
336  }
337  }
338  } // if the module is good!
339 
340  //Do the one where we don't exclude bad modules!
341  if (layer <= bounds::k_LayersAtTOBEnd) {
342  alllayerfound[layer] += num;
344  } else if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
345  if (tTopo_->tidSide(det) == 1) {
346  alllayerfound[layer] += num;
348  } else if (tTopo_->tidSide(det) == 2) {
349  alllayerfound[layer + 3] += num;
350  alllayertotal[layer + 3] += denom;
351  }
352  } else if (layer > bounds::k_LayersAtTIDEnd && layer <= bounds::k_LayersAtTECEnd) {
353  if (tTopo_->tecSide(det) == 1) {
354  alllayerfound[layer + 3] += num;
355  alllayertotal[layer + 3] += denom;
356  } else if (tTopo_->tecSide(det) == 2) {
359  }
360  }
361 
362  } // if denom
363  } // loop on DetIds
364 
365  if (autoIneffModTagging_) {
366  for (Long_t i = 1; i <= k_LayersAtTECEnd; i++) {
367  //Compute threshold to use for each layer
368  hEffInLayer[i]->getTH1()->GetXaxis()->SetRange(
369  3, hEffInLayer[i]->getNbinsX() + 1); // Remove from the avg modules below 1%
370  const double layer_min_eff = hEffInLayer[i]->getMean() - std::max(2.5 * hEffInLayer[i]->getRMS(), threshold_);
371  LOGPRINT << "Layer " << i << " threshold for bad modules: <" << layer_min_eff
372  << " (layer mean: " << hEffInLayer[i]->getMean() << " rms: " << hEffInLayer[i]->getRMS() << ")";
373 
374  hEffInLayer[i]->getTH1()->GetXaxis()->SetRange(1, hEffInLayer[i]->getNbinsX() + 1);
375 
376  for (auto det : stripDetIds_) {
377  if (stripQuality_->getBadApvs(det) == 0) {
378  const auto layer = ::checkLayer(det, tTopo_.get());
379  if (layer == i) {
380  const auto num = h_module_found->getValue(det);
381  const auto denom = h_module_total->getValue(det);
382  if (denom) {
383  const auto eff = num / denom;
384  const auto eff_up = TEfficiency::Bayesian(denom, num, .99, 1, 1, true);
385 
386  if ((denom >= nModsMin_) && (eff_up < layer_min_eff)) {
387  //We have a bad module, put it in the list!
388  badModules[det] = eff;
389  tkMapBad.fillc(det, 255, 0, 0);
390  if (!isAtPCL_ && doStoreOnTree_) {
391  t_isTaggedIneff = true;
392  }
393  } else {
394  //Fill the bad list with empty results for every module
395  tkMapBad.fillc(det, 255, 255, 255);
396  if (!isAtPCL_ && doStoreOnTree_) {
397  t_isTaggedIneff = false;
398  }
399  }
400  if (eff_up < layer_min_eff + 0.08) {
401  // printing message also for modules sligthly above (8%) the limit
402  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
403  << det.rawId() << " efficiency: " << eff << " , " << num << "/" << denom
404  << " , upper limit: " << eff_up;
405  }
406  if (denom < nModsMin_) {
407  LOGPRINT << "Layer " << layer << " (" << ::layerName(layer, showRings_, nTEClayers_) << ") module "
408  << det.rawId() << " layer " << layer << " is under occupancy at " << denom;
409  }
410 
411  if (!isAtPCL_ && doStoreOnTree_) {
412  t_DetId = det.rawId();
413  t_layer = layer;
414  t_found = num;
415  t_total = denom;
416  t_threshold = layer_min_eff;
417  tree->Fill();
418  } // if storing tree
419  } // if denom
420  } // layer = i
421  } // if there are no bad APVs
422  } // loop on detids
423  } // loop on layers
424  } // if auto tagging
425 
426  tkMap.save(true, 0, 0, "SiStripHitEffTKMap_NEW.png");
427  tkMapBad.save(true, 0, 0, "SiStripHitEffTKMapBad_NEW.png");
428  tkMapEff.save(true, tkMapMin_, 1., "SiStripHitEffTKMapEff_NEW.png");
429  tkMapNum.save(true, 0, 0, "SiStripHitEffTKMapNum_NEW.png");
430  tkMapDen.save(true, 0, 0, "SiStripHitEffTKMapDen_NEW.png");
431 
432  const auto detInfo =
434  SiStripQuality pQuality{detInfo};
435  //This is the list of the bad strips, use to mask out entire APVs
436  //Now simply go through the bad hit list and mask out things that
437  //are bad!
438  for (const auto it : badModules) {
439  const auto det = it.first;
440  std::vector<unsigned int> badStripList;
441  //We need to figure out how many strips are in this particular module
442  //To Mask correctly!
443  const auto nStrips = detInfo.getNumberOfApvsAndStripLength(det).first * sistrip::STRIPS_PER_APV;
444  LOGPRINT << "Number of strips module " << det << " is " << nStrips;
445  badStripList.push_back(pQuality.encode(0, nStrips, 0));
446  //Now compact into a single bad module
447  LOGPRINT << "ID1 shoudl match list of modules above " << det;
448  pQuality.compact(det, badStripList);
449  pQuality.put(det, SiStripQuality::Range(badStripList.begin(), badStripList.end()));
450  }
451  pQuality.fillBadComponents();
452  if (doStoreOnDB_) {
453  if (totalHits > 0u) {
454  writeBadStripPayload(pQuality);
455  } else {
456  edm::LogPrint("SiStripHitEfficiencyHarvester")
457  << __PRETTY_FUNCTION__ << " There are no SiStrip hits for a valid measurement, skipping!";
458  }
459  } else {
460  edm::LogInfo("SiStripHitEfficiencyHarvester") << "Will not produce payload!";
461  }
462 
463  printTotalStatistics(layerFound, layerTotal); // statistics by layer and subdetector
464  //LOGPRINT << "\n-----------------\nNew IOV starting from run " << e.id().run() << " event " << e.id().event()
465  // << " lumiBlock " << e.luminosityBlock() << " time " << e.time().value() << "\n-----------------\n";
466  printAndWriteBadModules(pQuality, detInfo); // TODO
467 
468  if (!isAtPCL_) {
470  makeSummary(getter, *fs); // TODO
471  makeSummaryVsBX(getter, *fs); // TODO
472  makeSummaryVsCM(getter, *fs); // TODO
473  }
474 
475  makeSummaryVsLumi(getter); // TODO
476 }
477 
479  const std::array<long, bounds::k_END_OF_LAYERS>& layerFound,
480  const std::array<long, bounds::k_END_OF_LAYERS>& layerTotal) const {
481  //Calculate the statistics by layer
482  int totalfound = 0;
483  int totaltotal = 0;
484  double layereff;
485  int subdetfound[5];
486  int subdettotal[5];
487 
488  for (Long_t i = 1; i < 5; i++) {
489  subdetfound[i] = 0;
490  subdettotal[i] = 0;
491  }
492 
493  for (Long_t i = 1; i <= bounds::k_LayersAtTECEnd; i++) {
494  layereff = double(layerFound[i]) / double(layerTotal[i]);
495  LOGPRINT << "Layer " << i << " (" << ::layerName(i, showRings_, nTEClayers_) << ") has total efficiency "
496  << layereff << " " << layerFound[i] << "/" << layerTotal[i];
497  totalfound += layerFound[i];
498  totaltotal += layerTotal[i];
499  if (i <= bounds::k_LayersAtTIBEnd) {
500  subdetfound[1] += layerFound[i];
501  subdettotal[1] += layerTotal[i];
502  }
503  if (i > bounds::k_LayersAtTIBEnd && i <= bounds::k_LayersAtTOBEnd) {
504  subdetfound[2] += layerFound[i];
505  subdettotal[2] += layerTotal[i];
506  }
507  if (i > bounds::k_LayersAtTOBEnd && i <= bounds::k_LayersAtTIDEnd) {
508  subdetfound[3] += layerFound[i];
509  subdettotal[3] += layerTotal[i];
510  }
511  if (i > bounds::k_LayersAtTIDEnd) {
512  subdetfound[4] += layerFound[i];
513  subdettotal[4] += layerTotal[i];
514  }
515  }
516 
517  LOGPRINT << "The total efficiency is " << double(totalfound) / double(totaltotal);
518  LOGPRINT << " TIB: " << double(subdetfound[1]) / subdettotal[1] << " " << subdetfound[1] << "/"
519  << subdettotal[1];
520  LOGPRINT << " TOB: " << double(subdetfound[2]) / subdettotal[2] << " " << subdetfound[2] << "/"
521  << subdettotal[2];
522  LOGPRINT << " TID: " << double(subdetfound[3]) / subdettotal[3] << " " << subdetfound[3] << "/"
523  << subdettotal[3];
524  LOGPRINT << " TEC: " << double(subdetfound[4]) / subdettotal[4] << " " << subdetfound[4] << "/"
525  << subdettotal[4];
526 }
527 
529  SiStripBadStrip pBadStrip{};
530  const auto pQdvBegin = quality.getDataVectorBegin();
531  for (auto rIt = quality.getRegistryVectorBegin(); rIt != quality.getRegistryVectorEnd(); ++rIt) {
532  const auto range = SiStripBadStrip::Range(pQdvBegin + rIt->ibegin, pQdvBegin + rIt->iend);
533  if (!pBadStrip.put(rIt->detid, range))
534  edm::LogError("SiStripHitEfficiencyHarvester") << "detid already exists in SiStripBadStrip";
535  }
537  if (poolDbService.isAvailable()) {
538  poolDbService->writeOneIOV(pBadStrip, poolDbService->currentTime(), record_);
539  } else {
540  throw cms::Exception("PoolDBService required");
541  }
542 }
543 
545  // use goodlayer_total/found and alllayer_total/found, collapse side and/or ring if needed
546 
547  unsigned int nLayers = 34;
548  if (showRings_)
549  nLayers = 30;
550  if (!showEndcapSides_) {
551  if (!showRings_)
552  nLayers = 22;
553  else
554  nLayers = 20;
555  }
556 
557  TH1F* found = fs.make<TH1F>("found", "found", nLayers + 1, 0, nLayers + 1);
558  TH1F* all = fs.make<TH1F>("all", "all", nLayers + 1, 0, nLayers + 1);
559  TH1F* found2 = fs.make<TH1F>("found2", "found2", nLayers + 1, 0, nLayers + 1);
560  TH1F* all2 = fs.make<TH1F>("all2", "all2", nLayers + 1, 0, nLayers + 1);
561  // first bin only to keep real data off the y axis so set to -1
562  found->SetBinContent(0, -1);
563  all->SetBinContent(0, 1);
564 
565  // new ROOT version: TGraph::Divide don't handle null or negative values
566  for (unsigned int i = 1; i < nLayers + 2; ++i) {
567  found->SetBinContent(i, 1e-6);
568  all->SetBinContent(i, 1);
569  found2->SetBinContent(i, 1e-6);
570  all2->SetBinContent(i, 1);
571  }
572 
573  TCanvas* c7 = new TCanvas("c7", " test ", 10, 10, 800, 600);
574  c7->SetFillColor(0);
575  c7->SetGrid();
576 
577  unsigned int nLayers_max = nLayers + 1; // barrel+endcap
578  if (!showEndcapSides_)
579  nLayers_max = 11; // barrel
580  for (unsigned int i = 1; i < nLayers_max; ++i) {
581  LOGPRINT << "Fill only good modules layer " << i << ": S = " << goodlayerfound[i]
582  << " B = " << goodlayertotal[i];
583  if (goodlayertotal[i] > 5) {
584  found->SetBinContent(i, goodlayerfound[i]);
585  all->SetBinContent(i, goodlayertotal[i]);
586  }
587 
588  LOGPRINT << "Filling all modules layer " << i << ": S = " << alllayerfound[i] << " B = " << alllayertotal[i];
589  if (alllayertotal[i] > 5) {
590  found2->SetBinContent(i, alllayerfound[i]);
591  all2->SetBinContent(i, alllayertotal[i]);
592  }
593  }
594 
595  // endcap - merging sides
596  if (!showEndcapSides_) {
597  for (unsigned int i = 11; i < 14; ++i) { // TID disks
598  LOGPRINT << "Fill only good modules layer " << i << ": S = " << goodlayerfound[i] + goodlayerfound[i + 3]
599  << " B = " << goodlayertotal[i] + goodlayertotal[i + 3];
600  if (goodlayertotal[i] + goodlayertotal[i + 3] > 5) {
601  found->SetBinContent(i, goodlayerfound[i] + goodlayerfound[i + 3]);
602  all->SetBinContent(i, goodlayertotal[i] + goodlayertotal[i + 3]);
603  }
604  LOGPRINT << "Filling all modules layer " << i << ": S = " << alllayerfound[i] + alllayerfound[i + 3]
605  << " B = " << alllayertotal[i] + alllayertotal[i + 3];
606  if (alllayertotal[i] + alllayertotal[i + 3] > 5) {
607  found2->SetBinContent(i, alllayerfound[i] + alllayerfound[i + 3]);
608  all2->SetBinContent(i, alllayertotal[i] + alllayertotal[i + 3]);
609  }
610  }
611  for (unsigned int i = 17; i < 17 + nTEClayers_; ++i) { // TEC disks
612  LOGPRINT << "Fill only good modules layer " << i - 3
613  << ": S = " << goodlayerfound[i] + goodlayerfound[i + nTEClayers_]
614  << " B = " << goodlayertotal[i] + goodlayertotal[i + nTEClayers_];
615  if (goodlayertotal[i] + goodlayertotal[i + nTEClayers_] > 5) {
616  found->SetBinContent(i - 3, goodlayerfound[i] + goodlayerfound[i + nTEClayers_]);
617  all->SetBinContent(i - 3, goodlayertotal[i] + goodlayertotal[i + nTEClayers_]);
618  }
619  LOGPRINT << "Filling all modules layer " << i - 3
620  << ": S = " << alllayerfound[i] + alllayerfound[i + nTEClayers_]
621  << " B = " << alllayertotal[i] + alllayertotal[i + nTEClayers_];
622  if (alllayertotal[i] + alllayertotal[i + nTEClayers_] > 5) {
623  found2->SetBinContent(i - 3, alllayerfound[i] + alllayerfound[i + nTEClayers_]);
624  all2->SetBinContent(i - 3, alllayertotal[i] + alllayertotal[i + nTEClayers_]);
625  }
626  }
627  }
628 
629  found->Sumw2();
630  all->Sumw2();
631 
632  found2->Sumw2();
633  all2->Sumw2();
634 
635  TGraphAsymmErrors* gr = fs.make<TGraphAsymmErrors>(nLayers + 1);
636  gr->SetName("eff_good");
637  gr->BayesDivide(found, all);
638 
639  TGraphAsymmErrors* gr2 = fs.make<TGraphAsymmErrors>(nLayers + 1);
640  gr2->SetName("eff_all");
641  gr2->BayesDivide(found2, all2);
642 
643  for (unsigned int j = 0; j < nLayers + 1; j++) {
644  gr->SetPointError(j, 0., 0., gr->GetErrorYlow(j), gr->GetErrorYhigh(j));
645  gr2->SetPointError(j, 0., 0., gr2->GetErrorYlow(j), gr2->GetErrorYhigh(j));
646  }
647 
648  gr->GetXaxis()->SetLimits(0, nLayers);
649  gr->SetMarkerColor(2);
650  gr->SetMarkerSize(1.2);
651  gr->SetLineColor(2);
652  gr->SetLineWidth(4);
653  gr->SetMarkerStyle(20);
654  gr->SetMinimum(effPlotMin_);
655  gr->SetMaximum(1.001);
656  gr->GetYaxis()->SetTitle("Efficiency");
657  gStyle->SetTitleFillColor(0);
658  gStyle->SetTitleBorderSize(0);
659  gr->SetTitle(title_.c_str());
660 
661  gr2->GetXaxis()->SetLimits(0, nLayers);
662  gr2->SetMarkerColor(1);
663  gr2->SetMarkerSize(1.2);
664  gr2->SetLineColor(1);
665  gr2->SetLineWidth(4);
666  gr2->SetMarkerStyle(21);
667  gr2->SetMinimum(effPlotMin_);
668  gr2->SetMaximum(1.001);
669  gr2->GetYaxis()->SetTitle("Efficiency");
670  gr2->SetTitle(title_.c_str());
671 
672  for (Long_t k = 1; k < nLayers + 1; k++) {
673  TString label;
674  if (showEndcapSides_)
675  label = ::layerSideName(k, showRings_, nTEClayers_);
676  else
678  if (!showTOB6TEC9_) {
679  if (k == 10)
680  label = "";
681  if (!showRings_ && k == nLayers)
682  label = "";
683  if (!showRings_ && showEndcapSides_ && k == 25)
684  label = "";
685  }
686  if (!showRings_) {
687  if (showEndcapSides_) {
688  gr->GetXaxis()->SetBinLabel(((k + 1) * 100 + 2) / (nLayers)-4, label);
689  gr2->GetXaxis()->SetBinLabel(((k + 1) * 100 + 2) / (nLayers)-4, label);
690  } else {
691  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-6, label);
692  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-6, label);
693  }
694  } else {
695  if (showEndcapSides_) {
696  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-4, label);
697  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-4, label);
698  } else {
699  gr->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-7, label);
700  gr2->GetXaxis()->SetBinLabel((k + 1) * 100 / (nLayers)-7, label);
701  }
702  }
703  }
704 
705  gr->Draw("AP");
706  gr->GetXaxis()->SetNdivisions(36);
707 
708  c7->cd();
709  TPad* overlay = new TPad("overlay", "", 0, 0, 1, 1);
710  overlay->SetFillStyle(4000);
711  overlay->SetFillColor(0);
712  overlay->SetFrameFillStyle(4000);
713  overlay->Draw("same");
714  overlay->cd();
716  gr2->Draw("AP");
717 
718  TLegend* leg = new TLegend(0.70, 0.27, 0.88, 0.40);
719  leg->AddEntry(gr, "Good Modules", "p");
721  leg->AddEntry(gr2, "All Modules", "p");
722  leg->SetTextSize(0.020);
723  leg->SetFillColor(0);
724  leg->Draw("same");
725 
726  c7->SaveAs("Summary.png");
727  c7->SaveAs("Summary.root");
728 }
729 
731  // use found/totalVsBx_layer%i [0,23)
732 }
733 
735  for (unsigned int iLayer = 1; iLayer != (showRings_ ? 20 : 22); ++iLayer) {
736  auto hfound = getter.get(fmt::format("{}/VsLumi/layerfound_vsLumi_layer_{}", inputFolder_, iLayer))->getTH1F();
737  auto htotal = getter.get(fmt::format("{}/VsLumi/layertotal_vsLumi_layer_{}", inputFolder_, iLayer))->getTH1F();
738 
739  if (hfound == nullptr or htotal == nullptr) {
740  if (hfound == nullptr)
741  edm::LogError("SiStripHitEfficiencyHarvester")
742  << fmt::format("{}/VsLumi/layerfound_vsLumi_layer_{}", inputFolder_, iLayer) << " was not found!";
743  if (htotal == nullptr)
744  edm::LogError("SiStripHitEfficiencyHarvester")
745  << fmt::format("{}/VsLumi/layertotal_vsLumi_layer_{}", inputFolder_, iLayer) << " was not found!";
746  // no input histograms -> continue in the loop
747  continue;
748  }
749 
750  if (!hfound->GetSumw2())
751  hfound->Sumw2();
752  if (!htotal->GetSumw2())
753  htotal->Sumw2();
754  for (Long_t i = 0; i != hfound->GetNbinsX() + 1; ++i) {
755  if (hfound->GetBinContent(i) == 0)
756  hfound->SetBinContent(i, 1e-6);
757  if (htotal->GetBinContent(i) == 0)
758  htotal->SetBinContent(i, 1);
759  }
760  LogDebug("SiStripHitEfficiencyHarvester")
761  << "Total hits for layer " << iLayer << " (vs lumi): " << htotal->GetEntries() << ", found "
762  << hfound->GetEntries();
763  }
764  // continue
765 }
766 
768 
769 namespace {
770  void setBadComponents(int i,
771  int comp,
773  std::stringstream ssV[4][19],
774  int nBad[4][19][4],
775  int nAPV) {
776  ssV[i][comp] << "\n\t\t " << bc.detid << " \t " << bc.BadModule << " \t " << ((bc.BadFibers) & 0x1) << " ";
777  if (nAPV == 4)
778  ssV[i][comp] << "x " << ((bc.BadFibers >> 1) & 0x1);
779 
780  if (nAPV == 6)
781  ssV[i][comp] << ((bc.BadFibers >> 1) & 0x1) << " " << ((bc.BadFibers >> 2) & 0x1);
782  ssV[i][comp] << " \t " << ((bc.BadApvs) & 0x1) << " " << ((bc.BadApvs >> 1) & 0x1) << " ";
783  if (nAPV == 4)
784  ssV[i][comp] << "x x " << ((bc.BadApvs >> 2) & 0x1) << " " << ((bc.BadApvs >> 3) & 0x1);
785  if (nAPV == 6)
786  ssV[i][comp] << ((bc.BadApvs >> 2) & 0x1) << " " << ((bc.BadApvs >> 3) & 0x1) << " " << ((bc.BadApvs >> 4) & 0x1)
787  << " " << ((bc.BadApvs >> 5) & 0x1) << " ";
788 
789  if (bc.BadApvs) {
790  nBad[i][0][2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
791  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
792  nBad[i][comp][2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
793  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
794  }
795  if (bc.BadFibers) {
796  nBad[i][0][1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
797  nBad[i][comp][1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
798  }
799  if (bc.BadModule) {
800  nBad[i][0][0]++;
801  nBad[i][comp][0]++;
802  }
803  }
804 } // namespace
805 
807  const SiStripDetInfo& detInfo) const {
809  //try to write out what's in the quality record
811  int nTkBadComp[4]; //k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
812  int nBadComp[4][19][4];
813  //legend: nBadComp[i][j][k]= SubSystem i, layer/disk/wheel j, BadModule/Fiber/Apv k
814  // i: 0=TIB, 1=TID, 2=TOB, 3=TEC
815  // k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
816  std::stringstream ssV[4][19];
817 
818  for (int i = 0; i < 4; ++i) {
819  nTkBadComp[i] = 0;
820  for (int j = 0; j < 19; ++j) {
821  ssV[i][j].str("");
822  for (int k = 0; k < 4; ++k)
823  nBadComp[i][j][k] = 0;
824  }
825  }
826 
827  for (const auto& bc : quality.getBadComponentList()) {
828  // Full Tk
829  if (bc.BadModule)
830  nTkBadComp[0]++;
831  if (bc.BadFibers)
832  nTkBadComp[1] += ((bc.BadFibers >> 2) & 0x1) + ((bc.BadFibers >> 1) & 0x1) + ((bc.BadFibers) & 0x1);
833  if (bc.BadApvs)
834  nTkBadComp[2] += ((bc.BadApvs >> 5) & 0x1) + ((bc.BadApvs >> 4) & 0x1) + ((bc.BadApvs >> 3) & 0x1) +
835  ((bc.BadApvs >> 2) & 0x1) + ((bc.BadApvs >> 1) & 0x1) + ((bc.BadApvs) & 0x1);
836  // single subsystem
837  DetId det(bc.detid);
838  if ((det.subdetId() >= SiStripSubdetector::TIB) && (det.subdetId() <= SiStripSubdetector::TEC)) {
839  const auto nAPV = detInfo.getNumberOfApvsAndStripLength(det).first;
840  switch (det.subdetId()) {
842  setBadComponents(0, tTopo_->tibLayer(det), bc, ssV, nBadComp, nAPV);
843  break;
846  (tTopo_->tidSide(det) == 2 ? tTopo_->tidWheel(det) : tTopo_->tidWheel(det) + 3),
847  bc,
848  ssV,
849  nBadComp,
850  nAPV);
851  break;
853  setBadComponents(2, tTopo_->tobLayer(det), bc, ssV, nBadComp, nAPV);
854  break;
857  (tTopo_->tecSide(det) == 2 ? tTopo_->tecWheel(det) : tTopo_->tecWheel(det) + 9),
858  bc,
859  ssV,
860  nBadComp,
861  nAPV);
862  break;
863  default:
864  break;
865  }
866  }
867  }
868  // single strip info
869  for (auto rp = quality.getRegistryVectorBegin(); rp != quality.getRegistryVectorEnd(); ++rp) {
870  DetId det{rp->detid};
871  int subdet = -999;
872  int component = -999;
873  switch (det.subdetId()) {
875  subdet = 0;
876  component = tTopo_->tibLayer(det);
877  break;
879  subdet = 1;
880  component = tTopo_->tidSide(det) == 2 ? tTopo_->tidWheel(det) : tTopo_->tidWheel(det) + 3;
881  break;
883  subdet = 2;
884  component = tTopo_->tobLayer(det);
885  break;
887  subdet = 3;
888  component = tTopo_->tecSide(det) == 2 ? tTopo_->tecWheel(det) : tTopo_->tecWheel(det) + 9;
889  break;
890  default:
891  break;
892  }
893 
894  const auto pQdvBegin = quality.getDataVectorBegin();
895  const auto sqrange = SiStripQuality::Range(pQdvBegin + rp->ibegin, pQdvBegin + rp->iend);
896  float percentage = 0;
897  for (int it = 0; it < sqrange.second - sqrange.first; it++) {
898  unsigned int range = quality.decode(*(sqrange.first + it)).range;
899  nTkBadComp[3] += range;
900  nBadComp[subdet][0][3] += range;
901  nBadComp[subdet][component][3] += range;
902  percentage += range;
903  }
904  if (percentage != 0)
905  percentage /= (sistrip::STRIPS_PER_APV * detInfo.getNumberOfApvsAndStripLength(det).first);
906  if (percentage > 1)
907  edm::LogError("SiStripHitEfficiencyHarvester") << "PROBLEM detid " << det.rawId() << " value " << percentage;
908  }
909 
910  // printout
911  std::ostringstream ss;
912  ss << "\n-----------------\nGlobal Info\n-----------------";
913  ss << "\nBadComp \t Modules \tFibers "
914  "\tApvs\tStrips\n----------------------------------------------------------------";
915  ss << "\nTracker:\t\t" << nTkBadComp[0] << "\t" << nTkBadComp[1] << "\t" << nTkBadComp[2] << "\t" << nTkBadComp[3];
916  ss << "\nTIB:\t\t\t" << nBadComp[0][0][0] << "\t" << nBadComp[0][0][1] << "\t" << nBadComp[0][0][2] << "\t"
917  << nBadComp[0][0][3];
918  ss << "\nTID:\t\t\t" << nBadComp[1][0][0] << "\t" << nBadComp[1][0][1] << "\t" << nBadComp[1][0][2] << "\t"
919  << nBadComp[1][0][3];
920  ss << "\nTOB:\t\t\t" << nBadComp[2][0][0] << "\t" << nBadComp[2][0][1] << "\t" << nBadComp[2][0][2] << "\t"
921  << nBadComp[2][0][3];
922  ss << "\nTEC:\t\t\t" << nBadComp[3][0][0] << "\t" << nBadComp[3][0][1] << "\t" << nBadComp[3][0][2] << "\t"
923  << nBadComp[3][0][3];
924  ss << "\n";
925 
926  for (int i = 1; i < 5; ++i)
927  ss << "\nTIB Layer " << i << " :\t\t" << nBadComp[0][i][0] << "\t" << nBadComp[0][i][1] << "\t" << nBadComp[0][i][2]
928  << "\t" << nBadComp[0][i][3];
929  ss << "\n";
930  for (int i = 1; i < 4; ++i)
931  ss << "\nTID+ Disk " << i << " :\t\t" << nBadComp[1][i][0] << "\t" << nBadComp[1][i][1] << "\t" << nBadComp[1][i][2]
932  << "\t" << nBadComp[1][i][3];
933  for (int i = 4; i < 7; ++i)
934  ss << "\nTID- Disk " << i - 3 << " :\t\t" << nBadComp[1][i][0] << "\t" << nBadComp[1][i][1] << "\t"
935  << nBadComp[1][i][2] << "\t" << nBadComp[1][i][3];
936  ss << "\n";
937  for (int i = 1; i < 7; ++i)
938  ss << "\nTOB Layer " << i << " :\t\t" << nBadComp[2][i][0] << "\t" << nBadComp[2][i][1] << "\t" << nBadComp[2][i][2]
939  << "\t" << nBadComp[2][i][3];
940  ss << "\n";
941  for (int i = 1; i < 10; ++i)
942  ss << "\nTEC+ Disk " << i << " :\t\t" << nBadComp[3][i][0] << "\t" << nBadComp[3][i][1] << "\t" << nBadComp[3][i][2]
943  << "\t" << nBadComp[3][i][3];
944  for (int i = 10; i < 19; ++i)
945  ss << "\nTEC- Disk " << i - 9 << " :\t\t" << nBadComp[3][i][0] << "\t" << nBadComp[3][i][1] << "\t"
946  << nBadComp[3][i][2] << "\t" << nBadComp[3][i][3];
947  ss << "\n";
948 
949  ss << "\n----------------------------------------------------------------\n\t\t Detid \tModules Fibers "
950  "Apvs\n----------------------------------------------------------------";
951  for (int i = 1; i < 5; ++i)
952  ss << "\nTIB Layer " << i << " :" << ssV[0][i].str();
953  ss << "\n";
954  for (int i = 1; i < 4; ++i)
955  ss << "\nTID+ Disk " << i << " :" << ssV[1][i].str();
956  for (int i = 4; i < 7; ++i)
957  ss << "\nTID- Disk " << i - 3 << " :" << ssV[1][i].str();
958  ss << "\n";
959  for (int i = 1; i < 7; ++i)
960  ss << "\nTOB Layer " << i << " :" << ssV[2][i].str();
961  ss << "\n";
962  for (int i = 1; i < 10; ++i)
963  ss << "\nTEC+ Disk " << i << " :" << ssV[3][i].str();
964  for (int i = 10; i < 19; ++i)
965  ss << "\nTEC- Disk " << i - 9 << " :" << ssV[3][i].str();
966 
967  LOGPRINT << ss.str();
968 
969  // store also bad modules in log file
970  std::ofstream badModules;
971  badModules.open("BadModules_NEW.log");
972  badModules << "\n----------------------------------------------------------------\n\t\t Detid \tModules Fibers "
973  "Apvs\n----------------------------------------------------------------";
974  for (int i = 1; i < 5; ++i)
975  badModules << "\nTIB Layer " << i << " :" << ssV[0][i].str();
976  badModules << "\n";
977  for (int i = 1; i < 4; ++i)
978  badModules << "\nTID+ Disk " << i << " :" << ssV[1][i].str();
979  for (int i = 4; i < 7; ++i)
980  badModules << "\nTID- Disk " << i - 3 << " :" << ssV[1][i].str();
981  badModules << "\n";
982  for (int i = 1; i < 7; ++i)
983  badModules << "\nTOB Layer " << i << " :" << ssV[2][i].str();
984  badModules << "\n";
985  for (int i = 1; i < 10; ++i)
986  badModules << "\nTEC+ Disk " << i << " :" << ssV[3][i].str();
987  for (int i = 10; i < 19; ++i)
988  badModules << "\nTEC- Disk " << i - 9 << " :" << ssV[3][i].str();
989  badModules.close();
990 }
991 
994  desc.add<std::string>("inputFolder", "AlCaReco/SiStripHitEfficiency");
995  desc.add<bool>("isAtPCL", false);
996  desc.add<bool>("doStoreOnDB", false);
997  desc.add<std::string>("Record", "SiStripBadStrip");
998  desc.add<double>("Threshold", 0.1);
999  desc.add<std::string>("Title", "Hit Efficiency");
1000  desc.add<int>("nModsMin", 5);
1001  desc.addUntracked<bool>("doStoreOnTree", false);
1002  desc.addUntracked<bool>("AutoIneffModTagging", false);
1003  desc.addUntracked<double>("TkMapMin", 0.9);
1004  desc.addUntracked<double>("EffPlotMin", 0.9);
1005  desc.addUntracked<bool>("ShowRings", false);
1006  desc.addUntracked<bool>("ShowEndcapSides", true);
1007  desc.addUntracked<bool>("ShowTOB6TEC9", false);
1008  desc.addUntracked<bool>("ShowOnlyGoodModules", false);
1009  descriptions.addWithDefaultLabel(desc);
1010 }
1011 
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
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...
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
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
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)
unsigned transform(const HcalDetId &id, unsigned transformCode)