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