File indexing completed on 2025-06-03 00:11:45
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010
0011 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0012 #include "CondCore/Utilities/interface/PayloadInspector.h"
0013 #include "CondCore/CondDB/interface/Time.h"
0014
0015
0016 #include "CondFormats/Alignment/interface/Alignments.h"
0017 #include "DataFormats/DetId/interface/DetId.h"
0018 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0019
0020
0021 #include "CLHEP/Vector/RotationInterfaces.h"
0022 #include "Alignment/CommonAlignment/interface/Utilities.h"
0023
0024
0025 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0026
0027
0028 #include "CondCore/AlignmentPlugins/interface/AlignmentPayloadInspectorHelper.h"
0029 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0030 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0031 #include "DQM/TrackerRemapper/interface/Phase1PixelSummaryMap.h"
0032 #include "DQM/TrackerRemapper/interface/Phase1PixelMaps.h"
0033
0034 #include <boost/range/adaptor/indexed.hpp>
0035 #include <iomanip> // std::setprecision
0036 #include <iostream>
0037 #include <memory>
0038 #include <sstream>
0039
0040
0041 #include "TH2F.h"
0042 #include "TGaxis.h"
0043 #include "TLegend.h"
0044 #include "TCanvas.h"
0045 #include "TLine.h"
0046 #include "TStyle.h"
0047 #include "TLatex.h"
0048 #include "TPave.h"
0049 #include "TMarker.h"
0050 #include "TPaveStats.h"
0051
0052 namespace {
0053
0054 using namespace cond::payloadInspector;
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065 const std::map<AlignmentPI::coordinate, float> hardcodeGPR = {
0066 {AlignmentPI::t_x, -9.00e-02}, {AlignmentPI::t_y, -1.10e-01}, {AlignmentPI::t_z, -1.70e-01}};
0067
0068
0069
0070
0071
0072
0073 enum RegionCategory { ALL = 0, INNER = 1, OUTER = 2 };
0074
0075 template <int ntags, IOVMultiplicity nIOVs, RegionCategory cat>
0076 class TrackerAlignmentCompareAll : public PlotImage<Alignments, nIOVs, ntags> {
0077 public:
0078 TrackerAlignmentCompareAll()
0079 : PlotImage<Alignments, nIOVs, ntags>("comparison of all coordinates between two geometries") {}
0080
0081 bool fill() override {
0082 TGaxis::SetExponentOffset(-0.12, 0.01, "y");
0083
0084
0085 auto theIOVs = PlotBase::getTag<0>().iovs;
0086 auto tagname1 = PlotBase::getTag<0>().name;
0087 std::string tagname2 = "";
0088 auto firstiov = theIOVs.front();
0089 std::tuple<cond::Time_t, cond::Hash> lastiov;
0090
0091
0092 assert(this->m_plotAnnotations.ntags < 3);
0093
0094 if (this->m_plotAnnotations.ntags == 2) {
0095 auto tag2iovs = PlotBase::getTag<1>().iovs;
0096 tagname2 = PlotBase::getTag<1>().name;
0097 lastiov = tag2iovs.front();
0098 } else {
0099 lastiov = theIOVs.back();
0100 }
0101
0102 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0103 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0104
0105 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0106 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0107
0108 std::vector<AlignTransform> ref_ali = first_payload->m_align;
0109 std::vector<AlignTransform> target_ali = last_payload->m_align;
0110
0111 const bool ph2 = (ref_ali.size() > AlignmentPI::phase1size);
0112
0113
0114 const char *path_toTopologyXML = nullptr;
0115 if (ph2) {
0116 if (AlignmentPI::isReorderedTFPXTEPX(ref_ali) && AlignmentPI::isReorderedTFPXTEPX(target_ali)) {
0117 edm::LogPrint("TrackerAlignment_PayloadInspector")
0118 << "Both reference and target alignments are reordered. Using the trackerParameters for the Reordered "
0119 "TFPX,TEPX.";
0120 path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/TFPXTEPXReordered/trackerParameters.xml";
0121 } else if (!AlignmentPI::isReorderedTFPXTEPX(ref_ali) && !AlignmentPI::isReorderedTFPXTEPX(target_ali)) {
0122 edm::LogPrint("TrackerAlignment_PayloadInspector")
0123 << "Neither reference nor target alignments are reordered. Using the standard trackerParameters.";
0124 path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
0125 } else {
0126 if (cat == RegionCategory::ALL || cat == RegionCategory::INNER) {
0127
0128 edm::LogWarning("TrackerAlignment_PayloadInspector")
0129 << "Mismatched alignments detected. One is reordered while the other is not. Unable to proceed.";
0130 return false;
0131 } else {
0132 edm::LogWarning("TrackerAlignment_PayloadInspector")
0133 << "Mismatched inner tracks alignments detected. One is reordered while the other is not. Ignoring as "
0134 "OT only comparison requested.";
0135 path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
0136 }
0137 }
0138 } else {
0139 path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0140 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0141 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0142 }
0143
0144
0145 if (cat != RegionCategory::ALL) {
0146 ref_ali.erase(std::remove_if(ref_ali.begin(),
0147 ref_ali.end(),
0148 [](const AlignTransform &transform) {
0149 int subid = DetId(transform.rawId()).subdetId();
0150 return (cat == RegionCategory::INNER) ? (subid > 2) : (subid <= 2);
0151 }),
0152 ref_ali.end());
0153
0154 target_ali.erase(std::remove_if(target_ali.begin(),
0155 target_ali.end(),
0156 [](const AlignTransform &transform) {
0157 int subid = DetId(transform.rawId()).subdetId();
0158 return (cat == RegionCategory::INNER) ? (subid > 2) : (subid <= 2);
0159 }),
0160 target_ali.end());
0161 }
0162 TCanvas canvas("Alignment Comparison", "Alignment Comparison", 2000, 1200);
0163 canvas.Divide(3, 2);
0164
0165 if (ref_ali.size() != target_ali.size()) {
0166 edm::LogError("TrackerAlignment_PayloadInspector")
0167 << "the size of the reference alignment (" << ref_ali.size()
0168 << ") is different from the one of the target (" << target_ali.size()
0169 << ")! You are probably trying to compare different underlying geometries. Exiting";
0170 return false;
0171 }
0172
0173 TrackerTopology tTopo =
0174 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0175
0176 for (const auto &ali : ref_ali) {
0177 auto mydetid = ali.rawId();
0178 if (DetId(mydetid).det() != DetId::Tracker) {
0179 edm::LogWarning("TrackerAlignment_PayloadInspector")
0180 << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0181 << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0182 << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0183 return false;
0184 }
0185 }
0186
0187 const std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
0188 AlignmentPI::t_y,
0189 AlignmentPI::t_z,
0190 AlignmentPI::rot_alpha,
0191 AlignmentPI::rot_beta,
0192 AlignmentPI::rot_gamma};
0193
0194 std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F>> diffs;
0195
0196
0197 for (const auto &coord : coords) {
0198 auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0199 std::string unit =
0200 (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0201
0202 diffs[coord] = std::make_unique<TH1F>(Form("comparison_%s", s_coord.c_str()),
0203 Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
0204 ref_ali.size(),
0205 -0.5,
0206 ref_ali.size() - 0.5);
0207 }
0208
0209
0210 std::map<int, AlignmentPI::partitions> boundaries;
0211 if (cat < RegionCategory::OUTER) {
0212 boundaries.insert({0, AlignmentPI::BPix});
0213 }
0214 AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs);
0215
0216 unsigned int subpad{1};
0217 TLegend legend = TLegend(0.17, 0.84, 0.95, 0.94);
0218 legend.SetTextSize(0.023);
0219 for (const auto &coord : coords) {
0220 canvas.cd(subpad);
0221 canvas.cd(subpad)->SetTopMargin(0.06);
0222 canvas.cd(subpad)->SetLeftMargin(0.17);
0223 canvas.cd(subpad)->SetRightMargin(0.05);
0224 canvas.cd(subpad)->SetBottomMargin(0.15);
0225 AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0226 auto max = diffs[coord]->GetMaximum();
0227 auto min = diffs[coord]->GetMinimum();
0228 auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
0229 if (range == 0.f)
0230 range = 0.1;
0231
0232
0233 diffs[coord]->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
0234 diffs[coord]->GetYaxis()->SetTitleOffset(1.5);
0235 diffs[coord]->SetMarkerStyle(20);
0236 diffs[coord]->SetMarkerSize(0.5);
0237 diffs[coord]->Draw("P");
0238
0239 if (subpad == 1) {
0240 if (this->m_plotAnnotations.ntags == 2) {
0241 legend.SetHeader("#bf{Two Tags Comparison}", "C");
0242 legend.AddEntry(
0243 diffs[coord].get(),
0244 ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}")
0245 .c_str(),
0246 "PL");
0247 } else {
0248 legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");
0249 legend.AddEntry(diffs[coord].get(),
0250 ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0251 "PL");
0252 }
0253 }
0254 subpad++;
0255 }
0256
0257 canvas.Update();
0258 canvas.cd();
0259 canvas.Modified();
0260
0261 bool doOnlyPixel = (cat == RegionCategory::INNER);
0262
0263 TLine l[6][boundaries.size()];
0264 TLatex tSubdet[6];
0265 for (unsigned int i = 0; i < 6; i++) {
0266 tSubdet[i].SetTextColor(kRed);
0267 tSubdet[i].SetNDC();
0268 tSubdet[i].SetTextAlign(21);
0269 tSubdet[i].SetTextSize(doOnlyPixel ? 0.05 : 0.03);
0270 tSubdet[i].SetTextAngle(90);
0271 }
0272
0273 subpad = 0;
0274 for (const auto &coord : coords) {
0275 auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0276 canvas.cd(subpad + 1);
0277 for (const auto &line : boundaries | boost::adaptors::indexed(0)) {
0278 const auto &index = line.index();
0279 const auto value = line.value();
0280 l[subpad][index] = TLine(diffs[coord]->GetBinLowEdge(value.first),
0281 canvas.cd(subpad + 1)->GetUymin(),
0282 diffs[coord]->GetBinLowEdge(value.first),
0283 canvas.cd(subpad + 1)->GetUymax() * 0.84);
0284 l[subpad][index].SetLineWidth(1);
0285 l[subpad][index].SetLineStyle(9);
0286 l[subpad][index].SetLineColor(2);
0287 l[subpad][index].Draw("same");
0288 }
0289
0290 for (const auto &elem : boundaries | boost::adaptors::indexed(0)) {
0291 const auto &lm = canvas.cd(subpad + 1)->GetLeftMargin();
0292 const auto &rm = 1 - canvas.cd(subpad + 1)->GetRightMargin();
0293 const auto &frac = float(elem.value().first) / ref_ali.size();
0294
0295 LogDebug("TrackerAlignmentCompareAll")
0296 << __PRETTY_FUNCTION__ << " left margin: " << lm << " right margin: " << rm << " fraction: " << frac;
0297
0298 float theX_ = lm + (rm - lm) * frac + ((elem.index() > 0 || doOnlyPixel) ? 0.025 : 0.01);
0299
0300 tSubdet[subpad].DrawLatex(
0301 theX_, 0.23, Form("%s", AlignmentPI::getStringFromPart(elem.value().second, ph2).c_str()));
0302 }
0303
0304 auto ltx = TLatex();
0305 ltx.SetTextFont(62);
0306 ltx.SetTextSize(0.042);
0307 ltx.SetTextAlign(11);
0308 ltx.DrawLatexNDC(canvas.cd(subpad + 1)->GetLeftMargin(),
0309 1 - canvas.cd(subpad + 1)->GetTopMargin() + 0.01,
0310 ("Tracker Alignment Compare : #color[4]{" + s_coord + "}").c_str());
0311 legend.Draw("same");
0312 subpad++;
0313 }
0314
0315 std::string fileName(this->m_imageFileName);
0316 canvas.SaveAs(fileName.c_str());
0317
0318
0319 return true;
0320 }
0321 };
0322
0323 typedef TrackerAlignmentCompareAll<1, MULTI_IOV, RegionCategory::ALL> TrackerAlignmentComparatorSingleTag;
0324 typedef TrackerAlignmentCompareAll<2, SINGLE_IOV, RegionCategory::ALL> TrackerAlignmentComparatorTwoTags;
0325
0326 typedef TrackerAlignmentCompareAll<1, MULTI_IOV, RegionCategory::INNER> PixelAlignmentComparatorSingleTag;
0327 typedef TrackerAlignmentCompareAll<2, SINGLE_IOV, RegionCategory::INNER> PixelAlignmentComparatorTwoTags;
0328
0329 typedef TrackerAlignmentCompareAll<1, MULTI_IOV, RegionCategory::OUTER> OTAlignmentComparatorSingleTag;
0330 typedef TrackerAlignmentCompareAll<2, SINGLE_IOV, RegionCategory::OUTER> OTAlignmentComparatorTwoTags;
0331
0332
0333
0334
0335
0336
0337 template <int ntags, IOVMultiplicity nIOVs, RegionCategory cat>
0338 class TrackerAlignmentCompareCylindricalBase : public PlotImage<Alignments, nIOVs, ntags> {
0339 enum coordinate {
0340 t_r = 1,
0341 t_phi = 2,
0342 t_z = 3,
0343 };
0344
0345 public:
0346 TrackerAlignmentCompareCylindricalBase()
0347 : PlotImage<Alignments, nIOVs, ntags>("comparison of cylindrical coordinates between two geometries") {}
0348
0349 bool fill() override {
0350 TGaxis::SetExponentOffset(-0.12, 0.01, "y");
0351
0352
0353 auto theIOVs = PlotBase::getTag<0>().iovs;
0354 auto tagname1 = PlotBase::getTag<0>().name;
0355 std::string tagname2 = "";
0356 auto firstiov = theIOVs.front();
0357 std::tuple<cond::Time_t, cond::Hash> lastiov;
0358
0359
0360 assert(this->m_plotAnnotations.ntags < 3);
0361
0362 if (this->m_plotAnnotations.ntags == 2) {
0363 auto tag2iovs = PlotBase::getTag<1>().iovs;
0364 tagname2 = PlotBase::getTag<1>().name;
0365 lastiov = tag2iovs.front();
0366 } else {
0367 lastiov = theIOVs.back();
0368 }
0369
0370 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0371 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0372
0373 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0374 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0375
0376 std::vector<AlignTransform> ref_ali = first_payload->m_align;
0377 std::vector<AlignTransform> target_ali = last_payload->m_align;
0378
0379 TCanvas canvas("Alignment Comparison", "", 2000, 600);
0380 canvas.Divide(3, 1);
0381
0382 const bool ph2 = (ref_ali.size() > AlignmentPI::phase1size);
0383
0384 const std::vector<coordinate> coords = {t_r, t_phi, t_z};
0385 std::unordered_map<coordinate, std::unique_ptr<TH1F>> diffs;
0386
0387
0388 if (cat != RegionCategory::ALL) {
0389 ref_ali.erase(std::remove_if(ref_ali.begin(),
0390 ref_ali.end(),
0391 [](const AlignTransform &transform) {
0392 int subid = DetId(transform.rawId()).subdetId();
0393 return (cat == RegionCategory::INNER) ? (subid > 2) : (subid <= 2);
0394 }),
0395 ref_ali.end());
0396
0397 target_ali.erase(std::remove_if(target_ali.begin(),
0398 target_ali.end(),
0399 [](const AlignTransform &transform) {
0400 int subid = DetId(transform.rawId()).subdetId();
0401 return (cat == RegionCategory::INNER) ? (subid > 2) : (subid <= 2);
0402 }),
0403 target_ali.end());
0404 }
0405
0406 auto h_deltaR = std::make_unique<TH1F>(
0407 "deltaR", Form(";Detector Id index; #DeltaR [#mum]"), ref_ali.size(), -0.5, ref_ali.size() - 0.5);
0408 auto h_deltaPhi = std::make_unique<TH1F>(
0409 "deltaPhi", Form(";Detector Id index; #Delta#phi [mrad]"), ref_ali.size(), -0.5, ref_ali.size() - 0.5);
0410 auto h_deltaZ = std::make_unique<TH1F>(
0411 "deltaZ", Form(";Detector Id index; #DeltaZ [#mum]"), ref_ali.size(), -0.5, ref_ali.size() - 0.5);
0412
0413 std::map<int, AlignmentPI::partitions> boundaries;
0414 if (cat < RegionCategory::OUTER) {
0415 boundaries.insert({0, AlignmentPI::BPix});
0416 }
0417
0418 int counter = 0;
0419 AlignmentPI::partitions currentPart = AlignmentPI::BPix;
0420 for (unsigned int i = 0; i < ref_ali.size(); i++) {
0421 if (ref_ali[i].rawId() == target_ali[i].rawId()) {
0422 counter++;
0423 int subid = DetId(ref_ali[i].rawId()).subdetId();
0424 auto thePart = static_cast<AlignmentPI::partitions>(subid);
0425
0426 if (thePart != currentPart) {
0427 currentPart = thePart;
0428 boundaries.insert({counter, thePart});
0429 }
0430
0431 const auto &deltaTrans = target_ali[i].translation() - ref_ali[i].translation();
0432 double dPhi = target_ali[i].translation().phi() - ref_ali[i].translation().phi();
0433 if (dPhi > M_PI) {
0434 dPhi -= 2.0 * M_PI;
0435 }
0436 if (dPhi < -M_PI) {
0437 dPhi += 2.0 * M_PI;
0438 }
0439
0440 h_deltaR->SetBinContent(i + 1, deltaTrans.perp() * AlignmentPI::cmToUm);
0441 h_deltaPhi->SetBinContent(i + 1, dPhi * AlignmentPI::tomRad);
0442 h_deltaZ->SetBinContent(i + 1, deltaTrans.z() * AlignmentPI::cmToUm);
0443 }
0444 }
0445
0446 diffs[t_r] = std::move(h_deltaR);
0447 diffs[t_phi] = std::move(h_deltaPhi);
0448 diffs[t_z] = std::move(h_deltaZ);
0449
0450 unsigned int subpad{1};
0451 TLegend legend = TLegend(0.17, 0.84, 0.95, 0.94);
0452 legend.SetTextSize(0.023);
0453 for (const auto &coord : coords) {
0454 canvas.cd(subpad);
0455 canvas.cd(subpad)->SetTopMargin(0.06);
0456 canvas.cd(subpad)->SetLeftMargin(0.17);
0457 canvas.cd(subpad)->SetRightMargin(0.05);
0458 canvas.cd(subpad)->SetBottomMargin(0.15);
0459 AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0460 auto max = diffs[coord]->GetMaximum();
0461 auto min = diffs[coord]->GetMinimum();
0462 auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
0463 if (range == 0.f) {
0464 range = 0.1;
0465 }
0466
0467
0468 if (coord != t_r) {
0469 diffs[coord]->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
0470 } else {
0471 diffs[coord]->GetYaxis()->SetRangeUser(0., range * 1.5);
0472 }
0473
0474 diffs[coord]->GetYaxis()->SetTitleOffset(1.5);
0475 diffs[coord]->SetMarkerStyle(20);
0476 diffs[coord]->SetMarkerSize(0.5);
0477 diffs[coord]->Draw("P");
0478
0479 if (subpad == 1) {
0480 if (this->m_plotAnnotations.ntags == 2) {
0481 legend.SetHeader("#bf{Two Tags Comparison}", "C");
0482 legend.AddEntry(
0483 diffs[coord].get(),
0484 ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}")
0485 .c_str(),
0486 "PL");
0487 } else {
0488 legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");
0489 legend.AddEntry(diffs[coord].get(),
0490 ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0491 "PL");
0492 }
0493 }
0494 subpad++;
0495 }
0496
0497 canvas.Update();
0498 canvas.cd();
0499 canvas.Modified();
0500
0501 bool doOnlyPixel = (cat == RegionCategory::INNER);
0502
0503 TLine l[6][boundaries.size()];
0504 TLatex tSubdet[6];
0505 for (unsigned int i = 0; i < 6; i++) {
0506 tSubdet[i].SetTextColor(kRed);
0507 tSubdet[i].SetNDC();
0508 tSubdet[i].SetTextAlign(21);
0509 tSubdet[i].SetTextSize(doOnlyPixel ? 0.05 : 0.03);
0510 tSubdet[i].SetTextAngle(90);
0511 }
0512
0513 subpad = 0;
0514 for (const auto &coord : coords) {
0515 auto s_coord = getStringFromCoordinate(coord);
0516 canvas.cd(subpad + 1);
0517 for (const auto &line : boundaries | boost::adaptors::indexed(0)) {
0518 const auto &index = line.index();
0519 const auto value = line.value();
0520 l[subpad][index] = TLine(diffs[coord]->GetBinLowEdge(value.first),
0521 canvas.cd(subpad + 1)->GetUymin(),
0522 diffs[coord]->GetBinLowEdge(value.first),
0523 canvas.cd(subpad + 1)->GetUymax() * 0.84);
0524 l[subpad][index].SetLineWidth(1);
0525 l[subpad][index].SetLineStyle(9);
0526 l[subpad][index].SetLineColor(2);
0527 l[subpad][index].Draw("same");
0528 }
0529
0530 for (const auto &elem : boundaries | boost::adaptors::indexed(0)) {
0531 const auto &lm = canvas.cd(subpad + 1)->GetLeftMargin();
0532 const auto &rm = 1 - canvas.cd(subpad + 1)->GetRightMargin();
0533 const auto &frac = float(elem.value().first) / ref_ali.size();
0534
0535 LogDebug("TrackerAlignmentCompareCylindricalBase")
0536 << __PRETTY_FUNCTION__ << " left margin: " << lm << " right margin: " << rm << " fraction: " << frac;
0537
0538 float theX_ = lm + (rm - lm) * frac + ((elem.index() > 0 || doOnlyPixel) ? 0.025 : 0.01);
0539
0540 tSubdet[subpad].DrawLatex(
0541 theX_, 0.23, Form("%s", AlignmentPI::getStringFromPart(elem.value().second, ph2).c_str()));
0542 }
0543
0544 auto ltx = TLatex();
0545 ltx.SetTextFont(62);
0546 ltx.SetTextSize(0.042);
0547 ltx.SetTextAlign(11);
0548 ltx.DrawLatexNDC(canvas.cd(subpad + 1)->GetLeftMargin(),
0549 1 - canvas.cd(subpad + 1)->GetTopMargin() + 0.01,
0550 ("Tracker Alignment Compare : #color[4]{" + s_coord + "}").c_str());
0551 legend.Draw("same");
0552 subpad++;
0553 }
0554
0555 std::string fileName(this->m_imageFileName);
0556 canvas.SaveAs(fileName.c_str());
0557
0558 return true;
0559 }
0560
0561 private:
0562
0563 inline std::string getStringFromCoordinate(coordinate coord)
0564
0565 {
0566 switch (coord) {
0567 case t_r:
0568 return "r-translation";
0569 case t_phi:
0570 return "#phi-rotation";
0571 case t_z:
0572 return "z-translation";
0573 default:
0574 return "should never be here!";
0575 }
0576 }
0577 };
0578
0579 typedef TrackerAlignmentCompareCylindricalBase<1, MULTI_IOV, RegionCategory::ALL>
0580 TrackerAlignmentCompareRPhiZSingleTag;
0581 typedef TrackerAlignmentCompareCylindricalBase<2, SINGLE_IOV, RegionCategory::ALL> TrackerAlignmentCompareRPhiZTwoTags;
0582
0583 typedef TrackerAlignmentCompareCylindricalBase<1, MULTI_IOV, RegionCategory::INNER>
0584 PixelAlignmentCompareRPhiZSingleTag;
0585 typedef TrackerAlignmentCompareCylindricalBase<2, SINGLE_IOV, RegionCategory::INNER> PixelAlignmentCompareRPhiZTwoTags;
0586
0587 typedef TrackerAlignmentCompareCylindricalBase<1, MULTI_IOV, RegionCategory::OUTER> OTAlignmentCompareRPhiZSingleTag;
0588 typedef TrackerAlignmentCompareCylindricalBase<2, SINGLE_IOV, RegionCategory::OUTER> OTAlignmentCompareRPhiZTwoTags;
0589
0590
0591
0592
0593
0594
0595 template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
0596 class TrackerAlignmentComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0597 public:
0598 TrackerAlignmentComparatorBase()
0599 : PlotImage<Alignments, nIOVs, ntags>("comparison of " + AlignmentPI::getStringFromCoordinate(coord) +
0600 " coordinate between two geometries") {}
0601
0602 bool fill() override {
0603 TGaxis::SetExponentOffset(-0.12, 0.01, "y");
0604
0605
0606 auto theIOVs = PlotBase::getTag<0>().iovs;
0607 auto tagname1 = PlotBase::getTag<0>().name;
0608 std::string tagname2 = "";
0609 auto firstiov = theIOVs.front();
0610 std::tuple<cond::Time_t, cond::Hash> lastiov;
0611
0612
0613 assert(this->m_plotAnnotations.ntags < 3);
0614
0615 if (this->m_plotAnnotations.ntags == 2) {
0616 auto tag2iovs = PlotBase::getTag<1>().iovs;
0617 tagname2 = PlotBase::getTag<1>().name;
0618 lastiov = tag2iovs.front();
0619 } else {
0620 lastiov = theIOVs.back();
0621 }
0622
0623 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0624 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0625
0626 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0627 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0628
0629 std::vector<AlignTransform> ref_ali = first_payload->m_align;
0630 std::vector<AlignTransform> target_ali = last_payload->m_align;
0631
0632 TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1200, 1200);
0633
0634 if (ref_ali.size() != target_ali.size()) {
0635 edm::LogError("TrackerAlignment_PayloadInspector")
0636 << "the size of the reference alignment (" << ref_ali.size()
0637 << ") is different from the one of the target (" << target_ali.size()
0638 << ")! You are probably trying to compare different underlying geometries. Exiting";
0639 return false;
0640 }
0641
0642 const bool ph2 = (ref_ali.size() > AlignmentPI::phase1size);
0643
0644
0645 const char *path_toTopologyXML = nullptr;
0646 if (ph2) {
0647 if (AlignmentPI::isReorderedTFPXTEPX(ref_ali) && AlignmentPI::isReorderedTFPXTEPX(target_ali)) {
0648 edm::LogPrint("TrackerAlignment_PayloadInspector")
0649 << "Both reference and target alignments are reordered. Using the trackerParameters for the Reordered "
0650 "TFPX,TEPX.";
0651 path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/TFPXTEPXReordered/trackerParameters.xml";
0652 } else if (!AlignmentPI::isReorderedTFPXTEPX(ref_ali) && !AlignmentPI::isReorderedTFPXTEPX(target_ali)) {
0653 edm::LogPrint("TrackerAlignment_PayloadInspector")
0654 << "Neither reference nor target alignments are reordered. Using the standard trackerParameters.";
0655 path_toTopologyXML = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
0656 } else {
0657
0658 edm::LogWarning("TrackerAlignment_PayloadInspector")
0659 << "Mismatched alignments detected. One is reordered while the other is not. Unable to proceed.";
0660 return false;
0661 }
0662 } else {
0663 path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0664 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0665 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0666 }
0667
0668 TrackerTopology tTopo =
0669 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0670
0671 for (const auto &ali : ref_ali) {
0672 auto mydetid = ali.rawId();
0673 if (DetId(mydetid).det() != DetId::Tracker) {
0674 edm::LogWarning("TrackerAlignment_PayloadInspector")
0675 << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0676 << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0677 << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0678 return false;
0679 }
0680 }
0681
0682 auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0683 std::string unit =
0684 (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0685
0686
0687 std::unique_ptr<TH1F> compare =
0688 std::make_unique<TH1F>("comparison",
0689 Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
0690 ref_ali.size(),
0691 -0.5,
0692 ref_ali.size() - 0.5);
0693
0694
0695 std::map<int, AlignmentPI::partitions> boundaries;
0696 boundaries.insert({0, AlignmentPI::BPix});
0697 AlignmentPI::fillComparisonHistogram(coord, boundaries, ref_ali, target_ali, compare);
0698
0699 canvas.cd();
0700
0701 canvas.SetTopMargin(0.06);
0702 canvas.SetLeftMargin(0.17);
0703 canvas.SetRightMargin(0.05);
0704 canvas.SetBottomMargin(0.15);
0705 AlignmentPI::makeNicePlotStyle(compare.get(), kBlack);
0706 auto max = compare->GetMaximum();
0707 auto min = compare->GetMinimum();
0708 auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
0709 if (range == 0.f)
0710 range = 0.1;
0711
0712
0713 compare->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
0714 compare->GetYaxis()->SetTitleOffset(1.5);
0715 compare->SetMarkerStyle(20);
0716 compare->SetMarkerSize(0.5);
0717 compare->Draw("P");
0718
0719 canvas.Update();
0720 canvas.cd();
0721
0722 TLine l[boundaries.size()];
0723 for (const auto &line : boundaries | boost::adaptors::indexed(0)) {
0724 const auto &index = line.index();
0725 const auto value = line.value();
0726 l[index] = TLine(compare->GetBinLowEdge(value.first),
0727 canvas.cd()->GetUymin(),
0728 compare->GetBinLowEdge(value.first),
0729 canvas.cd()->GetUymax());
0730 l[index].SetLineWidth(1);
0731 l[index].SetLineStyle(9);
0732 l[index].SetLineColor(2);
0733 l[index].Draw("same");
0734 }
0735
0736 TLatex tSubdet;
0737 tSubdet.SetNDC();
0738 tSubdet.SetTextAlign(21);
0739 tSubdet.SetTextSize(0.027);
0740 tSubdet.SetTextAngle(90);
0741
0742 for (const auto &elem : boundaries) {
0743 tSubdet.SetTextColor(kRed);
0744 auto myPair = AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(elem.first));
0745 float theX_ = elem.first != 0 ? myPair.first + 0.025 : myPair.first + 0.01;
0746 const bool isPhase2 = (ref_ali.size() > AlignmentPI::phase1size);
0747 tSubdet.DrawLatex(theX_, 0.20, Form("%s", AlignmentPI::getStringFromPart(elem.second, isPhase2).c_str()));
0748 }
0749
0750 TLegend legend = TLegend(0.17, 0.86, 0.95, 0.94);
0751 if (this->m_plotAnnotations.ntags == 2) {
0752 legend.SetHeader("#bf{Two Tags Comparison}", "C");
0753 legend.AddEntry(
0754 compare.get(),
0755 ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
0756 "PL");
0757 } else {
0758 legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");
0759 legend.AddEntry(compare.get(),
0760 ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0761 "PL");
0762 }
0763 legend.SetTextSize(0.020);
0764 legend.Draw("same");
0765
0766 auto ltx = TLatex();
0767 ltx.SetTextFont(62);
0768 ltx.SetTextSize(0.042);
0769 ltx.SetTextAlign(11);
0770 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0771 1 - gPad->GetTopMargin() + 0.01,
0772 ("Tracker Alignment Compare :#color[4]{" + s_coord + "}").c_str());
0773
0774 std::string fileName(this->m_imageFileName);
0775 canvas.SaveAs(fileName.c_str());
0776
0777 return true;
0778 }
0779 };
0780
0781 template <AlignmentPI::coordinate coord>
0782 using TrackerAlignmentCompare = TrackerAlignmentComparatorBase<coord, 1, MULTI_IOV>;
0783
0784 template <AlignmentPI::coordinate coord>
0785 using TrackerAlignmentCompareTwoTags = TrackerAlignmentComparatorBase<coord, 2, SINGLE_IOV>;
0786
0787 typedef TrackerAlignmentCompare<AlignmentPI::t_x> TrackerAlignmentCompareX;
0788 typedef TrackerAlignmentCompare<AlignmentPI::t_y> TrackerAlignmentCompareY;
0789 typedef TrackerAlignmentCompare<AlignmentPI::t_z> TrackerAlignmentCompareZ;
0790
0791 typedef TrackerAlignmentCompare<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlpha;
0792 typedef TrackerAlignmentCompare<AlignmentPI::rot_beta> TrackerAlignmentCompareBeta;
0793 typedef TrackerAlignmentCompare<AlignmentPI::rot_gamma> TrackerAlignmentCompareGamma;
0794
0795 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_x> TrackerAlignmentCompareXTwoTags;
0796 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_y> TrackerAlignmentCompareYTwoTags;
0797 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_z> TrackerAlignmentCompareZTwoTags;
0798
0799 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlphaTwoTags;
0800 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_beta> TrackerAlignmentCompareBetaTwoTags;
0801 typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_gamma> TrackerAlignmentCompareGammaTwoTags;
0802
0803
0804
0805
0806
0807 template <int ntags, IOVMultiplicity nIOVs, AlignmentPI::partitions q>
0808 class TrackerAlignmentSummaryBase : public PlotImage<Alignments, nIOVs, ntags> {
0809 public:
0810 TrackerAlignmentSummaryBase()
0811 : PlotImage<Alignments, nIOVs, ntags>("Comparison of all coordinates between two geometries for " +
0812 getStringFromPart(q)) {}
0813
0814 bool fill() override {
0815
0816 auto theIOVs = PlotBase::getTag<0>().iovs;
0817 auto tagname1 = PlotBase::getTag<0>().name;
0818 std::string tagname2 = "";
0819 auto firstiov = theIOVs.front();
0820 std::tuple<cond::Time_t, cond::Hash> lastiov;
0821
0822
0823 assert(this->m_plotAnnotations.ntags < 3);
0824
0825 if (this->m_plotAnnotations.ntags == 2) {
0826 auto tag2iovs = PlotBase::getTag<1>().iovs;
0827 tagname2 = PlotBase::getTag<1>().name;
0828 lastiov = tag2iovs.front();
0829 } else {
0830 lastiov = theIOVs.back();
0831 }
0832
0833 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0834 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0835
0836 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0837 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0838
0839 std::vector<AlignTransform> ref_ali = first_payload->m_align;
0840 std::vector<AlignTransform> target_ali = last_payload->m_align;
0841
0842 if (ref_ali.size() != target_ali.size()) {
0843 edm::LogError("TrackerAlignment_PayloadInspector")
0844 << "the size of the reference alignment (" << ref_ali.size()
0845 << ") is different from the one of the target (" << target_ali.size()
0846 << ")! You are probably trying to compare different underlying geometries. Exiting";
0847 return false;
0848 }
0849
0850
0851 const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0852 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0853 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0854 TrackerTopology tTopo =
0855 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0856
0857 for (const auto &ali : ref_ali) {
0858 auto mydetid = ali.rawId();
0859 if (DetId(mydetid).det() != DetId::Tracker) {
0860 edm::LogWarning("TrackerAlignment_PayloadInspector")
0861 << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0862 << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0863 << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0864 return false;
0865 }
0866 }
0867
0868 TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1800, 1200);
0869 canvas.Divide(3, 2);
0870
0871 std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F>> diffs;
0872 std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
0873 AlignmentPI::t_y,
0874 AlignmentPI::t_z,
0875 AlignmentPI::rot_alpha,
0876 AlignmentPI::rot_beta,
0877 AlignmentPI::rot_gamma};
0878
0879 for (const auto &coord : coords) {
0880 auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0881 std::string unit =
0882 (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0883
0884 diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s", s_coord.c_str()),
0885 Form(";#Delta%s %s;n. of modules", s_coord.c_str(), unit.c_str()),
0886 1001,
0887 -500.5,
0888 500.5);
0889 }
0890
0891
0892 std::map<int, AlignmentPI::partitions> boundaries;
0893 AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs, true, q);
0894
0895 int c_index = 1;
0896
0897
0898 auto legend = std::make_unique<TLegend>(0.14, 0.88, 0.96, 0.99);
0899 if (this->m_plotAnnotations.ntags == 2) {
0900 legend->SetHeader("#bf{Two Tags Comparison}", "C");
0901 legend->AddEntry(
0902 diffs[AlignmentPI::t_x].get(),
0903 ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
0904 "PL");
0905 } else {
0906 legend->SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");
0907 legend->AddEntry(diffs[AlignmentPI::t_x].get(),
0908 ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0909 "PL");
0910 }
0911 legend->SetTextSize(0.025);
0912
0913 for (const auto &coord : coords) {
0914 canvas.cd(c_index)->SetLogy();
0915 canvas.cd(c_index)->SetTopMargin(0.01);
0916 canvas.cd(c_index)->SetBottomMargin(0.15);
0917 canvas.cd(c_index)->SetLeftMargin(0.14);
0918 canvas.cd(c_index)->SetRightMargin(0.04);
0919 diffs[coord]->SetLineWidth(2);
0920 AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0921
0922
0923
0924
0925
0926
0927 int i_max = diffs[coord]->FindLastBinAbove(0.);
0928 int i_min = diffs[coord]->FindFirstBinAbove(0.);
0929 diffs[coord]->GetXaxis()->SetRange(std::max(1, i_min - 10), std::min(i_max + 10, diffs[coord]->GetNbinsX()));
0930 diffs[coord]->SetMaximum(diffs[coord]->GetMaximum() * 5);
0931 diffs[coord]->Draw("HIST");
0932 AlignmentPI::makeNiceStats(diffs[coord].get(), q, kBlack);
0933
0934 legend->Draw("same");
0935
0936 c_index++;
0937 }
0938
0939 std::string fileName(this->m_imageFileName);
0940 canvas.SaveAs(fileName.c_str());
0941
0942 return true;
0943 }
0944 };
0945
0946 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
0947 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
0948 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
0949
0950 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTID;
0951 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
0952 typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
0953
0954 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPixTwoTags;
0955 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPixTwoTags;
0956 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIBTwoTags;
0957
0958 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTIDTwoTags;
0959 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOBTwoTags;
0960 typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTECTwoTags;
0961
0962
0963
0964
0965 template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
0966 class PixelAlignmentComparisonMapBase : public PlotImage<Alignments, nIOVs, ntags> {
0967 public:
0968 PixelAlignmentComparisonMapBase()
0969 : PlotImage<Alignments, nIOVs, ntags>("SiPixel Comparison Map of " +
0970 AlignmentPI::getStringFromCoordinate(coord)) {
0971 label_ = "PixelAlignmentComparisonMap" + AlignmentPI::getStringFromCoordinate(coord);
0972 payloadString = "Tracker Alignment";
0973 }
0974
0975 bool fill() override {
0976 gStyle->SetPalette(1);
0977
0978
0979 auto theIOVs = PlotBase::getTag<0>().iovs;
0980 auto tagname1 = PlotBase::getTag<0>().name;
0981 std::string tagname2 = "";
0982 auto firstiov = theIOVs.front();
0983 std::tuple<cond::Time_t, cond::Hash> lastiov;
0984
0985
0986 assert(this->m_plotAnnotations.ntags < 3);
0987
0988 if (this->m_plotAnnotations.ntags == 2) {
0989 auto tag2iovs = PlotBase::getTag<1>().iovs;
0990 tagname2 = PlotBase::getTag<1>().name;
0991 lastiov = tag2iovs.front();
0992 } else {
0993 lastiov = theIOVs.back();
0994 }
0995
0996 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0997 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0998
0999 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1000 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1001
1002 const std::vector<AlignTransform> &ref_ali = first_payload->m_align;
1003 const std::vector<AlignTransform> &target_ali = last_payload->m_align;
1004
1005 if (last_payload.get() && first_payload.get()) {
1006 Phase1PixelSummaryMap fullMap(
1007 "",
1008 fmt::sprintf("%s %s", payloadString, AlignmentPI::getStringFromCoordinate(coord)),
1009 fmt::sprintf(
1010 "#Delta %s [%s]", AlignmentPI::getStringFromCoordinate(coord), (coord <= 3) ? "#mum" : "mrad"));
1011 fullMap.createTrackerBaseMap();
1012
1013 if (this->isPhase0(ref_ali) || this->isPhase0(target_ali)) {
1014 edm::LogError(label_) << "Pixel Tracker Alignment maps are not supported for non-Phase1 Pixel geometries !";
1015 TCanvas canvas("Canv", "Canv", 1200, 1000);
1016 AlignmentPI::displayNotSupported(canvas, 0);
1017 std::string fileName(this->m_imageFileName);
1018 canvas.SaveAs(fileName.c_str());
1019 return false;
1020 }
1021
1022
1023 std::map<uint32_t, double> diffPerDetid;
1024 this->fillPerDetIdDiff(coord, ref_ali, target_ali, diffPerDetid);
1025
1026
1027 for (const auto &elem : diffPerDetid) {
1028
1029 int subid = DetId(elem.first).subdetId();
1030 if (subid > 2) {
1031 continue;
1032 }
1033 fullMap.fillTrackerMap(elem.first, elem.second);
1034 }
1035
1036
1037
1038
1039 TCanvas canvas("Canv", "Canv", 3000, 2000);
1040 fullMap.printTrackerMap(canvas);
1041
1042 auto ltx = TLatex();
1043 ltx.SetTextFont(62);
1044 ltx.SetTextSize(0.025);
1045 ltx.SetTextAlign(11);
1046
1047 ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
1048 gPad->GetBottomMargin() + 0.01,
1049 ("#color[4]{" + tagname1 + "}, IOV: #color[4]{" + firstIOVsince + "} vs #color[4]{" +
1050 tagname2 + "}, IOV: #color[4]{" + lastIOVsince + "}")
1051 .c_str());
1052
1053 std::string fileName(this->m_imageFileName);
1054 canvas.SaveAs(fileName.c_str());
1055 }
1056 return true;
1057 }
1058
1059 protected:
1060 std::string payloadString;
1061 std::string label_;
1062
1063 private:
1064
1065 bool isPhase0(std::vector<AlignTransform> theAlis) {
1066 SiPixelDetInfoFileReader reader =
1067 SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
1068 const auto &p0detIds = reader.getAllDetIds();
1069
1070 std::vector<uint32_t> ownDetIds;
1071 std::transform(theAlis.begin(), theAlis.end(), std::back_inserter(ownDetIds), [](AlignTransform ali) -> uint32_t {
1072 return ali.rawId();
1073 });
1074
1075 for (const auto &det : ownDetIds) {
1076
1077 if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
1078 return true;
1079 }
1080 }
1081 return false;
1082 }
1083
1084
1085 void fillPerDetIdDiff(const AlignmentPI::coordinate &myCoord,
1086 const std::vector<AlignTransform> &ref_ali,
1087 const std::vector<AlignTransform> &target_ali,
1088 std::map<uint32_t, double> &diff)
1089
1090 {
1091 for (unsigned int i = 0; i < ref_ali.size(); i++) {
1092 uint32_t detid = ref_ali[i].rawId();
1093 if (ref_ali[i].rawId() == target_ali[i].rawId()) {
1094 CLHEP::HepRotation target_rot(target_ali[i].rotation());
1095 CLHEP::HepRotation ref_rot(ref_ali[i].rotation());
1096
1097 align::RotationType target_ROT(target_rot.xx(),
1098 target_rot.xy(),
1099 target_rot.xz(),
1100 target_rot.yx(),
1101 target_rot.yy(),
1102 target_rot.yz(),
1103 target_rot.zx(),
1104 target_rot.zy(),
1105 target_rot.zz());
1106
1107 align::RotationType ref_ROT(ref_rot.xx(),
1108 ref_rot.xy(),
1109 ref_rot.xz(),
1110 ref_rot.yx(),
1111 ref_rot.yy(),
1112 ref_rot.yz(),
1113 ref_rot.zx(),
1114 ref_rot.zy(),
1115 ref_rot.zz());
1116
1117 const std::vector<double> deltaRot = {
1118 ::deltaPhi(align::toAngles(target_ROT)[0], align::toAngles(ref_ROT)[0]),
1119 ::deltaPhi(align::toAngles(target_ROT)[1], align::toAngles(ref_ROT)[1]),
1120 ::deltaPhi(align::toAngles(target_ROT)[2], align::toAngles(ref_ROT)[2])};
1121
1122 const auto &deltaTrans = target_ali[i].translation() - ref_ali[i].translation();
1123
1124 switch (myCoord) {
1125 case AlignmentPI::t_x:
1126 diff.insert({detid, deltaTrans.x() * AlignmentPI::cmToUm});
1127 break;
1128 case AlignmentPI::t_y:
1129 diff.insert({detid, deltaTrans.y() * AlignmentPI::cmToUm});
1130 break;
1131 case AlignmentPI::t_z:
1132 diff.insert({detid, deltaTrans.z() * AlignmentPI::cmToUm});
1133 break;
1134 case AlignmentPI::rot_alpha:
1135 diff.insert({detid, deltaRot[0] * AlignmentPI::tomRad});
1136 break;
1137 case AlignmentPI::rot_beta:
1138 diff.insert({detid, deltaRot[1] * AlignmentPI::tomRad});
1139 break;
1140 case AlignmentPI::rot_gamma:
1141 diff.insert({detid, deltaRot[2] * AlignmentPI::tomRad});
1142 break;
1143 default:
1144 edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << myCoord << std::endl;
1145 break;
1146 }
1147 }
1148 }
1149 }
1150 };
1151
1152 template <AlignmentPI::coordinate coord>
1153 using PixelAlignmentCompareMap = PixelAlignmentComparisonMapBase<coord, 1, MULTI_IOV>;
1154
1155 template <AlignmentPI::coordinate coord>
1156 using PixelAlignmentCompareMapTwoTags = PixelAlignmentComparisonMapBase<coord, 2, SINGLE_IOV>;
1157
1158 typedef PixelAlignmentCompareMap<AlignmentPI::t_x> PixelAlignmentCompareMapX;
1159 typedef PixelAlignmentCompareMap<AlignmentPI::t_y> PixelAlignmentCompareMapY;
1160 typedef PixelAlignmentCompareMap<AlignmentPI::t_z> PixelAlignmentCompareMapZ;
1161
1162 typedef PixelAlignmentCompareMap<AlignmentPI::rot_alpha> PixelAlignmentCompareMapAlpha;
1163 typedef PixelAlignmentCompareMap<AlignmentPI::rot_beta> PixelAlignmentCompareMapBeta;
1164 typedef PixelAlignmentCompareMap<AlignmentPI::rot_gamma> PixelAlignmentCompareMapGamma;
1165
1166 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_x> PixelAlignmentCompareMapXTwoTags;
1167 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_y> PixelAlignmentCompareMapYTwoTags;
1168 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::t_z> PixelAlignmentCompareMapZTwoTags;
1169
1170 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_alpha> PixelAlignmentCompareMapAlphaTwoTags;
1171 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_beta> PixelAlignmentCompareMapBetaTwoTags;
1172 typedef PixelAlignmentCompareMapTwoTags<AlignmentPI::rot_gamma> PixelAlignmentCompareMapGammaTwoTags;
1173
1174
1175
1176
1177 template <AlignmentPI::partitions P, AlignmentPI::coordinate C, int NTAGS, IOVMultiplicity NIOVs>
1178 class PixelAlignmentCompareMapsBase : public PlotImage<Alignments, NIOVs, NTAGS> {
1179 public:
1180 PixelAlignmentCompareMapsBase()
1181 : PlotImage<Alignments, NIOVs, NTAGS>("SiPixel Comparison Map of " + AlignmentPI::getStringFromCoordinate(C)),
1182 label_("PixelAlignmentComparisonMap" + AlignmentPI::getStringFromCoordinate(C)),
1183 payloadStr_("Tracker Alignment") {}
1184
1185 bool fill() override {
1186 gStyle->SetPalette(1);
1187
1188 const auto &iovs = PlotBase::getTag<0>().iovs;
1189 const auto &tag1 = PlotBase::getTag<0>().name;
1190 const auto &[since1, hash1] = iovs.front();
1191 auto [since2, hash2] = (NTAGS == 2) ? PlotBase::getTag<1>().iovs.front() : iovs.back();
1192 const auto &tag2 = (NTAGS == 2) ? PlotBase::getTag<1>().name : tag1;
1193
1194 auto ref = this->fetchPayload(hash1);
1195 auto tgt = this->fetchPayload(hash2);
1196 if (!ref || !tgt)
1197 return false;
1198
1199 const auto &refAli = ref->m_align, tgtAli = tgt->m_align;
1200 if (isPhase0(refAli) || isPhase0(tgtAli)) {
1201 edm::LogError(label_) << "Not supported for non-Phase1 geometries.";
1202 TCanvas canv("Canv", "Canv", 1200, 1000);
1203 AlignmentPI::displayNotSupported(canv, 0);
1204 return false;
1205 }
1206
1207 Phase1PixelMaps maps("");
1208 maps.resetOption(P == AlignmentPI::INVALID ? "COLZA L" : "COLZL");
1209
1210 std::string coordStr = AlignmentPI::getStringFromCoordinate(C);
1211 std::string unit = (C <= 3) ? "#mum" : "mrad";
1212 std::string label = fmt::sprintf("#Delta %s [%s]", coordStr, unit);
1213
1214 std::string histoLabel = (P == AlignmentPI::BPix) ? "AlignmentBarrel"
1215 : (P == AlignmentPI::FPix) ? "AlignmentForward"
1216 : "Alignment";
1217
1218 if (P != AlignmentPI::FPix)
1219 maps.bookBarrelHistograms(histoLabel, (payloadStr_ + " " + coordStr).c_str(), label.c_str());
1220 if (P != AlignmentPI::BPix)
1221 maps.bookForwardHistograms(histoLabel, (payloadStr_ + " " + coordStr).c_str(), label.c_str());
1222
1223 std::map<uint32_t, double> diff;
1224 fillDiffs(refAli, tgtAli, diff);
1225
1226 for (const auto &[id, val] : diff) {
1227 int subid = DetId(id).subdetId();
1228 if (subid > 2)
1229 continue;
1230 if (P == AlignmentPI::INVALID || (P == AlignmentPI::BPix && subid == 1))
1231 maps.fillBarrelBin(histoLabel, id, val);
1232 if (P == AlignmentPI::INVALID || (P == AlignmentPI::FPix && subid == 2))
1233 maps.fillForwardBin(histoLabel, id, val);
1234 }
1235
1236 maps.beautifyAllHistograms();
1237
1238 TCanvas canvas("Canv", "Canv", (P == AlignmentPI::BPix ? 1200 : 1600), 1000);
1239 if (P == AlignmentPI::BPix)
1240 maps.drawBarrelMaps(histoLabel, canvas);
1241 else if (P == AlignmentPI::FPix)
1242 maps.drawForwardMaps(histoLabel, canvas);
1243 else
1244 maps.drawSummaryMaps(histoLabel, canvas);
1245
1246
1247 canvas.cd();
1248 TPaveText text(0, 0, 0.88, 0.04, "NDC");
1249 text.SetBorderSize(0);
1250 text.SetFillColor(0);
1251 text.AddText(
1252 fmt::sprintf(
1253 "#color[4]{%s}, IOV: #color[4]{%lu} vs #color[4]{%s}, IOV: #color[4]{%lu}", tag1, since1, tag2, since2)
1254 .c_str());
1255 text.Draw();
1256
1257 canvas.SaveAs(this->m_imageFileName.c_str());
1258 return true;
1259 }
1260
1261 protected:
1262 std::string label_, payloadStr_;
1263
1264 private:
1265 bool isPhase0(const std::vector<AlignTransform> &alis) const {
1266 static const auto phase0 =
1267 SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath()).getAllDetIds();
1268 return std::any_of(alis.begin(), alis.end(), [&](const auto &a) {
1269 return std::find(phase0.begin(), phase0.end(), a.rawId()) != phase0.end();
1270 });
1271 }
1272
1273 void fillDiffs(const std::vector<AlignTransform> &ref,
1274 const std::vector<AlignTransform> &tgt,
1275 std::map<uint32_t, double> &diff) const {
1276 for (size_t i = 0; i < ref.size(); ++i) {
1277 if (ref[i].rawId() != tgt[i].rawId())
1278 continue;
1279 auto id = ref[i].rawId();
1280 auto dt = tgt[i].translation() - ref[i].translation();
1281
1282 CLHEP::HepRotation dr = tgt[i].rotation();
1283 CLHEP::HepRotation rr = ref[i].rotation();
1284
1285 align::RotationType dr_align(dr.xx(), dr.xy(), dr.xz(), dr.yx(), dr.yy(), dr.yz(), dr.zx(), dr.zy(), dr.zz());
1286
1287 align::RotationType rr_align(rr.xx(), rr.xy(), rr.xz(), rr.yx(), rr.yy(), rr.yz(), rr.zx(), rr.zy(), rr.zz());
1288
1289 auto da = align::toAngles(dr_align);
1290 auto ra = align::toAngles(rr_align);
1291
1292 switch (C) {
1293 case AlignmentPI::t_x:
1294 diff[id] = dt.x() * AlignmentPI::cmToUm;
1295 break;
1296 case AlignmentPI::t_y:
1297 diff[id] = dt.y() * AlignmentPI::cmToUm;
1298 break;
1299 case AlignmentPI::t_z:
1300 diff[id] = dt.z() * AlignmentPI::cmToUm;
1301 break;
1302 case AlignmentPI::rot_alpha:
1303 diff[id] = deltaPhi(da[0], ra[0]) * AlignmentPI::tomRad;
1304 break;
1305 case AlignmentPI::rot_beta:
1306 diff[id] = deltaPhi(da[1], ra[1]) * AlignmentPI::tomRad;
1307 break;
1308 case AlignmentPI::rot_gamma:
1309 diff[id] = deltaPhi(da[2], ra[2]) * AlignmentPI::tomRad;
1310 break;
1311 default:
1312 break;
1313 }
1314 }
1315 }
1316 };
1317
1318
1319
1320 template <AlignmentPI::coordinate coord>
1321 using BPixAlignmentCompareMap = PixelAlignmentCompareMapsBase<AlignmentPI::BPix, coord, 1, MULTI_IOV>;
1322
1323 template <AlignmentPI::coordinate coord>
1324 using BPixAlignmentCompareMapTwoTags = PixelAlignmentCompareMapsBase<AlignmentPI::BPix, coord, 2, SINGLE_IOV>;
1325
1326 typedef BPixAlignmentCompareMap<AlignmentPI::t_x> BPixAlignmentCompareMapX;
1327 typedef BPixAlignmentCompareMap<AlignmentPI::t_y> BPixAlignmentCompareMapY;
1328 typedef BPixAlignmentCompareMap<AlignmentPI::t_z> BPixAlignmentCompareMapZ;
1329
1330 typedef BPixAlignmentCompareMap<AlignmentPI::rot_alpha> BPixAlignmentCompareMapAlpha;
1331 typedef BPixAlignmentCompareMap<AlignmentPI::rot_beta> BPixAlignmentCompareMapBeta;
1332 typedef BPixAlignmentCompareMap<AlignmentPI::rot_gamma> BPixAlignmentCompareMapGamma;
1333
1334 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::t_x> BPixAlignmentCompareMapXTwoTags;
1335 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::t_y> BPixAlignmentCompareMapYTwoTags;
1336 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::t_z> BPixAlignmentCompareMapZTwoTags;
1337
1338 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::rot_alpha> BPixAlignmentCompareMapAlphaTwoTags;
1339 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::rot_beta> BPixAlignmentCompareMapBetaTwoTags;
1340 typedef BPixAlignmentCompareMapTwoTags<AlignmentPI::rot_gamma> BPixAlignmentCompareMapGammaTwoTags;
1341
1342
1343
1344 template <AlignmentPI::coordinate coord>
1345 using FPixAlignmentCompareMap = PixelAlignmentCompareMapsBase<AlignmentPI::FPix, coord, 1, MULTI_IOV>;
1346
1347 template <AlignmentPI::coordinate coord>
1348 using FPixAlignmentCompareMapTwoTags = PixelAlignmentCompareMapsBase<AlignmentPI::FPix, coord, 2, SINGLE_IOV>;
1349
1350 typedef FPixAlignmentCompareMap<AlignmentPI::t_x> FPixAlignmentCompareMapX;
1351 typedef FPixAlignmentCompareMap<AlignmentPI::t_y> FPixAlignmentCompareMapY;
1352 typedef FPixAlignmentCompareMap<AlignmentPI::t_z> FPixAlignmentCompareMapZ;
1353
1354 typedef FPixAlignmentCompareMap<AlignmentPI::rot_alpha> FPixAlignmentCompareMapAlpha;
1355 typedef FPixAlignmentCompareMap<AlignmentPI::rot_beta> FPixAlignmentCompareMapBeta;
1356 typedef FPixAlignmentCompareMap<AlignmentPI::rot_gamma> FPixAlignmentCompareMapGamma;
1357
1358 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::t_x> FPixAlignmentCompareMapXTwoTags;
1359 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::t_y> FPixAlignmentCompareMapYTwoTags;
1360 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::t_z> FPixAlignmentCompareMapZTwoTags;
1361
1362 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::rot_alpha> FPixAlignmentCompareMapAlphaTwoTags;
1363 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::rot_beta> FPixAlignmentCompareMapBetaTwoTags;
1364 typedef FPixAlignmentCompareMapTwoTags<AlignmentPI::rot_gamma> FPixAlignmentCompareMapGammaTwoTags;
1365
1366
1367
1368
1369
1370 template <AlignmentPI::coordinate coord>
1371 class BPixBarycenterHistory : public HistoryPlot<Alignments, float> {
1372 public:
1373 BPixBarycenterHistory()
1374 : HistoryPlot<Alignments, float>(
1375 " Barrel Pixel " + AlignmentPI::getStringFromCoordinate(coord) + " positions vs time",
1376 AlignmentPI::getStringFromCoordinate(coord) + " position [cm]") {}
1377 ~BPixBarycenterHistory() override = default;
1378
1379 float getFromPayload(Alignments &payload) override {
1380 std::vector<AlignTransform> alignments = payload.m_align;
1381
1382 float barycenter = 0.;
1383 float nmodules(0.);
1384 for (const auto &ali : alignments) {
1385 if (DetId(ali.rawId()).det() != DetId::Tracker) {
1386 edm::LogWarning("TrackerAlignment_PayloadInspector")
1387 << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
1388 << " is different from " << DetId::Tracker << " - terminating ";
1389 return false;
1390 }
1391
1392 int subid = DetId(ali.rawId()).subdetId();
1393 if (subid != PixelSubdetector::PixelBarrel)
1394 continue;
1395
1396 nmodules++;
1397 switch (coord) {
1398 case AlignmentPI::t_x:
1399 barycenter += (ali.translation().x());
1400 break;
1401 case AlignmentPI::t_y:
1402 barycenter += (ali.translation().y());
1403 break;
1404 case AlignmentPI::t_z:
1405 barycenter += (ali.translation().z());
1406 break;
1407 default:
1408 edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
1409 break;
1410 }
1411 }
1412
1413 edm::LogInfo("TrackerAlignment_PayloadInspector") << "barycenter (" << barycenter << ")/n. modules (" << nmodules
1414 << ") = " << barycenter / nmodules << std::endl;
1415
1416
1417 barycenter /= nmodules;
1418
1419
1420 barycenter += hardcodeGPR.at(coord);
1421
1422 return barycenter;
1423
1424 }
1425 };
1426
1427 typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
1428 typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
1429 typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
1430
1431
1432
1433
1434 class TrackerAlignmentBarycenters : public PlotImage<Alignments, SINGLE_IOV> {
1435 public:
1436 TrackerAlignmentBarycenters() : PlotImage<Alignments, SINGLE_IOV>("Display of Tracker Alignment Barycenters") {}
1437
1438 bool fill() override {
1439 auto tag = PlotBase::getTag<0>();
1440 auto iov = tag.iovs.front();
1441 const auto &tagname = PlotBase::getTag<0>().name;
1442 std::shared_ptr<Alignments> payload = fetchPayload(std::get<1>(iov));
1443 unsigned int run = std::get<0>(iov);
1444
1445 TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1600, 1000);
1446 canvas.cd();
1447
1448 canvas.SetTopMargin(0.07);
1449 canvas.SetBottomMargin(0.06);
1450 canvas.SetLeftMargin(0.15);
1451 canvas.SetRightMargin(0.03);
1452 canvas.Modified();
1453 canvas.SetGrid();
1454
1455 auto h2_BarycenterParameters =
1456 std::make_unique<TH2F>("Parameters", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
1457
1458 auto h2_uncBarycenterParameters =
1459 std::make_unique<TH2F>("Parameters2", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
1460
1461 h2_BarycenterParameters->SetStats(false);
1462 h2_BarycenterParameters->SetTitle(nullptr);
1463 h2_uncBarycenterParameters->SetStats(false);
1464 h2_uncBarycenterParameters->SetTitle(nullptr);
1465
1466 std::vector<AlignTransform> alignments = payload->m_align;
1467
1468 isPhase0 = (alignments.size() == AlignmentPI::phase0size) ? true : false;
1469
1470
1471 const char *path_toTopologyXML = isPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1472 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1473
1474 TrackerTopology tTopo =
1475 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1476
1477 AlignmentPI::TkAlBarycenters barycenters;
1478
1479 barycenters.computeBarycenters(
1480 alignments, tTopo, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1481
1482 auto Xbarycenters = barycenters.getX();
1483 auto Ybarycenters = barycenters.getY();
1484 auto Zbarycenters = barycenters.getZ();
1485
1486
1487 barycenters.computeBarycenters(alignments, tTopo, hardcodeGPR);
1488
1489 auto c_Xbarycenters = barycenters.getX();
1490 auto c_Ybarycenters = barycenters.getY();
1491 auto c_Zbarycenters = barycenters.getZ();
1492
1493 h2_BarycenterParameters->GetXaxis()->SetBinLabel(1, "X [cm]");
1494 h2_BarycenterParameters->GetXaxis()->SetBinLabel(2, "Y [cm]");
1495 h2_BarycenterParameters->GetXaxis()->SetBinLabel(3, "Z [cm]");
1496 h2_BarycenterParameters->GetXaxis()->SetBinLabel(4, "X_{no GPR} [cm]");
1497 h2_BarycenterParameters->GetXaxis()->SetBinLabel(5, "Y_{no GPR} [cm]");
1498 h2_BarycenterParameters->GetXaxis()->SetBinLabel(6, "Z_{no GPR} [cm]");
1499
1500 bool isLikelyMC(false);
1501 int checkX =
1502 std::count_if(Xbarycenters.begin(), Xbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1503 int checkY =
1504 std::count_if(Ybarycenters.begin(), Ybarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1505 int checkZ =
1506 std::count_if(Zbarycenters.begin(), Zbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
1507
1508
1509
1510 if ((checkX + checkY + checkZ) == 0 && run == 1)
1511 isLikelyMC = true;
1512
1513 unsigned int yBin = 6;
1514 for (unsigned int i = 0; i < 6; i++) {
1515 auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
1516 std::string theLabel = getStringFromPart(thePart);
1517 h2_BarycenterParameters->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
1518 if (!isLikelyMC) {
1519 h2_BarycenterParameters->SetBinContent(1, yBin, c_Xbarycenters[i]);
1520 h2_BarycenterParameters->SetBinContent(2, yBin, c_Ybarycenters[i]);
1521 h2_BarycenterParameters->SetBinContent(3, yBin, c_Zbarycenters[i]);
1522 }
1523
1524 h2_uncBarycenterParameters->SetBinContent(4, yBin, Xbarycenters[i]);
1525 h2_uncBarycenterParameters->SetBinContent(5, yBin, Ybarycenters[i]);
1526 h2_uncBarycenterParameters->SetBinContent(6, yBin, Zbarycenters[i]);
1527 yBin--;
1528 }
1529
1530 h2_BarycenterParameters->GetXaxis()->LabelsOption("h");
1531 h2_BarycenterParameters->GetYaxis()->SetLabelSize(0.05);
1532 h2_BarycenterParameters->GetXaxis()->SetLabelSize(0.05);
1533 h2_BarycenterParameters->SetMarkerSize(1.5);
1534 h2_BarycenterParameters->Draw("TEXT");
1535
1536 h2_uncBarycenterParameters->SetMarkerSize(1.5);
1537 h2_uncBarycenterParameters->SetMarkerColor(kRed);
1538 h2_uncBarycenterParameters->Draw("TEXTsame");
1539
1540 TLatex t1;
1541 t1.SetNDC();
1542 t1.SetTextAlign(26);
1543 t1.SetTextSize(0.045);
1544 t1.DrawLatex(0.5, 0.96, Form("TkAl Barycenters, Tag: #color[4]{%s}, IOV #color[4]{%i}", tagname.c_str(), run));
1545 t1.SetTextSize(0.025);
1546
1547 std::string fileName(m_imageFileName);
1548 canvas.SaveAs(fileName.c_str());
1549
1550 return true;
1551 }
1552
1553 private:
1554 bool isPhase0;
1555 };
1556
1557
1558
1559
1560 template <int ntags, IOVMultiplicity nIOVs>
1561 class TrackerAlignmentBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
1562 public:
1563 TrackerAlignmentBarycentersComparatorBase()
1564 : PlotImage<Alignments, nIOVs, ntags>("Comparison of Tracker Alignment Barycenters") {}
1565
1566 bool fill() override {
1567
1568 auto theIOVs = PlotBase::getTag<0>().iovs;
1569 auto tagname1 = PlotBase::getTag<0>().name;
1570 std::string tagname2 = "";
1571 auto firstiov = theIOVs.front();
1572 std::tuple<cond::Time_t, cond::Hash> lastiov;
1573
1574
1575 assert(this->m_plotAnnotations.ntags < 3);
1576
1577 if (this->m_plotAnnotations.ntags == 2) {
1578 auto tag2iovs = PlotBase::getTag<1>().iovs;
1579 tagname2 = PlotBase::getTag<1>().name;
1580 lastiov = tag2iovs.front();
1581 } else {
1582 lastiov = theIOVs.back();
1583 }
1584
1585 unsigned int first_run = std::get<0>(firstiov);
1586 unsigned int last_run = std::get<0>(lastiov);
1587
1588 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
1589 std::vector<AlignTransform> last_alignments = last_payload->m_align;
1590
1591 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
1592 std::vector<AlignTransform> first_alignments = first_payload->m_align;
1593
1594 isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
1595 isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
1596
1597
1598 const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1599 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1600
1601 TrackerTopology tTopo_f =
1602 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1603
1604 path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1605 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1606
1607 TrackerTopology tTopo_l =
1608 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1609
1610 TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1200, 800);
1611 canvas.cd();
1612
1613 canvas.SetTopMargin(0.07);
1614 canvas.SetBottomMargin(0.06);
1615 canvas.SetLeftMargin(0.15);
1616 canvas.SetRightMargin(0.03);
1617 canvas.Modified();
1618 canvas.SetGrid();
1619
1620 auto h2_BarycenterDiff =
1621 std::make_unique<TH2F>("Parameters diff", "SubDetector Barycenter Difference", 3, 0.0, 3.0, 6, 0, 6.);
1622
1623 h2_BarycenterDiff->SetStats(false);
1624 h2_BarycenterDiff->SetTitle(nullptr);
1625 h2_BarycenterDiff->GetXaxis()->SetBinLabel(1, "X [#mum]");
1626 h2_BarycenterDiff->GetXaxis()->SetBinLabel(2, "Y [#mum]");
1627 h2_BarycenterDiff->GetXaxis()->SetBinLabel(3, "Z [#mum]");
1628
1629 AlignmentPI::TkAlBarycenters l_barycenters;
1630 l_barycenters.computeBarycenters(last_alignments, tTopo_l, hardcodeGPR);
1631
1632 AlignmentPI::TkAlBarycenters f_barycenters;
1633 f_barycenters.computeBarycenters(first_alignments, tTopo_f, hardcodeGPR);
1634
1635 unsigned int yBin = 6;
1636 for (unsigned int i = 0; i < 6; i++) {
1637 auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
1638 std::string theLabel = getStringFromPart(thePart);
1639 h2_BarycenterDiff->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
1640 h2_BarycenterDiff->SetBinContent(
1641 1, yBin, (l_barycenters.getX()[i] - f_barycenters.getX()[i]) * AlignmentPI::cmToUm);
1642 h2_BarycenterDiff->SetBinContent(
1643 2, yBin, (l_barycenters.getY()[i] - f_barycenters.getY()[i]) * AlignmentPI::cmToUm);
1644 h2_BarycenterDiff->SetBinContent(
1645 3, yBin, (l_barycenters.getZ()[i] - f_barycenters.getZ()[i]) * AlignmentPI::cmToUm);
1646 yBin--;
1647 }
1648
1649 h2_BarycenterDiff->GetXaxis()->LabelsOption("h");
1650 h2_BarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1651 h2_BarycenterDiff->GetXaxis()->SetLabelSize(0.05);
1652 h2_BarycenterDiff->SetMarkerSize(1.5);
1653 h2_BarycenterDiff->SetMarkerColor(kRed);
1654 h2_BarycenterDiff->Draw("TEXT");
1655
1656 TLatex t1;
1657 t1.SetNDC();
1658 t1.SetTextAlign(26);
1659 t1.SetTextSize(0.05);
1660 t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters Diff, IOV %i - IOV %i", last_run, first_run));
1661 t1.SetTextSize(0.025);
1662
1663 std::string fileName(this->m_imageFileName);
1664 canvas.SaveAs(fileName.c_str());
1665
1666 return true;
1667 }
1668
1669 private:
1670 bool isInitialPhase0;
1671 bool isFinalPhase0;
1672 };
1673
1674 using TrackerAlignmentBarycentersCompare = TrackerAlignmentBarycentersComparatorBase<1, MULTI_IOV>;
1675 using TrackerAlignmentBarycentersCompareTwoTags = TrackerAlignmentBarycentersComparatorBase<2, SINGLE_IOV>;
1676
1677
1678
1679
1680 template <int ntags, IOVMultiplicity nIOVs>
1681 class PixelBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
1682 public:
1683 PixelBarycentersComparatorBase() : PlotImage<Alignments, nIOVs, ntags>("Comparison of Pixel Barycenters") {}
1684
1685 bool fill() override {
1686
1687 auto theIOVs = PlotBase::getTag<0>().iovs;
1688 auto tagname1 = PlotBase::getTag<0>().name;
1689 std::string tagname2 = "";
1690 auto firstiov = theIOVs.front();
1691 std::tuple<cond::Time_t, cond::Hash> lastiov;
1692
1693
1694 assert(this->m_plotAnnotations.ntags < 3);
1695
1696 if (this->m_plotAnnotations.ntags == 2) {
1697 auto tag2iovs = PlotBase::getTag<1>().iovs;
1698 tagname2 = PlotBase::getTag<1>().name;
1699 lastiov = tag2iovs.front();
1700 } else {
1701 lastiov = theIOVs.back();
1702 }
1703
1704 unsigned int first_run = std::get<0>(firstiov);
1705 unsigned int last_run = std::get<0>(lastiov);
1706
1707 std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
1708 std::vector<AlignTransform> last_alignments = last_payload->m_align;
1709
1710 std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
1711 std::vector<AlignTransform> first_alignments = first_payload->m_align;
1712
1713 TCanvas canvas("Pixel Barycenter Summary", "Pixel Barycenter summary", 1200, 1200);
1714 canvas.Divide(2, 2);
1715 canvas.cd();
1716
1717 TLatex t1;
1718 t1.SetNDC();
1719 t1.SetTextAlign(26);
1720 t1.SetTextSize(0.03);
1721 t1.DrawLatex(0.5,
1722 0.97,
1723 ("Pixel Barycenters comparison, IOV: #color[2]{" + std::to_string(first_run) +
1724 "} vs IOV: #color[4]{" + std::to_string(last_run) + "}")
1725 .c_str());
1726 t1.SetTextSize(0.025);
1727
1728 for (unsigned int c = 1; c <= 4; c++) {
1729 canvas.cd(c)->SetTopMargin(0.07);
1730 canvas.cd(c)->SetBottomMargin(0.12);
1731 canvas.cd(c)->SetLeftMargin(0.15);
1732 canvas.cd(c)->SetRightMargin(0.03);
1733 canvas.cd(c)->Modified();
1734 canvas.cd(c)->SetGrid();
1735 }
1736
1737 std::array<std::string, 3> structures = {{"FPIX-", "BPIX", "FPIX+"}};
1738 std::array<std::unique_ptr<TH2F>, 3> histos;
1739
1740 isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
1741 isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
1742
1743
1744 const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1745 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1746
1747 TrackerTopology tTopo_f =
1748 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1749
1750 AlignmentPI::TkAlBarycenters myInitialBarycenters;
1751
1752 myInitialBarycenters.computeBarycenters(
1753 first_alignments, tTopo_f, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1754
1755 path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1756 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1757
1758 TrackerTopology tTopo_l =
1759 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1760
1761 AlignmentPI::TkAlBarycenters myFinalBarycenters;
1762
1763 myFinalBarycenters.computeBarycenters(
1764 last_alignments, tTopo_l, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1765
1766 if (isFinalPhase0 != isInitialPhase0) {
1767 edm::LogWarning("TrackerAlignment_PayloadInspector")
1768 << "the size of the reference alignment (" << first_alignments.size()
1769 << ") is different from the one of the target (" << last_alignments.size()
1770 << ")! You are probably trying to compare different underlying geometries.";
1771 }
1772
1773 unsigned int index(0);
1774 for (const auto &piece : structures) {
1775 const char *name = piece.c_str();
1776 histos[index] = std::make_unique<TH2F>(
1777 name,
1778 Form("%s x-y Barycenter Difference;x_{%s}-x_{TOB} [mm];y_{%s}-y_{TOB} [mm]", name, name, name),
1779 100,
1780 -3.,
1781 3.,
1782 100,
1783 -3.,
1784 3.);
1785
1786 histos[index]->SetStats(false);
1787 histos[index]->SetTitle(nullptr);
1788 histos[index]->GetYaxis()->SetLabelSize(0.05);
1789 histos[index]->GetXaxis()->SetLabelSize(0.05);
1790 histos[index]->GetYaxis()->SetTitleSize(0.06);
1791 histos[index]->GetXaxis()->SetTitleSize(0.06);
1792 histos[index]->GetYaxis()->CenterTitle();
1793 histos[index]->GetXaxis()->CenterTitle();
1794 histos[index]->GetXaxis()->SetTitleOffset(0.9);
1795 index++;
1796 }
1797
1798 auto h2_ZBarycenterDiff = std::make_unique<TH2F>(
1799 "Pixel_z_diff", "Pixel z-Barycenter Difference;; z_{Pixel-Ideal} -z_{TOB} [mm]", 3, -0.5, 2.5, 100, -10., 10.);
1800 h2_ZBarycenterDiff->SetStats(false);
1801 h2_ZBarycenterDiff->SetTitle(nullptr);
1802 h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(1, "FPIX -");
1803 h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(2, "BPIX");
1804 h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(3, "FPIX +");
1805 h2_ZBarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1806 h2_ZBarycenterDiff->GetXaxis()->SetLabelSize(0.07);
1807 h2_ZBarycenterDiff->GetYaxis()->SetTitleSize(0.06);
1808 h2_ZBarycenterDiff->GetXaxis()->SetTitleSize(0.06);
1809 h2_ZBarycenterDiff->GetYaxis()->CenterTitle();
1810 h2_ZBarycenterDiff->GetXaxis()->CenterTitle();
1811 h2_ZBarycenterDiff->GetYaxis()->SetTitleOffset(1.1);
1812
1813 std::function<GlobalPoint(int)> cutFunctorInitial = [&myInitialBarycenters](int index) {
1814 switch (index) {
1815 case 1:
1816 return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1817 case 2:
1818 return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1819 case 3:
1820 return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1821 default:
1822 return GlobalPoint(0, 0, 0);
1823 }
1824 };
1825
1826 std::function<GlobalPoint(int)> cutFunctorFinal = [&myFinalBarycenters](int index) {
1827 switch (index) {
1828 case 1:
1829 return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1830 case 2:
1831 return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1832 case 3:
1833 return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1834 default:
1835 return GlobalPoint(0, 0, 0);
1836 }
1837 };
1838
1839 float x0i, x0f, y0i, y0f;
1840
1841 t1.SetNDC(kFALSE);
1842 t1.SetTextSize(0.047);
1843 for (unsigned int c = 1; c <= 3; c++) {
1844 x0i = cutFunctorInitial(c).x() * 10;
1845 x0f = cutFunctorFinal(c).x() * 10;
1846 y0i = cutFunctorInitial(c).y() * 10;
1847 y0f = cutFunctorFinal(c).y() * 10;
1848
1849 canvas.cd(c);
1850 histos[c - 1]->Draw();
1851
1852 COUT << "initial x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0i << "," << y0i << ") mm"
1853 << std::endl;
1854 COUT << "final x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0f << "," << y0f << ") mm"
1855 << std::endl;
1856
1857 TMarker *initial = new TMarker(x0i, y0i, 21);
1858 TMarker *final = new TMarker(x0f, y0f, 20);
1859
1860 initial->SetMarkerColor(kRed);
1861 final->SetMarkerColor(kBlue);
1862 initial->SetMarkerSize(2.5);
1863 final->SetMarkerSize(2.5);
1864 t1.SetTextColor(kRed);
1865 initial->Draw();
1866 t1.DrawLatex(x0i, y0i + (y0i > y0f ? 0.3 : -0.5), Form("(%.2f,%.2f)", x0i, y0i));
1867 final->Draw("same");
1868 t1.SetTextColor(kBlue);
1869 t1.DrawLatex(x0f, y0f + (y0i > y0f ? -0.5 : 0.3), Form("(%.2f,%.2f)", x0f, y0f));
1870 }
1871
1872
1873 canvas.cd(4);
1874 h2_ZBarycenterDiff->Draw();
1875 float z0i, z0f;
1876
1877
1878
1879 std::array<double, 3> hardcodeIdealZPhase0 = {{-41.94909, 0., 41.94909}};
1880 std::array<double, 3> hardcodeIdealZPhase1 = {{-39.82911, 0., 39.82911}};
1881
1882 for (unsigned int c = 1; c <= 3; c++) {
1883
1884
1885 z0i =
1886 (cutFunctorInitial(c).z() - (isInitialPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) *
1887 10;
1888 z0f =
1889 (cutFunctorFinal(c).z() - (isFinalPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) * 10;
1890
1891 TMarker *initial = new TMarker(c - 1, z0i, 21);
1892 TMarker *final = new TMarker(c - 1, z0f, 20);
1893
1894 COUT << "initial z " << std::left << std::setw(7) << structures[c - 1] << " " << z0i << " mm" << std::endl;
1895 COUT << "final z " << std::left << std::setw(7) << structures[c - 1] << " " << z0f << " mm" << std::endl;
1896
1897 initial->SetMarkerColor(kRed);
1898 final->SetMarkerColor(kBlue);
1899 initial->SetMarkerSize(2.5);
1900 final->SetMarkerSize(2.5);
1901 initial->Draw();
1902 t1.SetTextColor(kRed);
1903 t1.DrawLatex(c - 1, z0i + (z0i > z0f ? 1. : -1.5), Form("(%.2f)", z0i));
1904 final->Draw("same");
1905 t1.SetTextColor(kBlue);
1906 t1.DrawLatex(c - 1, z0f + (z0i > z0f ? -1.5 : 1), Form("(%.2f)", z0f));
1907 }
1908
1909 std::string fileName(this->m_imageFileName);
1910 canvas.SaveAs(fileName.c_str());
1911
1912 return true;
1913 }
1914
1915 private:
1916 bool isInitialPhase0;
1917 bool isFinalPhase0;
1918 };
1919
1920 using PixelBarycentersCompare = PixelBarycentersComparatorBase<1, MULTI_IOV>;
1921 using PixelBarycentersCompareTwoTags = PixelBarycentersComparatorBase<2, SINGLE_IOV>;
1922
1923 }
1924
1925 PAYLOAD_INSPECTOR_MODULE(TrackerAlignment) {
1926 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentComparatorSingleTag);
1927 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentComparatorTwoTags);
1928 PAYLOAD_INSPECTOR_CLASS(OTAlignmentComparatorSingleTag);
1929 PAYLOAD_INSPECTOR_CLASS(OTAlignmentComparatorTwoTags);
1930 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorSingleTag);
1931 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorTwoTags);
1932 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareRPhiZSingleTag);
1933 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareRPhiZTwoTags);
1934 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareRPhiZSingleTag);
1935 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareRPhiZTwoTags);
1936 PAYLOAD_INSPECTOR_CLASS(OTAlignmentCompareRPhiZSingleTag);
1937 PAYLOAD_INSPECTOR_CLASS(OTAlignmentCompareRPhiZTwoTags);
1938 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
1939 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
1940 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
1941 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
1942 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
1943 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
1944 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareXTwoTags);
1945 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareYTwoTags);
1946 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZTwoTags);
1947 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlphaTwoTags);
1948 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBetaTwoTags);
1949 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGammaTwoTags);
1950 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapX);
1951 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapY);
1952 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapZ);
1953 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapAlpha);
1954 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapBeta);
1955 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapGamma);
1956 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapXTwoTags);
1957 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapYTwoTags);
1958 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapZTwoTags);
1959 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapAlphaTwoTags);
1960 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapBetaTwoTags);
1961 PAYLOAD_INSPECTOR_CLASS(PixelAlignmentCompareMapGammaTwoTags);
1962 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapX);
1963 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapY);
1964 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapZ);
1965 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapAlpha);
1966 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapBeta);
1967 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapGamma);
1968 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapXTwoTags);
1969 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapYTwoTags);
1970 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapZTwoTags);
1971 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapAlphaTwoTags);
1972 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapBetaTwoTags);
1973 PAYLOAD_INSPECTOR_CLASS(BPixAlignmentCompareMapGammaTwoTags);
1974 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapX);
1975 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapY);
1976 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapZ);
1977 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapAlpha);
1978 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapBeta);
1979 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapGamma);
1980 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapXTwoTags);
1981 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapYTwoTags);
1982 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapZTwoTags);
1983 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapAlphaTwoTags);
1984 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapBetaTwoTags);
1985 PAYLOAD_INSPECTOR_CLASS(FPixAlignmentCompareMapGammaTwoTags);
1986 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
1987 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
1988 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
1989 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
1990 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
1991 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
1992 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPixTwoTags);
1993 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPixTwoTags);
1994 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIBTwoTags);
1995 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIDTwoTags);
1996 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOBTwoTags);
1997 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTECTwoTags);
1998 PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
1999 PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
2000 PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
2001 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycenters);
2002 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompare);
2003 PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompareTwoTags);
2004 PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompare);
2005 PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompareTwoTags);
2006 }