CMS 3D CMS Logo

TrackerAlignment_PayloadInspector.cc
Go to the documentation of this file.
1 
10 
14 
15 // the data format of the condition to be inspected
19 
20 //#include "CondFormats/Alignment/interface/Definitions.h"
21 #include "CLHEP/Vector/RotationInterfaces.h"
23 
24 // needed for the tracker map
26 
27 // needed for mapping
32 
33 #include <boost/range/adaptor/indexed.hpp>
34 #include <iomanip> // std::setprecision
35 #include <iostream>
36 #include <memory>
37 #include <sstream>
38 
39 // include ROOT
40 #include "TH2F.h"
41 #include "TGaxis.h"
42 #include "TLegend.h"
43 #include "TCanvas.h"
44 #include "TLine.h"
45 #include "TStyle.h"
46 #include "TLatex.h"
47 #include "TPave.h"
48 #include "TMarker.h"
49 #include "TPaveStats.h"
50 
51 namespace {
52 
53  using namespace cond::payloadInspector;
54 
55  // M.M. 2017/09/29
56  // Hardcoded Tracker Global Position Record
57  // Without accessing the ES, it is not possible to access to the GPR with the PI technology,
58  // so this needs to be hardcoded.
59  // Anyway it is not likely to change until a new Tracker is installed.
60  // Details at:
61  // - https://indico.cern.ch/event/238026/contributions/513928/attachments/400000/556192/mm_TkAlMeeting_28_03_2013.pdf
62  // - https://twiki.cern.ch/twiki/bin/view/CMS/TkAlignmentPixelPosition
63 
64  const std::map<AlignmentPI::coordinate, float> hardcodeGPR = {
65  {AlignmentPI::t_x, -9.00e-02}, {AlignmentPI::t_y, -1.10e-01}, {AlignmentPI::t_z, -1.70e-01}};
66 
67  //*******************************************/
68  // Size of the movement over all partitions,
69  // one at a time
70  //******************************************//
71 
72  template <int ntags, IOVMultiplicity nIOVs, bool doOnlyPixel>
73  class TrackerAlignmentCompareAll : public PlotImage<Alignments, nIOVs, ntags> {
74  public:
75  TrackerAlignmentCompareAll()
76  : PlotImage<Alignments, nIOVs, ntags>("comparison of all coordinates between two geometries") {}
77 
78  bool fill() override {
79  TGaxis::SetExponentOffset(-0.12, 0.01, "y"); // Y offset
80 
81  // trick to deal with the multi-ioved tag and two tag case at the same time
82  auto theIOVs = PlotBase::getTag<0>().iovs;
83  auto tagname1 = PlotBase::getTag<0>().name;
84  std::string tagname2 = "";
85  auto firstiov = theIOVs.front();
86  std::tuple<cond::Time_t, cond::Hash> lastiov;
87 
88  // we don't support (yet) comparison with more than 2 tags
89  assert(this->m_plotAnnotations.ntags < 3);
90 
91  if (this->m_plotAnnotations.ntags == 2) {
92  auto tag2iovs = PlotBase::getTag<1>().iovs;
93  tagname2 = PlotBase::getTag<1>().name;
94  lastiov = tag2iovs.front();
95  } else {
96  lastiov = theIOVs.back();
97  }
98 
99  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
100  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
101 
102  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
103  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
104 
105  std::vector<AlignTransform> ref_ali = first_payload->m_align;
106  std::vector<AlignTransform> target_ali = last_payload->m_align;
107 
108  // Use remove_if along with a lambda expression to remove elements based on the condition (subid > 2)
109  if (doOnlyPixel) {
110  ref_ali.erase(std::remove_if(ref_ali.begin(),
111  ref_ali.end(),
112  [](const AlignTransform &transform) {
113  int subid = DetId(transform.rawId()).subdetId();
114  return subid > 2;
115  }),
116  ref_ali.end());
117 
118  target_ali.erase(std::remove_if(target_ali.begin(),
119  target_ali.end(),
120  [](const AlignTransform &transform) {
121  int subid = DetId(transform.rawId()).subdetId();
122  return subid > 2;
123  }),
124  target_ali.end());
125  }
126  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 2000, 1200);
127  canvas.Divide(3, 2);
128 
129  if (ref_ali.size() != target_ali.size()) {
130  edm::LogError("TrackerAlignment_PayloadInspector")
131  << "the size of the reference alignment (" << ref_ali.size()
132  << ") is different from the one of the target (" << target_ali.size()
133  << ")! You are probably trying to compare different underlying geometries. Exiting";
134  return false;
135  }
136 
137  const bool ph2 = (ref_ali.size() > AlignmentPI::phase1size);
138 
139  // check that the geomtery is a tracker one
140  const char *path_toTopologyXML = nullptr;
141  if (ph2) {
142  path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
143  } else {
144  path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
145  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
146  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
147  }
148 
149  TrackerTopology tTopo =
151 
152  for (const auto &ali : ref_ali) {
153  auto mydetid = ali.rawId();
154  if (DetId(mydetid).det() != DetId::Tracker) {
155  edm::LogWarning("TrackerAlignment_PayloadInspector")
156  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
157  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
158  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
159  return false;
160  }
161  }
162 
163  const std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
169 
170  std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
171 
172  // generate the map of histograms
173  for (const auto &coord : coords) {
174  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
175  std::string unit =
176  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
177 
178  diffs[coord] = std::make_unique<TH1F>(Form("comparison_%s", s_coord.c_str()),
179  Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
180  ref_ali.size(),
181  -0.5,
182  ref_ali.size() - 0.5);
183  }
184 
185  // fill all the histograms together
186  std::map<int, AlignmentPI::partitions> boundaries;
187  boundaries.insert({0, AlignmentPI::BPix}); // always start with BPix, not filled in the loop
188  AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs);
189 
190  unsigned int subpad{1};
191  TLegend legend = TLegend(0.17, 0.84, 0.95, 0.94);
192  legend.SetTextSize(0.023);
193  for (const auto &coord : coords) {
194  canvas.cd(subpad);
195  canvas.cd(subpad)->SetTopMargin(0.06);
196  canvas.cd(subpad)->SetLeftMargin(0.17);
197  canvas.cd(subpad)->SetRightMargin(0.05);
198  canvas.cd(subpad)->SetBottomMargin(0.15);
199  AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
200  auto max = diffs[coord]->GetMaximum();
201  auto min = diffs[coord]->GetMinimum();
202  auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
203  if (range == 0.f)
204  range = 0.1;
205  //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
206 
207  diffs[coord]->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
208  diffs[coord]->GetYaxis()->SetTitleOffset(1.5);
209  diffs[coord]->SetMarkerStyle(20);
210  diffs[coord]->SetMarkerSize(0.5);
211  diffs[coord]->Draw("P");
212 
213  if (subpad == 1) { /* fill the legend only at the first pass */
214  if (this->m_plotAnnotations.ntags == 2) {
215  legend.SetHeader("#bf{Two Tags Comparison}", "C"); // option "C" allows to center the header
216  legend.AddEntry(
217  diffs[coord].get(),
218  ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}")
219  .c_str(),
220  "PL");
221  } else {
222  legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C"); // option "C" allows to center the header
223  legend.AddEntry(diffs[coord].get(),
224  ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
225  "PL");
226  }
227  }
228  subpad++;
229  }
230 
231  canvas.Update();
232  canvas.cd();
233  canvas.Modified();
234 
235  TLine l[6][boundaries.size()];
236  TLatex tSubdet[6];
237  for (unsigned int i = 0; i < 6; i++) {
238  tSubdet[i].SetTextColor(kRed);
239  tSubdet[i].SetNDC();
240  tSubdet[i].SetTextAlign(21);
241  tSubdet[i].SetTextSize(doOnlyPixel ? 0.05 : 0.03);
242  tSubdet[i].SetTextAngle(90);
243  }
244 
245  subpad = 0;
246  for (const auto &coord : coords) {
247  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
248  canvas.cd(subpad + 1);
249  for (const auto &line : boundaries | boost::adaptors::indexed(0)) {
250  const auto &index = line.index();
251  const auto value = line.value();
252  l[subpad][index] = TLine(diffs[coord]->GetBinLowEdge(value.first),
253  canvas.cd(subpad + 1)->GetUymin(),
254  diffs[coord]->GetBinLowEdge(value.first),
255  canvas.cd(subpad + 1)->GetUymax() * 0.84);
256  l[subpad][index].SetLineWidth(1);
257  l[subpad][index].SetLineStyle(9);
258  l[subpad][index].SetLineColor(2);
259  l[subpad][index].Draw("same");
260  }
261 
262  for (const auto &elem : boundaries | boost::adaptors::indexed(0)) {
263  const auto &lm = canvas.cd(subpad + 1)->GetLeftMargin();
264  const auto &rm = 1 - canvas.cd(subpad + 1)->GetRightMargin();
265  const auto &frac = float(elem.value().first) / ref_ali.size();
266 
267  LogDebug("TrackerAlignmentCompareAll")
268  << __PRETTY_FUNCTION__ << " left margin: " << lm << " right margin: " << rm << " fraction: " << frac;
269 
270  float theX_ = lm + (rm - lm) * frac + ((elem.index() > 0 || doOnlyPixel) ? 0.025 : 0.01);
271 
272  tSubdet[subpad].DrawLatex(
273  theX_, 0.23, Form("%s", AlignmentPI::getStringFromPart(elem.value().second, /*is phase2?*/ ph2).c_str()));
274  }
275 
276  auto ltx = TLatex();
277  ltx.SetTextFont(62);
278  ltx.SetTextSize(0.042);
279  ltx.SetTextAlign(11);
280  ltx.DrawLatexNDC(canvas.cd(subpad + 1)->GetLeftMargin(),
281  1 - canvas.cd(subpad + 1)->GetTopMargin() + 0.01,
282  ("Tracker Alignment Compare : #color[4]{" + s_coord + "}").c_str());
283  legend.Draw("same");
284  subpad++;
285  } // loop on the coordinates
286 
287  std::string fileName(this->m_imageFileName);
288  canvas.SaveAs(fileName.c_str());
289  //canvas.SaveAs("out.root");
290 
291  return true;
292  }
293  };
294 
295  typedef TrackerAlignmentCompareAll<1, MULTI_IOV, false> TrackerAlignmentComparatorSingleTag;
296  typedef TrackerAlignmentCompareAll<2, SINGLE_IOV, false> TrackerAlignmentComparatorTwoTags;
297 
298  typedef TrackerAlignmentCompareAll<1, MULTI_IOV, true> PixelAlignmentComparatorSingleTag;
299  typedef TrackerAlignmentCompareAll<2, SINGLE_IOV, true> PixelAlignmentComparatorTwoTags;
300 
301  //*******************************************//
302  // Size of the movement over all partitions,
303  // one coordinate (x,y,z,...) at a time
304  //******************************************//
305 
306  template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
307  class TrackerAlignmentComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
308  public:
309  TrackerAlignmentComparatorBase()
310  : PlotImage<Alignments, nIOVs, ntags>("comparison of " + AlignmentPI::getStringFromCoordinate(coord) +
311  " coordinate between two geometries") {}
312 
313  bool fill() override {
314  TGaxis::SetExponentOffset(-0.12, 0.01, "y"); // Y offset
315 
316  // trick to deal with the multi-ioved tag and two tag case at the same time
317  auto theIOVs = PlotBase::getTag<0>().iovs;
318  auto tagname1 = PlotBase::getTag<0>().name;
319  std::string tagname2 = "";
320  auto firstiov = theIOVs.front();
321  std::tuple<cond::Time_t, cond::Hash> lastiov;
322 
323  // we don't support (yet) comparison with more than 2 tags
324  assert(this->m_plotAnnotations.ntags < 3);
325 
326  if (this->m_plotAnnotations.ntags == 2) {
327  auto tag2iovs = PlotBase::getTag<1>().iovs;
328  tagname2 = PlotBase::getTag<1>().name;
329  lastiov = tag2iovs.front();
330  } else {
331  lastiov = theIOVs.back();
332  }
333 
334  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
335  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
336 
337  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
338  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
339 
340  std::vector<AlignTransform> ref_ali = first_payload->m_align;
341  std::vector<AlignTransform> target_ali = last_payload->m_align;
342 
343  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1200, 1200);
344 
345  if (ref_ali.size() != target_ali.size()) {
346  edm::LogError("TrackerAlignment_PayloadInspector")
347  << "the size of the reference alignment (" << ref_ali.size()
348  << ") is different from the one of the target (" << target_ali.size()
349  << ")! You are probably trying to compare different underlying geometries. Exiting";
350  return false;
351  }
352 
353  // check that the geomtery is a tracker one
354  const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
355  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
356  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
357  TrackerTopology tTopo =
359 
360  for (const auto &ali : ref_ali) {
361  auto mydetid = ali.rawId();
362  if (DetId(mydetid).det() != DetId::Tracker) {
363  edm::LogWarning("TrackerAlignment_PayloadInspector")
364  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
365  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
366  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
367  return false;
368  }
369  }
370 
371  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
372  std::string unit =
373  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
374 
375  //std::unique_ptr<TH1F> compare = std::unique_ptr<TH1F>(new TH1F("comparison",Form("Comparison of %s;DetId index; #Delta%s %s",s_coord.c_str(),s_coord.c_str(),unit.c_str()),ref_ali.size(),-0.5,ref_ali.size()-0.5));
376  std::unique_ptr<TH1F> compare =
377  std::make_unique<TH1F>("comparison",
378  Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
379  ref_ali.size(),
380  -0.5,
381  ref_ali.size() - 0.5);
382 
383  // fill the histograms
384  std::map<int, AlignmentPI::partitions> boundaries;
385  boundaries.insert({0, AlignmentPI::BPix}); // always start with BPix, not filled in the loop
386  AlignmentPI::fillComparisonHistogram(coord, boundaries, ref_ali, target_ali, compare);
387 
388  canvas.cd();
389 
390  canvas.SetTopMargin(0.06);
391  canvas.SetLeftMargin(0.17);
392  canvas.SetRightMargin(0.05);
393  canvas.SetBottomMargin(0.15);
395  auto max = compare->GetMaximum();
396  auto min = compare->GetMinimum();
397  auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
398  if (range == 0.f)
399  range = 0.1;
400  //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
401 
402  compare->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
403  compare->GetYaxis()->SetTitleOffset(1.5);
404  compare->SetMarkerStyle(20);
405  compare->SetMarkerSize(0.5);
406  compare->Draw("P");
407 
408  canvas.Update();
409  canvas.cd();
410 
411  TLine l[boundaries.size()];
412  for (const auto &line : boundaries | boost::adaptors::indexed(0)) {
413  const auto &index = line.index();
414  const auto value = line.value();
415  l[index] = TLine(compare->GetBinLowEdge(value.first),
416  canvas.cd()->GetUymin(),
417  compare->GetBinLowEdge(value.first),
418  canvas.cd()->GetUymax());
419  l[index].SetLineWidth(1);
420  l[index].SetLineStyle(9);
421  l[index].SetLineColor(2);
422  l[index].Draw("same");
423  }
424 
425  TLatex tSubdet;
426  tSubdet.SetNDC();
427  tSubdet.SetTextAlign(21);
428  tSubdet.SetTextSize(0.027);
429  tSubdet.SetTextAngle(90);
430 
431  for (const auto &elem : boundaries) {
432  tSubdet.SetTextColor(kRed);
433  auto myPair = AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(elem.first));
434  float theX_ = elem.first != 0 ? myPair.first + 0.025 : myPair.first + 0.01;
435  const bool isPhase2 = (ref_ali.size() > AlignmentPI::phase1size);
436  tSubdet.DrawLatex(theX_, 0.20, Form("%s", AlignmentPI::getStringFromPart(elem.second, isPhase2).c_str()));
437  }
438 
439  TLegend legend = TLegend(0.17, 0.86, 0.95, 0.94);
440  if (this->m_plotAnnotations.ntags == 2) {
441  legend.SetHeader("#bf{Two Tags Comparison}", "C"); // option "C" allows to center the header
442  legend.AddEntry(
443  compare.get(),
444  ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
445  "PL");
446  } else {
447  legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C"); // option "C" allows to center the header
448  legend.AddEntry(compare.get(),
449  ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
450  "PL");
451  }
452  legend.SetTextSize(0.020);
453  legend.Draw("same");
454 
455  auto ltx = TLatex();
456  ltx.SetTextFont(62);
457  ltx.SetTextSize(0.042);
458  ltx.SetTextAlign(11);
459  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
460  1 - gPad->GetTopMargin() + 0.01,
461  ("Tracker Alignment Compare :#color[4]{" + s_coord + "}").c_str());
462 
463  std::string fileName(this->m_imageFileName);
464  canvas.SaveAs(fileName.c_str());
465 
466  return true;
467  }
468  };
469 
470  template <AlignmentPI::coordinate coord>
471  using TrackerAlignmentCompare = TrackerAlignmentComparatorBase<coord, 1, MULTI_IOV>;
472 
473  template <AlignmentPI::coordinate coord>
474  using TrackerAlignmentCompareTwoTags = TrackerAlignmentComparatorBase<coord, 2, SINGLE_IOV>;
475 
476  typedef TrackerAlignmentCompare<AlignmentPI::t_x> TrackerAlignmentCompareX;
477  typedef TrackerAlignmentCompare<AlignmentPI::t_y> TrackerAlignmentCompareY;
478  typedef TrackerAlignmentCompare<AlignmentPI::t_z> TrackerAlignmentCompareZ;
479 
480  typedef TrackerAlignmentCompare<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlpha;
481  typedef TrackerAlignmentCompare<AlignmentPI::rot_beta> TrackerAlignmentCompareBeta;
482  typedef TrackerAlignmentCompare<AlignmentPI::rot_gamma> TrackerAlignmentCompareGamma;
483 
484  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_x> TrackerAlignmentCompareXTwoTags;
485  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_y> TrackerAlignmentCompareYTwoTags;
486  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_z> TrackerAlignmentCompareZTwoTags;
487 
488  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlphaTwoTags;
489  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_beta> TrackerAlignmentCompareBetaTwoTags;
490  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_gamma> TrackerAlignmentCompareGammaTwoTags;
491 
492  //*******************************************//
493  // Summary canvas per subdetector
494  //******************************************//
495 
496  template <int ntags, IOVMultiplicity nIOVs, AlignmentPI::partitions q>
497  class TrackerAlignmentSummaryBase : public PlotImage<Alignments, nIOVs, ntags> {
498  public:
499  TrackerAlignmentSummaryBase()
500  : PlotImage<Alignments, nIOVs, ntags>("Comparison of all coordinates between two geometries for " +
501  getStringFromPart(q)) {}
502 
503  bool fill() override {
504  // trick to deal with the multi-ioved tag and two tag case at the same time
505  auto theIOVs = PlotBase::getTag<0>().iovs;
506  auto tagname1 = PlotBase::getTag<0>().name;
507  std::string tagname2 = "";
508  auto firstiov = theIOVs.front();
509  std::tuple<cond::Time_t, cond::Hash> lastiov;
510 
511  // we don't support (yet) comparison with more than 2 tags
512  assert(this->m_plotAnnotations.ntags < 3);
513 
514  if (this->m_plotAnnotations.ntags == 2) {
515  auto tag2iovs = PlotBase::getTag<1>().iovs;
516  tagname2 = PlotBase::getTag<1>().name;
517  lastiov = tag2iovs.front();
518  } else {
519  lastiov = theIOVs.back();
520  }
521 
522  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
523  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
524 
525  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
526  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
527 
528  std::vector<AlignTransform> ref_ali = first_payload->m_align;
529  std::vector<AlignTransform> target_ali = last_payload->m_align;
530 
531  if (ref_ali.size() != target_ali.size()) {
532  edm::LogError("TrackerAlignment_PayloadInspector")
533  << "the size of the reference alignment (" << ref_ali.size()
534  << ") is different from the one of the target (" << target_ali.size()
535  << ")! You are probably trying to compare different underlying geometries. Exiting";
536  return false;
537  }
538 
539  // check that the geomtery is a tracker one
540  const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
541  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
542  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
543  TrackerTopology tTopo =
545 
546  for (const auto &ali : ref_ali) {
547  auto mydetid = ali.rawId();
548  if (DetId(mydetid).det() != DetId::Tracker) {
549  edm::LogWarning("TrackerAlignment_PayloadInspector")
550  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
551  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
552  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
553  return false;
554  }
555  }
556 
557  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1800, 1200);
558  canvas.Divide(3, 2);
559 
560  std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
561  std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
567 
568  for (const auto &coord : coords) {
569  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
570  std::string unit =
571  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
572 
573  diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s", s_coord.c_str()),
574  Form(";#Delta%s %s;n. of modules", s_coord.c_str(), unit.c_str()),
575  1001,
576  -500.5,
577  500.5);
578  }
579 
580  // fill the comparison histograms
581  std::map<int, AlignmentPI::partitions> boundaries;
582  AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs, true, q);
583 
584  int c_index = 1;
585 
586  //TLegend (Double_t x1, Double_t y1, Double_t x2, Double_t y2, const char *header="", Option_t *option="brNDC")
587  auto legend = std::make_unique<TLegend>(0.14, 0.88, 0.96, 0.99);
588  if (this->m_plotAnnotations.ntags == 2) {
589  legend->SetHeader("#bf{Two Tags Comparison}", "C"); // option "C" allows to center the header
590  legend->AddEntry(
591  diffs[AlignmentPI::t_x].get(),
592  ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
593  "PL");
594  } else {
595  legend->SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C"); // option "C" allows to center the header
596  legend->AddEntry(diffs[AlignmentPI::t_x].get(),
597  ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
598  "PL");
599  }
600  legend->SetTextSize(0.025);
601 
602  for (const auto &coord : coords) {
603  canvas.cd(c_index)->SetLogy();
604  canvas.cd(c_index)->SetTopMargin(0.01);
605  canvas.cd(c_index)->SetBottomMargin(0.15);
606  canvas.cd(c_index)->SetLeftMargin(0.14);
607  canvas.cd(c_index)->SetRightMargin(0.04);
608  diffs[coord]->SetLineWidth(2);
609  AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
610 
611  //float x_max = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindLastBinAbove(0.));
612  //float x_min = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindFirstBinAbove(0.));
613  //float extremum = std::abs(x_max) > std::abs(x_min) ? std::abs(x_max) : std::abs(x_min);
614  //diffs[coord]->GetXaxis()->SetRangeUser(-extremum*2,extremum*2);
615 
616  int i_max = diffs[coord]->FindLastBinAbove(0.);
617  int i_min = diffs[coord]->FindFirstBinAbove(0.);
618  diffs[coord]->GetXaxis()->SetRange(std::max(1, i_min - 10), std::min(i_max + 10, diffs[coord]->GetNbinsX()));
619  diffs[coord]->SetMaximum(diffs[coord]->GetMaximum() * 5);
620  diffs[coord]->Draw("HIST");
621  AlignmentPI::makeNiceStats(diffs[coord].get(), q, kBlack);
622 
623  legend->Draw("same");
624 
625  c_index++;
626  }
627 
628  std::string fileName(this->m_imageFileName);
629  canvas.SaveAs(fileName.c_str());
630 
631  return true;
632  }
633  };
634 
635  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
636  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
637  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
638 
639  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTID;
640  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
641  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
642 
643  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPixTwoTags;
644  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPixTwoTags;
645  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIBTwoTags;
646 
647  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTIDTwoTags;
648  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOBTwoTags;
649  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTECTwoTags;
650 
651  /************************************************
652  Full Pixel Tracker Map Comparison of coordinates
653  *************************************************/
654  template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
655  class PixelAlignmentComparisonMapBase : public PlotImage<Alignments, nIOVs, ntags> {
656  public:
657  PixelAlignmentComparisonMapBase()
658  : PlotImage<Alignments, nIOVs, ntags>("SiPixel Comparison Map of " +
660  label_ = "PixelAlignmentComparisonMap" + AlignmentPI::getStringFromCoordinate(coord);
661  payloadString = "Tracker Alignment";
662  }
663 
664  bool fill() override {
665  gStyle->SetPalette(1);
666 
667  // trick to deal with the multi-ioved tag and two tag case at the same time
668  auto theIOVs = PlotBase::getTag<0>().iovs;
669  auto tagname1 = PlotBase::getTag<0>().name;
670  std::string tagname2 = "";
671  auto firstiov = theIOVs.front();
672  std::tuple<cond::Time_t, cond::Hash> lastiov;
673 
674  // we don't support (yet) comparison with more than 2 tags
675  assert(this->m_plotAnnotations.ntags < 3);
676 
677  if (this->m_plotAnnotations.ntags == 2) {
678  auto tag2iovs = PlotBase::getTag<1>().iovs;
679  tagname2 = PlotBase::getTag<1>().name;
680  lastiov = tag2iovs.front();
681  } else {
682  lastiov = theIOVs.back();
683  }
684 
685  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
686  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
687 
688  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
689  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
690 
691  const std::vector<AlignTransform> &ref_ali = first_payload->m_align;
692  const std::vector<AlignTransform> &target_ali = last_payload->m_align;
693 
694  if (last_payload.get() && first_payload.get()) {
695  Phase1PixelSummaryMap fullMap(
696  "",
697  fmt::sprintf("%s %s", payloadString, AlignmentPI::getStringFromCoordinate(coord)),
698  fmt::sprintf(
699  "#Delta %s [%s]", AlignmentPI::getStringFromCoordinate(coord), (coord <= 3) ? "#mum" : "mrad"));
700  fullMap.createTrackerBaseMap();
701 
702  if (this->isPhase0(ref_ali) || this->isPhase0(target_ali)) {
703  edm::LogError(label_) << "Pixel Tracker Alignment maps are not supported for non-Phase1 Pixel geometries !";
704  TCanvas canvas("Canv", "Canv", 1200, 1000);
706  std::string fileName(this->m_imageFileName);
707  canvas.SaveAs(fileName.c_str());
708  return false;
709  }
710 
711  // fill the map of differences
712  std::map<uint32_t, double> diffPerDetid;
713  this->fillPerDetIdDiff(coord, ref_ali, target_ali, diffPerDetid);
714 
715  // now fill the tracker map
716  for (const auto &elem : diffPerDetid) {
717  // reject Strips
718  int subid = DetId(elem.first).subdetId();
719  if (subid > 2) {
720  continue;
721  }
722  fullMap.fillTrackerMap(elem.first, elem.second);
723  }
724 
725  // limit the axis range (in case of need)
726  //fullMap.setZAxisRange(-50.f,50.f);
727 
728  TCanvas canvas("Canv", "Canv", 3000, 2000);
729  fullMap.printTrackerMap(canvas);
730 
731  auto ltx = TLatex();
732  ltx.SetTextFont(62);
733  ltx.SetTextSize(0.025);
734  ltx.SetTextAlign(11);
735 
736  ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
737  gPad->GetBottomMargin() + 0.01,
738  ("#color[4]{" + tagname1 + "}, IOV: #color[4]{" + firstIOVsince + "} vs #color[4]{" +
739  tagname2 + "}, IOV: #color[4]{" + lastIOVsince + "}")
740  .c_str());
741 
742  std::string fileName(this->m_imageFileName);
743  canvas.SaveAs(fileName.c_str());
744  }
745  return true;
746  }
747 
748  protected:
749  std::string payloadString;
750  std::string label_;
751 
752  private:
753  //_________________________________________________
754  bool isPhase0(std::vector<AlignTransform> theAlis) {
757  const auto &p0detIds = reader.getAllDetIds();
758 
759  std::vector<uint32_t> ownDetIds;
760  std::transform(theAlis.begin(), theAlis.end(), std::back_inserter(ownDetIds), [](AlignTransform ali) -> uint32_t {
761  return ali.rawId();
762  });
763 
764  for (const auto &det : ownDetIds) {
765  // if found at least one phase-0 detId early return
766  if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
767  return true;
768  }
769  }
770  return false;
771  }
772 
773  /*--------------------------------------------------------------------*/
774  void fillPerDetIdDiff(const AlignmentPI::coordinate &myCoord,
775  const std::vector<AlignTransform> &ref_ali,
776  const std::vector<AlignTransform> &target_ali,
777  std::map<uint32_t, double> &diff)
778  /*--------------------------------------------------------------------*/
779  {
780  for (unsigned int i = 0; i < ref_ali.size(); i++) {
781  uint32_t detid = ref_ali[i].rawId();
782  if (ref_ali[i].rawId() == target_ali[i].rawId()) {
783  CLHEP::HepRotation target_rot(target_ali[i].rotation());
784  CLHEP::HepRotation ref_rot(ref_ali[i].rotation());
785 
786  align::RotationType target_ROT(target_rot.xx(),
787  target_rot.xy(),
788  target_rot.xz(),
789  target_rot.yx(),
790  target_rot.yy(),
791  target_rot.yz(),
792  target_rot.zx(),
793  target_rot.zy(),
794  target_rot.zz());
795 
796  align::RotationType ref_ROT(ref_rot.xx(),
797  ref_rot.xy(),
798  ref_rot.xz(),
799  ref_rot.yx(),
800  ref_rot.yy(),
801  ref_rot.yz(),
802  ref_rot.zx(),
803  ref_rot.zy(),
804  ref_rot.zz());
805 
806  const std::vector<double> deltaRot = {
807  ::deltaPhi(align::toAngles(target_ROT)[0], align::toAngles(ref_ROT)[0]),
808  ::deltaPhi(align::toAngles(target_ROT)[1], align::toAngles(ref_ROT)[1]),
809  ::deltaPhi(align::toAngles(target_ROT)[2], align::toAngles(ref_ROT)[2])};
810 
811  const auto &deltaTrans = target_ali[i].translation() - ref_ali[i].translation();
812 
813  switch (myCoord) {
814  case AlignmentPI::t_x:
815  diff.insert({detid, deltaTrans.x() * AlignmentPI::cmToUm});
816  break;
817  case AlignmentPI::t_y:
818  diff.insert({detid, deltaTrans.y() * AlignmentPI::cmToUm});
819  break;
820  case AlignmentPI::t_z:
821  diff.insert({detid, deltaTrans.z() * AlignmentPI::cmToUm});
822  break;
824  diff.insert({detid, deltaRot[0] * AlignmentPI::tomRad});
825  break;
827  diff.insert({detid, deltaRot[1] * AlignmentPI::tomRad});
828  break;
830  diff.insert({detid, deltaRot[2] * AlignmentPI::tomRad});
831  break;
832  default:
833  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << myCoord << std::endl;
834  break;
835  } // switch on the coordinate
836  } // check on the same detID
837  } // loop on the components
838  }
839  };
840 
841  template <AlignmentPI::coordinate coord>
842  using PixelAlignmentCompareMap = PixelAlignmentComparisonMapBase<coord, 1, MULTI_IOV>;
843 
844  template <AlignmentPI::coordinate coord>
845  using PixelAlignmentCompareMapTwoTags = PixelAlignmentComparisonMapBase<coord, 2, SINGLE_IOV>;
846 
847  typedef PixelAlignmentCompareMap<AlignmentPI::t_x> PixelAlignmentCompareMapX;
848  typedef PixelAlignmentCompareMap<AlignmentPI::t_y> PixelAlignmentCompareMapY;
849  typedef PixelAlignmentCompareMap<AlignmentPI::t_z> PixelAlignmentCompareMapZ;
850 
851  typedef PixelAlignmentCompareMap<AlignmentPI::rot_alpha> PixelAlignmentCompareMapAlpha;
852  typedef PixelAlignmentCompareMap<AlignmentPI::rot_beta> PixelAlignmentCompareMapBeta;
853  typedef PixelAlignmentCompareMap<AlignmentPI::rot_gamma> PixelAlignmentCompareMapGamma;
854 
855  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_x> PixelAlignmentCompareMapXTwoTags;
856  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_y> PixelAlignmentCompareMapYTwoTags;
857  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_z> PixelAlignmentCompareMapZTwoTags;
858 
859  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_alpha> PixelAlignmentCompareMapAlphaTwoTags;
860  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_beta> PixelAlignmentCompareMapBetaTwoTags;
861  typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_gamma> PixelAlignmentCompareMapGammaTwoTags;
862 
863  //*******************************************//
864  // History of the position of the BPix Barycenter
865  //******************************************//
866 
867  template <AlignmentPI::coordinate coord>
868  class BPixBarycenterHistory : public HistoryPlot<Alignments, float> {
869  public:
870  BPixBarycenterHistory()
872  " Barrel Pixel " + AlignmentPI::getStringFromCoordinate(coord) + " positions vs time",
873  AlignmentPI::getStringFromCoordinate(coord) + " position [cm]") {}
874  ~BPixBarycenterHistory() override = default;
875 
876  float getFromPayload(Alignments &payload) override {
877  std::vector<AlignTransform> alignments = payload.m_align;
878 
879  float barycenter = 0.;
880  float nmodules(0.);
881  for (const auto &ali : alignments) {
882  if (DetId(ali.rawId()).det() != DetId::Tracker) {
883  edm::LogWarning("TrackerAlignment_PayloadInspector")
884  << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
885  << " is different from " << DetId::Tracker << " - terminating ";
886  return false;
887  }
888 
889  int subid = DetId(ali.rawId()).subdetId();
890  if (subid != PixelSubdetector::PixelBarrel)
891  continue;
892 
893  nmodules++;
894  switch (coord) {
895  case AlignmentPI::t_x:
896  barycenter += (ali.translation().x());
897  break;
898  case AlignmentPI::t_y:
899  barycenter += (ali.translation().y());
900  break;
901  case AlignmentPI::t_z:
902  barycenter += (ali.translation().z());
903  break;
904  default:
905  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
906  break;
907  } // switch on the coordinate (only X,Y,Z are interesting)
908  } // ends loop on the alignments
909 
910  edm::LogInfo("TrackerAlignment_PayloadInspector") << "barycenter (" << barycenter << ")/n. modules (" << nmodules
911  << ") = " << barycenter / nmodules << std::endl;
912 
913  // take the mean
914  barycenter /= nmodules;
915 
916  // applied GPR correction to move barycenter to global CMS coordinates
917  barycenter += hardcodeGPR.at(coord);
918 
919  return barycenter;
920 
921  } // payload
922  };
923 
924  typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
925  typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
926  typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
927 
928  /************************************************
929  Display of Tracker Detector barycenters
930  *************************************************/
931  class TrackerAlignmentBarycenters : public PlotImage<Alignments, SINGLE_IOV> {
932  public:
933  TrackerAlignmentBarycenters() : PlotImage<Alignments, SINGLE_IOV>("Display of Tracker Alignment Barycenters") {}
934 
935  bool fill() override {
936  auto tag = PlotBase::getTag<0>();
937  auto iov = tag.iovs.front();
938  const auto &tagname = PlotBase::getTag<0>().name;
939  std::shared_ptr<Alignments> payload = fetchPayload(std::get<1>(iov));
940  unsigned int run = std::get<0>(iov);
941 
942  TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1600, 1000);
943  canvas.cd();
944 
945  canvas.SetTopMargin(0.07);
946  canvas.SetBottomMargin(0.06);
947  canvas.SetLeftMargin(0.15);
948  canvas.SetRightMargin(0.03);
949  canvas.Modified();
950  canvas.SetGrid();
951 
952  auto h2_BarycenterParameters =
953  std::make_unique<TH2F>("Parameters", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
954 
955  auto h2_uncBarycenterParameters =
956  std::make_unique<TH2F>("Parameters2", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
957 
958  h2_BarycenterParameters->SetStats(false);
959  h2_BarycenterParameters->SetTitle(nullptr);
960  h2_uncBarycenterParameters->SetStats(false);
961  h2_uncBarycenterParameters->SetTitle(nullptr);
962 
963  std::vector<AlignTransform> alignments = payload->m_align;
964 
965  isPhase0 = (alignments.size() == AlignmentPI::phase0size) ? true : false;
966 
967  // check that the geomtery is a tracker one
968  const char *path_toTopologyXML = isPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
969  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
970 
971  TrackerTopology tTopo =
973 
974  AlignmentPI::TkAlBarycenters barycenters;
975  // compute uncorrected barycenter
976  barycenters.computeBarycenters(
977  alignments, tTopo, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
978 
979  auto Xbarycenters = barycenters.getX();
980  auto Ybarycenters = barycenters.getY();
981  auto Zbarycenters = barycenters.getZ();
982 
983  // compute barycenter corrected for the GPR
984  barycenters.computeBarycenters(alignments, tTopo, hardcodeGPR);
985 
986  auto c_Xbarycenters = barycenters.getX();
987  auto c_Ybarycenters = barycenters.getY();
988  auto c_Zbarycenters = barycenters.getZ();
989 
990  h2_BarycenterParameters->GetXaxis()->SetBinLabel(1, "X [cm]");
991  h2_BarycenterParameters->GetXaxis()->SetBinLabel(2, "Y [cm]");
992  h2_BarycenterParameters->GetXaxis()->SetBinLabel(3, "Z [cm]");
993  h2_BarycenterParameters->GetXaxis()->SetBinLabel(4, "X_{no GPR} [cm]");
994  h2_BarycenterParameters->GetXaxis()->SetBinLabel(5, "Y_{no GPR} [cm]");
995  h2_BarycenterParameters->GetXaxis()->SetBinLabel(6, "Z_{no GPR} [cm]");
996 
997  bool isLikelyMC(false);
998  int checkX =
999  std::count_if(Xbarycenters.begin(), Xbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1000  int checkY =
1001  std::count_if(Ybarycenters.begin(), Ybarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1002  int checkZ =
1003  std::count_if(Zbarycenters.begin(), Zbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1004 
1005  // if all the coordinate barycenters for both BPix and FPix are below 10um
1006  // this is very likely a MC payload
1007  if ((checkX + checkY + checkZ) == 0 && run == 1)
1008  isLikelyMC = true;
1009 
1010  unsigned int yBin = 6;
1011  for (unsigned int i = 0; i < 6; i++) {
1012  auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
1013  std::string theLabel = getStringFromPart(thePart);
1014  h2_BarycenterParameters->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
1015  if (!isLikelyMC) {
1016  h2_BarycenterParameters->SetBinContent(1, yBin, c_Xbarycenters[i]);
1017  h2_BarycenterParameters->SetBinContent(2, yBin, c_Ybarycenters[i]);
1018  h2_BarycenterParameters->SetBinContent(3, yBin, c_Zbarycenters[i]);
1019  }
1020 
1021  h2_uncBarycenterParameters->SetBinContent(4, yBin, Xbarycenters[i]);
1022  h2_uncBarycenterParameters->SetBinContent(5, yBin, Ybarycenters[i]);
1023  h2_uncBarycenterParameters->SetBinContent(6, yBin, Zbarycenters[i]);
1024  yBin--;
1025  }
1026 
1027  h2_BarycenterParameters->GetXaxis()->LabelsOption("h");
1028  h2_BarycenterParameters->GetYaxis()->SetLabelSize(0.05);
1029  h2_BarycenterParameters->GetXaxis()->SetLabelSize(0.05);
1030  h2_BarycenterParameters->SetMarkerSize(1.5);
1031  h2_BarycenterParameters->Draw("TEXT");
1032 
1033  h2_uncBarycenterParameters->SetMarkerSize(1.5);
1034  h2_uncBarycenterParameters->SetMarkerColor(kRed);
1035  h2_uncBarycenterParameters->Draw("TEXTsame");
1036 
1037  TLatex t1;
1038  t1.SetNDC();
1039  t1.SetTextAlign(26);
1040  t1.SetTextSize(0.045);
1041  t1.DrawLatex(0.5, 0.96, Form("TkAl Barycenters, Tag: #color[4]{%s}, IOV #color[4]{%i}", tagname.c_str(), run));
1042  t1.SetTextSize(0.025);
1043 
1044  std::string fileName(m_imageFileName);
1045  canvas.SaveAs(fileName.c_str());
1046 
1047  return true;
1048  }
1049 
1050  private:
1051  bool isPhase0;
1052  };
1053 
1054  /************************************************
1055  Comparator of Tracker Detector barycenters
1056  *************************************************/
1057  template <int ntags, IOVMultiplicity nIOVs>
1058  class TrackerAlignmentBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
1059  public:
1060  TrackerAlignmentBarycentersComparatorBase()
1061  : PlotImage<Alignments, nIOVs, ntags>("Comparison of Tracker Alignment Barycenters") {}
1062 
1063  bool fill() override {
1064  // trick to deal with the multi-ioved tag and two tag case at the same time
1065  auto theIOVs = PlotBase::getTag<0>().iovs;
1066  auto tagname1 = PlotBase::getTag<0>().name;
1067  std::string tagname2 = "";
1068  auto firstiov = theIOVs.front();
1069  std::tuple<cond::Time_t, cond::Hash> lastiov;
1070 
1071  // we don't support (yet) comparison with more than 2 tags
1072  assert(this->m_plotAnnotations.ntags < 3);
1073 
1074  if (this->m_plotAnnotations.ntags == 2) {
1075  auto tag2iovs = PlotBase::getTag<1>().iovs;
1076  tagname2 = PlotBase::getTag<1>().name;
1077  lastiov = tag2iovs.front();
1078  } else {
1079  lastiov = theIOVs.back();
1080  }
1081 
1082  unsigned int first_run = std::get<0>(firstiov);
1083  unsigned int last_run = std::get<0>(lastiov);
1084 
1085  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
1086  std::vector<AlignTransform> last_alignments = last_payload->m_align;
1087 
1088  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
1089  std::vector<AlignTransform> first_alignments = first_payload->m_align;
1090 
1091  isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
1092  isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
1093 
1094  // check that the geomtery is a tracker one
1095  const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1096  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1097 
1098  TrackerTopology tTopo_f =
1100 
1101  path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1102  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1103 
1104  TrackerTopology tTopo_l =
1106 
1107  TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1200, 800);
1108  canvas.cd();
1109 
1110  canvas.SetTopMargin(0.07);
1111  canvas.SetBottomMargin(0.06);
1112  canvas.SetLeftMargin(0.15);
1113  canvas.SetRightMargin(0.03);
1114  canvas.Modified();
1115  canvas.SetGrid();
1116 
1117  auto h2_BarycenterDiff =
1118  std::make_unique<TH2F>("Parameters diff", "SubDetector Barycenter Difference", 3, 0.0, 3.0, 6, 0, 6.);
1119 
1120  h2_BarycenterDiff->SetStats(false);
1121  h2_BarycenterDiff->SetTitle(nullptr);
1122  h2_BarycenterDiff->GetXaxis()->SetBinLabel(1, "X [#mum]");
1123  h2_BarycenterDiff->GetXaxis()->SetBinLabel(2, "Y [#mum]");
1124  h2_BarycenterDiff->GetXaxis()->SetBinLabel(3, "Z [#mum]");
1125 
1126  AlignmentPI::TkAlBarycenters l_barycenters;
1127  l_barycenters.computeBarycenters(last_alignments, tTopo_l, hardcodeGPR);
1128 
1129  AlignmentPI::TkAlBarycenters f_barycenters;
1130  f_barycenters.computeBarycenters(first_alignments, tTopo_f, hardcodeGPR);
1131 
1132  unsigned int yBin = 6;
1133  for (unsigned int i = 0; i < 6; i++) {
1134  auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
1135  std::string theLabel = getStringFromPart(thePart);
1136  h2_BarycenterDiff->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
1137  h2_BarycenterDiff->SetBinContent(
1138  1, yBin, (l_barycenters.getX()[i] - f_barycenters.getX()[i]) * AlignmentPI::cmToUm);
1139  h2_BarycenterDiff->SetBinContent(
1140  2, yBin, (l_barycenters.getY()[i] - f_barycenters.getY()[i]) * AlignmentPI::cmToUm);
1141  h2_BarycenterDiff->SetBinContent(
1142  3, yBin, (l_barycenters.getZ()[i] - f_barycenters.getZ()[i]) * AlignmentPI::cmToUm);
1143  yBin--;
1144  }
1145 
1146  h2_BarycenterDiff->GetXaxis()->LabelsOption("h");
1147  h2_BarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1148  h2_BarycenterDiff->GetXaxis()->SetLabelSize(0.05);
1149  h2_BarycenterDiff->SetMarkerSize(1.5);
1150  h2_BarycenterDiff->SetMarkerColor(kRed);
1151  h2_BarycenterDiff->Draw("TEXT");
1152 
1153  TLatex t1;
1154  t1.SetNDC();
1155  t1.SetTextAlign(26);
1156  t1.SetTextSize(0.05);
1157  t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters Diff, IOV %i - IOV %i", last_run, first_run));
1158  t1.SetTextSize(0.025);
1159 
1160  std::string fileName(this->m_imageFileName);
1161  canvas.SaveAs(fileName.c_str());
1162 
1163  return true;
1164  }
1165 
1166  private:
1167  bool isInitialPhase0;
1168  bool isFinalPhase0;
1169  };
1170 
1171  using TrackerAlignmentBarycentersCompare = TrackerAlignmentBarycentersComparatorBase<1, MULTI_IOV>;
1172  using TrackerAlignmentBarycentersCompareTwoTags = TrackerAlignmentBarycentersComparatorBase<2, SINGLE_IOV>;
1173 
1174  /************************************************
1175  Comparator of Pixel Tracker Detector barycenters
1176  *************************************************/
1177  template <int ntags, IOVMultiplicity nIOVs>
1178  class PixelBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
1179  public:
1180  PixelBarycentersComparatorBase() : PlotImage<Alignments, nIOVs, ntags>("Comparison of Pixel Barycenters") {}
1181 
1182  bool fill() override {
1183  // trick to deal with the multi-ioved tag and two tag case at the same time
1184  auto theIOVs = PlotBase::getTag<0>().iovs;
1185  auto tagname1 = PlotBase::getTag<0>().name;
1186  std::string tagname2 = "";
1187  auto firstiov = theIOVs.front();
1188  std::tuple<cond::Time_t, cond::Hash> lastiov;
1189 
1190  // we don't support (yet) comparison with more than 2 tags
1191  assert(this->m_plotAnnotations.ntags < 3);
1192 
1193  if (this->m_plotAnnotations.ntags == 2) {
1194  auto tag2iovs = PlotBase::getTag<1>().iovs;
1195  tagname2 = PlotBase::getTag<1>().name;
1196  lastiov = tag2iovs.front();
1197  } else {
1198  lastiov = theIOVs.back();
1199  }
1200 
1201  unsigned int first_run = std::get<0>(firstiov);
1202  unsigned int last_run = std::get<0>(lastiov);
1203 
1204  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
1205  std::vector<AlignTransform> last_alignments = last_payload->m_align;
1206 
1207  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
1208  std::vector<AlignTransform> first_alignments = first_payload->m_align;
1209 
1210  TCanvas canvas("Pixel Barycenter Summary", "Pixel Barycenter summary", 1200, 1200);
1211  canvas.Divide(2, 2);
1212  canvas.cd();
1213 
1214  TLatex t1;
1215  t1.SetNDC();
1216  t1.SetTextAlign(26);
1217  t1.SetTextSize(0.03);
1218  t1.DrawLatex(0.5,
1219  0.97,
1220  ("Pixel Barycenters comparison, IOV: #color[2]{" + std::to_string(first_run) +
1221  "} vs IOV: #color[4]{" + std::to_string(last_run) + "}")
1222  .c_str());
1223  t1.SetTextSize(0.025);
1224 
1225  for (unsigned int c = 1; c <= 4; c++) {
1226  canvas.cd(c)->SetTopMargin(0.07);
1227  canvas.cd(c)->SetBottomMargin(0.12);
1228  canvas.cd(c)->SetLeftMargin(0.15);
1229  canvas.cd(c)->SetRightMargin(0.03);
1230  canvas.cd(c)->Modified();
1231  canvas.cd(c)->SetGrid();
1232  }
1233 
1234  std::array<std::string, 3> structures = {{"FPIX-", "BPIX", "FPIX+"}};
1235  std::array<std::unique_ptr<TH2F>, 3> histos;
1236 
1237  isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
1238  isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
1239 
1240  // check that the geomtery is a tracker one
1241  const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1242  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1243 
1244  TrackerTopology tTopo_f =
1246 
1247  AlignmentPI::TkAlBarycenters myInitialBarycenters;
1248  //myInitialBarycenters.computeBarycenters(first_alignments,tTopo_f,hardcodeGPR);
1249  myInitialBarycenters.computeBarycenters(
1250  first_alignments, tTopo_f, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1251 
1252  path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1253  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1254 
1255  TrackerTopology tTopo_l =
1257 
1258  AlignmentPI::TkAlBarycenters myFinalBarycenters;
1259  //myFinalBarycenters.computeBarycenters(last_alignments,tTopo_l,hardcodeGPR);
1260  myFinalBarycenters.computeBarycenters(
1261  last_alignments, tTopo_l, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1262 
1263  if (isFinalPhase0 != isInitialPhase0) {
1264  edm::LogWarning("TrackerAlignment_PayloadInspector")
1265  << "the size of the reference alignment (" << first_alignments.size()
1266  << ") is different from the one of the target (" << last_alignments.size()
1267  << ")! You are probably trying to compare different underlying geometries.";
1268  }
1269 
1270  unsigned int index(0);
1271  for (const auto &piece : structures) {
1272  const char *name = piece.c_str();
1273  histos[index] = std::make_unique<TH2F>(
1274  name,
1275  Form("%s x-y Barycenter Difference;x_{%s}-x_{TOB} [mm];y_{%s}-y_{TOB} [mm]", name, name, name),
1276  100,
1277  -3.,
1278  3.,
1279  100,
1280  -3.,
1281  3.);
1282 
1283  histos[index]->SetStats(false);
1284  histos[index]->SetTitle(nullptr);
1285  histos[index]->GetYaxis()->SetLabelSize(0.05);
1286  histos[index]->GetXaxis()->SetLabelSize(0.05);
1287  histos[index]->GetYaxis()->SetTitleSize(0.06);
1288  histos[index]->GetXaxis()->SetTitleSize(0.06);
1289  histos[index]->GetYaxis()->CenterTitle();
1290  histos[index]->GetXaxis()->CenterTitle();
1291  histos[index]->GetXaxis()->SetTitleOffset(0.9);
1292  index++;
1293  }
1294 
1295  auto h2_ZBarycenterDiff = std::make_unique<TH2F>(
1296  "Pixel_z_diff", "Pixel z-Barycenter Difference;; z_{Pixel-Ideal} -z_{TOB} [mm]", 3, -0.5, 2.5, 100, -10., 10.);
1297  h2_ZBarycenterDiff->SetStats(false);
1298  h2_ZBarycenterDiff->SetTitle(nullptr);
1299  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(1, "FPIX -");
1300  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(2, "BPIX");
1301  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(3, "FPIX +");
1302  h2_ZBarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1303  h2_ZBarycenterDiff->GetXaxis()->SetLabelSize(0.07);
1304  h2_ZBarycenterDiff->GetYaxis()->SetTitleSize(0.06);
1305  h2_ZBarycenterDiff->GetXaxis()->SetTitleSize(0.06);
1306  h2_ZBarycenterDiff->GetYaxis()->CenterTitle();
1307  h2_ZBarycenterDiff->GetXaxis()->CenterTitle();
1308  h2_ZBarycenterDiff->GetYaxis()->SetTitleOffset(1.1);
1309 
1310  std::function<GlobalPoint(int)> cutFunctorInitial = [&myInitialBarycenters](int index) {
1311  switch (index) {
1312  case 1:
1313  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1314  case 2:
1315  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1316  case 3:
1317  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1318  default:
1319  return GlobalPoint(0, 0, 0);
1320  }
1321  };
1322 
1323  std::function<GlobalPoint(int)> cutFunctorFinal = [&myFinalBarycenters](int index) {
1324  switch (index) {
1325  case 1:
1326  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1327  case 2:
1328  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1329  case 3:
1330  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1331  default:
1332  return GlobalPoint(0, 0, 0);
1333  }
1334  };
1335 
1336  float x0i, x0f, y0i, y0f;
1337 
1338  t1.SetNDC(kFALSE);
1339  t1.SetTextSize(0.047);
1340  for (unsigned int c = 1; c <= 3; c++) {
1341  x0i = cutFunctorInitial(c).x() * 10; // transform cm to mm (x10)
1342  x0f = cutFunctorFinal(c).x() * 10;
1343  y0i = cutFunctorInitial(c).y() * 10;
1344  y0f = cutFunctorFinal(c).y() * 10;
1345 
1346  canvas.cd(c);
1347  histos[c - 1]->Draw();
1348 
1349  COUT << "initial x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0i << "," << y0i << ") mm"
1350  << std::endl;
1351  COUT << "final x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0f << "," << y0f << ") mm"
1352  << std::endl;
1353 
1354  TMarker *initial = new TMarker(x0i, y0i, 21);
1355  TMarker *final = new TMarker(x0f, y0f, 20);
1356 
1357  initial->SetMarkerColor(kRed);
1358  final->SetMarkerColor(kBlue);
1359  initial->SetMarkerSize(2.5);
1360  final->SetMarkerSize(2.5);
1361  t1.SetTextColor(kRed);
1362  initial->Draw();
1363  t1.DrawLatex(x0i, y0i + (y0i > y0f ? 0.3 : -0.5), Form("(%.2f,%.2f)", x0i, y0i));
1364  final->Draw("same");
1365  t1.SetTextColor(kBlue);
1366  t1.DrawLatex(x0f, y0f + (y0i > y0f ? -0.5 : 0.3), Form("(%.2f,%.2f)", x0f, y0f));
1367  }
1368 
1369  // fourth pad is a special case for the z coordinate
1370  canvas.cd(4);
1371  h2_ZBarycenterDiff->Draw();
1372  float z0i, z0f;
1373 
1374  // numbers do agree with https://twiki.cern.ch/twiki/bin/view/CMSPublic/TkAlignmentPerformancePhaseIStartUp17#Pixel_Barycentre_Positions
1375 
1376  std::array<double, 3> hardcodeIdealZPhase0 = {{-41.94909, 0., 41.94909}}; // units are cm
1377  std::array<double, 3> hardcodeIdealZPhase1 = {{-39.82911, 0., 39.82911}}; // units are cm
1378 
1379  for (unsigned int c = 1; c <= 3; c++) {
1380  // less than pretty but needed to remove the z position of the FPix barycenters != 0
1381 
1382  z0i =
1383  (cutFunctorInitial(c).z() - (isInitialPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) *
1384  10; // convert to mm
1385  z0f =
1386  (cutFunctorFinal(c).z() - (isFinalPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) * 10;
1387 
1388  TMarker *initial = new TMarker(c - 1, z0i, 21);
1389  TMarker *final = new TMarker(c - 1, z0f, 20);
1390 
1391  COUT << "initial z " << std::left << std::setw(7) << structures[c - 1] << " " << z0i << " mm" << std::endl;
1392  COUT << "final z " << std::left << std::setw(7) << structures[c - 1] << " " << z0f << " mm" << std::endl;
1393 
1394  initial->SetMarkerColor(kRed);
1395  final->SetMarkerColor(kBlue);
1396  initial->SetMarkerSize(2.5);
1397  final->SetMarkerSize(2.5);
1398  initial->Draw();
1399  t1.SetTextColor(kRed);
1400  t1.DrawLatex(c - 1, z0i + (z0i > z0f ? 1. : -1.5), Form("(%.2f)", z0i));
1401  final->Draw("same");
1402  t1.SetTextColor(kBlue);
1403  t1.DrawLatex(c - 1, z0f + (z0i > z0f ? -1.5 : 1), Form("(%.2f)", z0f));
1404  }
1405 
1406  std::string fileName(this->m_imageFileName);
1407  canvas.SaveAs(fileName.c_str());
1408 
1409  return true;
1410  }
1411 
1412  private:
1413  bool isInitialPhase0;
1414  bool isFinalPhase0;
1415  };
1416 
1417  using PixelBarycentersCompare = PixelBarycentersComparatorBase<1, MULTI_IOV>;
1418  using PixelBarycentersCompareTwoTags = PixelBarycentersComparatorBase<2, SINGLE_IOV>;
1419 
1420 } // namespace
1421 
1423  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentComparatorSingleTag);
1424  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentComparatorTwoTags);
1425  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorSingleTag);
1426  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorTwoTags);
1427  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
1428  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
1429  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
1430  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
1431  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
1432  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
1433  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareXTwoTags);
1434  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareYTwoTags);
1435  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZTwoTags);
1436  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlphaTwoTags);
1437  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBetaTwoTags);
1438  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGammaTwoTags);
1439  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapX);
1440  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapY);
1441  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapZ);
1442  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapAlpha);
1443  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapBeta);
1444  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapGamma);
1445  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapXTwoTags);
1446  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapYTwoTags);
1447  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapZTwoTags);
1448  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapAlphaTwoTags);
1449  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapBetaTwoTags);
1450  PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapGammaTwoTags);
1451  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
1452  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
1453  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
1454  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
1455  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
1456  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
1457  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPixTwoTags);
1458  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPixTwoTags);
1459  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIBTwoTags);
1460  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIDTwoTags);
1461  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOBTwoTags);
1462  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTECTwoTags);
1463  PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
1464  PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
1465  PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
1466  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycenters);
1467  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompare);
1468  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompareTwoTags);
1469  PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompare);
1470  PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompareTwoTags);
1471 }
const std::array< double, 6 > getY()
def rm(path, rec=False)
Definition: eostools.py:363
bool tidIsDoubleSide(const DetId &id) const
void makeNicePlotStyle(TH1 *hist, int color)
void computeBarycenters(const std::vector< AlignTransform > &input, const TrackerTopology &tTopo, const std::map< AlignmentPI::coordinate, float > &GPR)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
reader
Definition: DQM.py:105
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
static std::string to_string(const XMLCh *ch)
static const float cmToUm
void makeNiceStats(TH1F *hist, AlignmentPI::partitions part, int color)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
const std::array< double, 6 > getX()
static const unsigned int phase0size
static const float tomRad
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static const unsigned int phase1size
Definition: value.py:1
EulerAngles toAngles(const RotationType &)
Convert rotation matrix to angles about x-, y-, z-axes (frame rotation).
Definition: Utilities.cc:8
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
#define COUT
Basic3DVector unit() const
Log< level::Info, false > LogInfo
void fillComparisonHistograms(std::map< int, AlignmentPI::partitions > &boundaries, const std::vector< AlignTransform > &ref_ali, const std::vector< AlignTransform > &target_ali, std::unordered_map< AlignmentPI::coordinate, std::unique_ptr< TH1F > > &compare, bool diff=false, AlignmentPI::partitions checkPart=AlignmentPI::INVALID)
GlobalPoint getPartitionAvg(AlignmentPI::PARTITION p)
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::string getStringFromCoordinate(AlignmentPI::coordinate coord)
histos
Definition: combine.py:4
std::pair< double, double > calculatePosition(TVirtualPad *myPad, int boundary)
double a
Definition: hdecay.h:121
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
static constexpr char const *const kPh0DefaultFile
void displayNotSupported(TCanvas &canv, const unsigned int size)
void fillComparisonHistogram(const AlignmentPI::coordinate &coord, std::map< int, AlignmentPI::partitions > &boundaries, const std::vector< AlignTransform > &ref_ali, const std::vector< AlignTransform > &target_ali, std::unique_ptr< TH1F > &compare)
Log< level::Warning, false > LogWarning
const std::array< double, 6 > getZ()
std::string getStringFromPart(AlignmentPI::partitions i, bool isPhase2=false)
#define LogDebug(id)
unsigned transform(const HcalDetId &id, unsigned transformCode)