Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-10 23:04:48

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