CMS 3D CMS Logo

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