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
30 
31 #include <memory>
32 #include <sstream>
33 #include <iostream>
34 #include <iomanip> // std::setprecision
35 
36 // include ROOT
37 #include "TH2F.h"
38 #include "TGaxis.h"
39 #include "TLegend.h"
40 #include "TCanvas.h"
41 #include "TLine.h"
42 #include "TStyle.h"
43 #include "TLatex.h"
44 #include "TPave.h"
45 #include "TMarker.h"
46 #include "TPaveStats.h"
47 
48 namespace {
49 
50  using namespace cond::payloadInspector;
51 
52  // M.M. 2017/09/29
53  // Hardcoded Tracker Global Position Record
54  // Without accessing the ES, it is not possible to access to the GPR with the PI technology,
55  // so this needs to be hardcoded.
56  // Anyway it is not likely to change until a new Tracker is installed.
57  // Details at:
58  // - https://indico.cern.ch/event/238026/contributions/513928/attachments/400000/556192/mm_TkAlMeeting_28_03_2013.pdf
59  // - https://twiki.cern.ch/twiki/bin/view/CMS/TkAlignmentPixelPosition
60 
61  const std::map<AlignmentPI::coordinate, float> hardcodeGPR = {
62  {AlignmentPI::t_x, -9.00e-02}, {AlignmentPI::t_y, -1.10e-01}, {AlignmentPI::t_z, -1.70e-01}};
63 
64  //*******************************************/
65  // Size of the movement over all partitions,
66  // one at a time
67  //******************************************//
68 
69  template <int ntags, IOVMultiplicity nIOVs>
70  class TrackerAlignmentCompareAll : public PlotImage<Alignments, nIOVs, ntags> {
71  public:
72  TrackerAlignmentCompareAll()
73  : PlotImage<Alignments, nIOVs, ntags>("comparison of all coordinates between two geometries") {}
74 
75  bool fill() override {
76  TGaxis::SetExponentOffset(-0.12, 0.01, "y"); // Y offset
77 
78  // trick to deal with the multi-ioved tag and two tag case at the same time
79  auto theIOVs = PlotBase::getTag<0>().iovs;
80  auto tagname1 = PlotBase::getTag<0>().name;
81  std::string tagname2 = "";
82  auto firstiov = theIOVs.front();
83  std::tuple<cond::Time_t, cond::Hash> lastiov;
84 
85  // we don't support (yet) comparison with more than 2 tags
86  assert(this->m_plotAnnotations.ntags < 3);
87 
88  if (this->m_plotAnnotations.ntags == 2) {
89  auto tag2iovs = PlotBase::getTag<1>().iovs;
90  tagname2 = PlotBase::getTag<1>().name;
91  lastiov = tag2iovs.front();
92  } else {
93  lastiov = theIOVs.back();
94  }
95 
96  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
97  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
98 
99  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
100  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
101 
102  std::vector<AlignTransform> ref_ali = first_payload->m_align;
103  std::vector<AlignTransform> target_ali = last_payload->m_align;
104 
105  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 2000, 1200);
106  canvas.Divide(3, 2);
107 
108  if (ref_ali.size() != target_ali.size()) {
109  edm::LogError("TrackerAlignment_PayloadInspector")
110  << "the size of the reference alignment (" << ref_ali.size()
111  << ") is different from the one of the target (" << target_ali.size()
112  << ")! You are probably trying to compare different underlying geometries. Exiting";
113  return false;
114  }
115 
116  // check that the geomtery is a tracker one
117  const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
118  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
119  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
120  TrackerTopology tTopo =
122 
123  for (const auto &ali : ref_ali) {
124  auto mydetid = ali.rawId();
125  if (DetId(mydetid).det() != DetId::Tracker) {
126  edm::LogWarning("TrackerAlignment_PayloadInspector")
127  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
128  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
129  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
130  return false;
131  }
132  }
133 
134  const std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
140 
141  std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
142 
143  // generate the map of histograms
144  for (const auto &coord : coords) {
145  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
146  std::string unit =
147  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
148 
149  diffs[coord] = std::make_unique<TH1F>(Form("comparison_%s", s_coord.c_str()),
150  Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
151  ref_ali.size(),
152  -0.5,
153  ref_ali.size() - 0.5);
154  }
155 
156  // fill all the histograms together
157  std::vector<int> boundaries;
158  AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs);
159 
160  unsigned int subpad{1};
161  TLegend legend = TLegend(0.17, 0.84, 0.95, 0.94);
162  legend.SetTextSize(0.023);
163  for (const auto &coord : coords) {
164  canvas.cd(subpad);
165  canvas.cd(subpad)->SetTopMargin(0.06);
166  canvas.cd(subpad)->SetLeftMargin(0.17);
167  canvas.cd(subpad)->SetRightMargin(0.05);
168  canvas.cd(subpad)->SetBottomMargin(0.15);
169  AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
170  auto max = diffs[coord]->GetMaximum();
171  auto min = diffs[coord]->GetMinimum();
172  auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
173  if (range == 0.f)
174  range = 0.1;
175  //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
176 
177  diffs[coord]->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
178  diffs[coord]->GetYaxis()->SetTitleOffset(1.5);
179  diffs[coord]->SetMarkerStyle(20);
180  diffs[coord]->SetMarkerSize(0.5);
181  diffs[coord]->Draw("P");
182 
183  if (subpad == 1) { /* fill the legend only at the first pass */
184  if (this->m_plotAnnotations.ntags == 2) {
185  legend.SetHeader("#bf{Two Tags Comparison}", "C"); // option "C" allows to center the header
186  legend.AddEntry(
187  diffs[coord].get(),
188  ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}")
189  .c_str(),
190  "PL");
191  } else {
192  legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C"); // option "C" allows to center the header
193  legend.AddEntry(diffs[coord].get(),
194  ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
195  "PL");
196  }
197  }
198  subpad++;
199  }
200 
201  canvas.Update();
202  canvas.cd();
203  canvas.Modified();
204 
205  TLine l[6][boundaries.size()];
206  TLatex tSubdet[6];
207  for (unsigned int i = 0; i < 6; i++) {
208  tSubdet[i].SetTextColor(kRed);
209  tSubdet[i].SetNDC();
210  tSubdet[i].SetTextAlign(21);
211  tSubdet[i].SetTextSize(0.03);
212  tSubdet[i].SetTextAngle(90);
213  }
214 
215  subpad = 0;
216  for (const auto &coord : coords) {
217  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
218  canvas.cd(subpad + 1);
219 
220  unsigned int i = 0;
221  for (const auto &line : boundaries) {
222  l[subpad][i] = TLine(diffs[coord]->GetBinLowEdge(line),
223  canvas.cd(subpad + 1)->GetUymin(),
224  diffs[coord]->GetBinLowEdge(line),
225  canvas.cd(subpad + 1)->GetUymax() * 0.84);
226  l[subpad][i].SetLineWidth(1);
227  l[subpad][i].SetLineStyle(9);
228  l[subpad][i].SetLineColor(2);
229  l[subpad][i].Draw("same");
230  i++;
231  }
232 
233  const auto &theX_ = {0.2, 0.24, 0.31, 0.4, 0.55, 0.8};
234  for (unsigned int j = 1; j <= 7; j++) {
235  auto thePart = static_cast<AlignmentPI::partitions>(j);
236  tSubdet[subpad].DrawLatex(
237  theX_.begin()[j - 1], 0.20, Form("%s", (AlignmentPI::getStringFromPart(thePart)).c_str()));
238  }
239 
240  auto ltx = TLatex();
241  ltx.SetTextFont(62);
242  ltx.SetTextSize(0.042);
243  ltx.SetTextAlign(11);
244  ltx.DrawLatexNDC(canvas.cd(subpad + 1)->GetLeftMargin(),
245  1 - canvas.cd(subpad + 1)->GetTopMargin() + 0.01,
246  ("Tracker Alignment Compare : #color[4]{" + s_coord + "}").c_str());
247  legend.Draw("same");
248  subpad++;
249  } // loop on the coordinates
250 
251  std::string fileName(this->m_imageFileName);
252  canvas.SaveAs(fileName.c_str());
253  //canvas.SaveAs("out.root");
254 
255  return true;
256  }
257  };
258 
259  typedef TrackerAlignmentCompareAll<1, MULTI_IOV> TrackerAlignmentComparatorSingleTag;
260  typedef TrackerAlignmentCompareAll<2, SINGLE_IOV> TrackerAlignmentComparatorTwoTags;
261 
262  //*******************************************//
263  // Size of the movement over all partitions,
264  // one coordinate (x,y,z,...) at a time
265  //******************************************//
266 
267  template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
268  class TrackerAlignmentComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
269  public:
270  TrackerAlignmentComparatorBase()
271  : PlotImage<Alignments, nIOVs, ntags>("comparison of " + AlignmentPI::getStringFromCoordinate(coord) +
272  " coordinate between two geometries") {}
273 
274  bool fill() override {
275  TGaxis::SetExponentOffset(-0.12, 0.01, "y"); // Y offset
276 
277  // trick to deal with the multi-ioved tag and two tag case at the same time
278  auto theIOVs = PlotBase::getTag<0>().iovs;
279  auto tagname1 = PlotBase::getTag<0>().name;
280  std::string tagname2 = "";
281  auto firstiov = theIOVs.front();
282  std::tuple<cond::Time_t, cond::Hash> lastiov;
283 
284  // we don't support (yet) comparison with more than 2 tags
285  assert(this->m_plotAnnotations.ntags < 3);
286 
287  if (this->m_plotAnnotations.ntags == 2) {
288  auto tag2iovs = PlotBase::getTag<1>().iovs;
289  tagname2 = PlotBase::getTag<1>().name;
290  lastiov = tag2iovs.front();
291  } else {
292  lastiov = theIOVs.back();
293  }
294 
295  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
296  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
297 
298  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
299  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
300 
301  std::vector<AlignTransform> ref_ali = first_payload->m_align;
302  std::vector<AlignTransform> target_ali = last_payload->m_align;
303 
304  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1200, 1200);
305 
306  if (ref_ali.size() != target_ali.size()) {
307  edm::LogError("TrackerAlignment_PayloadInspector")
308  << "the size of the reference alignment (" << ref_ali.size()
309  << ") is different from the one of the target (" << target_ali.size()
310  << ")! You are probably trying to compare different underlying geometries. Exiting";
311  return false;
312  }
313 
314  // check that the geomtery is a tracker one
315  const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
316  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
317  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
318  TrackerTopology tTopo =
320 
321  for (const auto &ali : ref_ali) {
322  auto mydetid = ali.rawId();
323  if (DetId(mydetid).det() != DetId::Tracker) {
324  edm::LogWarning("TrackerAlignment_PayloadInspector")
325  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
326  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
327  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
328  return false;
329  }
330  }
331 
332  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
333  std::string unit =
334  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
335 
336  //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));
337  std::unique_ptr<TH1F> compare =
338  std::make_unique<TH1F>("comparison",
339  Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
340  ref_ali.size(),
341  -0.5,
342  ref_ali.size() - 0.5);
343 
344  // fill the histograms
345  std::vector<int> boundaries;
346  AlignmentPI::fillComparisonHistogram(coord, boundaries, ref_ali, target_ali, compare);
347 
348  canvas.cd();
349 
350  canvas.SetTopMargin(0.06);
351  canvas.SetLeftMargin(0.17);
352  canvas.SetRightMargin(0.05);
353  canvas.SetBottomMargin(0.15);
355  auto max = compare->GetMaximum();
356  auto min = compare->GetMinimum();
357  auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
358  if (range == 0.f)
359  range = 0.1;
360  //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
361 
362  compare->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
363  compare->GetYaxis()->SetTitleOffset(1.5);
364  compare->SetMarkerStyle(20);
365  compare->SetMarkerSize(0.5);
366  compare->Draw("P");
367 
368  canvas.Update();
369  canvas.cd();
370 
371  TLine l[boundaries.size()];
372  unsigned int i = 0;
373  for (const auto &line : boundaries) {
374  l[i] = TLine(compare->GetBinLowEdge(line),
375  canvas.cd()->GetUymin(),
376  compare->GetBinLowEdge(line),
377  canvas.cd()->GetUymax());
378  l[i].SetLineWidth(1);
379  l[i].SetLineStyle(9);
380  l[i].SetLineColor(2);
381  l[i].Draw("same");
382  i++;
383  }
384 
385  TLatex tSubdet;
386  tSubdet.SetNDC();
387  tSubdet.SetTextAlign(21);
388  tSubdet.SetTextSize(0.027);
389  tSubdet.SetTextAngle(90);
390  for (unsigned int j = 1; j <= 6; j++) {
391  auto thePart = static_cast<AlignmentPI::partitions>(j);
392  tSubdet.SetTextColor(kRed);
393  auto myPair = (j > 1) ? AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(boundaries[j - 2]))
394  : AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(0));
395  float theX_ = myPair.first + 0.025;
396  tSubdet.DrawLatex(theX_, 0.20, Form("%s", (AlignmentPI::getStringFromPart(thePart)).c_str()));
397  }
398 
399  TLegend legend = TLegend(0.17, 0.86, 0.95, 0.94);
400  if (this->m_plotAnnotations.ntags == 2) {
401  legend.SetHeader("#bf{Two Tags Comparison}", "C"); // option "C" allows to center the header
402  legend.AddEntry(
403  compare.get(),
404  ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
405  "PL");
406  } else {
407  legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C"); // option "C" allows to center the header
408  legend.AddEntry(compare.get(),
409  ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
410  "PL");
411  }
412  legend.SetTextSize(0.020);
413  legend.Draw("same");
414 
415  auto ltx = TLatex();
416  ltx.SetTextFont(62);
417  ltx.SetTextSize(0.042);
418  ltx.SetTextAlign(11);
419  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
420  1 - gPad->GetTopMargin() + 0.01,
421  ("Tracker Alignment Compare :#color[4]{" + s_coord + "}").c_str());
422 
423  std::string fileName(this->m_imageFileName);
424  canvas.SaveAs(fileName.c_str());
425 
426  return true;
427  }
428  };
429 
430  template <AlignmentPI::coordinate coord>
431  using TrackerAlignmentCompare = TrackerAlignmentComparatorBase<coord, 1, MULTI_IOV>;
432 
433  template <AlignmentPI::coordinate coord>
434  using TrackerAlignmentCompareTwoTags = TrackerAlignmentComparatorBase<coord, 2, SINGLE_IOV>;
435 
436  typedef TrackerAlignmentCompare<AlignmentPI::t_x> TrackerAlignmentCompareX;
437  typedef TrackerAlignmentCompare<AlignmentPI::t_y> TrackerAlignmentCompareY;
438  typedef TrackerAlignmentCompare<AlignmentPI::t_z> TrackerAlignmentCompareZ;
439 
440  typedef TrackerAlignmentCompare<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlpha;
441  typedef TrackerAlignmentCompare<AlignmentPI::rot_beta> TrackerAlignmentCompareBeta;
442  typedef TrackerAlignmentCompare<AlignmentPI::rot_gamma> TrackerAlignmentCompareGamma;
443 
444  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_x> TrackerAlignmentCompareXTwoTags;
445  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_y> TrackerAlignmentCompareYTwoTags;
446  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_z> TrackerAlignmentCompareZTwoTags;
447 
448  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlphaTwoTags;
449  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_beta> TrackerAlignmentCompareBetaTwoTags;
450  typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_gamma> TrackerAlignmentCompareGammaTwoTags;
451 
452  //*******************************************//
453  // Summary canvas per subdetector
454  //******************************************//
455 
456  template <int ntags, IOVMultiplicity nIOVs, AlignmentPI::partitions q>
457  class TrackerAlignmentSummaryBase : public PlotImage<Alignments, nIOVs, ntags> {
458  public:
459  TrackerAlignmentSummaryBase()
460  : PlotImage<Alignments, nIOVs, ntags>("Comparison of all coordinates between two geometries for " +
461  getStringFromPart(q)) {}
462 
463  bool fill() override {
464  // trick to deal with the multi-ioved tag and two tag case at the same time
465  auto theIOVs = PlotBase::getTag<0>().iovs;
466  auto tagname1 = PlotBase::getTag<0>().name;
467  std::string tagname2 = "";
468  auto firstiov = theIOVs.front();
469  std::tuple<cond::Time_t, cond::Hash> lastiov;
470 
471  // we don't support (yet) comparison with more than 2 tags
472  assert(this->m_plotAnnotations.ntags < 3);
473 
474  if (this->m_plotAnnotations.ntags == 2) {
475  auto tag2iovs = PlotBase::getTag<1>().iovs;
476  tagname2 = PlotBase::getTag<1>().name;
477  lastiov = tag2iovs.front();
478  } else {
479  lastiov = theIOVs.back();
480  }
481 
482  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
483  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
484 
485  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
486  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
487 
488  std::vector<AlignTransform> ref_ali = first_payload->m_align;
489  std::vector<AlignTransform> target_ali = last_payload->m_align;
490 
491  if (ref_ali.size() != target_ali.size()) {
492  edm::LogError("TrackerAlignment_PayloadInspector")
493  << "the size of the reference alignment (" << ref_ali.size()
494  << ") is different from the one of the target (" << target_ali.size()
495  << ")! You are probably trying to compare different underlying geometries. Exiting";
496  return false;
497  }
498 
499  // check that the geomtery is a tracker one
500  const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
501  ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
502  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
503  TrackerTopology tTopo =
505 
506  for (const auto &ali : ref_ali) {
507  auto mydetid = ali.rawId();
508  if (DetId(mydetid).det() != DetId::Tracker) {
509  edm::LogWarning("TrackerAlignment_PayloadInspector")
510  << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
511  << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
512  << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
513  return false;
514  }
515  }
516 
517  TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1800, 1200);
518  canvas.Divide(3, 2);
519 
520  std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
521  std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
527 
528  for (const auto &coord : coords) {
529  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
530  std::string unit =
531  (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
532 
533  diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s", s_coord.c_str()),
534  Form(";#Delta%s %s;n. of modules", s_coord.c_str(), unit.c_str()),
535  1000,
536  -500.,
537  500.);
538  }
539 
540  // fill the comparison histograms
541  std::vector<int> boundaries{};
542  AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs, true, q);
543 
544  int c_index = 1;
545 
546  auto legend = std::make_unique<TLegend>(0.14, 0.93, 0.55, 0.98);
547  legend->AddEntry(
548  diffs[AlignmentPI::t_x].get(),
549  ("#DeltaIOV: " + std::to_string(std::get<0>(lastiov)) + "-" + std::to_string(std::get<0>(firstiov))).c_str(),
550  "L");
551  legend->SetTextSize(0.03);
552 
553  for (const auto &coord : coords) {
554  canvas.cd(c_index)->SetLogy();
555  canvas.cd(c_index)->SetTopMargin(0.02);
556  canvas.cd(c_index)->SetBottomMargin(0.15);
557  canvas.cd(c_index)->SetLeftMargin(0.14);
558  canvas.cd(c_index)->SetRightMargin(0.04);
559  diffs[coord]->SetLineWidth(2);
560  AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
561 
562  //float x_max = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindLastBinAbove(0.));
563  //float x_min = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindFirstBinAbove(0.));
564  //float extremum = std::abs(x_max) > std::abs(x_min) ? std::abs(x_max) : std::abs(x_min);
565  //diffs[coord]->GetXaxis()->SetRangeUser(-extremum*2,extremum*2);
566 
567  int i_max = diffs[coord]->FindLastBinAbove(0.);
568  int i_min = diffs[coord]->FindFirstBinAbove(0.);
569  diffs[coord]->GetXaxis()->SetRange(std::max(1, i_min - 10), std::min(i_max + 10, diffs[coord]->GetNbinsX()));
570  diffs[coord]->Draw("HIST");
571  AlignmentPI::makeNiceStats(diffs[coord].get(), q, kBlack);
572 
573  legend->Draw("same");
574 
575  c_index++;
576  }
577 
578  std::string fileName(this->m_imageFileName);
579  canvas.SaveAs(fileName.c_str());
580 
581  return true;
582  }
583  };
584 
585  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
586  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
587  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
588 
589  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTID;
590  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
591  typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
592 
593  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPixTwoTags;
594  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPixTwoTags;
595  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIBTwoTags;
596 
597  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTIDTwoTags;
598  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOBTwoTags;
599  typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTECTwoTags;
600 
601  //*******************************************//
602  // History of the position of the BPix Barycenter
603  //******************************************//
604 
605  template <AlignmentPI::coordinate coord>
606  class BPixBarycenterHistory : public HistoryPlot<Alignments, float> {
607  public:
608  BPixBarycenterHistory()
610  " Barrel Pixel " + AlignmentPI::getStringFromCoordinate(coord) + " positions vs time",
611  AlignmentPI::getStringFromCoordinate(coord) + " position [cm]") {}
612  ~BPixBarycenterHistory() override = default;
613 
614  float getFromPayload(Alignments &payload) override {
615  std::vector<AlignTransform> alignments = payload.m_align;
616 
617  float barycenter = 0.;
618  float nmodules(0.);
619  for (const auto &ali : alignments) {
620  if (DetId(ali.rawId()).det() != DetId::Tracker) {
621  edm::LogWarning("TrackerAlignment_PayloadInspector")
622  << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
623  << " is different from " << DetId::Tracker << " - terminating ";
624  return false;
625  }
626 
627  int subid = DetId(ali.rawId()).subdetId();
628  if (subid != PixelSubdetector::PixelBarrel)
629  continue;
630 
631  nmodules++;
632  switch (coord) {
633  case AlignmentPI::t_x:
634  barycenter += (ali.translation().x());
635  break;
636  case AlignmentPI::t_y:
637  barycenter += (ali.translation().y());
638  break;
639  case AlignmentPI::t_z:
640  barycenter += (ali.translation().z());
641  break;
642  default:
643  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
644  break;
645  } // switch on the coordinate (only X,Y,Z are interesting)
646  } // ends loop on the alignments
647 
648  edm::LogInfo("TrackerAlignment_PayloadInspector") << "barycenter (" << barycenter << ")/n. modules (" << nmodules
649  << ") = " << barycenter / nmodules << std::endl;
650 
651  // take the mean
652  barycenter /= nmodules;
653 
654  // applied GPR correction to move barycenter to global CMS coordinates
655  barycenter += hardcodeGPR.at(coord);
656 
657  return barycenter;
658 
659  } // payload
660  };
661 
662  typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
663  typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
664  typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
665 
666  /************************************************
667  Display of Tracker Detector barycenters
668  *************************************************/
669  class TrackerAlignmentBarycenters : public PlotImage<Alignments, SINGLE_IOV> {
670  public:
671  TrackerAlignmentBarycenters() : PlotImage<Alignments, SINGLE_IOV>("Display of Tracker Alignment Barycenters") {}
672 
673  bool fill() override {
674  auto tag = PlotBase::getTag<0>();
675  auto iov = tag.iovs.front();
676  std::shared_ptr<Alignments> payload = fetchPayload(std::get<1>(iov));
677  unsigned int run = std::get<0>(iov);
678 
679  TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1600, 1000);
680  canvas.cd();
681 
682  canvas.SetTopMargin(0.07);
683  canvas.SetBottomMargin(0.06);
684  canvas.SetLeftMargin(0.15);
685  canvas.SetRightMargin(0.03);
686  canvas.Modified();
687  canvas.SetGrid();
688 
689  auto h2_BarycenterParameters =
690  std::make_unique<TH2F>("Parameters", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
691 
692  auto h2_uncBarycenterParameters =
693  std::make_unique<TH2F>("Parameters2", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
694 
695  h2_BarycenterParameters->SetStats(false);
696  h2_BarycenterParameters->SetTitle(nullptr);
697  h2_uncBarycenterParameters->SetStats(false);
698  h2_uncBarycenterParameters->SetTitle(nullptr);
699 
700  std::vector<AlignTransform> alignments = payload->m_align;
701 
702  isPhase0 = (alignments.size() == AlignmentPI::phase0size) ? true : false;
703 
704  // check that the geomtery is a tracker one
705  const char *path_toTopologyXML = isPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
706  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
707 
708  TrackerTopology tTopo =
710 
711  AlignmentPI::TkAlBarycenters barycenters;
712  // compute uncorrected barycenter
713  barycenters.computeBarycenters(
714  alignments, tTopo, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
715 
716  auto Xbarycenters = barycenters.getX();
717  auto Ybarycenters = barycenters.getY();
718  auto Zbarycenters = barycenters.getZ();
719 
720  // compute barycenter corrected for the GPR
721  barycenters.computeBarycenters(alignments, tTopo, hardcodeGPR);
722 
723  auto c_Xbarycenters = barycenters.getX();
724  auto c_Ybarycenters = barycenters.getY();
725  auto c_Zbarycenters = barycenters.getZ();
726 
727  h2_BarycenterParameters->GetXaxis()->SetBinLabel(1, "X [cm]");
728  h2_BarycenterParameters->GetXaxis()->SetBinLabel(2, "Y [cm]");
729  h2_BarycenterParameters->GetXaxis()->SetBinLabel(3, "Z [cm]");
730  h2_BarycenterParameters->GetXaxis()->SetBinLabel(4, "X_{no GPR} [cm]");
731  h2_BarycenterParameters->GetXaxis()->SetBinLabel(5, "Y_{no GPR} [cm]");
732  h2_BarycenterParameters->GetXaxis()->SetBinLabel(6, "Z_{no GPR} [cm]");
733 
734  bool isLikelyMC(false);
735  int checkX =
736  std::count_if(Xbarycenters.begin(), Xbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
737  int checkY =
738  std::count_if(Ybarycenters.begin(), Ybarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
739  int checkZ =
740  std::count_if(Zbarycenters.begin(), Zbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
741 
742  // if all the coordinate barycenters for both BPix and FPix are below 10um
743  // this is very likely a MC payload
744  if ((checkX + checkY + checkZ) == 0 && run == 1)
745  isLikelyMC = true;
746 
747  unsigned int yBin = 6;
748  for (unsigned int i = 0; i < 6; i++) {
749  auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
750  std::string theLabel = getStringFromPart(thePart);
751  h2_BarycenterParameters->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
752  if (!isLikelyMC) {
753  h2_BarycenterParameters->SetBinContent(1, yBin, c_Xbarycenters[i]);
754  h2_BarycenterParameters->SetBinContent(2, yBin, c_Ybarycenters[i]);
755  h2_BarycenterParameters->SetBinContent(3, yBin, c_Zbarycenters[i]);
756  }
757 
758  h2_uncBarycenterParameters->SetBinContent(4, yBin, Xbarycenters[i]);
759  h2_uncBarycenterParameters->SetBinContent(5, yBin, Ybarycenters[i]);
760  h2_uncBarycenterParameters->SetBinContent(6, yBin, Zbarycenters[i]);
761  yBin--;
762  }
763 
764  h2_BarycenterParameters->GetXaxis()->LabelsOption("h");
765  h2_BarycenterParameters->GetYaxis()->SetLabelSize(0.05);
766  h2_BarycenterParameters->GetXaxis()->SetLabelSize(0.05);
767  h2_BarycenterParameters->SetMarkerSize(1.5);
768  h2_BarycenterParameters->Draw("TEXT");
769 
770  h2_uncBarycenterParameters->SetMarkerSize(1.5);
771  h2_uncBarycenterParameters->SetMarkerColor(kRed);
772  h2_uncBarycenterParameters->Draw("TEXTsame");
773 
774  TLatex t1;
775  t1.SetNDC();
776  t1.SetTextAlign(26);
777  t1.SetTextSize(0.05);
778  t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters, IOV %i", run));
779  t1.SetTextSize(0.025);
780 
781  std::string fileName(m_imageFileName);
782  canvas.SaveAs(fileName.c_str());
783 
784  return true;
785  }
786 
787  private:
788  bool isPhase0;
789  };
790 
791  /************************************************
792  Comparator of Tracker Detector barycenters
793  *************************************************/
794  template <int ntags, IOVMultiplicity nIOVs>
795  class TrackerAlignmentBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
796  public:
797  TrackerAlignmentBarycentersComparatorBase()
798  : PlotImage<Alignments, nIOVs, ntags>("Comparison of Tracker Alignment Barycenters") {}
799 
800  bool fill() override {
801  // trick to deal with the multi-ioved tag and two tag case at the same time
802  auto theIOVs = PlotBase::getTag<0>().iovs;
803  auto tagname1 = PlotBase::getTag<0>().name;
804  std::string tagname2 = "";
805  auto firstiov = theIOVs.front();
806  std::tuple<cond::Time_t, cond::Hash> lastiov;
807 
808  // we don't support (yet) comparison with more than 2 tags
809  assert(this->m_plotAnnotations.ntags < 3);
810 
811  if (this->m_plotAnnotations.ntags == 2) {
812  auto tag2iovs = PlotBase::getTag<1>().iovs;
813  tagname2 = PlotBase::getTag<1>().name;
814  lastiov = tag2iovs.front();
815  } else {
816  lastiov = theIOVs.back();
817  }
818 
819  unsigned int first_run = std::get<0>(firstiov);
820  unsigned int last_run = std::get<0>(lastiov);
821 
822  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
823  std::vector<AlignTransform> last_alignments = last_payload->m_align;
824 
825  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
826  std::vector<AlignTransform> first_alignments = first_payload->m_align;
827 
828  isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
829  isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
830 
831  // check that the geomtery is a tracker one
832  const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
833  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
834 
835  TrackerTopology tTopo_f =
837 
838  path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
839  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
840 
841  TrackerTopology tTopo_l =
843 
844  TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1200, 800);
845  canvas.cd();
846 
847  canvas.SetTopMargin(0.07);
848  canvas.SetBottomMargin(0.06);
849  canvas.SetLeftMargin(0.15);
850  canvas.SetRightMargin(0.03);
851  canvas.Modified();
852  canvas.SetGrid();
853 
854  auto h2_BarycenterDiff =
855  std::make_unique<TH2F>("Parameters diff", "SubDetector Barycenter Difference", 3, 0.0, 3.0, 6, 0, 6.);
856 
857  h2_BarycenterDiff->SetStats(false);
858  h2_BarycenterDiff->SetTitle(nullptr);
859  h2_BarycenterDiff->GetXaxis()->SetBinLabel(1, "X [#mum]");
860  h2_BarycenterDiff->GetXaxis()->SetBinLabel(2, "Y [#mum]");
861  h2_BarycenterDiff->GetXaxis()->SetBinLabel(3, "Z [#mum]");
862 
863  AlignmentPI::TkAlBarycenters l_barycenters;
864  l_barycenters.computeBarycenters(last_alignments, tTopo_l, hardcodeGPR);
865 
866  AlignmentPI::TkAlBarycenters f_barycenters;
867  f_barycenters.computeBarycenters(first_alignments, tTopo_f, hardcodeGPR);
868 
869  unsigned int yBin = 6;
870  for (unsigned int i = 0; i < 6; i++) {
871  auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
872  std::string theLabel = getStringFromPart(thePart);
873  h2_BarycenterDiff->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
874  h2_BarycenterDiff->SetBinContent(
875  1, yBin, (l_barycenters.getX()[i] - f_barycenters.getX()[i]) * AlignmentPI::cmToUm);
876  h2_BarycenterDiff->SetBinContent(
877  2, yBin, (l_barycenters.getY()[i] - f_barycenters.getY()[i]) * AlignmentPI::cmToUm);
878  h2_BarycenterDiff->SetBinContent(
879  3, yBin, (l_barycenters.getZ()[i] - f_barycenters.getZ()[i]) * AlignmentPI::cmToUm);
880  yBin--;
881  }
882 
883  h2_BarycenterDiff->GetXaxis()->LabelsOption("h");
884  h2_BarycenterDiff->GetYaxis()->SetLabelSize(0.05);
885  h2_BarycenterDiff->GetXaxis()->SetLabelSize(0.05);
886  h2_BarycenterDiff->SetMarkerSize(1.5);
887  h2_BarycenterDiff->SetMarkerColor(kRed);
888  h2_BarycenterDiff->Draw("TEXT");
889 
890  TLatex t1;
891  t1.SetNDC();
892  t1.SetTextAlign(26);
893  t1.SetTextSize(0.05);
894  t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters Diff, IOV %i - IOV %i", last_run, first_run));
895  t1.SetTextSize(0.025);
896 
897  std::string fileName(this->m_imageFileName);
898  canvas.SaveAs(fileName.c_str());
899 
900  return true;
901  }
902 
903  private:
904  bool isInitialPhase0;
905  bool isFinalPhase0;
906  };
907 
908  using TrackerAlignmentBarycentersCompare = TrackerAlignmentBarycentersComparatorBase<1, MULTI_IOV>;
909  using TrackerAlignmentBarycentersCompareTwoTags = TrackerAlignmentBarycentersComparatorBase<2, SINGLE_IOV>;
910 
911  /************************************************
912  Comparator of Pixel Tracker Detector barycenters
913  *************************************************/
914  template <int ntags, IOVMultiplicity nIOVs>
915  class PixelBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
916  public:
917  PixelBarycentersComparatorBase() : PlotImage<Alignments, nIOVs, ntags>("Comparison of Pixel Barycenters") {}
918 
919  bool fill() override {
920  // trick to deal with the multi-ioved tag and two tag case at the same time
921  auto theIOVs = PlotBase::getTag<0>().iovs;
922  auto tagname1 = PlotBase::getTag<0>().name;
923  std::string tagname2 = "";
924  auto firstiov = theIOVs.front();
925  std::tuple<cond::Time_t, cond::Hash> lastiov;
926 
927  // we don't support (yet) comparison with more than 2 tags
928  assert(this->m_plotAnnotations.ntags < 3);
929 
930  if (this->m_plotAnnotations.ntags == 2) {
931  auto tag2iovs = PlotBase::getTag<1>().iovs;
932  tagname2 = PlotBase::getTag<1>().name;
933  lastiov = tag2iovs.front();
934  } else {
935  lastiov = theIOVs.back();
936  }
937 
938  unsigned int first_run = std::get<0>(firstiov);
939  unsigned int last_run = std::get<0>(lastiov);
940 
941  std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
942  std::vector<AlignTransform> last_alignments = last_payload->m_align;
943 
944  std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
945  std::vector<AlignTransform> first_alignments = first_payload->m_align;
946 
947  TCanvas canvas("Pixel Barycenter Summary", "Pixel Barycenter summary", 1200, 1200);
948  canvas.Divide(2, 2);
949  canvas.cd();
950 
951  TLatex t1;
952  t1.SetNDC();
953  t1.SetTextAlign(26);
954  t1.SetTextSize(0.03);
955  t1.DrawLatex(0.5,
956  0.97,
957  ("Pixel Barycenters comparison, IOV: #color[2]{" + std::to_string(first_run) +
958  "} vs IOV: #color[4]{" + std::to_string(last_run) + "}")
959  .c_str());
960  t1.SetTextSize(0.025);
961 
962  for (unsigned int c = 1; c <= 4; c++) {
963  canvas.cd(c)->SetTopMargin(0.07);
964  canvas.cd(c)->SetBottomMargin(0.12);
965  canvas.cd(c)->SetLeftMargin(0.15);
966  canvas.cd(c)->SetRightMargin(0.03);
967  canvas.cd(c)->Modified();
968  canvas.cd(c)->SetGrid();
969  }
970 
971  std::array<std::string, 3> structures = {{"FPIX-", "BPIX", "FPIX+"}};
972  std::array<std::unique_ptr<TH2F>, 3> histos;
973 
974  isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
975  isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
976 
977  // check that the geomtery is a tracker one
978  const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
979  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
980 
981  TrackerTopology tTopo_f =
983 
984  AlignmentPI::TkAlBarycenters myInitialBarycenters;
985  //myInitialBarycenters.computeBarycenters(first_alignments,tTopo_f,hardcodeGPR);
986  myInitialBarycenters.computeBarycenters(
987  first_alignments, tTopo_f, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
988 
989  path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
990  : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
991 
992  TrackerTopology tTopo_l =
994 
995  AlignmentPI::TkAlBarycenters myFinalBarycenters;
996  //myFinalBarycenters.computeBarycenters(last_alignments,tTopo_l,hardcodeGPR);
997  myFinalBarycenters.computeBarycenters(
998  last_alignments, tTopo_l, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
999 
1000  if (isFinalPhase0 != isInitialPhase0) {
1001  edm::LogWarning("TrackerAlignment_PayloadInspector")
1002  << "the size of the reference alignment (" << first_alignments.size()
1003  << ") is different from the one of the target (" << last_alignments.size()
1004  << ")! You are probably trying to compare different underlying geometries.";
1005  }
1006 
1007  unsigned int index(0);
1008  for (const auto &piece : structures) {
1009  const char *name = piece.c_str();
1010  histos[index] = std::make_unique<TH2F>(
1011  name,
1012  Form("%s x-y Barycenter Difference;x_{%s}-x_{TOB} [mm];y_{%s}-y_{TOB} [mm]", name, name, name),
1013  100,
1014  -3.,
1015  3.,
1016  100,
1017  -3.,
1018  3.);
1019 
1020  histos[index]->SetStats(false);
1021  histos[index]->SetTitle(nullptr);
1022  histos[index]->GetYaxis()->SetLabelSize(0.05);
1023  histos[index]->GetXaxis()->SetLabelSize(0.05);
1024  histos[index]->GetYaxis()->SetTitleSize(0.06);
1025  histos[index]->GetXaxis()->SetTitleSize(0.06);
1026  histos[index]->GetYaxis()->CenterTitle();
1027  histos[index]->GetXaxis()->CenterTitle();
1028  histos[index]->GetXaxis()->SetTitleOffset(0.9);
1029  index++;
1030  }
1031 
1032  auto h2_ZBarycenterDiff = std::make_unique<TH2F>(
1033  "Pixel_z_diff", "Pixel z-Barycenter Difference;; z_{Pixel-Ideal} -z_{TOB} [mm]", 3, -0.5, 2.5, 100, -10., 10.);
1034  h2_ZBarycenterDiff->SetStats(false);
1035  h2_ZBarycenterDiff->SetTitle(nullptr);
1036  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(1, "FPIX -");
1037  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(2, "BPIX");
1038  h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(3, "FPIX +");
1039  h2_ZBarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1040  h2_ZBarycenterDiff->GetXaxis()->SetLabelSize(0.07);
1041  h2_ZBarycenterDiff->GetYaxis()->SetTitleSize(0.06);
1042  h2_ZBarycenterDiff->GetXaxis()->SetTitleSize(0.06);
1043  h2_ZBarycenterDiff->GetYaxis()->CenterTitle();
1044  h2_ZBarycenterDiff->GetXaxis()->CenterTitle();
1045  h2_ZBarycenterDiff->GetYaxis()->SetTitleOffset(1.1);
1046 
1047  std::function<GlobalPoint(int)> cutFunctorInitial = [&myInitialBarycenters](int index) {
1048  switch (index) {
1049  case 1:
1050  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1051  case 2:
1052  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1053  case 3:
1054  return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1055  default:
1056  return GlobalPoint(0, 0, 0);
1057  }
1058  };
1059 
1060  std::function<GlobalPoint(int)> cutFunctorFinal = [&myFinalBarycenters](int index) {
1061  switch (index) {
1062  case 1:
1063  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1064  case 2:
1065  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1066  case 3:
1067  return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1068  default:
1069  return GlobalPoint(0, 0, 0);
1070  }
1071  };
1072 
1073  float x0i, x0f, y0i, y0f;
1074 
1075  t1.SetNDC(kFALSE);
1076  t1.SetTextSize(0.047);
1077  for (unsigned int c = 1; c <= 3; c++) {
1078  x0i = cutFunctorInitial(c).x() * 10; // transform cm to mm (x10)
1079  x0f = cutFunctorFinal(c).x() * 10;
1080  y0i = cutFunctorInitial(c).y() * 10;
1081  y0f = cutFunctorFinal(c).y() * 10;
1082 
1083  canvas.cd(c);
1084  histos[c - 1]->Draw();
1085 
1086  COUT << "initial x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0i << "," << y0i << ") mm"
1087  << std::endl;
1088  COUT << "final x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0f << "," << y0f << ") mm"
1089  << std::endl;
1090 
1091  TMarker *initial = new TMarker(x0i, y0i, 21);
1092  TMarker *final = new TMarker(x0f, y0f, 20);
1093 
1094  initial->SetMarkerColor(kRed);
1095  final->SetMarkerColor(kBlue);
1096  initial->SetMarkerSize(2.5);
1097  final->SetMarkerSize(2.5);
1098  t1.SetTextColor(kRed);
1099  initial->Draw();
1100  t1.DrawLatex(x0i, y0i + (y0i > y0f ? 0.3 : -0.5), Form("(%.2f,%.2f)", x0i, y0i));
1101  final->Draw("same");
1102  t1.SetTextColor(kBlue);
1103  t1.DrawLatex(x0f, y0f + (y0i > y0f ? -0.5 : 0.3), Form("(%.2f,%.2f)", x0f, y0f));
1104  }
1105 
1106  // fourth pad is a special case for the z coordinate
1107  canvas.cd(4);
1108  h2_ZBarycenterDiff->Draw();
1109  float z0i, z0f;
1110 
1111  // numbers do agree with https://twiki.cern.ch/twiki/bin/view/CMSPublic/TkAlignmentPerformancePhaseIStartUp17#Pixel_Barycentre_Positions
1112 
1113  std::array<double, 3> hardcodeIdealZPhase0 = {{-41.94909, 0., 41.94909}}; // units are cm
1114  std::array<double, 3> hardcodeIdealZPhase1 = {{-39.82911, 0., 39.82911}}; // units are cm
1115 
1116  for (unsigned int c = 1; c <= 3; c++) {
1117  // less than pretty but needed to remove the z position of the FPix barycenters != 0
1118 
1119  z0i =
1120  (cutFunctorInitial(c).z() - (isInitialPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) *
1121  10; // convert to mm
1122  z0f =
1123  (cutFunctorFinal(c).z() - (isFinalPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) * 10;
1124 
1125  TMarker *initial = new TMarker(c - 1, z0i, 21);
1126  TMarker *final = new TMarker(c - 1, z0f, 20);
1127 
1128  COUT << "initial z " << std::left << std::setw(7) << structures[c - 1] << " " << z0i << " mm" << std::endl;
1129  COUT << "final z " << std::left << std::setw(7) << structures[c - 1] << " " << z0f << " mm" << std::endl;
1130 
1131  initial->SetMarkerColor(kRed);
1132  final->SetMarkerColor(kBlue);
1133  initial->SetMarkerSize(2.5);
1134  final->SetMarkerSize(2.5);
1135  initial->Draw();
1136  t1.SetTextColor(kRed);
1137  t1.DrawLatex(c - 1, z0i + (z0i > z0f ? 1. : -1.5), Form("(%.2f)", z0i));
1138  final->Draw("same");
1139  t1.SetTextColor(kBlue);
1140  t1.DrawLatex(c - 1, z0f + (z0i > z0f ? -1.5 : 1), Form("(%.2f)", z0f));
1141  }
1142 
1143  std::string fileName(this->m_imageFileName);
1144  canvas.SaveAs(fileName.c_str());
1145 
1146  return true;
1147  }
1148 
1149  private:
1150  bool isInitialPhase0;
1151  bool isFinalPhase0;
1152  };
1153 
1154  using PixelBarycentersCompare = PixelBarycentersComparatorBase<1, MULTI_IOV>;
1155  using PixelBarycentersCompareTwoTags = PixelBarycentersComparatorBase<2, SINGLE_IOV>;
1156 
1157 } // namespace
1158 
1160  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorSingleTag);
1161  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorTwoTags);
1162  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
1163  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
1164  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
1165  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
1166  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
1167  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
1168  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareXTwoTags);
1169  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareYTwoTags);
1170  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZTwoTags);
1171  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlphaTwoTags);
1172  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBetaTwoTags);
1173  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGammaTwoTags);
1174  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
1175  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
1176  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
1177  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
1178  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
1179  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
1180  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPixTwoTags);
1181  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPixTwoTags);
1182  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIBTwoTags);
1183  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIDTwoTags);
1184  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOBTwoTags);
1185  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTECTwoTags);
1186  PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
1187  PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
1188  PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
1189  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycenters);
1190  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompare);
1191  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompareTwoTags);
1192  PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompare);
1193  PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompareTwoTags);
1194 }
const std::array< double, 6 > getY()
bool tidIsDoubleSide(const DetId &id) const
void fillComparisonHistogram(const AlignmentPI::coordinate &coord, std::vector< int > &boundaries, const std::vector< AlignTransform > &ref_ali, const std::vector< AlignTransform > &target_ali, std::unique_ptr< TH1F > &compare)
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
std::string getStringFromPart(AlignmentPI::partitions i)
std::string to_string(const V &value)
Definition: OMSAccess.h:77
Log< level::Error, false > LogError
assert(be >=bs)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
static const float cmToUm
void makeNiceStats(TH1F *hist, AlignmentPI::partitions part, int color)
void fillComparisonHistograms(std::vector< int > &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)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
const std::array< double, 6 > getX()
static const unsigned int phase0size
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
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
#define COUT
Basic3DVector unit() const
Log< level::Info, false > LogInfo
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:119
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
Log< level::Warning, false > LogWarning
const std::array< double, 6 > getZ()