Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-10-20 23:30:36

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 
0031 #include <memory>
0032 #include <sstream>
0033 #include <iostream>
0034 #include <iomanip>  // std::setprecision
0035 
0036 // include ROOT
0037 #include "TH2F.h"
0038 #include "TGaxis.h"
0039 #include "TLegend.h"
0040 #include "TCanvas.h"
0041 #include "TLine.h"
0042 #include "TStyle.h"
0043 #include "TLatex.h"
0044 #include "TPave.h"
0045 #include "TMarker.h"
0046 #include "TPaveStats.h"
0047 
0048 namespace {
0049 
0050   using namespace cond::payloadInspector;
0051 
0052   // M.M. 2017/09/29
0053   // Hardcoded Tracker Global Position Record
0054   // Without accessing the ES, it is not possible to access to the GPR with the PI technology,
0055   // so this needs to be hardcoded.
0056   // Anyway it is not likely to change until a new Tracker is installed.
0057   // Details at:
0058   // - https://indico.cern.ch/event/238026/contributions/513928/attachments/400000/556192/mm_TkAlMeeting_28_03_2013.pdf
0059   // - https://twiki.cern.ch/twiki/bin/view/CMS/TkAlignmentPixelPosition
0060 
0061   const std::map<AlignmentPI::coordinate, float> hardcodeGPR = {
0062       {AlignmentPI::t_x, -9.00e-02}, {AlignmentPI::t_y, -1.10e-01}, {AlignmentPI::t_z, -1.70e-01}};
0063 
0064   //*******************************************/
0065   // Size of the movement over all partitions,
0066   // one at a time
0067   //******************************************//
0068 
0069   template <int ntags, IOVMultiplicity nIOVs>
0070   class TrackerAlignmentCompareAll : public PlotImage<Alignments, nIOVs, ntags> {
0071   public:
0072     TrackerAlignmentCompareAll()
0073         : PlotImage<Alignments, nIOVs, ntags>("comparison of all coordinates between two geometries") {}
0074 
0075     bool fill() override {
0076       TGaxis::SetExponentOffset(-0.12, 0.01, "y");  // Y offset
0077 
0078       // trick to deal with the multi-ioved tag and two tag case at the same time
0079       auto theIOVs = PlotBase::getTag<0>().iovs;
0080       auto tagname1 = PlotBase::getTag<0>().name;
0081       std::string tagname2 = "";
0082       auto firstiov = theIOVs.front();
0083       std::tuple<cond::Time_t, cond::Hash> lastiov;
0084 
0085       // we don't support (yet) comparison with more than 2 tags
0086       assert(this->m_plotAnnotations.ntags < 3);
0087 
0088       if (this->m_plotAnnotations.ntags == 2) {
0089         auto tag2iovs = PlotBase::getTag<1>().iovs;
0090         tagname2 = PlotBase::getTag<1>().name;
0091         lastiov = tag2iovs.front();
0092       } else {
0093         lastiov = theIOVs.back();
0094       }
0095 
0096       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0097       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0098 
0099       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0100       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0101 
0102       std::vector<AlignTransform> ref_ali = first_payload->m_align;
0103       std::vector<AlignTransform> target_ali = last_payload->m_align;
0104 
0105       TCanvas canvas("Alignment Comparison", "Alignment Comparison", 2000, 1200);
0106       canvas.Divide(3, 2);
0107 
0108       if (ref_ali.size() != target_ali.size()) {
0109         edm::LogError("TrackerAlignment_PayloadInspector")
0110             << "the size of the reference alignment (" << ref_ali.size()
0111             << ") is different from the one of the target (" << target_ali.size()
0112             << ")! You are probably trying to compare different underlying geometries. Exiting";
0113         return false;
0114       }
0115 
0116       // check that the geomtery is a tracker one
0117       const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0118                                            ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0119                                            : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0120       TrackerTopology tTopo =
0121           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0122 
0123       for (const auto &ali : ref_ali) {
0124         auto mydetid = ali.rawId();
0125         if (DetId(mydetid).det() != DetId::Tracker) {
0126           edm::LogWarning("TrackerAlignment_PayloadInspector")
0127               << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0128               << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0129               << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0130           return false;
0131         }
0132       }
0133 
0134       const std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
0135                                                            AlignmentPI::t_y,
0136                                                            AlignmentPI::t_z,
0137                                                            AlignmentPI::rot_alpha,
0138                                                            AlignmentPI::rot_beta,
0139                                                            AlignmentPI::rot_gamma};
0140 
0141       std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
0142 
0143       // generate the map of histograms
0144       for (const auto &coord : coords) {
0145         auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0146         std::string unit =
0147             (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0148 
0149         diffs[coord] = std::make_unique<TH1F>(Form("comparison_%s", s_coord.c_str()),
0150                                               Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
0151                                               ref_ali.size(),
0152                                               -0.5,
0153                                               ref_ali.size() - 0.5);
0154       }
0155 
0156       // fill all the histograms together
0157       std::vector<int> boundaries;
0158       AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs);
0159 
0160       unsigned int subpad{1};
0161       TLegend legend = TLegend(0.17, 0.84, 0.95, 0.94);
0162       legend.SetTextSize(0.023);
0163       for (const auto &coord : coords) {
0164         canvas.cd(subpad);
0165         canvas.cd(subpad)->SetTopMargin(0.06);
0166         canvas.cd(subpad)->SetLeftMargin(0.17);
0167         canvas.cd(subpad)->SetRightMargin(0.05);
0168         canvas.cd(subpad)->SetBottomMargin(0.15);
0169         AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0170         auto max = diffs[coord]->GetMaximum();
0171         auto min = diffs[coord]->GetMinimum();
0172         auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
0173         if (range == 0.f)
0174           range = 0.1;
0175         //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
0176 
0177         diffs[coord]->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
0178         diffs[coord]->GetYaxis()->SetTitleOffset(1.5);
0179         diffs[coord]->SetMarkerStyle(20);
0180         diffs[coord]->SetMarkerSize(0.5);
0181         diffs[coord]->Draw("P");
0182 
0183         if (subpad == 1) { /* fill the legend only at the first pass */
0184           if (this->m_plotAnnotations.ntags == 2) {
0185             legend.SetHeader("#bf{Two Tags Comparison}", "C");  // option "C" allows to center the header
0186             legend.AddEntry(
0187                 diffs[coord].get(),
0188                 ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}")
0189                     .c_str(),
0190                 "PL");
0191           } else {
0192             legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");  // option "C" allows to center the header
0193             legend.AddEntry(diffs[coord].get(),
0194                             ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0195                             "PL");
0196           }
0197         }
0198         subpad++;
0199       }
0200 
0201       canvas.Update();
0202       canvas.cd();
0203       canvas.Modified();
0204 
0205       TLine l[6][boundaries.size()];
0206       TLatex tSubdet[6];
0207       for (unsigned int i = 0; i < 6; i++) {
0208         tSubdet[i].SetTextColor(kRed);
0209         tSubdet[i].SetNDC();
0210         tSubdet[i].SetTextAlign(21);
0211         tSubdet[i].SetTextSize(0.03);
0212         tSubdet[i].SetTextAngle(90);
0213       }
0214 
0215       subpad = 0;
0216       for (const auto &coord : coords) {
0217         auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0218         canvas.cd(subpad + 1);
0219 
0220         unsigned int i = 0;
0221         for (const auto &line : boundaries) {
0222           l[subpad][i] = TLine(diffs[coord]->GetBinLowEdge(line),
0223                                canvas.cd(subpad + 1)->GetUymin(),
0224                                diffs[coord]->GetBinLowEdge(line),
0225                                canvas.cd(subpad + 1)->GetUymax() * 0.84);
0226           l[subpad][i].SetLineWidth(1);
0227           l[subpad][i].SetLineStyle(9);
0228           l[subpad][i].SetLineColor(2);
0229           l[subpad][i].Draw("same");
0230           i++;
0231         }
0232 
0233         const auto &theX_ = {0.2, 0.24, 0.31, 0.4, 0.55, 0.8};
0234         for (unsigned int j = 1; j <= 7; j++) {
0235           auto thePart = static_cast<AlignmentPI::partitions>(j);
0236           tSubdet[subpad].DrawLatex(
0237               theX_.begin()[j - 1], 0.20, Form("%s", (AlignmentPI::getStringFromPart(thePart)).c_str()));
0238         }
0239 
0240         auto ltx = TLatex();
0241         ltx.SetTextFont(62);
0242         ltx.SetTextSize(0.042);
0243         ltx.SetTextAlign(11);
0244         ltx.DrawLatexNDC(canvas.cd(subpad + 1)->GetLeftMargin(),
0245                          1 - canvas.cd(subpad + 1)->GetTopMargin() + 0.01,
0246                          ("Tracker Alignment Compare : #color[4]{" + s_coord + "}").c_str());
0247         legend.Draw("same");
0248         subpad++;
0249       }  // loop on the coordinates
0250 
0251       std::string fileName(this->m_imageFileName);
0252       canvas.SaveAs(fileName.c_str());
0253       //canvas.SaveAs("out.root");
0254 
0255       return true;
0256     }
0257   };
0258 
0259   typedef TrackerAlignmentCompareAll<1, MULTI_IOV> TrackerAlignmentComparatorSingleTag;
0260   typedef TrackerAlignmentCompareAll<2, SINGLE_IOV> TrackerAlignmentComparatorTwoTags;
0261 
0262   //*******************************************//
0263   // Size of the movement over all partitions,
0264   // one coordinate (x,y,z,...) at a time
0265   //******************************************//
0266 
0267   template <AlignmentPI::coordinate coord, int ntags, IOVMultiplicity nIOVs>
0268   class TrackerAlignmentComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0269   public:
0270     TrackerAlignmentComparatorBase()
0271         : PlotImage<Alignments, nIOVs, ntags>("comparison of " + AlignmentPI::getStringFromCoordinate(coord) +
0272                                               " coordinate between two geometries") {}
0273 
0274     bool fill() override {
0275       TGaxis::SetExponentOffset(-0.12, 0.01, "y");  // Y offset
0276 
0277       // trick to deal with the multi-ioved tag and two tag case at the same time
0278       auto theIOVs = PlotBase::getTag<0>().iovs;
0279       auto tagname1 = PlotBase::getTag<0>().name;
0280       std::string tagname2 = "";
0281       auto firstiov = theIOVs.front();
0282       std::tuple<cond::Time_t, cond::Hash> lastiov;
0283 
0284       // we don't support (yet) comparison with more than 2 tags
0285       assert(this->m_plotAnnotations.ntags < 3);
0286 
0287       if (this->m_plotAnnotations.ntags == 2) {
0288         auto tag2iovs = PlotBase::getTag<1>().iovs;
0289         tagname2 = PlotBase::getTag<1>().name;
0290         lastiov = tag2iovs.front();
0291       } else {
0292         lastiov = theIOVs.back();
0293       }
0294 
0295       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0296       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0297 
0298       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0299       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0300 
0301       std::vector<AlignTransform> ref_ali = first_payload->m_align;
0302       std::vector<AlignTransform> target_ali = last_payload->m_align;
0303 
0304       TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1200, 1200);
0305 
0306       if (ref_ali.size() != target_ali.size()) {
0307         edm::LogError("TrackerAlignment_PayloadInspector")
0308             << "the size of the reference alignment (" << ref_ali.size()
0309             << ") is different from the one of the target (" << target_ali.size()
0310             << ")! You are probably trying to compare different underlying geometries. Exiting";
0311         return false;
0312       }
0313 
0314       // check that the geomtery is a tracker one
0315       const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0316                                            ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0317                                            : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0318       TrackerTopology tTopo =
0319           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0320 
0321       for (const auto &ali : ref_ali) {
0322         auto mydetid = ali.rawId();
0323         if (DetId(mydetid).det() != DetId::Tracker) {
0324           edm::LogWarning("TrackerAlignment_PayloadInspector")
0325               << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0326               << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0327               << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0328           return false;
0329         }
0330       }
0331 
0332       auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0333       std::string unit =
0334           (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0335 
0336       //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));
0337       std::unique_ptr<TH1F> compare =
0338           std::make_unique<TH1F>("comparison",
0339                                  Form(";Detector Id index; #Delta%s %s", s_coord.c_str(), unit.c_str()),
0340                                  ref_ali.size(),
0341                                  -0.5,
0342                                  ref_ali.size() - 0.5);
0343 
0344       // fill the histograms
0345       std::vector<int> boundaries;
0346       AlignmentPI::fillComparisonHistogram(coord, boundaries, ref_ali, target_ali, compare);
0347 
0348       canvas.cd();
0349 
0350       canvas.SetTopMargin(0.06);
0351       canvas.SetLeftMargin(0.17);
0352       canvas.SetRightMargin(0.05);
0353       canvas.SetBottomMargin(0.15);
0354       AlignmentPI::makeNicePlotStyle(compare.get(), kBlack);
0355       auto max = compare->GetMaximum();
0356       auto min = compare->GetMinimum();
0357       auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
0358       if (range == 0.f)
0359         range = 0.1;
0360       //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
0361 
0362       compare->GetYaxis()->SetRangeUser(-range * 1.5, range * 1.5);
0363       compare->GetYaxis()->SetTitleOffset(1.5);
0364       compare->SetMarkerStyle(20);
0365       compare->SetMarkerSize(0.5);
0366       compare->Draw("P");
0367 
0368       canvas.Update();
0369       canvas.cd();
0370 
0371       TLine l[boundaries.size()];
0372       unsigned int i = 0;
0373       for (const auto &line : boundaries) {
0374         l[i] = TLine(compare->GetBinLowEdge(line),
0375                      canvas.cd()->GetUymin(),
0376                      compare->GetBinLowEdge(line),
0377                      canvas.cd()->GetUymax());
0378         l[i].SetLineWidth(1);
0379         l[i].SetLineStyle(9);
0380         l[i].SetLineColor(2);
0381         l[i].Draw("same");
0382         i++;
0383       }
0384 
0385       TLatex tSubdet;
0386       tSubdet.SetNDC();
0387       tSubdet.SetTextAlign(21);
0388       tSubdet.SetTextSize(0.027);
0389       tSubdet.SetTextAngle(90);
0390       for (unsigned int j = 1; j <= 6; j++) {
0391         auto thePart = static_cast<AlignmentPI::partitions>(j);
0392         tSubdet.SetTextColor(kRed);
0393         auto myPair = (j > 1) ? AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(boundaries[j - 2]))
0394                               : AlignmentPI::calculatePosition(gPad, compare->GetBinLowEdge(0));
0395         float theX_ = myPair.first + 0.025;
0396         tSubdet.DrawLatex(theX_, 0.20, Form("%s", (AlignmentPI::getStringFromPart(thePart)).c_str()));
0397       }
0398 
0399       TLegend legend = TLegend(0.17, 0.86, 0.95, 0.94);
0400       if (this->m_plotAnnotations.ntags == 2) {
0401         legend.SetHeader("#bf{Two Tags Comparison}", "C");  // option "C" allows to center the header
0402         legend.AddEntry(
0403             compare.get(),
0404             ("#splitline{" + tagname1 + " : " + firstIOVsince + "}{" + tagname2 + " : " + lastIOVsince + "}").c_str(),
0405             "PL");
0406       } else {
0407         legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");  // option "C" allows to center the header
0408         legend.AddEntry(compare.get(),
0409                         ("#splitline{IOV since: " + firstIOVsince + "}{IOV since: " + lastIOVsince + "}").c_str(),
0410                         "PL");
0411       }
0412       legend.SetTextSize(0.020);
0413       legend.Draw("same");
0414 
0415       auto ltx = TLatex();
0416       ltx.SetTextFont(62);
0417       ltx.SetTextSize(0.042);
0418       ltx.SetTextAlign(11);
0419       ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0420                        1 - gPad->GetTopMargin() + 0.01,
0421                        ("Tracker Alignment Compare :#color[4]{" + s_coord + "}").c_str());
0422 
0423       std::string fileName(this->m_imageFileName);
0424       canvas.SaveAs(fileName.c_str());
0425 
0426       return true;
0427     }
0428   };
0429 
0430   template <AlignmentPI::coordinate coord>
0431   using TrackerAlignmentCompare = TrackerAlignmentComparatorBase<coord, 1, MULTI_IOV>;
0432 
0433   template <AlignmentPI::coordinate coord>
0434   using TrackerAlignmentCompareTwoTags = TrackerAlignmentComparatorBase<coord, 2, SINGLE_IOV>;
0435 
0436   typedef TrackerAlignmentCompare<AlignmentPI::t_x> TrackerAlignmentCompareX;
0437   typedef TrackerAlignmentCompare<AlignmentPI::t_y> TrackerAlignmentCompareY;
0438   typedef TrackerAlignmentCompare<AlignmentPI::t_z> TrackerAlignmentCompareZ;
0439 
0440   typedef TrackerAlignmentCompare<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlpha;
0441   typedef TrackerAlignmentCompare<AlignmentPI::rot_beta> TrackerAlignmentCompareBeta;
0442   typedef TrackerAlignmentCompare<AlignmentPI::rot_gamma> TrackerAlignmentCompareGamma;
0443 
0444   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_x> TrackerAlignmentCompareXTwoTags;
0445   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_y> TrackerAlignmentCompareYTwoTags;
0446   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::t_z> TrackerAlignmentCompareZTwoTags;
0447 
0448   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlphaTwoTags;
0449   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_beta> TrackerAlignmentCompareBetaTwoTags;
0450   typedef TrackerAlignmentCompareTwoTags<AlignmentPI::rot_gamma> TrackerAlignmentCompareGammaTwoTags;
0451 
0452   //*******************************************//
0453   // Summary canvas per subdetector
0454   //******************************************//
0455 
0456   template <int ntags, IOVMultiplicity nIOVs, AlignmentPI::partitions q>
0457   class TrackerAlignmentSummaryBase : public PlotImage<Alignments, nIOVs, ntags> {
0458   public:
0459     TrackerAlignmentSummaryBase()
0460         : PlotImage<Alignments, nIOVs, ntags>("Comparison of all coordinates between two geometries for " +
0461                                               getStringFromPart(q)) {}
0462 
0463     bool fill() override {
0464       // trick to deal with the multi-ioved tag and two tag case at the same time
0465       auto theIOVs = PlotBase::getTag<0>().iovs;
0466       auto tagname1 = PlotBase::getTag<0>().name;
0467       std::string tagname2 = "";
0468       auto firstiov = theIOVs.front();
0469       std::tuple<cond::Time_t, cond::Hash> lastiov;
0470 
0471       // we don't support (yet) comparison with more than 2 tags
0472       assert(this->m_plotAnnotations.ntags < 3);
0473 
0474       if (this->m_plotAnnotations.ntags == 2) {
0475         auto tag2iovs = PlotBase::getTag<1>().iovs;
0476         tagname2 = PlotBase::getTag<1>().name;
0477         lastiov = tag2iovs.front();
0478       } else {
0479         lastiov = theIOVs.back();
0480       }
0481 
0482       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0483       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0484 
0485       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0486       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0487 
0488       std::vector<AlignTransform> ref_ali = first_payload->m_align;
0489       std::vector<AlignTransform> target_ali = last_payload->m_align;
0490 
0491       if (ref_ali.size() != target_ali.size()) {
0492         edm::LogError("TrackerAlignment_PayloadInspector")
0493             << "the size of the reference alignment (" << ref_ali.size()
0494             << ") is different from the one of the target (" << target_ali.size()
0495             << ")! You are probably trying to compare different underlying geometries. Exiting";
0496         return false;
0497       }
0498 
0499       // check that the geomtery is a tracker one
0500       const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0501                                            ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0502                                            : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0503       TrackerTopology tTopo =
0504           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0505 
0506       for (const auto &ali : ref_ali) {
0507         auto mydetid = ali.rawId();
0508         if (DetId(mydetid).det() != DetId::Tracker) {
0509           edm::LogWarning("TrackerAlignment_PayloadInspector")
0510               << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0511               << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0512               << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0513           return false;
0514         }
0515       }
0516 
0517       TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1800, 1200);
0518       canvas.Divide(3, 2);
0519 
0520       std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
0521       std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
0522                                                      AlignmentPI::t_y,
0523                                                      AlignmentPI::t_z,
0524                                                      AlignmentPI::rot_alpha,
0525                                                      AlignmentPI::rot_beta,
0526                                                      AlignmentPI::rot_gamma};
0527 
0528       for (const auto &coord : coords) {
0529         auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0530         std::string unit =
0531             (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0532 
0533         diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s", s_coord.c_str()),
0534                                               Form(";#Delta%s %s;n. of modules", s_coord.c_str(), unit.c_str()),
0535                                               1000,
0536                                               -500.,
0537                                               500.);
0538       }
0539 
0540       // fill the comparison histograms
0541       std::vector<int> boundaries{};
0542       AlignmentPI::fillComparisonHistograms(boundaries, ref_ali, target_ali, diffs, true, q);
0543 
0544       int c_index = 1;
0545 
0546       auto legend = std::make_unique<TLegend>(0.14, 0.93, 0.55, 0.98);
0547       legend->AddEntry(
0548           diffs[AlignmentPI::t_x].get(),
0549           ("#DeltaIOV: " + std::to_string(std::get<0>(lastiov)) + "-" + std::to_string(std::get<0>(firstiov))).c_str(),
0550           "L");
0551       legend->SetTextSize(0.03);
0552 
0553       for (const auto &coord : coords) {
0554         canvas.cd(c_index)->SetLogy();
0555         canvas.cd(c_index)->SetTopMargin(0.02);
0556         canvas.cd(c_index)->SetBottomMargin(0.15);
0557         canvas.cd(c_index)->SetLeftMargin(0.14);
0558         canvas.cd(c_index)->SetRightMargin(0.04);
0559         diffs[coord]->SetLineWidth(2);
0560         AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0561 
0562         //float x_max = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindLastBinAbove(0.));
0563         //float x_min = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindFirstBinAbove(0.));
0564         //float extremum = std::abs(x_max) > std::abs(x_min) ? std::abs(x_max) : std::abs(x_min);
0565         //diffs[coord]->GetXaxis()->SetRangeUser(-extremum*2,extremum*2);
0566 
0567         int i_max = diffs[coord]->FindLastBinAbove(0.);
0568         int i_min = diffs[coord]->FindFirstBinAbove(0.);
0569         diffs[coord]->GetXaxis()->SetRange(std::max(1, i_min - 10), std::min(i_max + 10, diffs[coord]->GetNbinsX()));
0570         diffs[coord]->Draw("HIST");
0571         AlignmentPI::makeNiceStats(diffs[coord].get(), q, kBlack);
0572 
0573         legend->Draw("same");
0574 
0575         c_index++;
0576       }
0577 
0578       std::string fileName(this->m_imageFileName);
0579       canvas.SaveAs(fileName.c_str());
0580 
0581       return true;
0582     }
0583   };
0584 
0585   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
0586   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
0587   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
0588 
0589   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTID;
0590   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
0591   typedef TrackerAlignmentSummaryBase<1, MULTI_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
0592 
0593   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::BPix> TrackerAlignmentSummaryBPixTwoTags;
0594   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::FPix> TrackerAlignmentSummaryFPixTwoTags;
0595   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TIB> TrackerAlignmentSummaryTIBTwoTags;
0596 
0597   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TID> TrackerAlignmentSummaryTIDTwoTags;
0598   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TOB> TrackerAlignmentSummaryTOBTwoTags;
0599   typedef TrackerAlignmentSummaryBase<2, SINGLE_IOV, AlignmentPI::TEC> TrackerAlignmentSummaryTECTwoTags;
0600 
0601   //*******************************************//
0602   // History of the position of the BPix Barycenter
0603   //******************************************//
0604 
0605   template <AlignmentPI::coordinate coord>
0606   class BPixBarycenterHistory : public HistoryPlot<Alignments, float> {
0607   public:
0608     BPixBarycenterHistory()
0609         : HistoryPlot<Alignments, float>(
0610               " Barrel Pixel " + AlignmentPI::getStringFromCoordinate(coord) + " positions vs time",
0611               AlignmentPI::getStringFromCoordinate(coord) + " position [cm]") {}
0612     ~BPixBarycenterHistory() override = default;
0613 
0614     float getFromPayload(Alignments &payload) override {
0615       std::vector<AlignTransform> alignments = payload.m_align;
0616 
0617       float barycenter = 0.;
0618       float nmodules(0.);
0619       for (const auto &ali : alignments) {
0620         if (DetId(ali.rawId()).det() != DetId::Tracker) {
0621           edm::LogWarning("TrackerAlignment_PayloadInspector")
0622               << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
0623               << " is different from " << DetId::Tracker << "  - terminating ";
0624           return false;
0625         }
0626 
0627         int subid = DetId(ali.rawId()).subdetId();
0628         if (subid != PixelSubdetector::PixelBarrel)
0629           continue;
0630 
0631         nmodules++;
0632         switch (coord) {
0633           case AlignmentPI::t_x:
0634             barycenter += (ali.translation().x());
0635             break;
0636           case AlignmentPI::t_y:
0637             barycenter += (ali.translation().y());
0638             break;
0639           case AlignmentPI::t_z:
0640             barycenter += (ali.translation().z());
0641             break;
0642           default:
0643             edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
0644             break;
0645         }  // switch on the coordinate (only X,Y,Z are interesting)
0646       }    // ends loop on the alignments
0647 
0648       edm::LogInfo("TrackerAlignment_PayloadInspector") << "barycenter (" << barycenter << ")/n. modules (" << nmodules
0649                                                         << ") =  " << barycenter / nmodules << std::endl;
0650 
0651       // take the mean
0652       barycenter /= nmodules;
0653 
0654       // applied GPR correction to move barycenter to global CMS coordinates
0655       barycenter += hardcodeGPR.at(coord);
0656 
0657       return barycenter;
0658 
0659     }  // payload
0660   };
0661 
0662   typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
0663   typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
0664   typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
0665 
0666   /************************************************
0667     Display of Tracker Detector barycenters
0668   *************************************************/
0669   class TrackerAlignmentBarycenters : public PlotImage<Alignments, SINGLE_IOV> {
0670   public:
0671     TrackerAlignmentBarycenters() : PlotImage<Alignments, SINGLE_IOV>("Display of Tracker Alignment Barycenters") {}
0672 
0673     bool fill() override {
0674       auto tag = PlotBase::getTag<0>();
0675       auto iov = tag.iovs.front();
0676       std::shared_ptr<Alignments> payload = fetchPayload(std::get<1>(iov));
0677       unsigned int run = std::get<0>(iov);
0678 
0679       TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1600, 1000);
0680       canvas.cd();
0681 
0682       canvas.SetTopMargin(0.07);
0683       canvas.SetBottomMargin(0.06);
0684       canvas.SetLeftMargin(0.15);
0685       canvas.SetRightMargin(0.03);
0686       canvas.Modified();
0687       canvas.SetGrid();
0688 
0689       auto h2_BarycenterParameters =
0690           std::make_unique<TH2F>("Parameters", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
0691 
0692       auto h2_uncBarycenterParameters =
0693           std::make_unique<TH2F>("Parameters2", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
0694 
0695       h2_BarycenterParameters->SetStats(false);
0696       h2_BarycenterParameters->SetTitle(nullptr);
0697       h2_uncBarycenterParameters->SetStats(false);
0698       h2_uncBarycenterParameters->SetTitle(nullptr);
0699 
0700       std::vector<AlignTransform> alignments = payload->m_align;
0701 
0702       isPhase0 = (alignments.size() == AlignmentPI::phase0size) ? true : false;
0703 
0704       // check that the geomtery is a tracker one
0705       const char *path_toTopologyXML = isPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0706                                                 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0707 
0708       TrackerTopology tTopo =
0709           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0710 
0711       AlignmentPI::TkAlBarycenters barycenters;
0712       // compute uncorrected barycenter
0713       barycenters.computeBarycenters(
0714           alignments, tTopo, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
0715 
0716       auto Xbarycenters = barycenters.getX();
0717       auto Ybarycenters = barycenters.getY();
0718       auto Zbarycenters = barycenters.getZ();
0719 
0720       // compute barycenter corrected for the GPR
0721       barycenters.computeBarycenters(alignments, tTopo, hardcodeGPR);
0722 
0723       auto c_Xbarycenters = barycenters.getX();
0724       auto c_Ybarycenters = barycenters.getY();
0725       auto c_Zbarycenters = barycenters.getZ();
0726 
0727       h2_BarycenterParameters->GetXaxis()->SetBinLabel(1, "X [cm]");
0728       h2_BarycenterParameters->GetXaxis()->SetBinLabel(2, "Y [cm]");
0729       h2_BarycenterParameters->GetXaxis()->SetBinLabel(3, "Z [cm]");
0730       h2_BarycenterParameters->GetXaxis()->SetBinLabel(4, "X_{no GPR} [cm]");
0731       h2_BarycenterParameters->GetXaxis()->SetBinLabel(5, "Y_{no GPR} [cm]");
0732       h2_BarycenterParameters->GetXaxis()->SetBinLabel(6, "Z_{no GPR} [cm]");
0733 
0734       bool isLikelyMC(false);
0735       int checkX =
0736           std::count_if(Xbarycenters.begin(), Xbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0737       int checkY =
0738           std::count_if(Ybarycenters.begin(), Ybarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0739       int checkZ =
0740           std::count_if(Zbarycenters.begin(), Zbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0741 
0742       // if all the coordinate barycenters for both BPix and FPix are below 10um
0743       // this is very likely a MC payload
0744       if ((checkX + checkY + checkZ) == 0 && run == 1)
0745         isLikelyMC = true;
0746 
0747       unsigned int yBin = 6;
0748       for (unsigned int i = 0; i < 6; i++) {
0749         auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
0750         std::string theLabel = getStringFromPart(thePart);
0751         h2_BarycenterParameters->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0752         if (!isLikelyMC) {
0753           h2_BarycenterParameters->SetBinContent(1, yBin, c_Xbarycenters[i]);
0754           h2_BarycenterParameters->SetBinContent(2, yBin, c_Ybarycenters[i]);
0755           h2_BarycenterParameters->SetBinContent(3, yBin, c_Zbarycenters[i]);
0756         }
0757 
0758         h2_uncBarycenterParameters->SetBinContent(4, yBin, Xbarycenters[i]);
0759         h2_uncBarycenterParameters->SetBinContent(5, yBin, Ybarycenters[i]);
0760         h2_uncBarycenterParameters->SetBinContent(6, yBin, Zbarycenters[i]);
0761         yBin--;
0762       }
0763 
0764       h2_BarycenterParameters->GetXaxis()->LabelsOption("h");
0765       h2_BarycenterParameters->GetYaxis()->SetLabelSize(0.05);
0766       h2_BarycenterParameters->GetXaxis()->SetLabelSize(0.05);
0767       h2_BarycenterParameters->SetMarkerSize(1.5);
0768       h2_BarycenterParameters->Draw("TEXT");
0769 
0770       h2_uncBarycenterParameters->SetMarkerSize(1.5);
0771       h2_uncBarycenterParameters->SetMarkerColor(kRed);
0772       h2_uncBarycenterParameters->Draw("TEXTsame");
0773 
0774       TLatex t1;
0775       t1.SetNDC();
0776       t1.SetTextAlign(26);
0777       t1.SetTextSize(0.05);
0778       t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters, IOV %i", run));
0779       t1.SetTextSize(0.025);
0780 
0781       std::string fileName(m_imageFileName);
0782       canvas.SaveAs(fileName.c_str());
0783 
0784       return true;
0785     }
0786 
0787   private:
0788     bool isPhase0;
0789   };
0790 
0791   /************************************************
0792     Comparator of Tracker Detector barycenters
0793   *************************************************/
0794   template <int ntags, IOVMultiplicity nIOVs>
0795   class TrackerAlignmentBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0796   public:
0797     TrackerAlignmentBarycentersComparatorBase()
0798         : PlotImage<Alignments, nIOVs, ntags>("Comparison of Tracker Alignment Barycenters") {}
0799 
0800     bool fill() override {
0801       // trick to deal with the multi-ioved tag and two tag case at the same time
0802       auto theIOVs = PlotBase::getTag<0>().iovs;
0803       auto tagname1 = PlotBase::getTag<0>().name;
0804       std::string tagname2 = "";
0805       auto firstiov = theIOVs.front();
0806       std::tuple<cond::Time_t, cond::Hash> lastiov;
0807 
0808       // we don't support (yet) comparison with more than 2 tags
0809       assert(this->m_plotAnnotations.ntags < 3);
0810 
0811       if (this->m_plotAnnotations.ntags == 2) {
0812         auto tag2iovs = PlotBase::getTag<1>().iovs;
0813         tagname2 = PlotBase::getTag<1>().name;
0814         lastiov = tag2iovs.front();
0815       } else {
0816         lastiov = theIOVs.back();
0817       }
0818 
0819       unsigned int first_run = std::get<0>(firstiov);
0820       unsigned int last_run = std::get<0>(lastiov);
0821 
0822       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0823       std::vector<AlignTransform> last_alignments = last_payload->m_align;
0824 
0825       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0826       std::vector<AlignTransform> first_alignments = first_payload->m_align;
0827 
0828       isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
0829       isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
0830 
0831       // check that the geomtery is a tracker one
0832       const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0833                                                        : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0834 
0835       TrackerTopology tTopo_f =
0836           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0837 
0838       path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0839                                          : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0840 
0841       TrackerTopology tTopo_l =
0842           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0843 
0844       TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1200, 800);
0845       canvas.cd();
0846 
0847       canvas.SetTopMargin(0.07);
0848       canvas.SetBottomMargin(0.06);
0849       canvas.SetLeftMargin(0.15);
0850       canvas.SetRightMargin(0.03);
0851       canvas.Modified();
0852       canvas.SetGrid();
0853 
0854       auto h2_BarycenterDiff =
0855           std::make_unique<TH2F>("Parameters diff", "SubDetector Barycenter Difference", 3, 0.0, 3.0, 6, 0, 6.);
0856 
0857       h2_BarycenterDiff->SetStats(false);
0858       h2_BarycenterDiff->SetTitle(nullptr);
0859       h2_BarycenterDiff->GetXaxis()->SetBinLabel(1, "X [#mum]");
0860       h2_BarycenterDiff->GetXaxis()->SetBinLabel(2, "Y [#mum]");
0861       h2_BarycenterDiff->GetXaxis()->SetBinLabel(3, "Z [#mum]");
0862 
0863       AlignmentPI::TkAlBarycenters l_barycenters;
0864       l_barycenters.computeBarycenters(last_alignments, tTopo_l, hardcodeGPR);
0865 
0866       AlignmentPI::TkAlBarycenters f_barycenters;
0867       f_barycenters.computeBarycenters(first_alignments, tTopo_f, hardcodeGPR);
0868 
0869       unsigned int yBin = 6;
0870       for (unsigned int i = 0; i < 6; i++) {
0871         auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
0872         std::string theLabel = getStringFromPart(thePart);
0873         h2_BarycenterDiff->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0874         h2_BarycenterDiff->SetBinContent(
0875             1, yBin, (l_barycenters.getX()[i] - f_barycenters.getX()[i]) * AlignmentPI::cmToUm);
0876         h2_BarycenterDiff->SetBinContent(
0877             2, yBin, (l_barycenters.getY()[i] - f_barycenters.getY()[i]) * AlignmentPI::cmToUm);
0878         h2_BarycenterDiff->SetBinContent(
0879             3, yBin, (l_barycenters.getZ()[i] - f_barycenters.getZ()[i]) * AlignmentPI::cmToUm);
0880         yBin--;
0881       }
0882 
0883       h2_BarycenterDiff->GetXaxis()->LabelsOption("h");
0884       h2_BarycenterDiff->GetYaxis()->SetLabelSize(0.05);
0885       h2_BarycenterDiff->GetXaxis()->SetLabelSize(0.05);
0886       h2_BarycenterDiff->SetMarkerSize(1.5);
0887       h2_BarycenterDiff->SetMarkerColor(kRed);
0888       h2_BarycenterDiff->Draw("TEXT");
0889 
0890       TLatex t1;
0891       t1.SetNDC();
0892       t1.SetTextAlign(26);
0893       t1.SetTextSize(0.05);
0894       t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters Diff, IOV %i - IOV %i", last_run, first_run));
0895       t1.SetTextSize(0.025);
0896 
0897       std::string fileName(this->m_imageFileName);
0898       canvas.SaveAs(fileName.c_str());
0899 
0900       return true;
0901     }
0902 
0903   private:
0904     bool isInitialPhase0;
0905     bool isFinalPhase0;
0906   };
0907 
0908   using TrackerAlignmentBarycentersCompare = TrackerAlignmentBarycentersComparatorBase<1, MULTI_IOV>;
0909   using TrackerAlignmentBarycentersCompareTwoTags = TrackerAlignmentBarycentersComparatorBase<2, SINGLE_IOV>;
0910 
0911   /************************************************
0912     Comparator of Pixel Tracker Detector barycenters
0913   *************************************************/
0914   template <int ntags, IOVMultiplicity nIOVs>
0915   class PixelBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0916   public:
0917     PixelBarycentersComparatorBase() : PlotImage<Alignments, nIOVs, ntags>("Comparison of Pixel Barycenters") {}
0918 
0919     bool fill() override {
0920       // trick to deal with the multi-ioved tag and two tag case at the same time
0921       auto theIOVs = PlotBase::getTag<0>().iovs;
0922       auto tagname1 = PlotBase::getTag<0>().name;
0923       std::string tagname2 = "";
0924       auto firstiov = theIOVs.front();
0925       std::tuple<cond::Time_t, cond::Hash> lastiov;
0926 
0927       // we don't support (yet) comparison with more than 2 tags
0928       assert(this->m_plotAnnotations.ntags < 3);
0929 
0930       if (this->m_plotAnnotations.ntags == 2) {
0931         auto tag2iovs = PlotBase::getTag<1>().iovs;
0932         tagname2 = PlotBase::getTag<1>().name;
0933         lastiov = tag2iovs.front();
0934       } else {
0935         lastiov = theIOVs.back();
0936       }
0937 
0938       unsigned int first_run = std::get<0>(firstiov);
0939       unsigned int last_run = std::get<0>(lastiov);
0940 
0941       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0942       std::vector<AlignTransform> last_alignments = last_payload->m_align;
0943 
0944       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0945       std::vector<AlignTransform> first_alignments = first_payload->m_align;
0946 
0947       TCanvas canvas("Pixel Barycenter Summary", "Pixel Barycenter summary", 1200, 1200);
0948       canvas.Divide(2, 2);
0949       canvas.cd();
0950 
0951       TLatex t1;
0952       t1.SetNDC();
0953       t1.SetTextAlign(26);
0954       t1.SetTextSize(0.03);
0955       t1.DrawLatex(0.5,
0956                    0.97,
0957                    ("Pixel Barycenters comparison, IOV: #color[2]{" + std::to_string(first_run) +
0958                     "} vs IOV: #color[4]{" + std::to_string(last_run) + "}")
0959                        .c_str());
0960       t1.SetTextSize(0.025);
0961 
0962       for (unsigned int c = 1; c <= 4; c++) {
0963         canvas.cd(c)->SetTopMargin(0.07);
0964         canvas.cd(c)->SetBottomMargin(0.12);
0965         canvas.cd(c)->SetLeftMargin(0.15);
0966         canvas.cd(c)->SetRightMargin(0.03);
0967         canvas.cd(c)->Modified();
0968         canvas.cd(c)->SetGrid();
0969       }
0970 
0971       std::array<std::string, 3> structures = {{"FPIX-", "BPIX", "FPIX+"}};
0972       std::array<std::unique_ptr<TH2F>, 3> histos;
0973 
0974       isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
0975       isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
0976 
0977       // check that the geomtery is a tracker one
0978       const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0979                                                        : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0980 
0981       TrackerTopology tTopo_f =
0982           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0983 
0984       AlignmentPI::TkAlBarycenters myInitialBarycenters;
0985       //myInitialBarycenters.computeBarycenters(first_alignments,tTopo_f,hardcodeGPR);
0986       myInitialBarycenters.computeBarycenters(
0987           first_alignments, tTopo_f, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
0988 
0989       path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0990                                          : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0991 
0992       TrackerTopology tTopo_l =
0993           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0994 
0995       AlignmentPI::TkAlBarycenters myFinalBarycenters;
0996       //myFinalBarycenters.computeBarycenters(last_alignments,tTopo_l,hardcodeGPR);
0997       myFinalBarycenters.computeBarycenters(
0998           last_alignments, tTopo_l, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
0999 
1000       if (isFinalPhase0 != isInitialPhase0) {
1001         edm::LogWarning("TrackerAlignment_PayloadInspector")
1002             << "the size of the reference alignment (" << first_alignments.size()
1003             << ") is different from the one of the target (" << last_alignments.size()
1004             << ")! You are probably trying to compare different underlying geometries.";
1005       }
1006 
1007       unsigned int index(0);
1008       for (const auto &piece : structures) {
1009         const char *name = piece.c_str();
1010         histos[index] = std::make_unique<TH2F>(
1011             name,
1012             Form("%s x-y Barycenter Difference;x_{%s}-x_{TOB} [mm];y_{%s}-y_{TOB} [mm]", name, name, name),
1013             100,
1014             -3.,
1015             3.,
1016             100,
1017             -3.,
1018             3.);
1019 
1020         histos[index]->SetStats(false);
1021         histos[index]->SetTitle(nullptr);
1022         histos[index]->GetYaxis()->SetLabelSize(0.05);
1023         histos[index]->GetXaxis()->SetLabelSize(0.05);
1024         histos[index]->GetYaxis()->SetTitleSize(0.06);
1025         histos[index]->GetXaxis()->SetTitleSize(0.06);
1026         histos[index]->GetYaxis()->CenterTitle();
1027         histos[index]->GetXaxis()->CenterTitle();
1028         histos[index]->GetXaxis()->SetTitleOffset(0.9);
1029         index++;
1030       }
1031 
1032       auto h2_ZBarycenterDiff = std::make_unique<TH2F>(
1033           "Pixel_z_diff", "Pixel z-Barycenter Difference;; z_{Pixel-Ideal} -z_{TOB} [mm]", 3, -0.5, 2.5, 100, -10., 10.);
1034       h2_ZBarycenterDiff->SetStats(false);
1035       h2_ZBarycenterDiff->SetTitle(nullptr);
1036       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(1, "FPIX -");
1037       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(2, "BPIX");
1038       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(3, "FPIX +");
1039       h2_ZBarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1040       h2_ZBarycenterDiff->GetXaxis()->SetLabelSize(0.07);
1041       h2_ZBarycenterDiff->GetYaxis()->SetTitleSize(0.06);
1042       h2_ZBarycenterDiff->GetXaxis()->SetTitleSize(0.06);
1043       h2_ZBarycenterDiff->GetYaxis()->CenterTitle();
1044       h2_ZBarycenterDiff->GetXaxis()->CenterTitle();
1045       h2_ZBarycenterDiff->GetYaxis()->SetTitleOffset(1.1);
1046 
1047       std::function<GlobalPoint(int)> cutFunctorInitial = [&myInitialBarycenters](int index) {
1048         switch (index) {
1049           case 1:
1050             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1051           case 2:
1052             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1053           case 3:
1054             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1055           default:
1056             return GlobalPoint(0, 0, 0);
1057         }
1058       };
1059 
1060       std::function<GlobalPoint(int)> cutFunctorFinal = [&myFinalBarycenters](int index) {
1061         switch (index) {
1062           case 1:
1063             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1064           case 2:
1065             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1066           case 3:
1067             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1068           default:
1069             return GlobalPoint(0, 0, 0);
1070         }
1071       };
1072 
1073       float x0i, x0f, y0i, y0f;
1074 
1075       t1.SetNDC(kFALSE);
1076       t1.SetTextSize(0.047);
1077       for (unsigned int c = 1; c <= 3; c++) {
1078         x0i = cutFunctorInitial(c).x() * 10;  // transform cm to mm (x10)
1079         x0f = cutFunctorFinal(c).x() * 10;
1080         y0i = cutFunctorInitial(c).y() * 10;
1081         y0f = cutFunctorFinal(c).y() * 10;
1082 
1083         canvas.cd(c);
1084         histos[c - 1]->Draw();
1085 
1086         COUT << "initial x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0i << "," << y0i << ") mm"
1087              << std::endl;
1088         COUT << "final   x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0f << "," << y0f << ") mm"
1089              << std::endl;
1090 
1091         TMarker *initial = new TMarker(x0i, y0i, 21);
1092         TMarker *final = new TMarker(x0f, y0f, 20);
1093 
1094         initial->SetMarkerColor(kRed);
1095         final->SetMarkerColor(kBlue);
1096         initial->SetMarkerSize(2.5);
1097         final->SetMarkerSize(2.5);
1098         t1.SetTextColor(kRed);
1099         initial->Draw();
1100         t1.DrawLatex(x0i, y0i + (y0i > y0f ? 0.3 : -0.5), Form("(%.2f,%.2f)", x0i, y0i));
1101         final->Draw("same");
1102         t1.SetTextColor(kBlue);
1103         t1.DrawLatex(x0f, y0f + (y0i > y0f ? -0.5 : 0.3), Form("(%.2f,%.2f)", x0f, y0f));
1104       }
1105 
1106       // fourth pad is a special case for the z coordinate
1107       canvas.cd(4);
1108       h2_ZBarycenterDiff->Draw();
1109       float z0i, z0f;
1110 
1111       // numbers do agree with https://twiki.cern.ch/twiki/bin/view/CMSPublic/TkAlignmentPerformancePhaseIStartUp17#Pixel_Barycentre_Positions
1112 
1113       std::array<double, 3> hardcodeIdealZPhase0 = {{-41.94909, 0., 41.94909}};  // units are cm
1114       std::array<double, 3> hardcodeIdealZPhase1 = {{-39.82911, 0., 39.82911}};  // units are cm
1115 
1116       for (unsigned int c = 1; c <= 3; c++) {
1117         // less than pretty but needed to remove the z position of the FPix barycenters != 0
1118 
1119         z0i =
1120             (cutFunctorInitial(c).z() - (isInitialPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) *
1121             10;  // convert to mm
1122         z0f =
1123             (cutFunctorFinal(c).z() - (isFinalPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) * 10;
1124 
1125         TMarker *initial = new TMarker(c - 1, z0i, 21);
1126         TMarker *final = new TMarker(c - 1, z0f, 20);
1127 
1128         COUT << "initial   z " << std::left << std::setw(7) << structures[c - 1] << " " << z0i << " mm" << std::endl;
1129         COUT << "final     z " << std::left << std::setw(7) << structures[c - 1] << " " << z0f << " mm" << std::endl;
1130 
1131         initial->SetMarkerColor(kRed);
1132         final->SetMarkerColor(kBlue);
1133         initial->SetMarkerSize(2.5);
1134         final->SetMarkerSize(2.5);
1135         initial->Draw();
1136         t1.SetTextColor(kRed);
1137         t1.DrawLatex(c - 1, z0i + (z0i > z0f ? 1. : -1.5), Form("(%.2f)", z0i));
1138         final->Draw("same");
1139         t1.SetTextColor(kBlue);
1140         t1.DrawLatex(c - 1, z0f + (z0i > z0f ? -1.5 : 1), Form("(%.2f)", z0f));
1141       }
1142 
1143       std::string fileName(this->m_imageFileName);
1144       canvas.SaveAs(fileName.c_str());
1145 
1146       return true;
1147     }
1148 
1149   private:
1150     bool isInitialPhase0;
1151     bool isFinalPhase0;
1152   };
1153 
1154   using PixelBarycentersCompare = PixelBarycentersComparatorBase<1, MULTI_IOV>;
1155   using PixelBarycentersCompareTwoTags = PixelBarycentersComparatorBase<2, SINGLE_IOV>;
1156 
1157 }  // namespace
1158 
1159 PAYLOAD_INSPECTOR_MODULE(TrackerAlignment) {
1160   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorSingleTag);
1161   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorTwoTags);
1162   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
1163   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
1164   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
1165   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
1166   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
1167   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
1168   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareXTwoTags);
1169   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareYTwoTags);
1170   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZTwoTags);
1171   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlphaTwoTags);
1172   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBetaTwoTags);
1173   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGammaTwoTags);
1174   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
1175   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
1176   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
1177   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
1178   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
1179   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
1180   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPixTwoTags);
1181   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPixTwoTags);
1182   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIBTwoTags);
1183   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIDTwoTags);
1184   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOBTwoTags);
1185   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTECTwoTags);
1186   PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
1187   PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
1188   PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
1189   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycenters);
1190   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompare);
1191   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompareTwoTags);
1192   PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompare);
1193   PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompareTwoTags);
1194 }