Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-11 23:27:41

0001 #ifndef CONDCORE_PCLCONFIGPLUGINS_SIPIXELALIPCLTHRESHOLDSPAYLOADINSPECTORHELPER_H
0002 #define CONDCORE_PCLCONFIGPLUGINS_SIPIXELALIPCLTHRESHOLDSPAYLOADINSPECTORHELPER_H
0003 
0004 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0005 #include "CondCore/Utilities/interface/PayloadInspector.h"
0006 #include "CondCore/CondDB/interface/Time.h"
0007 
0008 // the data format of the condition to be inspected
0009 #include "CondFormats/PCLConfig/interface/AlignPCLThresholds.h"
0010 #include "CondFormats/PCLConfig/interface/AlignPCLThresholdsHG.h"
0011 
0012 #include <memory>
0013 #include <sstream>
0014 #include <iostream>
0015 #include <functional>
0016 #include <fmt/printf.h>
0017 
0018 // include ROOT
0019 #include "TH2F.h"
0020 #include "TLegend.h"
0021 #include "TCanvas.h"
0022 #include "TLatex.h"
0023 #include "TLine.h"
0024 #include "TStyle.h"
0025 
0026 namespace AlignPCLThresholdPlotHelper {
0027 
0028   enum types { DELTA = 0, SIG = 1, MAXMOVE = 2, MAXERR = 3, FRACTION_CUT = 4, END_OF_TYPES = 5 };
0029 
0030   /************************************************/
0031   inline const std::string getStringFromCoordEnum(const AlignPCLThresholds::coordType& coord) {
0032     switch (coord) {
0033       case AlignPCLThresholds::X:
0034         return "X";
0035       case AlignPCLThresholds::Y:
0036         return "Y";
0037       case AlignPCLThresholds::Z:
0038         return "Z";
0039       case AlignPCLThresholds::theta_X:
0040         return "#theta_{X}";
0041       case AlignPCLThresholds::theta_Y:
0042         return "#theta_{Y}";
0043       case AlignPCLThresholds::theta_Z:
0044         return "#theta_{Z}";
0045       default:
0046         return "should never be here";
0047     }
0048   }
0049 
0050   /************************************************/
0051   inline const std::string getStringFromTypeEnum(const types& type) {
0052     switch (type) {
0053       case types::DELTA:
0054         return "#Delta";
0055       case types::SIG:
0056         return "#Delta/#sigma ";
0057       case types::MAXMOVE:
0058         return "max. move ";
0059       case types::MAXERR:
0060         return "max. err ";
0061       case types::FRACTION_CUT:
0062         return "fraction cut ";
0063       default:
0064         return "should never be here";
0065     }
0066   }
0067 
0068   /************************************************/
0069   inline std::string replaceAll(const std::string& str, const std::string& from, const std::string& to) {
0070     std::string out(str);
0071 
0072     if (from.empty())
0073       return out;
0074     size_t start_pos = 0;
0075     while ((start_pos = out.find(from, start_pos)) != std::string::npos) {
0076       out.replace(start_pos, from.length(), to);
0077       start_pos += to.length();  // In case 'to' contains 'from', like replacing 'x' with 'yx'
0078     }
0079     return out;
0080   }
0081 
0082   /************************************************
0083     Display of AlignPCLThreshold
0084   *************************************************/
0085   template <class PayloadType>
0086   class AlignPCLThresholds_DisplayBase
0087       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
0088   public:
0089     AlignPCLThresholds_DisplayBase()
0090         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
0091               "Display of threshold parameters for SiPixelAli PCL") {
0092       if constexpr (std::is_same_v<PayloadType, AlignPCLThresholdsHG>) {
0093         isHighGranularity_ = true;
0094         label_ = "AlignPCLThresholdsHG_PayloadInspector";
0095       } else {
0096         isHighGranularity_ = false;
0097         label_ = "AlignPCLThresholds_PayloadInspector";
0098       }
0099     }
0100 
0101     bool fill() override {
0102       // to be able to see zeros
0103       gStyle->SetHistMinimumZero(kTRUE);
0104 
0105       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
0106       auto iov = tag.iovs.front();
0107       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
0108       auto alignables = payload->getAlignableList();
0109 
0110       TCanvas canvas("Alignment PCL thresholds summary", "Alignment PCL thresholds summary", 1500, 800);
0111       canvas.cd();
0112 
0113       canvas.SetTopMargin(0.07);
0114       canvas.SetBottomMargin(isHighGranularity_ ? 0.20 : 0.06);
0115       canvas.SetLeftMargin(0.11);
0116       canvas.SetRightMargin(0.05);
0117       canvas.Modified();
0118       canvas.SetGrid();
0119 
0120       // needed for the internal loop
0121       const int local_end_of_types = static_cast<int>(isHighGranularity_ ? END_OF_TYPES : FRACTION_CUT);
0122       const int N_Y_BINS = static_cast<int>(AlignPCLThresholds::extra_DOF) * local_end_of_types;
0123 
0124       auto Thresholds =
0125           std::make_unique<TH2F>("Thresholds", "", alignables.size(), 0, alignables.size(), N_Y_BINS, 0, N_Y_BINS);
0126       Thresholds->SetStats(false);
0127       Thresholds->GetXaxis()->SetLabelSize(0.028);
0128 
0129       std::function<float(types, std::string, AlignPCLThresholds::coordType)> cutFunctor =
0130           [&payload](types my_type, std::string alignable, AlignPCLThresholds::coordType coord) {
0131             float ret(-999.);
0132             switch (my_type) {
0133               case DELTA:
0134                 return payload->getCut(alignable, coord);
0135               case SIG:
0136                 return payload->getSigCut(alignable, coord);
0137               case MAXMOVE:
0138                 return payload->getMaxMoveCut(alignable, coord);
0139               case MAXERR:
0140                 return payload->getMaxErrorCut(alignable, coord);
0141               case FRACTION_CUT: {
0142                 if constexpr (std::is_same_v<PayloadType, AlignPCLThresholdsHG>) {
0143                   const AlignPCLThresholdsHG::param_map& floatMap = payload->getFloatMap();
0144                   if (floatMap.find(alignable) != floatMap.end()) {
0145                     return payload->getFractionCut(alignable, coord);
0146                   } else {
0147                     return 0.f;
0148                   }
0149                 } else {
0150                   assert(false); /* cannot be here */
0151                 }
0152               }
0153               case END_OF_TYPES:
0154                 return ret;
0155               default:
0156                 return ret;
0157             }
0158           };
0159 
0160       unsigned int xBin = 0;
0161       for (const auto& alignable : alignables) {
0162         xBin++;
0163         auto xLabel = replaceAll(replaceAll(alignable, "minus", "(-)"), "plus", "(+)");
0164         Thresholds->GetXaxis()->SetBinLabel(xBin, (xLabel).c_str());
0165         unsigned int yBin = N_Y_BINS;
0166         for (int foo = PayloadType::X; foo != PayloadType::extra_DOF; foo++) {
0167           AlignPCLThresholds::coordType coord = static_cast<AlignPCLThresholds::coordType>(foo);
0168           for (int bar = DELTA; bar != local_end_of_types; bar++) {
0169             types type = static_cast<types>(bar);
0170             std::string theLabel = getStringFromTypeEnum(type) + getStringFromCoordEnum(coord);
0171             if (xBin == 1) {
0172               Thresholds->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0173             }
0174 
0175             Thresholds->SetBinContent(xBin, yBin, cutFunctor(type, alignable, coord));
0176 
0177             yBin--;
0178           }  // loop on types
0179         }    // loop on coordinates
0180       }      // loop on alignables
0181 
0182       Thresholds->GetXaxis()->LabelsOption(isHighGranularity_ ? "v" : "h");
0183       Thresholds->Draw("TEXT");
0184 
0185       auto ltx = TLatex();
0186       ltx.SetTextFont(62);
0187       //ltx.SetTextColor(kBlue);
0188       ltx.SetTextSize(0.047);
0189       ltx.SetTextAlign(11);
0190       std::string ltxText =
0191           fmt::sprintf("#color[4]{%s} IOV: #color[4]{%s}", tag.name, std::to_string(std::get<0>(iov)));
0192       ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
0193 
0194       std::string fileName(this->m_imageFileName);
0195       canvas.SaveAs(fileName.c_str());
0196 
0197       return true;
0198     }
0199 
0200   private:
0201     bool isHighGranularity_;
0202     std::string label_;
0203   };
0204 
0205   /************************************************
0206     Compare AlignPCLThresholdsHG mapping
0207   *************************************************/
0208   template <class PayloadType, cond::payloadInspector::IOVMultiplicity nIOVs, int ntags>
0209   class AlignPCLThresholds_CompareBase : public cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags> {
0210   public:
0211     AlignPCLThresholds_CompareBase()
0212         : cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags>("Table of AlignPCLThresholdsHG comparison") {
0213       if constexpr (std::is_same_v<PayloadType, AlignPCLThresholdsHG>) {
0214         isHighGranularity_ = true;
0215         label_ = "AlignPCLThresholdsHG_PayloadInspector";
0216       } else {
0217         isHighGranularity_ = false;
0218         label_ = "AlignPCLThresholds_PayloadInspector";
0219       }
0220     }
0221 
0222     bool fill() override {
0223       gStyle->SetPalette(kTemperatureMap);
0224 
0225       // to be able to see zeros
0226       gStyle->SetHistMinimumZero(kTRUE);
0227 
0228       // trick to deal with the multi-ioved tag and two tag case at the same time
0229       auto theIOVs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
0230       auto f_tagname = cond::payloadInspector::PlotBase::getTag<0>().name;
0231       std::string l_tagname = "";
0232       auto firstiov = theIOVs.front();
0233       std::tuple<cond::Time_t, cond::Hash> lastiov;
0234 
0235       // we don't support (yet) comparison with more than 2 tags
0236       assert(this->m_plotAnnotations.ntags < 3);
0237 
0238       if (this->m_plotAnnotations.ntags == 2) {
0239         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
0240         l_tagname = cond::payloadInspector::PlotBase::getTag<1>().name;
0241         lastiov = tag2iovs.front();
0242       } else {
0243         lastiov = theIOVs.back();
0244       }
0245 
0246       std::shared_ptr<PayloadType> l_payload = this->fetchPayload(std::get<1>(lastiov));
0247       std::shared_ptr<PayloadType> f_payload = this->fetchPayload(std::get<1>(firstiov));
0248 
0249       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0250       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0251 
0252       const auto& alignables = l_payload->getAlignableList();
0253       const auto& alignables2 = f_payload->getAlignableList();
0254 
0255       std::vector<std::string> v_intersection;
0256 
0257       if (!isEqual(alignables, alignables2)) {
0258         edm::LogWarning(label_)
0259             << "Cannot compare directly the two AlignPCLThresholds objects, as the list of alignables differs";
0260         std::set_intersection(alignables.begin(),
0261                               alignables.end(),
0262                               alignables2.begin(),
0263                               alignables2.end(),
0264                               std::back_inserter(v_intersection));
0265 
0266         std::vector<std::string> not_in_first_keys, not_in_last_keys;
0267 
0268         // find the elements not in common
0269         std::set_difference(alignables.begin(),
0270                             alignables.end(),
0271                             alignables2.begin(),
0272                             alignables2.end(),
0273                             std::inserter(not_in_last_keys, not_in_last_keys.begin()));
0274 
0275         std::stringstream ss;
0276         ss << "the following keys are not in the last IoV: ";
0277         for (const auto& key : not_in_last_keys) {
0278           ss << key << ",";
0279         }
0280         ss << std::endl;
0281         edm::LogWarning(label_) << ss.str();
0282         ss.str(std::string()); /* clear the stream */
0283 
0284         std::set_difference(alignables2.begin(),
0285                             alignables2.end(),
0286                             alignables.begin(),
0287                             alignables.end(),
0288                             std::inserter(not_in_first_keys, not_in_first_keys.begin()));
0289 
0290         ss << "the following keys are not in the first IoV: ";
0291         for (const auto& key : not_in_first_keys) {
0292           ss << key << ",";
0293         }
0294         ss << std::endl;
0295         edm::LogWarning(label_) << ss.str();
0296       } else {
0297         // vectors are the same, just copy one
0298         v_intersection = alignables;
0299       }
0300 
0301       TCanvas canvas("Alignment PCL thresholds summary", "Alignment PCL thresholds summary", 1500, 800);
0302       canvas.cd();
0303 
0304       canvas.SetTopMargin(0.07);
0305       canvas.SetBottomMargin(isHighGranularity_ ? 0.20 : 0.06);
0306       canvas.SetLeftMargin(0.11);
0307       canvas.SetRightMargin(0.12);
0308       canvas.Modified();
0309       canvas.SetGrid();
0310 
0311       // needed for the internal loop
0312       const int local_end_of_types = static_cast<int>(isHighGranularity_ ? END_OF_TYPES : FRACTION_CUT);
0313       const int N_Y_BINS = static_cast<int>(AlignPCLThresholds::extra_DOF) * local_end_of_types;
0314 
0315       auto Thresholds = std::make_unique<TH2F>(
0316           "Thresholds", "", v_intersection.size(), 0, v_intersection.size(), N_Y_BINS, 0, N_Y_BINS);
0317       Thresholds->SetStats(false);
0318       Thresholds->GetXaxis()->SetLabelSize(0.028);
0319 
0320       auto ThresholdsColor = std::make_unique<TH2F>(
0321           "ThresholdsC", "", v_intersection.size(), 0, v_intersection.size(), N_Y_BINS, 0, N_Y_BINS);
0322       ThresholdsColor->SetStats(false);
0323       ThresholdsColor->GetXaxis()->SetLabelSize(0.028);
0324 
0325       std::function<float(types, std::string, AlignPCLThresholds::coordType)> cutFunctor =
0326           [&f_payload, &l_payload](types my_type, std::string alignable, AlignPCLThresholds::coordType coord) {
0327             float ret(-999.);
0328             switch (my_type) {
0329               case DELTA:
0330                 return l_payload->getCut(alignable, coord) - f_payload->getCut(alignable, coord);
0331               case SIG:
0332                 return l_payload->getSigCut(alignable, coord) - f_payload->getSigCut(alignable, coord);
0333               case MAXMOVE:
0334                 return l_payload->getMaxMoveCut(alignable, coord) - f_payload->getMaxMoveCut(alignable, coord);
0335               case MAXERR:
0336                 return l_payload->getMaxErrorCut(alignable, coord) - f_payload->getMaxErrorCut(alignable, coord);
0337               case FRACTION_CUT: {
0338                 if constexpr (std::is_same_v<PayloadType, AlignPCLThresholdsHG>) {
0339                   const AlignPCLThresholdsHG::param_map& f_floatMap = f_payload->getFloatMap();
0340                   const AlignPCLThresholdsHG::param_map& l_floatMap = l_payload->getFloatMap();
0341                   if (f_floatMap.find(alignable) != f_floatMap.end() &&
0342                       l_floatMap.find(alignable) != l_floatMap.end()) {
0343                     return l_payload->getFractionCut(alignable, coord) - f_payload->getFractionCut(alignable, coord);
0344                   } else {
0345                     return +999.f;
0346                   }
0347                 } else {
0348                   assert(false); /* cannot be here */
0349                 }
0350               }
0351               case END_OF_TYPES:
0352                 return ret;
0353               default:
0354                 return ret;
0355             }
0356           };
0357 
0358       unsigned int xBin = 0;
0359       for (const auto& alignable : v_intersection) {
0360         xBin++;
0361         auto xLabel = replaceAll(replaceAll(alignable, "minus", "(-)"), "plus", "(+)");
0362         Thresholds->GetXaxis()->SetBinLabel(xBin, (xLabel).c_str());
0363         ThresholdsColor->GetXaxis()->SetBinLabel(xBin, (xLabel).c_str());
0364         unsigned int yBin = N_Y_BINS;
0365         for (int foo = PayloadType::X; foo != PayloadType::extra_DOF; foo++) {
0366           AlignPCLThresholds::coordType coord = static_cast<AlignPCLThresholds::coordType>(foo);
0367           for (int bar = DELTA; bar != local_end_of_types; bar++) {
0368             types type = static_cast<types>(bar);
0369             std::string theLabel = getStringFromTypeEnum(type) + getStringFromCoordEnum(coord);
0370             if (xBin == 1) {
0371               Thresholds->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0372               ThresholdsColor->GetYaxis()->SetBinLabel(yBin, theLabel.c_str());
0373             }
0374 
0375             const auto& value = cutFunctor(type, alignable, coord);
0376             Thresholds->SetBinContent(xBin, yBin, value);
0377             ThresholdsColor->SetBinContent(xBin, yBin, value);
0378 
0379             yBin--;
0380           }  // loop on types
0381         }    // loop on coordinates
0382       }      // loop on alignables
0383 
0384       ThresholdsColor->Draw("COLZ0");
0385       ThresholdsColor->GetXaxis()->LabelsOption(isHighGranularity_ ? "v" : "h");
0386       Thresholds->Draw("TEXTsame");
0387       Thresholds->GetXaxis()->LabelsOption(isHighGranularity_ ? "v" : "h");
0388 
0389       auto ltx = TLatex();
0390       ltx.SetTextFont(62);
0391       //ltx.SetTextColor(kBlue);
0392       ltx.SetTextSize(0.047);
0393       ltx.SetTextAlign(11);
0394       std::string ltxText;
0395       if (this->m_plotAnnotations.ntags == 2) {
0396         ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
0397                                f_tagname,
0398                                std::to_string(std::get<0>(firstiov)),
0399                                l_tagname,
0400                                std::to_string(std::get<0>(lastiov)));
0401       } else {
0402         ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
0403                                f_tagname,
0404                                std::to_string(std::get<0>(firstiov)),
0405                                std::to_string(std::get<0>(lastiov)));
0406       }
0407       ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
0408 
0409       std::string fileName(this->m_imageFileName);
0410       canvas.SaveAs(fileName.c_str());
0411 
0412       return true;
0413     }
0414 
0415   private:
0416     template <typename T>
0417     bool isEqual(std::vector<T> const& v1, std::vector<T> const& v2) {
0418       return (v1.size() == v2.size() && std::equal(v1.begin(), v1.end(), v2.begin()));
0419     }
0420     bool isHighGranularity_;
0421     std::string label_;
0422   };
0423 }  // namespace AlignPCLThresholdPlotHelper
0424 
0425 #endif