Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:24

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