CMS 3D CMS Logo

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