Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-22 04:02:31

0001 /**
0002  * @package   Alignment/MillePedeAlignmentAlgorithm
0003  * @file      MillePedeDQMModule.cc
0004  *
0005  * @author    Max Stark (max.stark@cern.ch)
0006  * @date      Feb 19, 2016
0007  */
0008 
0009 /*** header-file ***/
0010 #include "Alignment/MillePedeAlignmentAlgorithm/plugins/MillePedeDQMModule.h"
0011 
0012 /*** ROOT objects ***/
0013 #include "TH1F.h"
0014 
0015 /*** Core framework functionality ***/
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 
0019 /*** Alignment ***/
0020 #include "Alignment/MillePedeAlignmentAlgorithm/interface/PedeLabelerBase.h"
0021 #include "Alignment/MillePedeAlignmentAlgorithm/interface/PedeLabelerPluginFactory.h"
0022 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
0023 
0024 /*** Necessary Framework infrastructure ***/
0025 #include "FWCore/Framework/interface/ProcessBlock.h"
0026 #include "DataFormats/Alignment/interface/AlignmentToken.h"
0027 
0028 MillePedeDQMModule ::MillePedeDQMModule(const edm::ParameterSet& config)
0029     : tTopoToken_(esConsumes<edm::Transition::BeginRun>()),
0030       gDetToken_(esConsumes<edm::Transition::BeginRun>()),
0031       ptpToken_(esConsumes<edm::Transition::BeginRun>()),
0032       ptitpToken_(esConsumes<edm::Transition::BeginRun>()),
0033       aliThrToken_(esConsumes<edm::Transition::BeginRun>()),
0034       geomToken_(esConsumes<edm::Transition::BeginRun>()),
0035       outputFolder_(config.getParameter<std::string>("outputFolder")),
0036       mpReaderConfig_(config.getParameter<edm::ParameterSet>("MillePedeFileReader")),
0037       isHG_(mpReaderConfig_.getParameter<bool>("isHG")) {
0038   consumes<AlignmentToken, edm::InProcess>(config.getParameter<edm::InputTag>("alignmentTokenSrc"));
0039 }
0040 
0041 MillePedeDQMModule ::~MillePedeDQMModule() {}
0042 
0043 //=============================================================================
0044 //===   INTERFACE IMPLEMENTATION                                            ===
0045 //=============================================================================
0046 
0047 void MillePedeDQMModule ::bookHistograms(DQMStore::IBooker& booker) {
0048   edm::LogInfo("MillePedeDQMModule") << "Booking histograms";
0049 
0050   booker.cd();
0051   if (!isHG_) {
0052     if (outputFolder_.find("HG") != std::string::npos) {
0053       throw cms::Exception("LogicError")
0054           << "MillePedeDQMModule is configured as Low Granularity but the outputfolder is for High Granularity";
0055     }
0056 
0057     booker.setCurrentFolder(outputFolder_);
0058     h_xPos = booker.book1D("Xpos", "Alignment fit #DeltaX;;#mum", 36, 0., 36.);
0059     h_xRot = booker.book1D("Xrot", "Alignment fit #Delta#theta_{X};;#murad", 36, 0., 36.);
0060     h_yPos = booker.book1D("Ypos", "Alignment fit #DeltaY;;#mum", 36, 0., 36.);
0061     h_yRot = booker.book1D("Yrot", "Alignment fit #Delta#theta_{Y};;#murad", 36, 0., 36.);
0062     h_zPos = booker.book1D("Zpos", "Alignment fit #DeltaZ;;#mum", 36, 0., 36.);
0063     h_zRot = booker.book1D("Zrot", "Alignment fit #Delta#theta_{Z};;#murad", 36, 0., 36.);
0064     statusResults = booker.book2D("statusResults", "Status of SiPixelAli PCL workflow;;", 6, 0., 6., 1, 0., 1.);
0065   } else {
0066     if (outputFolder_.find("HG") == std::string::npos) {
0067       throw cms::Exception("LogicError")
0068           << "MillePedeDQMModule is configured as High Granularity but the outputfolder is for Low Granularity";
0069     }
0070 
0071     booker.setCurrentFolder(outputFolder_);
0072 
0073     layerVec = {{"Layer1", pixelTopologyMap_->getPXBLadders(1)},
0074                 {"Layer2", pixelTopologyMap_->getPXBLadders(2)},
0075                 {"Layer3", pixelTopologyMap_->getPXBLadders(3)},
0076                 {"Layer4", pixelTopologyMap_->getPXBLadders(4)},
0077                 {"Disk-3", pixelTopologyMap_->getPXFBlades(-3) * 2},
0078                 {"Disk-2", pixelTopologyMap_->getPXFBlades(-2) * 2},
0079                 {"Disk-1", pixelTopologyMap_->getPXFBlades(-1) * 2},
0080                 {"Disk1", pixelTopologyMap_->getPXFBlades(1) * 2},
0081                 {"Disk2", pixelTopologyMap_->getPXFBlades(2) * 2},
0082                 {"Disk3", pixelTopologyMap_->getPXFBlades(3) * 2}};
0083 
0084     for (const auto& layer : layerVec) {
0085       h_xPos_HG[layer.first] = booker.book1D("Xpos_HG_" + layer.first,
0086                                              "Alignment fit #DeltaX for " + layer.first + ";;#mum",
0087                                              layer.second + 5,
0088                                              0.,
0089                                              layer.second + 5);
0090       h_xRot_HG[layer.first] = booker.book1D("Xrot_HG_" + layer.first,
0091                                              "Alignment fit #Delta#theta_{X} for " + layer.first + ";;#murad",
0092                                              layer.second + 5,
0093                                              0.,
0094                                              layer.second + 5);
0095       h_yPos_HG[layer.first] = booker.book1D("Ypos_HG_" + layer.first,
0096                                              "Alignment fit #DeltaY for " + layer.first + ";;#mum",
0097                                              layer.second + 5,
0098                                              0.,
0099                                              layer.second + 5);
0100       h_yRot_HG[layer.first] = booker.book1D("Yrot_HG_" + layer.first,
0101                                              "Alignment fit #Delta#theta_{Y} for " + layer.first + ";;#murad",
0102                                              layer.second + 5,
0103                                              0.,
0104                                              layer.second + 5);
0105       h_zPos_HG[layer.first] = booker.book1D("Zpos_HG_" + layer.first,
0106                                              "Alignment fit #DeltaZ for " + layer.first + ";;#mum",
0107                                              layer.second + 5,
0108                                              0.,
0109                                              layer.second + 5);
0110       h_zRot_HG[layer.first] = booker.book1D("Zrot_HG_" + layer.first,
0111                                              "Alignment fit #Delta#theta_{Z} for " + layer.first + ";;#murad",
0112                                              layer.second + 5,
0113                                              0.,
0114                                              layer.second + 5);
0115     }
0116 
0117     statusResults =
0118         booker.book2D("statusResults", "Fraction threshold check for SiPixelAliHG PCL;;", 6, 0., 6., 10, 0., 10.);
0119   }
0120 
0121   binariesAvalaible = booker.bookInt("BinariesFound");
0122   exitCode = booker.bookString("PedeExitCode", "");
0123   isVetoed = booker.bookString("IsVetoed", "");
0124 
0125   booker.cd();
0126 }
0127 
0128 void MillePedeDQMModule ::dqmEndJob(DQMStore::IBooker& booker, DQMStore::IGetter&) {
0129   bookHistograms(booker);
0130   if (mpReader_) {
0131     mpReader_->read();
0132   } else {
0133     throw cms::Exception("LogicError") << "@SUB=MillePedeDQMModule::dqmEndJob\n"
0134                                        << "Try to read MillePede results before initializing MillePedeFileReader";
0135   }
0136   if (!isHG_) {
0137     fillExpertHistos();
0138     fillStatusHisto(statusResults);
0139   } else {
0140     fillExpertHistos_HG();
0141     fillStatusHistoHG(statusResults);
0142   }
0143   binariesAvalaible->Fill(mpReader_->binariesAmount());
0144   auto theResults = mpReader_->getResults();
0145   std::string exitCodeStr = theResults.getExitMessage();
0146 
0147   std::string vetoStr{};
0148   if (mpReader_->storeAlignments()) {
0149     vetoStr = "DB Updated!"; /* easy peasy, fait accompli an alignment is there */
0150   } else {
0151     if (theResults.isHighGranularity()) { /* HG case */
0152       if (theResults.getDBVetoed() && theResults.getDBUpdated()) {
0153         vetoStr = "DB Update Vetoed"; /* this can happen in the HG PCL case */
0154       } else {
0155         vetoStr = "N/A";
0156       }
0157     } else { /* LG case */
0158       if (theResults.exceedsCutoffs()) {
0159         vetoStr = "DB Update Vetoed"; /* this can happen in the LG PCL case */
0160       } else {
0161         vetoStr = "N/A";
0162       }  // if the alignment exceeds the cutoffs
0163     }    // LG case
0164   }      // if the alignment was not stored
0165 
0166   exitCode->Fill(exitCodeStr);
0167   isVetoed->Fill(vetoStr);
0168 }
0169 
0170 //=============================================================================
0171 //===   PRIVATE METHOD IMPLEMENTATION                                       ===
0172 //=============================================================================
0173 
0174 void MillePedeDQMModule ::beginRun(const edm::Run&, const edm::EventSetup& setup) {
0175   if (!setupChanged(setup))
0176     return;
0177 
0178   const TrackerTopology* const tTopo = &setup.getData(tTopoToken_);
0179   const GeometricDet* geometricDet = &setup.getData(gDetToken_);
0180   const PTrackerParameters* ptp = &setup.getData(ptpToken_);
0181   const PTrackerAdditionalParametersPerDet* ptitp = &setup.getData(ptitpToken_);
0182   const TrackerGeometry* geom = &setup.getData(geomToken_);
0183 
0184   pixelTopologyMap_ = std::make_shared<PixelTopologyMap>(geom, tTopo);
0185 
0186   // take the thresholds from DB
0187   const auto& thresholds_ = &setup.getData(aliThrToken_);
0188 
0189   auto myThresholds = std::make_shared<AlignPCLThresholdsHG>();
0190   myThresholds->setAlignPCLThresholds(thresholds_->getNrecords(), thresholds_->getThreshold_Map());
0191   myThresholds->setFloatMap(thresholds_->getFloatMap());
0192 
0193   TrackerGeomBuilderFromGeometricDet builder;
0194 
0195   const auto trackerGeometry = builder.build(geometricDet, ptitp, *ptp, tTopo);
0196   tracker_ = std::make_unique<AlignableTracker>(trackerGeometry, tTopo);
0197 
0198   const std::string labelerPlugin{"PedeLabeler"};
0199   edm::ParameterSet labelerConfig{};
0200   labelerConfig.addUntrackedParameter("plugin", labelerPlugin);
0201   labelerConfig.addUntrackedParameter("RunRangeSelection", edm::VParameterSet{});
0202 
0203   std::shared_ptr<PedeLabelerBase> pedeLabeler{PedeLabelerPluginFactory::get()->create(
0204       labelerPlugin, PedeLabelerBase::TopLevelAlignables(tracker_.get(), nullptr, nullptr), labelerConfig)};
0205 
0206   mpReader_ = std::make_unique<MillePedeFileReader>(
0207       mpReaderConfig_, pedeLabeler, std::shared_ptr<const AlignPCLThresholdsHG>(myThresholds), pixelTopologyMap_);
0208 }
0209 
0210 void MillePedeDQMModule ::fillStatusHisto(MonitorElement* statusHisto) {
0211   TH2F* histo_status = statusHisto->getTH2F();
0212   auto theResults = mpReader_->getResults();
0213   theResults.print();
0214   histo_status->SetBinContent(1, 1, theResults.getDBUpdated());
0215   histo_status->GetXaxis()->SetBinLabel(1, "DB updated");
0216   histo_status->SetBinContent(2, 1, theResults.exceedsCutoffs());
0217   histo_status->GetXaxis()->SetBinLabel(2, "significant movement");
0218   histo_status->SetBinContent(3, 1, theResults.getDBVetoed());
0219   histo_status->GetXaxis()->SetBinLabel(3, "DB update vetoed");
0220   histo_status->SetBinContent(4, 1, !theResults.exceedsThresholds());
0221   histo_status->GetXaxis()->SetBinLabel(4, "within max movement");
0222   histo_status->SetBinContent(5, 1, !theResults.exceedsMaxError());
0223   histo_status->GetXaxis()->SetBinLabel(5, "within max error");
0224   histo_status->SetBinContent(6, 1, !theResults.belowSignificance());
0225   histo_status->GetXaxis()->SetBinLabel(6, "above significance");
0226 }
0227 
0228 void MillePedeDQMModule ::fillStatusHistoHG(MonitorElement* statusHisto) {
0229   TH2F* histo_status = statusHisto->getTH2F();
0230   auto& theResults = mpReader_->getResultsHG();
0231   histo_status->GetXaxis()->SetBinLabel(1, "#DeltaX");
0232   histo_status->GetXaxis()->SetBinLabel(2, "#DeltaY");
0233   histo_status->GetXaxis()->SetBinLabel(3, "#DeltaZ");
0234   histo_status->GetXaxis()->SetBinLabel(4, "#Delta#theta_{X}");
0235   histo_status->GetXaxis()->SetBinLabel(5, "#Delta#theta_{Y}");
0236   histo_status->GetXaxis()->SetBinLabel(6, "#Delta#theta_{Z}");
0237 
0238   int i = 0;
0239   for (const auto& result : theResults) {
0240     histo_status->GetYaxis()->SetBinLabel(i + 1, result.first.data());
0241     for (std::size_t j = 0; j < result.second.size(); ++j) {
0242       histo_status->SetBinContent(j + 1, i + 1, result.second[j]);
0243     }
0244     i++;
0245   }
0246 }
0247 
0248 void MillePedeDQMModule ::fillExpertHistos() {
0249   std::array<double, SIZE_INDEX> Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_;
0250   std::array<double, SIZE_INDEX> tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_;
0251 
0252   std::array<double, SIZE_INDEX> Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_;
0253   std::array<double, SIZE_INDEX> tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_;
0254 
0255   std::array<double, SIZE_INDEX> Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_;
0256   std::array<double, SIZE_INDEX> tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_;
0257 
0258   auto myMap = mpReader_->getThresholdMap();
0259 
0260   std::vector<std::string> alignablesList;
0261   for (auto it = myMap.begin(); it != myMap.end(); ++it) {
0262     alignablesList.push_back(it->first);
0263   }
0264 
0265   for (auto& alignable : alignablesList) {
0266     int detIndex = getIndexFromString(alignable);
0267 
0268     Xcut_[detIndex] = myMap[alignable].getXcut();
0269     sigXcut_[detIndex] = myMap[alignable].getSigXcut();
0270     maxMoveXcut_[detIndex] = myMap[alignable].getMaxMoveXcut();
0271     maxErrorXcut_[detIndex] = myMap[alignable].getErrorXcut();
0272 
0273     Ycut_[detIndex] = myMap[alignable].getYcut();
0274     sigYcut_[detIndex] = myMap[alignable].getSigYcut();
0275     maxMoveYcut_[detIndex] = myMap[alignable].getMaxMoveYcut();
0276     maxErrorYcut_[detIndex] = myMap[alignable].getErrorYcut();
0277 
0278     Zcut_[detIndex] = myMap[alignable].getZcut();
0279     sigZcut_[detIndex] = myMap[alignable].getSigZcut();
0280     maxMoveZcut_[detIndex] = myMap[alignable].getMaxMoveZcut();
0281     maxErrorZcut_[detIndex] = myMap[alignable].getErrorZcut();
0282 
0283     tXcut_[detIndex] = myMap[alignable].getThetaXcut();
0284     sigtXcut_[detIndex] = myMap[alignable].getSigThetaXcut();
0285     maxMovetXcut_[detIndex] = myMap[alignable].getMaxMoveThetaXcut();
0286     maxErrortXcut_[detIndex] = myMap[alignable].getErrorThetaXcut();
0287 
0288     tYcut_[detIndex] = myMap[alignable].getThetaYcut();
0289     sigtYcut_[detIndex] = myMap[alignable].getSigThetaYcut();
0290     maxMovetYcut_[detIndex] = myMap[alignable].getMaxMoveThetaYcut();
0291     maxErrortYcut_[detIndex] = myMap[alignable].getErrorThetaYcut();
0292 
0293     tZcut_[detIndex] = myMap[alignable].getThetaZcut();
0294     sigtZcut_[detIndex] = myMap[alignable].getSigThetaZcut();
0295     maxMovetZcut_[detIndex] = myMap[alignable].getMaxMoveThetaZcut();
0296     maxErrortZcut_[detIndex] = myMap[alignable].getErrorThetaZcut();
0297   }
0298 
0299   fillExpertHisto(h_xPos, Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_, mpReader_->getXobs(), mpReader_->getXobsErr());
0300   fillExpertHisto(
0301       h_xRot, tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_, mpReader_->getTXobs(), mpReader_->getTXobsErr());
0302 
0303   fillExpertHisto(h_yPos, Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_, mpReader_->getYobs(), mpReader_->getYobsErr());
0304   fillExpertHisto(
0305       h_yRot, tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_, mpReader_->getTYobs(), mpReader_->getTYobsErr());
0306 
0307   fillExpertHisto(h_zPos, Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_, mpReader_->getZobs(), mpReader_->getZobsErr());
0308   fillExpertHisto(
0309       h_zRot, tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_, mpReader_->getTZobs(), mpReader_->getTZobsErr());
0310 }
0311 
0312 void MillePedeDQMModule ::fillExpertHisto(MonitorElement* histo,
0313                                           const std::array<double, SIZE_INDEX>& cut,
0314                                           const std::array<double, SIZE_INDEX>& sigCut,
0315                                           const std::array<double, SIZE_INDEX>& maxMoveCut,
0316                                           const std::array<double, SIZE_INDEX>& maxErrorCut,
0317                                           const std::array<double, SIZE_LG_STRUCTS>& obs,
0318                                           const std::array<double, SIZE_LG_STRUCTS>& obsErr) {
0319   TH1F* histo_0 = histo->getTH1F();
0320 
0321   double max_ = *std::max_element(maxMoveCut.begin(), maxMoveCut.end());
0322 
0323   histo_0->SetMinimum(-(max_));
0324   histo_0->SetMaximum(max_);
0325 
0326   //  Schematics of the bin contents
0327   //
0328   //  XX XX XX XX XX XX    OO OO OO OO    II II II II
0329   // |--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|--|
0330   // | 1| 2| 3| 4| 5| 6| 7| 8| 9|10|11|12|13|14|15|16|17|  ...
0331   //
0332   // |-----------------|  |-----------|  |-----------|
0333   // |observed movement|  |thresholds1|  |thresholds2|
0334 
0335   for (size_t i = 0; i < obs.size(); ++i) {
0336     // fist obs.size() bins for observed movements
0337     histo_0->SetBinContent(i + 1, obs[i]);
0338     histo_0->SetBinError(i + 1, obsErr[i]);
0339 
0340     // then at bin 8,8+5,8+10,... for cutoffs
0341     // 5 bins is the space allocated for the 4 other thresholds + 1 empty separation bin
0342     histo_0->SetBinContent(8 + i * 5, cut[i]);
0343 
0344     // then at bin 9,9+5,9+10,... for significances
0345     histo_0->SetBinContent(9 + i * 5, sigCut[i]);
0346 
0347     // then at bin 10,10+5,10+10,... for maximum movements
0348     histo_0->SetBinContent(10 + i * 5, maxMoveCut[i]);
0349 
0350     // then at bin 11,11+5,11+10,... for maximum errors
0351     histo_0->SetBinContent(11 + i * 5, maxErrorCut[i]);
0352   }
0353 }
0354 
0355 void MillePedeDQMModule ::fillExpertHistos_HG() {
0356   std::array<double, SIZE_INDEX> Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_;
0357   std::array<double, SIZE_INDEX> tXcut_, sigtXcut_, maxMovetXcut_, maxErrortXcut_;
0358 
0359   std::array<double, SIZE_INDEX> Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_;
0360   std::array<double, SIZE_INDEX> tYcut_, sigtYcut_, maxMovetYcut_, maxErrortYcut_;
0361 
0362   std::array<double, SIZE_INDEX> Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_;
0363   std::array<double, SIZE_INDEX> tZcut_, sigtZcut_, maxMovetZcut_, maxErrortZcut_;
0364 
0365   auto myMap = mpReader_->getThresholdMap();
0366 
0367   std::vector<std::string> alignablesList;
0368   for (auto it = myMap.begin(); it != myMap.end(); ++it) {
0369     alignablesList.push_back(it->first);
0370   }
0371 
0372   for (auto& alignable : alignablesList) {
0373     int detIndex = getIndexFromString(alignable);
0374 
0375     Xcut_[detIndex] = myMap[alignable].getXcut();
0376     sigXcut_[detIndex] = myMap[alignable].getSigXcut();
0377     maxMoveXcut_[detIndex] = myMap[alignable].getMaxMoveXcut();
0378     maxErrorXcut_[detIndex] = myMap[alignable].getErrorXcut();
0379 
0380     Ycut_[detIndex] = myMap[alignable].getYcut();
0381     sigYcut_[detIndex] = myMap[alignable].getSigYcut();
0382     maxMoveYcut_[detIndex] = myMap[alignable].getMaxMoveYcut();
0383     maxErrorYcut_[detIndex] = myMap[alignable].getErrorYcut();
0384 
0385     Zcut_[detIndex] = myMap[alignable].getZcut();
0386     sigZcut_[detIndex] = myMap[alignable].getSigZcut();
0387     maxMoveZcut_[detIndex] = myMap[alignable].getMaxMoveZcut();
0388     maxErrorZcut_[detIndex] = myMap[alignable].getErrorZcut();
0389 
0390     tXcut_[detIndex] = myMap[alignable].getThetaXcut();
0391     sigtXcut_[detIndex] = myMap[alignable].getSigThetaXcut();
0392     maxMovetXcut_[detIndex] = myMap[alignable].getMaxMoveThetaXcut();
0393     maxErrortXcut_[detIndex] = myMap[alignable].getErrorThetaXcut();
0394 
0395     tYcut_[detIndex] = myMap[alignable].getThetaYcut();
0396     sigtYcut_[detIndex] = myMap[alignable].getSigThetaYcut();
0397     maxMovetYcut_[detIndex] = myMap[alignable].getMaxMoveThetaYcut();
0398     maxErrortYcut_[detIndex] = myMap[alignable].getErrorThetaYcut();
0399 
0400     tZcut_[detIndex] = myMap[alignable].getThetaZcut();
0401     sigtZcut_[detIndex] = myMap[alignable].getSigThetaZcut();
0402     maxMovetZcut_[detIndex] = myMap[alignable].getMaxMoveThetaZcut();
0403     maxErrortZcut_[detIndex] = myMap[alignable].getErrorThetaZcut();
0404   }
0405 
0406   fillExpertHisto_HG(
0407       h_xPos_HG, Xcut_, sigXcut_, maxMoveXcut_, maxErrorXcut_, mpReader_->getXobs_HG(), mpReader_->getXobsErr_HG());
0408   fillExpertHisto_HG(h_xRot_HG,
0409                      tXcut_,
0410                      sigtXcut_,
0411                      maxMovetXcut_,
0412                      maxErrortXcut_,
0413                      mpReader_->getTXobs_HG(),
0414                      mpReader_->getTXobsErr_HG());
0415 
0416   fillExpertHisto_HG(
0417       h_yPos_HG, Ycut_, sigYcut_, maxMoveYcut_, maxErrorYcut_, mpReader_->getYobs_HG(), mpReader_->getYobsErr_HG());
0418   fillExpertHisto_HG(h_yRot_HG,
0419                      tYcut_,
0420                      sigtYcut_,
0421                      maxMovetYcut_,
0422                      maxErrortYcut_,
0423                      mpReader_->getTYobs_HG(),
0424                      mpReader_->getTYobsErr_HG());
0425 
0426   fillExpertHisto_HG(
0427       h_zPos_HG, Zcut_, sigZcut_, maxMoveZcut_, maxErrorZcut_, mpReader_->getZobs_HG(), mpReader_->getZobsErr_HG());
0428   fillExpertHisto_HG(h_zRot_HG,
0429                      tZcut_,
0430                      sigtZcut_,
0431                      maxMovetZcut_,
0432                      maxErrortZcut_,
0433                      mpReader_->getTZobs_HG(),
0434                      mpReader_->getTZobsErr_HG());
0435 }
0436 
0437 void MillePedeDQMModule ::fillExpertHisto_HG(std::map<std::string, MonitorElement*>& histo_map,
0438                                              const std::array<double, SIZE_INDEX>& cut,
0439                                              const std::array<double, SIZE_INDEX>& sigCut,
0440                                              const std::array<double, SIZE_INDEX>& maxMoveCut,
0441                                              const std::array<double, SIZE_INDEX>& maxErrorCut,
0442                                              const std::array<double, SIZE_HG_STRUCTS>& obs,
0443                                              const std::array<double, SIZE_HG_STRUCTS>& obsErr) {
0444   int currentStart = 0;
0445   int bin = 0;
0446   double max_ = 0;
0447 
0448   for (const auto& layer : layerVec) {
0449     TH1F* histo_0 = histo_map[layer.first]->getTH1F();
0450 
0451     max_ = -1;
0452     for (int i = currentStart; i < (currentStart + layer.second); ++i) {
0453       // first obs.size() bins for observed movements
0454       bin = i - currentStart + 1;
0455 
0456       // fill observed values
0457       histo_0->SetBinContent(bin, obs[i]);
0458       histo_0->SetBinError(bin, obsErr[i]);
0459 
0460       if (std::abs(obs[i]) > max_) {
0461         max_ = std::abs(obs[i]);
0462       }
0463     }
0464 
0465     // five extra bins at the end, one empty, one with threshold, one with sigCut, one with maxMoveCut, one with MaxErrorCut
0466     histo_0->SetBinContent(bin + 1, 0);
0467     histo_0->SetBinError(bin + 1, 0);
0468 
0469     int detIndex;
0470     if (layer.first.find("Disk") != std::string::npos) {
0471       // 7 is the detId for panels, see getIndexFromString
0472       detIndex = 7;
0473       histo_0->GetXaxis()->SetTitle("Panel");
0474     } else {
0475       // 6 is the detId for ladders, see getIndexFromString
0476       detIndex = 6;
0477       histo_0->GetXaxis()->SetTitle("Ladder");
0478     }
0479 
0480     histo_0->SetBinContent(bin + 2, cut[detIndex]);
0481     histo_0->SetBinError(bin + 2, 0);
0482     histo_0->SetBinContent(bin + 3, sigCut[detIndex]);
0483     histo_0->SetBinError(bin + 3, 0);
0484     histo_0->SetBinContent(bin + 4, maxMoveCut[detIndex]);
0485     histo_0->SetBinError(bin + 4, 0);
0486     histo_0->SetBinContent(bin + 5, maxErrorCut[detIndex]);
0487     histo_0->SetBinError(bin + 5, 0);
0488 
0489     // always scale so the cutoff is visible
0490     max_ = std::max(cut[detIndex] * 1.2, max_);
0491 
0492     histo_0->SetMinimum(-(max_) * 1.2);
0493     histo_0->SetMaximum(max_ * 1.2);
0494 
0495     currentStart += layer.second;
0496   }
0497 }
0498 
0499 bool MillePedeDQMModule ::setupChanged(const edm::EventSetup& setup) {
0500   bool changed{false};
0501 
0502   if (watchIdealGeometryRcd_.check(setup))
0503     changed = true;
0504   if (watchTrackerTopologyRcd_.check(setup))
0505     changed = true;
0506   if (watchPTrackerParametersRcd_.check(setup))
0507     changed = true;
0508 
0509   return changed;
0510 }
0511 
0512 int MillePedeDQMModule ::getIndexFromString(const std::string& alignableId) {
0513   if (alignableId == "TPBHalfBarrelXminus") {
0514     return 3;
0515   } else if (alignableId == "TPBHalfBarrelXplus") {
0516     return 2;
0517   } else if (alignableId == "TPEHalfCylinderXminusZminus") {
0518     return 1;
0519   } else if (alignableId == "TPEHalfCylinderXplusZminus") {
0520     return 0;
0521   } else if (alignableId == "TPEHalfCylinderXminusZplus") {
0522     return 5;
0523   } else if (alignableId == "TPEHalfCylinderXplusZplus") {
0524     return 4;
0525   } else if (alignableId.rfind("TPBLadder", 0) == 0) {
0526     return 6;
0527   } else if (alignableId.rfind("TPEPanel", 0) == 0) {
0528     return 7;
0529   } else {
0530     throw cms::Exception("LogicError") << "@SUB=MillePedeDQMModule::getIndexFromString\n"
0531                                        << "Retrieving conversion for not supported Alignable partition" << alignableId;
0532   }
0533 }