Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-09-18 22:46:39

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 <AlignmentPI::partitions q>
0457   class TrackerAlignmentSummary : public PlotImage<Alignments, MULTI_IOV> {
0458   public:
0459     TrackerAlignmentSummary()
0460         : PlotImage<Alignments, MULTI_IOV>("Comparison of all coordinates between two geometries for " +
0461                                            getStringFromPart(q)) {}
0462 
0463     bool fill() override {
0464       auto tag = PlotBase::getTag<0>();
0465       auto sorted_iovs = tag.iovs;
0466 
0467       // make absolute sure the IOVs are sortd by since
0468       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const &t1, auto const &t2) {
0469         return std::get<0>(t1) < std::get<0>(t2);
0470       });
0471 
0472       auto firstiov = sorted_iovs.front();
0473       auto lastiov = sorted_iovs.back();
0474 
0475       std::shared_ptr<Alignments> last_payload = fetchPayload(std::get<1>(lastiov));
0476       std::shared_ptr<Alignments> first_payload = fetchPayload(std::get<1>(firstiov));
0477 
0478       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0479       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0480 
0481       std::vector<AlignTransform> ref_ali = first_payload->m_align;
0482       std::vector<AlignTransform> target_ali = last_payload->m_align;
0483 
0484       if (ref_ali.size() != target_ali.size()) {
0485         edm::LogError("TrackerAlignment_PayloadInspector")
0486             << "the size of the reference alignment (" << ref_ali.size()
0487             << ") is different from the one of the target (" << target_ali.size()
0488             << ")! You are probably trying to compare different underlying geometries. Exiting";
0489         return false;
0490       }
0491 
0492       // check that the geomtery is a tracker one
0493       const char *path_toTopologyXML = (ref_ali.size() == AlignmentPI::phase0size)
0494                                            ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0495                                            : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0496       TrackerTopology tTopo =
0497           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0498 
0499       for (const auto &ali : ref_ali) {
0500         auto mydetid = ali.rawId();
0501         if (DetId(mydetid).det() != DetId::Tracker) {
0502           edm::LogWarning("TrackerAlignment_PayloadInspector")
0503               << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() << " (" << DetId(mydetid).det()
0504               << ") is different from " << DetId::Tracker << " (is DoubleSide: " << tTopo.tidIsDoubleSide(mydetid)
0505               << "); subdetId " << DetId(mydetid).subdetId() << " - terminating ";
0506           return false;
0507         }
0508       }
0509 
0510       TCanvas canvas("Alignment Comparison", "Alignment Comparison", 1800, 1200);
0511       canvas.Divide(3, 2);
0512 
0513       std::unordered_map<AlignmentPI::coordinate, std::unique_ptr<TH1F> > diffs;
0514       std::vector<AlignmentPI::coordinate> coords = {AlignmentPI::t_x,
0515                                                      AlignmentPI::t_y,
0516                                                      AlignmentPI::t_z,
0517                                                      AlignmentPI::rot_alpha,
0518                                                      AlignmentPI::rot_beta,
0519                                                      AlignmentPI::rot_gamma};
0520 
0521       for (const auto &coord : coords) {
0522         auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
0523         std::string unit =
0524             (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z) ? "[#mum]" : "[mrad]";
0525 
0526         diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s", s_coord.c_str()),
0527                                               Form(";#Delta%s %s;n. of modules", s_coord.c_str(), unit.c_str()),
0528                                               1000,
0529                                               -500.,
0530                                               500.);
0531       }
0532 
0533       int loopedComponents(0);
0534       for (unsigned int i = 0; i < ref_ali.size(); i++) {
0535         if (ref_ali[i].rawId() == target_ali[i].rawId()) {
0536           loopedComponents++;
0537           int subid = DetId(ref_ali[i].rawId()).subdetId();
0538           auto thePart = static_cast<AlignmentPI::partitions>(subid);
0539           if (thePart != q)
0540             continue;
0541 
0542           CLHEP::HepRotation target_rot(target_ali[i].rotation());
0543           CLHEP::HepRotation ref_rot(ref_ali[i].rotation());
0544 
0545           align::RotationType target_rotation(target_rot.xx(),
0546                                               target_rot.xy(),
0547                                               target_rot.xz(),
0548                                               target_rot.yx(),
0549                                               target_rot.yy(),
0550                                               target_rot.yz(),
0551                                               target_rot.zx(),
0552                                               target_rot.zy(),
0553                                               target_rot.zz());
0554 
0555           align::RotationType ref_rotation(ref_rot.xx(),
0556                                            ref_rot.xy(),
0557                                            ref_rot.xz(),
0558                                            ref_rot.yx(),
0559                                            ref_rot.yy(),
0560                                            ref_rot.yz(),
0561                                            ref_rot.zx(),
0562                                            ref_rot.zy(),
0563                                            ref_rot.zz());
0564 
0565           align::EulerAngles target_eulerAngles = align::toAngles(target_rotation);
0566           align::EulerAngles ref_eulerAngles = align::toAngles(ref_rotation);
0567 
0568           for (const auto &coord : coords) {
0569             switch (coord) {
0570               case AlignmentPI::t_x:
0571                 diffs[coord]->Fill((target_ali[i].translation().x() - ref_ali[i].translation().x()) *
0572                                    AlignmentPI::cmToUm);
0573                 break;
0574               case AlignmentPI::t_y:
0575                 diffs[coord]->Fill((target_ali[i].translation().y() - ref_ali[i].translation().y()) *
0576                                    AlignmentPI::cmToUm);
0577                 break;
0578               case AlignmentPI::t_z:
0579                 diffs[coord]->Fill((target_ali[i].translation().z() - ref_ali[i].translation().z()) *
0580                                    AlignmentPI::cmToUm);
0581                 break;
0582               case AlignmentPI::rot_alpha: {
0583                 auto deltaRot = target_eulerAngles[0] - ref_eulerAngles[0];
0584                 diffs[coord]->Fill(AlignmentPI::returnZeroIfNear2PI(deltaRot) * AlignmentPI::tomRad);
0585                 break;
0586               }
0587               case AlignmentPI::rot_beta: {
0588                 auto deltaRot = target_eulerAngles[1] - ref_eulerAngles[1];
0589                 diffs[coord]->Fill(AlignmentPI::returnZeroIfNear2PI(deltaRot) * AlignmentPI::tomRad);
0590                 break;
0591               }
0592               case AlignmentPI::rot_gamma: {
0593                 auto deltaRot = target_eulerAngles[2] - ref_eulerAngles[2];
0594                 diffs[coord]->Fill(AlignmentPI::returnZeroIfNear2PI(deltaRot) * AlignmentPI::tomRad);
0595                 break;
0596               }
0597               default:
0598                 edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
0599                 break;
0600             }  // switch on the coordinate
0601           }
0602         }  // check on the same detID
0603       }    // loop on the components
0604 
0605       int c_index = 1;
0606 
0607       auto legend = std::make_unique<TLegend>(0.14, 0.93, 0.55, 0.98);
0608       legend->AddEntry(
0609           diffs[AlignmentPI::t_x].get(),
0610           ("#DeltaIOV: " + std::to_string(std::get<0>(lastiov)) + "-" + std::to_string(std::get<0>(firstiov))).c_str(),
0611           "L");
0612       legend->SetTextSize(0.03);
0613 
0614       for (const auto &coord : coords) {
0615         canvas.cd(c_index)->SetLogy();
0616         canvas.cd(c_index)->SetTopMargin(0.02);
0617         canvas.cd(c_index)->SetBottomMargin(0.15);
0618         canvas.cd(c_index)->SetLeftMargin(0.14);
0619         canvas.cd(c_index)->SetRightMargin(0.04);
0620         diffs[coord]->SetLineWidth(2);
0621         AlignmentPI::makeNicePlotStyle(diffs[coord].get(), kBlack);
0622 
0623         //float x_max = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindLastBinAbove(0.));
0624         //float x_min = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindFirstBinAbove(0.));
0625         //float extremum = std::abs(x_max) > std::abs(x_min) ? std::abs(x_max) : std::abs(x_min);
0626         //diffs[coord]->GetXaxis()->SetRangeUser(-extremum*2,extremum*2);
0627 
0628         int i_max = diffs[coord]->FindLastBinAbove(0.);
0629         int i_min = diffs[coord]->FindFirstBinAbove(0.);
0630         diffs[coord]->GetXaxis()->SetRange(std::max(1, i_min - 10), std::min(i_max + 10, diffs[coord]->GetNbinsX()));
0631         diffs[coord]->Draw("HIST");
0632         AlignmentPI::makeNiceStats(diffs[coord].get(), q, kBlack);
0633 
0634         legend->Draw("same");
0635 
0636         c_index++;
0637       }
0638 
0639       std::string fileName(m_imageFileName);
0640       canvas.SaveAs(fileName.c_str());
0641 
0642       return true;
0643     }
0644   };
0645 
0646   typedef TrackerAlignmentSummary<AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
0647   typedef TrackerAlignmentSummary<AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
0648   typedef TrackerAlignmentSummary<AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
0649 
0650   typedef TrackerAlignmentSummary<AlignmentPI::TID> TrackerAlignmentSummaryTID;
0651   typedef TrackerAlignmentSummary<AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
0652   typedef TrackerAlignmentSummary<AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
0653 
0654   //*******************************************//
0655   // History of the position of the BPix Barycenter
0656   //******************************************//
0657 
0658   template <AlignmentPI::coordinate coord>
0659   class BPixBarycenterHistory : public HistoryPlot<Alignments, float> {
0660   public:
0661     BPixBarycenterHistory()
0662         : HistoryPlot<Alignments, float>(
0663               " Barrel Pixel " + AlignmentPI::getStringFromCoordinate(coord) + " positions vs time",
0664               AlignmentPI::getStringFromCoordinate(coord) + " position [cm]") {}
0665     ~BPixBarycenterHistory() override = default;
0666 
0667     float getFromPayload(Alignments &payload) override {
0668       std::vector<AlignTransform> alignments = payload.m_align;
0669 
0670       float barycenter = 0.;
0671       float nmodules(0.);
0672       for (const auto &ali : alignments) {
0673         if (DetId(ali.rawId()).det() != DetId::Tracker) {
0674           edm::LogWarning("TrackerAlignment_PayloadInspector")
0675               << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
0676               << " is different from " << DetId::Tracker << "  - terminating ";
0677           return false;
0678         }
0679 
0680         int subid = DetId(ali.rawId()).subdetId();
0681         if (subid != PixelSubdetector::PixelBarrel)
0682           continue;
0683 
0684         nmodules++;
0685         switch (coord) {
0686           case AlignmentPI::t_x:
0687             barycenter += (ali.translation().x());
0688             break;
0689           case AlignmentPI::t_y:
0690             barycenter += (ali.translation().y());
0691             break;
0692           case AlignmentPI::t_z:
0693             barycenter += (ali.translation().z());
0694             break;
0695           default:
0696             edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate " << coord << std::endl;
0697             break;
0698         }  // switch on the coordinate (only X,Y,Z are interesting)
0699       }    // ends loop on the alignments
0700 
0701       edm::LogInfo("TrackerAlignment_PayloadInspector") << "barycenter (" << barycenter << ")/n. modules (" << nmodules
0702                                                         << ") =  " << barycenter / nmodules << std::endl;
0703 
0704       // take the mean
0705       barycenter /= nmodules;
0706 
0707       // applied GPR correction to move barycenter to global CMS coordinates
0708       barycenter += hardcodeGPR.at(coord);
0709 
0710       return barycenter;
0711 
0712     }  // payload
0713   };
0714 
0715   typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
0716   typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
0717   typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
0718 
0719   /************************************************
0720     Display of Tracker Detector barycenters
0721   *************************************************/
0722   class TrackerAlignmentBarycenters : public PlotImage<Alignments, SINGLE_IOV> {
0723   public:
0724     TrackerAlignmentBarycenters() : PlotImage<Alignments, SINGLE_IOV>("Display of Tracker Alignment Barycenters") {}
0725 
0726     bool fill() override {
0727       auto tag = PlotBase::getTag<0>();
0728       auto iov = tag.iovs.front();
0729       std::shared_ptr<Alignments> payload = fetchPayload(std::get<1>(iov));
0730       unsigned int run = std::get<0>(iov);
0731 
0732       TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1600, 1000);
0733       canvas.cd();
0734 
0735       canvas.SetTopMargin(0.07);
0736       canvas.SetBottomMargin(0.06);
0737       canvas.SetLeftMargin(0.15);
0738       canvas.SetRightMargin(0.03);
0739       canvas.Modified();
0740       canvas.SetGrid();
0741 
0742       auto h2_BarycenterParameters =
0743           std::make_unique<TH2F>("Parameters", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
0744 
0745       auto h2_uncBarycenterParameters =
0746           std::make_unique<TH2F>("Parameters2", "SubDetector Barycenter summary", 6, 0.0, 6.0, 6, 0, 6.);
0747 
0748       h2_BarycenterParameters->SetStats(false);
0749       h2_BarycenterParameters->SetTitle(nullptr);
0750       h2_uncBarycenterParameters->SetStats(false);
0751       h2_uncBarycenterParameters->SetTitle(nullptr);
0752 
0753       std::vector<AlignTransform> alignments = payload->m_align;
0754 
0755       isPhase0 = (alignments.size() == AlignmentPI::phase0size) ? true : false;
0756 
0757       // check that the geomtery is a tracker one
0758       const char *path_toTopologyXML = isPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0759                                                 : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0760 
0761       TrackerTopology tTopo =
0762           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0763 
0764       AlignmentPI::TkAlBarycenters barycenters;
0765       // compute uncorrected barycenter
0766       barycenters.computeBarycenters(
0767           alignments, tTopo, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
0768 
0769       auto Xbarycenters = barycenters.getX();
0770       auto Ybarycenters = barycenters.getY();
0771       auto Zbarycenters = barycenters.getZ();
0772 
0773       // compute barycenter corrected for the GPR
0774       barycenters.computeBarycenters(alignments, tTopo, hardcodeGPR);
0775 
0776       auto c_Xbarycenters = barycenters.getX();
0777       auto c_Ybarycenters = barycenters.getY();
0778       auto c_Zbarycenters = barycenters.getZ();
0779 
0780       h2_BarycenterParameters->GetXaxis()->SetBinLabel(1, "X [cm]");
0781       h2_BarycenterParameters->GetXaxis()->SetBinLabel(2, "Y [cm]");
0782       h2_BarycenterParameters->GetXaxis()->SetBinLabel(3, "Z [cm]");
0783       h2_BarycenterParameters->GetXaxis()->SetBinLabel(4, "X_{no GPR} [cm]");
0784       h2_BarycenterParameters->GetXaxis()->SetBinLabel(5, "Y_{no GPR} [cm]");
0785       h2_BarycenterParameters->GetXaxis()->SetBinLabel(6, "Z_{no GPR} [cm]");
0786 
0787       bool isLikelyMC(false);
0788       int checkX =
0789           std::count_if(Xbarycenters.begin(), Xbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0790       int checkY =
0791           std::count_if(Ybarycenters.begin(), Ybarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0792       int checkZ =
0793           std::count_if(Zbarycenters.begin(), Zbarycenters.begin() + 2, [](float a) { return (std::abs(a) >= 1.e-4); });
0794 
0795       // if all the coordinate barycenters for both BPix and FPix are below 10um
0796       // this is very likely a MC payload
0797       if ((checkX + checkY + checkZ) == 0 && run == 1)
0798         isLikelyMC = true;
0799 
0800       unsigned int yBin = 6;
0801       for (unsigned int i = 0; i < 6; i++) {
0802         auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
0803         std::string theLabel = getStringFromPart(thePart);
0804         h2_BarycenterParameters->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0805         if (!isLikelyMC) {
0806           h2_BarycenterParameters->SetBinContent(1, yBin, c_Xbarycenters[i]);
0807           h2_BarycenterParameters->SetBinContent(2, yBin, c_Ybarycenters[i]);
0808           h2_BarycenterParameters->SetBinContent(3, yBin, c_Zbarycenters[i]);
0809         }
0810 
0811         h2_uncBarycenterParameters->SetBinContent(4, yBin, Xbarycenters[i]);
0812         h2_uncBarycenterParameters->SetBinContent(5, yBin, Ybarycenters[i]);
0813         h2_uncBarycenterParameters->SetBinContent(6, yBin, Zbarycenters[i]);
0814         yBin--;
0815       }
0816 
0817       h2_BarycenterParameters->GetXaxis()->LabelsOption("h");
0818       h2_BarycenterParameters->GetYaxis()->SetLabelSize(0.05);
0819       h2_BarycenterParameters->GetXaxis()->SetLabelSize(0.05);
0820       h2_BarycenterParameters->SetMarkerSize(1.5);
0821       h2_BarycenterParameters->Draw("TEXT");
0822 
0823       h2_uncBarycenterParameters->SetMarkerSize(1.5);
0824       h2_uncBarycenterParameters->SetMarkerColor(kRed);
0825       h2_uncBarycenterParameters->Draw("TEXTsame");
0826 
0827       TLatex t1;
0828       t1.SetNDC();
0829       t1.SetTextAlign(26);
0830       t1.SetTextSize(0.05);
0831       t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters, IOV %i", run));
0832       t1.SetTextSize(0.025);
0833 
0834       std::string fileName(m_imageFileName);
0835       canvas.SaveAs(fileName.c_str());
0836 
0837       return true;
0838     }
0839 
0840   private:
0841     bool isPhase0;
0842   };
0843 
0844   /************************************************
0845     Comparator of Tracker Detector barycenters
0846   *************************************************/
0847   template <int ntags, IOVMultiplicity nIOVs>
0848   class TrackerAlignmentBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0849   public:
0850     TrackerAlignmentBarycentersComparatorBase()
0851         : PlotImage<Alignments, nIOVs, ntags>("Comparison of Tracker Alignment Barycenters") {}
0852 
0853     bool fill() override {
0854       // trick to deal with the multi-ioved tag and two tag case at the same time
0855       auto theIOVs = PlotBase::getTag<0>().iovs;
0856       auto tagname1 = PlotBase::getTag<0>().name;
0857       std::string tagname2 = "";
0858       auto firstiov = theIOVs.front();
0859       std::tuple<cond::Time_t, cond::Hash> lastiov;
0860 
0861       // we don't support (yet) comparison with more than 2 tags
0862       assert(this->m_plotAnnotations.ntags < 3);
0863 
0864       if (this->m_plotAnnotations.ntags == 2) {
0865         auto tag2iovs = PlotBase::getTag<1>().iovs;
0866         tagname2 = PlotBase::getTag<1>().name;
0867         lastiov = tag2iovs.front();
0868       } else {
0869         lastiov = theIOVs.back();
0870       }
0871 
0872       unsigned int first_run = std::get<0>(firstiov);
0873       unsigned int last_run = std::get<0>(lastiov);
0874 
0875       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0876       std::vector<AlignTransform> last_alignments = last_payload->m_align;
0877 
0878       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0879       std::vector<AlignTransform> first_alignments = first_payload->m_align;
0880 
0881       isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
0882       isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
0883 
0884       // check that the geomtery is a tracker one
0885       const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0886                                                        : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0887 
0888       TrackerTopology tTopo_f =
0889           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0890 
0891       path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
0892                                          : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0893 
0894       TrackerTopology tTopo_l =
0895           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0896 
0897       TCanvas canvas("Tracker Alignment Barycenter Summary", "Tracker Alignment Barycenter summary", 1200, 800);
0898       canvas.cd();
0899 
0900       canvas.SetTopMargin(0.07);
0901       canvas.SetBottomMargin(0.06);
0902       canvas.SetLeftMargin(0.15);
0903       canvas.SetRightMargin(0.03);
0904       canvas.Modified();
0905       canvas.SetGrid();
0906 
0907       auto h2_BarycenterDiff =
0908           std::make_unique<TH2F>("Parameters diff", "SubDetector Barycenter Difference", 3, 0.0, 3.0, 6, 0, 6.);
0909 
0910       h2_BarycenterDiff->SetStats(false);
0911       h2_BarycenterDiff->SetTitle(nullptr);
0912       h2_BarycenterDiff->GetXaxis()->SetBinLabel(1, "X [#mum]");
0913       h2_BarycenterDiff->GetXaxis()->SetBinLabel(2, "Y [#mum]");
0914       h2_BarycenterDiff->GetXaxis()->SetBinLabel(3, "Z [#mum]");
0915 
0916       AlignmentPI::TkAlBarycenters l_barycenters;
0917       l_barycenters.computeBarycenters(last_alignments, tTopo_l, hardcodeGPR);
0918 
0919       AlignmentPI::TkAlBarycenters f_barycenters;
0920       f_barycenters.computeBarycenters(first_alignments, tTopo_f, hardcodeGPR);
0921 
0922       unsigned int yBin = 6;
0923       for (unsigned int i = 0; i < 6; i++) {
0924         auto thePart = static_cast<AlignmentPI::partitions>(i + 1);
0925         std::string theLabel = getStringFromPart(thePart);
0926         h2_BarycenterDiff->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0927         h2_BarycenterDiff->SetBinContent(
0928             1, yBin, (l_barycenters.getX()[i] - f_barycenters.getX()[i]) * AlignmentPI::cmToUm);
0929         h2_BarycenterDiff->SetBinContent(
0930             2, yBin, (l_barycenters.getY()[i] - f_barycenters.getY()[i]) * AlignmentPI::cmToUm);
0931         h2_BarycenterDiff->SetBinContent(
0932             3, yBin, (l_barycenters.getZ()[i] - f_barycenters.getZ()[i]) * AlignmentPI::cmToUm);
0933         yBin--;
0934       }
0935 
0936       h2_BarycenterDiff->GetXaxis()->LabelsOption("h");
0937       h2_BarycenterDiff->GetYaxis()->SetLabelSize(0.05);
0938       h2_BarycenterDiff->GetXaxis()->SetLabelSize(0.05);
0939       h2_BarycenterDiff->SetMarkerSize(1.5);
0940       h2_BarycenterDiff->SetMarkerColor(kRed);
0941       h2_BarycenterDiff->Draw("TEXT");
0942 
0943       TLatex t1;
0944       t1.SetNDC();
0945       t1.SetTextAlign(26);
0946       t1.SetTextSize(0.05);
0947       t1.DrawLatex(0.5, 0.96, Form("Tracker Alignment Barycenters Diff, IOV %i - IOV %i", last_run, first_run));
0948       t1.SetTextSize(0.025);
0949 
0950       std::string fileName(this->m_imageFileName);
0951       canvas.SaveAs(fileName.c_str());
0952 
0953       return true;
0954     }
0955 
0956   private:
0957     bool isInitialPhase0;
0958     bool isFinalPhase0;
0959   };
0960 
0961   using TrackerAlignmentBarycentersCompare = TrackerAlignmentBarycentersComparatorBase<1, MULTI_IOV>;
0962   using TrackerAlignmentBarycentersCompareTwoTags = TrackerAlignmentBarycentersComparatorBase<2, SINGLE_IOV>;
0963 
0964   /************************************************
0965     Comparator of Pixel Tracker Detector barycenters
0966   *************************************************/
0967   template <int ntags, IOVMultiplicity nIOVs>
0968   class PixelBarycentersComparatorBase : public PlotImage<Alignments, nIOVs, ntags> {
0969   public:
0970     PixelBarycentersComparatorBase() : PlotImage<Alignments, nIOVs, ntags>("Comparison of Pixel Barycenters") {}
0971 
0972     bool fill() override {
0973       // trick to deal with the multi-ioved tag and two tag case at the same time
0974       auto theIOVs = PlotBase::getTag<0>().iovs;
0975       auto tagname1 = PlotBase::getTag<0>().name;
0976       std::string tagname2 = "";
0977       auto firstiov = theIOVs.front();
0978       std::tuple<cond::Time_t, cond::Hash> lastiov;
0979 
0980       // we don't support (yet) comparison with more than 2 tags
0981       assert(this->m_plotAnnotations.ntags < 3);
0982 
0983       if (this->m_plotAnnotations.ntags == 2) {
0984         auto tag2iovs = PlotBase::getTag<1>().iovs;
0985         tagname2 = PlotBase::getTag<1>().name;
0986         lastiov = tag2iovs.front();
0987       } else {
0988         lastiov = theIOVs.back();
0989       }
0990 
0991       unsigned int first_run = std::get<0>(firstiov);
0992       unsigned int last_run = std::get<0>(lastiov);
0993 
0994       std::shared_ptr<Alignments> last_payload = this->fetchPayload(std::get<1>(lastiov));
0995       std::vector<AlignTransform> last_alignments = last_payload->m_align;
0996 
0997       std::shared_ptr<Alignments> first_payload = this->fetchPayload(std::get<1>(firstiov));
0998       std::vector<AlignTransform> first_alignments = first_payload->m_align;
0999 
1000       TCanvas canvas("Pixel Barycenter Summary", "Pixel Barycenter summary", 1200, 1200);
1001       canvas.Divide(2, 2);
1002       canvas.cd();
1003 
1004       TLatex t1;
1005       t1.SetNDC();
1006       t1.SetTextAlign(26);
1007       t1.SetTextSize(0.03);
1008       t1.DrawLatex(0.5,
1009                    0.97,
1010                    ("Pixel Barycenters comparison, IOV: #color[2]{" + std::to_string(first_run) +
1011                     "} vs IOV: #color[4]{" + std::to_string(last_run) + "}")
1012                        .c_str());
1013       t1.SetTextSize(0.025);
1014 
1015       for (unsigned int c = 1; c <= 4; c++) {
1016         canvas.cd(c)->SetTopMargin(0.07);
1017         canvas.cd(c)->SetBottomMargin(0.12);
1018         canvas.cd(c)->SetLeftMargin(0.15);
1019         canvas.cd(c)->SetRightMargin(0.03);
1020         canvas.cd(c)->Modified();
1021         canvas.cd(c)->SetGrid();
1022       }
1023 
1024       std::array<std::string, 3> structures = {{"FPIX-", "BPIX", "FPIX+"}};
1025       std::array<std::unique_ptr<TH2F>, 3> histos;
1026 
1027       isInitialPhase0 = (first_alignments.size() == AlignmentPI::phase0size) ? true : false;
1028       isFinalPhase0 = (last_alignments.size() == AlignmentPI::phase0size) ? true : false;
1029 
1030       // check that the geomtery is a tracker one
1031       const char *path_toTopologyXML = isInitialPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1032                                                        : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1033 
1034       TrackerTopology tTopo_f =
1035           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1036 
1037       AlignmentPI::TkAlBarycenters myInitialBarycenters;
1038       //myInitialBarycenters.computeBarycenters(first_alignments,tTopo_f,hardcodeGPR);
1039       myInitialBarycenters.computeBarycenters(
1040           first_alignments, tTopo_f, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1041 
1042       path_toTopologyXML = isFinalPhase0 ? "Geometry/TrackerCommonData/data/trackerParameters.xml"
1043                                          : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
1044 
1045       TrackerTopology tTopo_l =
1046           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1047 
1048       AlignmentPI::TkAlBarycenters myFinalBarycenters;
1049       //myFinalBarycenters.computeBarycenters(last_alignments,tTopo_l,hardcodeGPR);
1050       myFinalBarycenters.computeBarycenters(
1051           last_alignments, tTopo_l, {{AlignmentPI::t_x, 0.0}, {AlignmentPI::t_y, 0.0}, {AlignmentPI::t_z, 0.0}});
1052 
1053       if (isFinalPhase0 != isInitialPhase0) {
1054         edm::LogWarning("TrackerAlignment_PayloadInspector")
1055             << "the size of the reference alignment (" << first_alignments.size()
1056             << ") is different from the one of the target (" << last_alignments.size()
1057             << ")! You are probably trying to compare different underlying geometries.";
1058       }
1059 
1060       unsigned int index(0);
1061       for (const auto &piece : structures) {
1062         const char *name = piece.c_str();
1063         histos[index] = std::make_unique<TH2F>(
1064             name,
1065             Form("%s x-y Barycenter Difference;x_{%s}-x_{TOB} [mm];y_{%s}-y_{TOB} [mm]", name, name, name),
1066             100,
1067             -3.,
1068             3.,
1069             100,
1070             -3.,
1071             3.);
1072 
1073         histos[index]->SetStats(false);
1074         histos[index]->SetTitle(nullptr);
1075         histos[index]->GetYaxis()->SetLabelSize(0.05);
1076         histos[index]->GetXaxis()->SetLabelSize(0.05);
1077         histos[index]->GetYaxis()->SetTitleSize(0.06);
1078         histos[index]->GetXaxis()->SetTitleSize(0.06);
1079         histos[index]->GetYaxis()->CenterTitle();
1080         histos[index]->GetXaxis()->CenterTitle();
1081         histos[index]->GetXaxis()->SetTitleOffset(0.9);
1082         index++;
1083       }
1084 
1085       auto h2_ZBarycenterDiff = std::make_unique<TH2F>(
1086           "Pixel_z_diff", "Pixel z-Barycenter Difference;; z_{Pixel-Ideal} -z_{TOB} [mm]", 3, -0.5, 2.5, 100, -10., 10.);
1087       h2_ZBarycenterDiff->SetStats(false);
1088       h2_ZBarycenterDiff->SetTitle(nullptr);
1089       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(1, "FPIX -");
1090       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(2, "BPIX");
1091       h2_ZBarycenterDiff->GetXaxis()->SetBinLabel(3, "FPIX +");
1092       h2_ZBarycenterDiff->GetYaxis()->SetLabelSize(0.05);
1093       h2_ZBarycenterDiff->GetXaxis()->SetLabelSize(0.07);
1094       h2_ZBarycenterDiff->GetYaxis()->SetTitleSize(0.06);
1095       h2_ZBarycenterDiff->GetXaxis()->SetTitleSize(0.06);
1096       h2_ZBarycenterDiff->GetYaxis()->CenterTitle();
1097       h2_ZBarycenterDiff->GetXaxis()->CenterTitle();
1098       h2_ZBarycenterDiff->GetYaxis()->SetTitleOffset(1.1);
1099 
1100       std::function<GlobalPoint(int)> cutFunctorInitial = [&myInitialBarycenters](int index) {
1101         switch (index) {
1102           case 1:
1103             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1104           case 2:
1105             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1106           case 3:
1107             return myInitialBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1108           default:
1109             return GlobalPoint(0, 0, 0);
1110         }
1111       };
1112 
1113       std::function<GlobalPoint(int)> cutFunctorFinal = [&myFinalBarycenters](int index) {
1114         switch (index) {
1115           case 1:
1116             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXm);
1117           case 2:
1118             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::BPIX);
1119           case 3:
1120             return myFinalBarycenters.getPartitionAvg(AlignmentPI::PARTITION::FPIXp);
1121           default:
1122             return GlobalPoint(0, 0, 0);
1123         }
1124       };
1125 
1126       float x0i, x0f, y0i, y0f;
1127 
1128       t1.SetNDC(kFALSE);
1129       t1.SetTextSize(0.047);
1130       for (unsigned int c = 1; c <= 3; c++) {
1131         x0i = cutFunctorInitial(c).x() * 10;  // transform cm to mm (x10)
1132         x0f = cutFunctorFinal(c).x() * 10;
1133         y0i = cutFunctorInitial(c).y() * 10;
1134         y0f = cutFunctorFinal(c).y() * 10;
1135 
1136         canvas.cd(c);
1137         histos[c - 1]->Draw();
1138 
1139         COUT << "initial x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0i << "," << y0i << ") mm"
1140              << std::endl;
1141         COUT << "final   x,y " << std::left << std::setw(7) << structures[c - 1] << " (" << x0f << "," << y0f << ") mm"
1142              << std::endl;
1143 
1144         TMarker *initial = new TMarker(x0i, y0i, 21);
1145         TMarker *final = new TMarker(x0f, y0f, 20);
1146 
1147         initial->SetMarkerColor(kRed);
1148         final->SetMarkerColor(kBlue);
1149         initial->SetMarkerSize(2.5);
1150         final->SetMarkerSize(2.5);
1151         t1.SetTextColor(kRed);
1152         initial->Draw();
1153         t1.DrawLatex(x0i, y0i + (y0i > y0f ? 0.3 : -0.5), Form("(%.2f,%.2f)", x0i, y0i));
1154         final->Draw("same");
1155         t1.SetTextColor(kBlue);
1156         t1.DrawLatex(x0f, y0f + (y0i > y0f ? -0.5 : 0.3), Form("(%.2f,%.2f)", x0f, y0f));
1157       }
1158 
1159       // fourth pad is a special case for the z coordinate
1160       canvas.cd(4);
1161       h2_ZBarycenterDiff->Draw();
1162       float z0i, z0f;
1163 
1164       // numbers do agree with https://twiki.cern.ch/twiki/bin/view/CMSPublic/TkAlignmentPerformancePhaseIStartUp17#Pixel_Barycentre_Positions
1165 
1166       std::array<double, 3> hardcodeIdealZPhase0 = {{-41.94909, 0., 41.94909}};  // units are cm
1167       std::array<double, 3> hardcodeIdealZPhase1 = {{-39.82911, 0., 39.82911}};  // units are cm
1168 
1169       for (unsigned int c = 1; c <= 3; c++) {
1170         // less than pretty but needed to remove the z position of the FPix barycenters != 0
1171 
1172         z0i =
1173             (cutFunctorInitial(c).z() - (isInitialPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) *
1174             10;  // convert to mm
1175         z0f =
1176             (cutFunctorFinal(c).z() - (isFinalPhase0 ? hardcodeIdealZPhase0[c - 1] : hardcodeIdealZPhase1[c - 1])) * 10;
1177 
1178         TMarker *initial = new TMarker(c - 1, z0i, 21);
1179         TMarker *final = new TMarker(c - 1, z0f, 20);
1180 
1181         COUT << "initial   z " << std::left << std::setw(7) << structures[c - 1] << " " << z0i << " mm" << std::endl;
1182         COUT << "final     z " << std::left << std::setw(7) << structures[c - 1] << " " << z0f << " mm" << std::endl;
1183 
1184         initial->SetMarkerColor(kRed);
1185         final->SetMarkerColor(kBlue);
1186         initial->SetMarkerSize(2.5);
1187         final->SetMarkerSize(2.5);
1188         initial->Draw();
1189         t1.SetTextColor(kRed);
1190         t1.DrawLatex(c - 1, z0i + (z0i > z0f ? 1. : -1.5), Form("(%.2f)", z0i));
1191         final->Draw("same");
1192         t1.SetTextColor(kBlue);
1193         t1.DrawLatex(c - 1, z0f + (z0i > z0f ? -1.5 : 1), Form("(%.2f)", z0f));
1194       }
1195 
1196       std::string fileName(this->m_imageFileName);
1197       canvas.SaveAs(fileName.c_str());
1198 
1199       return true;
1200     }
1201 
1202   private:
1203     bool isInitialPhase0;
1204     bool isFinalPhase0;
1205   };
1206 
1207   using PixelBarycentersCompare = PixelBarycentersComparatorBase<1, MULTI_IOV>;
1208   using PixelBarycentersCompareTwoTags = PixelBarycentersComparatorBase<2, SINGLE_IOV>;
1209 
1210 }  // namespace
1211 
1212 PAYLOAD_INSPECTOR_MODULE(TrackerAlignment) {
1213   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorSingleTag);
1214   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentComparatorTwoTags);
1215   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
1216   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
1217   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
1218   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
1219   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
1220   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
1221   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareXTwoTags);
1222   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareYTwoTags);
1223   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZTwoTags);
1224   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlphaTwoTags);
1225   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBetaTwoTags);
1226   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGammaTwoTags);
1227   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
1228   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
1229   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
1230   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
1231   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
1232   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
1233   PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
1234   PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
1235   PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
1236   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycenters);
1237   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompare);
1238   PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentBarycentersCompareTwoTags);
1239   PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompare);
1240   PAYLOAD_INSPECTOR_CLASS(PixelBarycentersCompareTwoTags);
1241 }