CMS 3D CMS Logo

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