CMS 3D CMS Logo

MillePedeDQMModule.cc
Go to the documentation of this file.
1 
9 /*** header-file ***/
11 
12 /*** ROOT objects ***/
13 #include "TH1F.h"
14 
15 /*** Core framework functionality ***/
18 
19 /*** Alignment ***/
23 
24 /*** Necessary Framework infrastructure ***/
27 
29  : tTopoToken_(esConsumes<edm::Transition::BeginRun>()),
30  gDetToken_(esConsumes<edm::Transition::BeginRun>()),
31  ptpToken_(esConsumes<edm::Transition::BeginRun>()),
32  aliThrToken_(esConsumes<edm::Transition::BeginRun>()),
33  geomToken_(esConsumes<edm::Transition::BeginRun>()),
34  outputFolder_(config.getParameter<std::string>("outputFolder")),
35  mpReaderConfig_(config.getParameter<edm::ParameterSet>("MillePedeFileReader")),
36  isHG_(mpReaderConfig_.getParameter<bool>("isHG")) {
37  consumes<AlignmentToken, edm::InProcess>(config.getParameter<edm::InputTag>("alignmentTokenSrc"));
38 }
39 
41 
42 //=============================================================================
43 //=== INTERFACE IMPLEMENTATION ===
44 //=============================================================================
45 
47  edm::LogInfo("MillePedeDQMModule") << "Booking histograms";
48 
49  booker.cd();
50  if (!isHG_) {
51  if (outputFolder_.find("HG") != std::string::npos) {
52  throw cms::Exception("LogicError")
53  << "MillePedeDQMModule is configured as Low Granularity but the outputfolder is for High Granularity";
54  }
55 
57  h_xPos = booker.book1D("Xpos", "Alignment fit #DeltaX;;#mum", 36, 0., 36.);
58  h_xRot = booker.book1D("Xrot", "Alignment fit #Delta#theta_{X};;#murad", 36, 0., 36.);
59  h_yPos = booker.book1D("Ypos", "Alignment fit #DeltaY;;#mum", 36, 0., 36.);
60  h_yRot = booker.book1D("Yrot", "Alignment fit #Delta#theta_{Y};;#murad", 36, 0., 36.);
61  h_zPos = booker.book1D("Zpos", "Alignment fit #DeltaZ;;#mum", 36, 0., 36.);
62  h_zRot = booker.book1D("Zrot", "Alignment fit #Delta#theta_{Z};;#murad", 36, 0., 36.);
63  statusResults = booker.book2D("statusResults", "Status of SiPixelAli PCL workflow;;", 6, 0., 6., 1, 0., 1.);
64  } else {
65  if (outputFolder_.find("HG") == std::string::npos) {
66  throw cms::Exception("LogicError")
67  << "MillePedeDQMModule is configured as High Granularity but the outputfolder is for Low Granularity";
68  }
69 
71 
72  layerVec = {{"Layer1", pixelTopologyMap_->getPXBLadders(1)},
73  {"Layer2", pixelTopologyMap_->getPXBLadders(2)},
74  {"Layer3", pixelTopologyMap_->getPXBLadders(3)},
75  {"Layer4", pixelTopologyMap_->getPXBLadders(4)},
76  {"Disk-3", pixelTopologyMap_->getPXFBlades(-3) * 2},
77  {"Disk-2", pixelTopologyMap_->getPXFBlades(-2) * 2},
78  {"Disk-1", pixelTopologyMap_->getPXFBlades(-1) * 2},
79  {"Disk1", pixelTopologyMap_->getPXFBlades(1) * 2},
80  {"Disk2", pixelTopologyMap_->getPXFBlades(2) * 2},
81  {"Disk3", pixelTopologyMap_->getPXFBlades(3) * 2}};
82 
83  for (const auto& layer : layerVec) {
84  h_xPos_HG[layer.first] = booker.book1D("Xpos_HG_" + layer.first,
85  "Alignment fit #DeltaX for " + layer.first + ";;#mum",
86  layer.second + 5,
87  0.,
88  layer.second + 5);
89  h_xRot_HG[layer.first] = booker.book1D("Xrot_HG_" + layer.first,
90  "Alignment fit #Delta#theta_{X} for " + layer.first + ";;#murad",
91  layer.second + 5,
92  0.,
93  layer.second + 5);
94  h_yPos_HG[layer.first] = booker.book1D("Ypos_HG_" + layer.first,
95  "Alignment fit #DeltaY for " + layer.first + ";;#mum",
96  layer.second + 5,
97  0.,
98  layer.second + 5);
99  h_yRot_HG[layer.first] = booker.book1D("Yrot_HG_" + layer.first,
100  "Alignment fit #Delta#theta_{Y} for " + layer.first + ";;#murad",
101  layer.second + 5,
102  0.,
103  layer.second + 5);
104  h_zPos_HG[layer.first] = booker.book1D("Zpos_HG_" + layer.first,
105  "Alignment fit #DeltaZ for " + layer.first + ";;#mum",
106  layer.second + 5,
107  0.,
108  layer.second + 5);
109  h_zRot_HG[layer.first] = booker.book1D("Zrot_HG_" + layer.first,
110  "Alignment fit #Delta#theta_{Z} for " + layer.first + ";;#murad",
111  layer.second + 5,
112  0.,
113  layer.second + 5);
114  }
115 
116  statusResults =
117  booker.book2D("statusResults", "Fraction threshold check for SiPixelAliHG PCL;;", 6, 0., 6., 10, 0., 10.);
118  }
119 
120  binariesAvalaible = booker.bookInt("BinariesFound");
121  exitCode = booker.bookString("PedeExitCode", "");
122  isVetoed = booker.bookString("IsVetoed", "");
123 
124  booker.cd();
125 }
126 
128  bookHistograms(booker);
129  if (mpReader_) {
130  mpReader_->read();
131  } else {
132  throw cms::Exception("LogicError") << "@SUB=MillePedeDQMModule::dqmEndJob\n"
133  << "Try to read MillePede results before initializing MillePedeFileReader";
134  }
135  if (!isHG_) {
138  } else {
141  }
142  binariesAvalaible->Fill(mpReader_->binariesAmount());
143  auto theResults = mpReader_->getResults();
144  std::string exitCodeStr = theResults.getExitMessage();
145 
146  std::string vetoStr{};
147  if (mpReader_->storeAlignments()) {
148  vetoStr = "DB Updated!"; /* easy peasy, fait accompli an alignment is there */
149  } else {
150  if (theResults.isHighGranularity()) { /* HG case */
151  if (theResults.getDBVetoed() && theResults.getDBUpdated()) {
152  vetoStr = "DB Update Vetoed"; /* this can happen in the HG PCL case */
153  } else {
154  vetoStr = "N/A";
155  }
156  } else { /* LG case */
157  if (theResults.exceedsCutoffs()) {
158  vetoStr = "DB Update Vetoed"; /* this can happen in the LG PCL case */
159  } else {
160  vetoStr = "N/A";
161  } // if the alignment exceeds the cutoffs
162  } // LG case
163  } // if the alignment was not stored
164 
165  exitCode->Fill(exitCodeStr);
166  isVetoed->Fill(vetoStr);
167 }
168 
169 //=============================================================================
170 //=== PRIVATE METHOD IMPLEMENTATION ===
171 //=============================================================================
172 
174  if (!setupChanged(setup))
175  return;
176 
177  const TrackerTopology* const tTopo = &setup.getData(tTopoToken_);
178  const GeometricDet* geometricDet = &setup.getData(gDetToken_);
179  const PTrackerParameters* ptp = &setup.getData(ptpToken_);
180  const TrackerGeometry* geom = &setup.getData(geomToken_);
181 
182  pixelTopologyMap_ = std::make_shared<PixelTopologyMap>(geom, tTopo);
183 
184  // take the thresholds from DB
185  const auto& thresholds_ = &setup.getData(aliThrToken_);
186 
187  auto myThresholds = std::make_shared<AlignPCLThresholdsHG>();
188  myThresholds->setAlignPCLThresholds(thresholds_->getNrecords(), thresholds_->getThreshold_Map());
189  myThresholds->setFloatMap(thresholds_->getFloatMap());
190 
192 
193  const auto trackerGeometry = builder.build(geometricDet, *ptp, tTopo);
194  tracker_ = std::make_unique<AlignableTracker>(trackerGeometry, tTopo);
195 
196  const std::string labelerPlugin{"PedeLabeler"};
197  edm::ParameterSet labelerConfig{};
198  labelerConfig.addUntrackedParameter("plugin", labelerPlugin);
199  labelerConfig.addUntrackedParameter("RunRangeSelection", edm::VParameterSet{});
200 
201  std::shared_ptr<PedeLabelerBase> pedeLabeler{PedeLabelerPluginFactory::get()->create(
202  labelerPlugin, PedeLabelerBase::TopLevelAlignables(tracker_.get(), nullptr, nullptr), labelerConfig)};
203 
204  mpReader_ = std::make_unique<MillePedeFileReader>(
205  mpReaderConfig_, pedeLabeler, std::shared_ptr<const AlignPCLThresholdsHG>(myThresholds), pixelTopologyMap_);
206 }
207 
209  TH2F* histo_status = statusHisto->getTH2F();
210  auto theResults = mpReader_->getResults();
211  theResults.print();
212  histo_status->SetBinContent(1, 1, theResults.getDBUpdated());
213  histo_status->GetXaxis()->SetBinLabel(1, "DB updated");
214  histo_status->SetBinContent(2, 1, theResults.exceedsCutoffs());
215  histo_status->GetXaxis()->SetBinLabel(2, "significant movement");
216  histo_status->SetBinContent(3, 1, theResults.getDBVetoed());
217  histo_status->GetXaxis()->SetBinLabel(3, "DB update vetoed");
218  histo_status->SetBinContent(4, 1, !theResults.exceedsThresholds());
219  histo_status->GetXaxis()->SetBinLabel(4, "within max movement");
220  histo_status->SetBinContent(5, 1, !theResults.exceedsMaxError());
221  histo_status->GetXaxis()->SetBinLabel(5, "within max error");
222  histo_status->SetBinContent(6, 1, !theResults.belowSignificance());
223  histo_status->GetXaxis()->SetBinLabel(6, "above significance");
224 }
225 
227  TH2F* histo_status = statusHisto->getTH2F();
228  auto& theResults = mpReader_->getResultsHG();
229  histo_status->GetXaxis()->SetBinLabel(1, "#DeltaX");
230  histo_status->GetXaxis()->SetBinLabel(2, "#DeltaY");
231  histo_status->GetXaxis()->SetBinLabel(3, "#DeltaZ");
232  histo_status->GetXaxis()->SetBinLabel(4, "#Delta#theta_{X}");
233  histo_status->GetXaxis()->SetBinLabel(5, "#Delta#theta_{Y}");
234  histo_status->GetXaxis()->SetBinLabel(6, "#Delta#theta_{Z}");
235 
236  int i = 0;
237  for (const auto& result : theResults) {
238  histo_status->GetYaxis()->SetBinLabel(i + 1, result.first.data());
239  for (std::size_t j = 0; j < result.second.size(); ++j) {
240  histo_status->SetBinContent(j + 1, i + 1, result.second[j]);
241  }
242  i++;
243  }
244 }
245 
247  std::array<double, SIZE_INDEX> Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_;
248  std::array<double, SIZE_INDEX> tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_;
249 
250  std::array<double, SIZE_INDEX> Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_;
251  std::array<double, SIZE_INDEX> tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_;
252 
253  std::array<double, SIZE_INDEX> Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_;
254  std::array<double, SIZE_INDEX> tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_;
255 
256  auto myMap = mpReader_->getThresholdMap();
257 
258  std::vector<std::string> alignablesList;
259  for (auto it = myMap.begin(); it != myMap.end(); ++it) {
260  alignablesList.push_back(it->first);
261  }
262 
263  for (auto& alignable : alignablesList) {
264  int detIndex = getIndexFromString(alignable);
265 
266  Xcut_[detIndex] = myMap[alignable].getXcut();
267  sigXcut_[detIndex] = myMap[alignable].getSigXcut();
268  maxMoveXcut_[detIndex] = myMap[alignable].getMaxMoveXcut();
269  maxErrorXcut_[detIndex] = myMap[alignable].getErrorXcut();
270 
271  Ycut_[detIndex] = myMap[alignable].getYcut();
272  sigYcut_[detIndex] = myMap[alignable].getSigYcut();
273  maxMoveYcut_[detIndex] = myMap[alignable].getMaxMoveYcut();
274  maxErrorYcut_[detIndex] = myMap[alignable].getErrorYcut();
275 
276  Zcut_[detIndex] = myMap[alignable].getZcut();
277  sigZcut_[detIndex] = myMap[alignable].getSigZcut();
278  maxMoveZcut_[detIndex] = myMap[alignable].getMaxMoveZcut();
279  maxErrorZcut_[detIndex] = myMap[alignable].getErrorZcut();
280 
281  tXcut_[detIndex] = myMap[alignable].getThetaXcut();
282  sigtXcut_[detIndex] = myMap[alignable].getSigThetaXcut();
283  maxMovetXcut_[detIndex] = myMap[alignable].getMaxMoveThetaXcut();
284  maxErrortXcut_[detIndex] = myMap[alignable].getErrorThetaXcut();
285 
286  tYcut_[detIndex] = myMap[alignable].getThetaYcut();
287  sigtYcut_[detIndex] = myMap[alignable].getSigThetaYcut();
288  maxMovetYcut_[detIndex] = myMap[alignable].getMaxMoveThetaYcut();
289  maxErrortYcut_[detIndex] = myMap[alignable].getErrorThetaYcut();
290 
291  tZcut_[detIndex] = myMap[alignable].getThetaZcut();
292  sigtZcut_[detIndex] = myMap[alignable].getSigThetaZcut();
293  maxMovetZcut_[detIndex] = myMap[alignable].getMaxMoveThetaZcut();
294  maxErrortZcut_[detIndex] = myMap[alignable].getErrorThetaZcut();
295  }
296 
297  fillExpertHisto(h_xPos, Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_, mpReader_->getXobs(), mpReader_->getXobsErr());
299  h_xRot, tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_, mpReader_->getTXobs(), mpReader_->getTXobsErr());
300 
301  fillExpertHisto(h_yPos, Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_, mpReader_->getYobs(), mpReader_->getYobsErr());
303  h_yRot, tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_, mpReader_->getTYobs(), mpReader_->getTYobsErr());
304 
305  fillExpertHisto(h_zPos, Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_, mpReader_->getZobs(), mpReader_->getZobsErr());
307  h_zRot, tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_, mpReader_->getTZobs(), mpReader_->getTZobsErr());
308 }
309 
311  const std::array<double, SIZE_INDEX>& cut,
312  const std::array<double, SIZE_INDEX>& sigCut,
313  const std::array<double, SIZE_INDEX>& maxMoveCut,
314  const std::array<double, SIZE_INDEX>& maxErrorCut,
315  const std::array<double, SIZE_LG_STRUCTS>& obs,
316  const std::array<double, SIZE_LG_STRUCTS>& obsErr) {
317  TH1F* histo_0 = histo->getTH1F();
318 
319  double max_ = *std::max_element(maxMoveCut.begin(), maxMoveCut.end());
320 
321  histo_0->SetMinimum(-(max_));
322  histo_0->SetMaximum(max_);
323 
324  // Schematics of the bin contents
325  //
326  // XX XX XX XX XX XX OO OO OO OO II II II II
327  // |--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|
328  // | 1| 2| 3| 4| 5| 6| 7| 8| 9|10|11|12|13|14|15|16|17| ...
329  //
330  // |-----------------| |-----------| |-----------|
331  // |observed movement| |thresholds1| |thresholds2|
332 
333  for (size_t i = 0; i < obs.size(); ++i) {
334  // fist obs.size() bins for observed movements
335  histo_0->SetBinContent(i + 1, obs[i]);
336  histo_0->SetBinError(i + 1, obsErr[i]);
337 
338  // then at bin 8,8+5,8+10,... for cutoffs
339  // 5 bins is the space allocated for the 4 other thresholds + 1 empty separation bin
340  histo_0->SetBinContent(8 + i * 5, cut[i]);
341 
342  // then at bin 9,9+5,9+10,... for significances
343  histo_0->SetBinContent(9 + i * 5, sigCut[i]);
344 
345  // then at bin 10,10+5,10+10,... for maximum movements
346  histo_0->SetBinContent(10 + i * 5, maxMoveCut[i]);
347 
348  // then at bin 11,11+5,11+10,... for maximum errors
349  histo_0->SetBinContent(11 + i * 5, maxErrorCut[i]);
350  }
351 }
352 
354  std::array<double, SIZE_INDEX> Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_;
355  std::array<double, SIZE_INDEX> tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_;
356 
357  std::array<double, SIZE_INDEX> Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_;
358  std::array<double, SIZE_INDEX> tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_;
359 
360  std::array<double, SIZE_INDEX> Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_;
361  std::array<double, SIZE_INDEX> tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_;
362 
363  auto myMap = mpReader_->getThresholdMap();
364 
365  std::vector<std::string> alignablesList;
366  for (auto it = myMap.begin(); it != myMap.end(); ++it) {
367  alignablesList.push_back(it->first);
368  }
369 
370  for (auto& alignable : alignablesList) {
371  int detIndex = getIndexFromString(alignable);
372 
373  Xcut_[detIndex] = myMap[alignable].getXcut();
374  sigXcut_[detIndex] = myMap[alignable].getSigXcut();
375  maxMoveXcut_[detIndex] = myMap[alignable].getMaxMoveXcut();
376  maxErrorXcut_[detIndex] = myMap[alignable].getErrorXcut();
377 
378  Ycut_[detIndex] = myMap[alignable].getYcut();
379  sigYcut_[detIndex] = myMap[alignable].getSigYcut();
380  maxMoveYcut_[detIndex] = myMap[alignable].getMaxMoveYcut();
381  maxErrorYcut_[detIndex] = myMap[alignable].getErrorYcut();
382 
383  Zcut_[detIndex] = myMap[alignable].getZcut();
384  sigZcut_[detIndex] = myMap[alignable].getSigZcut();
385  maxMoveZcut_[detIndex] = myMap[alignable].getMaxMoveZcut();
386  maxErrorZcut_[detIndex] = myMap[alignable].getErrorZcut();
387 
388  tXcut_[detIndex] = myMap[alignable].getThetaXcut();
389  sigtXcut_[detIndex] = myMap[alignable].getSigThetaXcut();
390  maxMovetXcut_[detIndex] = myMap[alignable].getMaxMoveThetaXcut();
391  maxErrortXcut_[detIndex] = myMap[alignable].getErrorThetaXcut();
392 
393  tYcut_[detIndex] = myMap[alignable].getThetaYcut();
394  sigtYcut_[detIndex] = myMap[alignable].getSigThetaYcut();
395  maxMovetYcut_[detIndex] = myMap[alignable].getMaxMoveThetaYcut();
396  maxErrortYcut_[detIndex] = myMap[alignable].getErrorThetaYcut();
397 
398  tZcut_[detIndex] = myMap[alignable].getThetaZcut();
399  sigtZcut_[detIndex] = myMap[alignable].getSigThetaZcut();
400  maxMovetZcut_[detIndex] = myMap[alignable].getMaxMoveThetaZcut();
401  maxErrortZcut_[detIndex] = myMap[alignable].getErrorThetaZcut();
402  }
403 
405  h_xPos_HG, Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_, mpReader_->getXobs_HG(), mpReader_->getXobsErr_HG());
407  tXcut_,
408  sigtXcut_,
409  maxMovetXcut_,
410  maxErrortXcut_,
411  mpReader_->getTXobs_HG(),
412  mpReader_->getTXobsErr_HG());
413 
415  h_yPos_HG, Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_, mpReader_->getYobs_HG(), mpReader_->getYobsErr_HG());
417  tYcut_,
418  sigtYcut_,
419  maxMovetYcut_,
420  maxErrortYcut_,
421  mpReader_->getTYobs_HG(),
422  mpReader_->getTYobsErr_HG());
423 
425  h_zPos_HG, Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_, mpReader_->getZobs_HG(), mpReader_->getZobsErr_HG());
427  tZcut_,
428  sigtZcut_,
429  maxMovetZcut_,
430  maxErrortZcut_,
431  mpReader_->getTZobs_HG(),
432  mpReader_->getTZobsErr_HG());
433 }
434 
435 void MillePedeDQMModule ::fillExpertHisto_HG(std::map<std::string, MonitorElement*>& histo_map,
436  const std::array<double, SIZE_INDEX>& cut,
437  const std::array<double, SIZE_INDEX>& sigCut,
438  const std::array<double, SIZE_INDEX>& maxMoveCut,
439  const std::array<double, SIZE_INDEX>& maxErrorCut,
440  const std::array<double, SIZE_HG_STRUCTS>& obs,
441  const std::array<double, SIZE_HG_STRUCTS>& obsErr) {
442  int currentStart = 0;
443  int bin = 0;
444  double max_ = 0;
445 
446  for (const auto& layer : layerVec) {
447  TH1F* histo_0 = histo_map[layer.first]->getTH1F();
448 
449  max_ = -1;
450  for (int i = currentStart; i < (currentStart + layer.second); ++i) {
451  // first obs.size() bins for observed movements
452  bin = i - currentStart + 1;
453 
454  // fill observed values
455  histo_0->SetBinContent(bin, obs[i]);
456  histo_0->SetBinError(bin, obsErr[i]);
457 
458  if (std::abs(obs[i]) > max_) {
459  max_ = std::abs(obs[i]);
460  }
461  }
462 
463  // five extra bins at the end, one empty, one with threshold, one with sigCut, one with maxMoveCut, one with MaxErrorCut
464  histo_0->SetBinContent(bin + 1, 0);
465  histo_0->SetBinError(bin + 1, 0);
466 
467  int detIndex;
468  if (layer.first.find("Disk") != std::string::npos) {
469  // 7 is the detId for panels, see getIndexFromString
470  detIndex = 7;
471  histo_0->GetXaxis()->SetTitle("Panel");
472  } else {
473  // 6 is the detId for ladders, see getIndexFromString
474  detIndex = 6;
475  histo_0->GetXaxis()->SetTitle("Ladder");
476  }
477 
478  histo_0->SetBinContent(bin + 2, cut[detIndex]);
479  histo_0->SetBinError(bin + 2, 0);
480  histo_0->SetBinContent(bin + 3, sigCut[detIndex]);
481  histo_0->SetBinError(bin + 3, 0);
482  histo_0->SetBinContent(bin + 4, maxMoveCut[detIndex]);
483  histo_0->SetBinError(bin + 4, 0);
484  histo_0->SetBinContent(bin + 5, maxErrorCut[detIndex]);
485  histo_0->SetBinError(bin + 5, 0);
486 
487  // always scale so the cutoff is visible
488  max_ = std::max(cut[detIndex] * 1.2, max_);
489 
490  histo_0->SetMinimum(-(max_) * 1.2);
491  histo_0->SetMaximum(max_ * 1.2);
492 
493  currentStart += layer.second;
494  }
495 }
496 
498  bool changed{false};
499 
501  changed = true;
503  changed = true;
505  changed = true;
506 
507  return changed;
508 }
509 
511  if (alignableId == "TPBHalfBarrelXminus") {
512  return 3;
513  } else if (alignableId == "TPBHalfBarrelXplus") {
514  return 2;
515  } else if (alignableId == "TPEHalfCylinderXminusZminus") {
516  return 1;
517  } else if (alignableId == "TPEHalfCylinderXplusZminus") {
518  return 0;
519  } else if (alignableId == "TPEHalfCylinderXminusZplus") {
520  return 5;
521  } else if (alignableId == "TPEHalfCylinderXplusZplus") {
522  return 4;
523  } else if (alignableId.rfind("TPBLadder", 0) == 0) {
524  return 6;
525  } else if (alignableId.rfind("TPEPanel", 0) == 0) {
526  return 7;
527  } else {
528  throw cms::Exception("LogicError") << "@SUB=MillePedeDQMModule::getIndexFromString\n"
529  << "Retrieving conversion for not supported Alignable partition" << alignableId;
530  }
531 }
void fillExpertHisto(MonitorElement *histo, const std::array< double, SIZE_INDEX > &cut, const std::array< double, SIZE_INDEX > &sigCut, const std::array< double, SIZE_INDEX > &maxMoveCut, const std::array< double, SIZE_INDEX > &maxErrorCut, const std::array< double, SIZE_LG_STRUCTS > &obs, const std::array< double, SIZE_LG_STRUCTS > &obsErr)
const edm::ParameterSet mpReaderConfig_
void fillStatusHisto(MonitorElement *statusHisto)
std::map< std::string, MonitorElement * > h_xRot_HG
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< std::pair< std::string, int > > layerVec
std::map< std::string, MonitorElement * > h_yPos_HG
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:35
bool setupChanged(const edm::EventSetup &)
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
MonitorElement * h_zPos
DQM Plotter for PCL-Alignment.
TrackerGeometry * build(const GeometricDet *gd, const PTrackerParameters &ptp, const TrackerTopology *tTopo)
const edm::ESGetToken< GeometricDet, IdealGeometryRecord > gDetToken_
Definition: config.py:1
std::map< std::string, MonitorElement * > h_zRot_HG
MonitorElement * binariesAvalaible
MonitorElement * bookString(TString const &name, TString const &value, FUNC onbooking=NOOP())
Definition: DQMStore.h:87
void Fill(long long x)
virtual TH2F * getTH2F() const
edm::ESWatcher< IdealGeometryRecord > watchIdealGeometryRcd_
MonitorElement * isVetoed
MillePedeDQMModule(const edm::ParameterSet &)
std::map< std::string, MonitorElement * > h_yRot_HG
void bookHistograms(DQMStore::IBooker &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
std::map< std::string, MonitorElement * > h_xPos_HG
MonitorElement * h_xPos
void beginRun(const edm::Run &, const edm::EventSetup &) override
MonitorElement * h_yPos
Log< level::Info, false > LogInfo
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
std::shared_ptr< PixelTopologyMap > pixelTopologyMap_
void addUntrackedParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:193
const edm::ESGetToken< AlignPCLThresholdsHG, AlignPCLThresholdsHGRcd > aliThrToken_
MonitorElement * statusResults
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
MonitorElement * bookInt(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:73
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
edm::ESWatcher< TrackerTopologyRcd > watchTrackerTopologyRcd_
int currentStart
Definition: mps_splice.py:66
edm::ESWatcher< PTrackerParametersRcd > watchPTrackerParametersRcd_
MonitorElement * h_xRot
HLT enums.
const edm::ESGetToken< PTrackerParameters, PTrackerParametersRcd > ptpToken_
std::unique_ptr< AlignableTracker > tracker_
std::unique_ptr< MillePedeFileReader > mpReader_
int getIndexFromString(const std::string &alignableId)
#define get
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * h_zRot
MonitorElement * exitCode
std::map< std::string, MonitorElement * > h_zPos_HG
~MillePedeDQMModule() override
const std::string outputFolder_
void fillStatusHistoHG(MonitorElement *statusHisto)
MonitorElement * h_yRot
Definition: Run.h:45
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
void fillExpertHisto_HG(std::map< std::string, MonitorElement *> &histo_map, const std::array< double, SIZE_INDEX > &cut, const std::array< double, SIZE_INDEX > &sigCut, const std::array< double, SIZE_INDEX > &maxMoveCut, const std::array< double, SIZE_INDEX > &maxErrorCut, const std::array< double, SIZE_HG_STRUCTS > &obs, const std::array< double, SIZE_HG_STRUCTS > &obsErr)