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