Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:06

0001 
0002 
0003 /*
0004  *  See header file for a description of this class.
0005  *
0006  *  \author G. Mila - INFN Torino
0007  *
0008  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
0009  *
0010  */
0011 
0012 #include "DQM/DTMonitorClient/src/DTResolutionTest.h"
0013 
0014 // Framework
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "FWCore/Framework/interface/ESHandle.h"
0018 #include "FWCore/Framework/interface/EventSetup.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 
0021 // Geometry
0022 #include "Geometry/DTGeometry/interface/DTChamber.h"
0023 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0024 
0025 #include "DQMServices/Core/interface/DQMStore.h"
0026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0027 
0028 #include <iostream>
0029 #include <cstdio>
0030 #include <string>
0031 #include <sstream>
0032 #include <cmath>
0033 #include "TF1.h"
0034 
0035 using namespace edm;
0036 using namespace std;
0037 
0038 DTResolutionTest::DTResolutionTest(const edm::ParameterSet& ps)
0039     : muonGeomToken_(esConsumes<edm::Transition::EndLuminosityBlock>()) {
0040   edm::LogVerbatim("resolution") << "[DTResolutionTest]: Constructor";
0041   parameters = ps;
0042 
0043   //FR: no idea if this input file needs to be used! comment it for now
0044   //-  if(ps.getUntrackedParameter<bool>("readFile", false))
0045   //-   dbe->open(ps.getUntrackedParameter<string>("inputFile", "residuals.root"));
0046 
0047   prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
0048 
0049   percentual = parameters.getUntrackedParameter<int>("BadSLpercentual", 10);
0050 
0051   nevents = 0;
0052 
0053   bookingdone = false;
0054 }
0055 
0056 DTResolutionTest::~DTResolutionTest() {
0057   edm::LogVerbatim("resolution") << "DTResolutionTest: analyzed " << nevents << " events";
0058 }
0059 
0060 void DTResolutionTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0061                                              DQMStore::IGetter& igetter,
0062                                              edm::LuminosityBlock const& lumiSeg,
0063                                              edm::EventSetup const& context) {
0064   // counts number of updates (online mode) or number of events (standalone mode)
0065   //nevents++;
0066   // if running in standalone perform diagnostic only after a reasonalbe amount of events
0067   //if ( parameters.getUntrackedParameter<bool>("runningStandalone", false) &&
0068   //   nevents%parameters.getUntrackedParameter<int>("diagnosticPrescale", 1000) != 0 ) return;
0069   //edm::LogVerbatim ("resolution") << "[DTResolutionTest]: "<<nevents<<" updates";
0070 
0071   if (!bookingdone) {
0072     // Get the geometry
0073     muonGeom = &context.getData(muonGeomToken_);
0074     // book the histos
0075     for (int wheel = -2; wheel < 3; wheel++) {
0076       bookHistos(ibooker, wheel);
0077     }
0078     vector<const DTChamber*> chambers = muonGeom->chambers();
0079     for (vector<const DTChamber*>::const_iterator chamber = chambers.begin(); chamber != chambers.end(); ++chamber) {
0080       bookHistos(ibooker, (*chamber)->id());
0081     }
0082   }
0083   bookingdone = true;
0084 
0085   edm::LogVerbatim("resolution") << "[DTResolutionTest]: End of LS transition, performing the DQM client operation";
0086 
0087   // counts number of lumiSegs
0088   nLumiSegs = lumiSeg.id().luminosityBlock();
0089 
0090   // prescale factor
0091   if (nLumiSegs % prescaleFactor != 0)
0092     return;
0093 
0094   for (map<int, MonitorElement*>::const_iterator histo = wheelMeanHistos.begin(); histo != wheelMeanHistos.end();
0095        histo++) {
0096     (*histo).second->Reset();
0097   }
0098   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0099     for (map<int, MonitorElement*>::const_iterator histo = wheelSigmaHistos.begin(); histo != wheelSigmaHistos.end();
0100          histo++) {
0101       (*histo).second->Reset();
0102     }
0103   }
0104   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0105     for (map<int, MonitorElement*>::const_iterator histo = wheelSlopeHistos.begin(); histo != wheelSlopeHistos.end();
0106          histo++) {
0107       (*histo).second->Reset();
0108     }
0109   }
0110 
0111   cmsMeanHistos.clear();
0112   for (int i = -2; i < 3; i++) {
0113     for (int j = 1; j < 15; j++) {
0114       MeanFilled[make_pair(i, j)] = false;
0115     }
0116   }
0117   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0118     cmsSigmaHistos.clear();
0119     for (int i = -2; i < 3; i++) {
0120       for (int j = 1; j < 15; j++) {
0121         SigmaFilled[make_pair(i, j)] = false;
0122       }
0123     }
0124   }
0125   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0126     cmsSlopeHistos.clear();
0127     for (int i = -2; i < 3; i++) {
0128       for (int j = 1; j < 15; j++) {
0129         SlopeFilled[make_pair(i, j)] = false;
0130       }
0131     }
0132   }
0133 
0134   // summary histos initialization
0135   for (int wh = -2; wh <= 3; wh++) {
0136     if (wh != 3) {
0137       for (int xBin = 0; xBin < 14; xBin++) {
0138         for (int yBin = 0; yBin < 11; yBin++) {
0139           wheelMeanHistos[wh]->setBinContent(xBin, yBin, 0);
0140           if (parameters.getUntrackedParameter<bool>("sigmaTest"))
0141             wheelSigmaHistos[wh]->setBinContent(xBin, yBin, 0);
0142           if (parameters.getUntrackedParameter<bool>("slopeTest"))
0143             wheelSlopeHistos[wh]->setBinContent(xBin, yBin, 0);
0144         }
0145       }
0146     } else {
0147       for (int xBin = 0; xBin < 14; xBin++) {
0148         for (int yBin = -2; yBin < 3; yBin++) {
0149           wheelMeanHistos[wh]->setBinContent(xBin, yBin, 0);
0150           if (parameters.getUntrackedParameter<bool>("sigmaTest"))
0151             wheelSigmaHistos[wh]->setBinContent(xBin, yBin, 0);
0152           if (parameters.getUntrackedParameter<bool>("slopeTest"))
0153             wheelSlopeHistos[wh]->setBinContent(xBin, yBin, 0);
0154         }
0155       }
0156     }
0157   }
0158 
0159   edm::LogVerbatim("resolution") << "[DTResolutionTest]: " << nLumiSegs << " updates";
0160 
0161   vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0162   vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0163 
0164   edm::LogVerbatim("resolution") << "[DTResolutionTest]: Residual Distribution tests results";
0165 
0166   for (; ch_it != ch_end; ++ch_it) {
0167     DTChamberId chID = (*ch_it)->id();
0168 
0169     // Fill the test histos
0170     int entry = -1;
0171     if (chID.station() == 1)
0172       entry = 0;
0173     if (chID.station() == 2)
0174       entry = 3;
0175     if (chID.station() == 3)
0176       entry = 6;
0177     if (chID.station() == 4)
0178       entry = 9;
0179 
0180     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
0181     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
0182 
0183     for (; sl_it != sl_end; ++sl_it) {
0184       DTSuperLayerId slID = (*sl_it)->id();
0185 
0186       edm::LogVerbatim("resolution") << "[DTResolutionTest]: Superlayer: " << slID;
0187 
0188       stringstream wheel;
0189       wheel << slID.wheel();
0190       stringstream station;
0191       station << slID.station();
0192       stringstream sector;
0193       sector << slID.sector();
0194       stringstream superLayer;
0195       superLayer << slID.superlayer();
0196 
0197       string HistoName = "W" + wheel.str() + "_Sec" + sector.str();
0198       string supLayer = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str();
0199 
0200       MonitorElement* res_histo = igetter.get(getMEName(slID));
0201       if (res_histo) {
0202         // gaussian test
0203         string GaussianCriterionName =
0204             parameters.getUntrackedParameter<string>("resDistributionTestName", "ResidualsDistributionGaussianTest");
0205         const QReport* GaussianReport = res_histo->getQReport(GaussianCriterionName);
0206         if (GaussianReport) {
0207           // FIXE ME: if the quality test fails this cout returns a null pointer
0208           //edm::LogWarning ("resolution") << "-------- SuperLayer : "<<supLayer<<"  "
0209           //                               <<GaussianReport->getMessage()<<" ------- "<<GaussianReport->getStatus();
0210         }
0211         int BinNumber = entry + slID.superLayer();
0212         if (BinNumber == 12)
0213           BinNumber = 11;
0214         float mean = (*res_histo).getMean(1);
0215         float sigma = (*res_histo).getRMS(1);
0216         MeanHistos.find(make_pair(slID.wheel(), slID.sector()))->second->setBinContent(BinNumber, mean);
0217         if (parameters.getUntrackedParameter<bool>("sigmaTest"))
0218           SigmaHistos.find(make_pair(slID.wheel(), slID.sector()))->second->setBinContent(BinNumber, sigma);
0219       }
0220 
0221       if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0222         MonitorElement* res_histo_2D = igetter.get(getMEName2D(slID));
0223         if (res_histo_2D) {
0224           TH2F* res_histo_2D_root = res_histo_2D->getTH2F();
0225           int BinNumber = entry + slID.superLayer();
0226           if (BinNumber == 12)
0227             BinNumber = 11;
0228           TProfile* prof = res_histo_2D_root->ProfileX();
0229           prof->GetXaxis()->SetRangeUser(0, 2);
0230           //prof->Fit("pol1","Q0");
0231           //need our own copy to avoid threading problems
0232           TF1 fitting("mypol1", "pol1");
0233           try {
0234             prof->Fit(&fitting, "Q0");
0235           } catch (cms::Exception& iException) {
0236             edm::LogError("resolution") << "[DTResolutionTest]: Exception when fitting..."
0237                                         << "SuperLayer : " << slID << "\n"
0238                                         << "                    STEP : "
0239                                         << parameters.getUntrackedParameter<string>("STEP", "STEP3") << "\n"
0240                                         << "Filling slope histogram with standard value -99. for bin " << BinNumber;
0241             SlopeHistos.find(make_pair(slID.wheel(), slID.sector()))->second->setBinContent(BinNumber, -99.);
0242             continue;
0243           }
0244           double slope = fitting.GetParameter(1);
0245           SlopeHistos.find(make_pair(slID.wheel(), slID.sector()))->second->setBinContent(BinNumber, slope);
0246         }
0247       }
0248     }
0249   }
0250 
0251   // Mean test
0252   string MeanCriterionName = parameters.getUntrackedParameter<string>("meanTestName", "ResidualsMeanInRange");
0253   for (map<pair<int, int>, MonitorElement*>::const_iterator hMean = MeanHistos.begin(); hMean != MeanHistos.end();
0254        hMean++) {
0255     const QReport* theMeanQReport = (*hMean).second->getQReport(MeanCriterionName);
0256     stringstream wheel;
0257     wheel << (*hMean).first.first;
0258     stringstream sector;
0259     sector << (*hMean).first.second;
0260     // Report the channels failing the test on the mean
0261     if (theMeanQReport) {
0262       vector<dqm::me_util::Channel> badChannels = theMeanQReport->getBadChannels();
0263       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0264            channel++) {
0265         edm::LogError("resolution") << "Bad mean channel: wh: " << wheel.str()
0266                                     << " st: " << stationFromBin((*channel).getBin()) << " sect: " << sector.str()
0267                                     << " sl: " << slFromBin((*channel).getBin())
0268                                     << " mean (cm): " << (*channel).getContents();
0269         string HistoName = "W" + wheel.str() + "_Sec" + sector.str();
0270         if (parameters.getUntrackedParameter<bool>("meanWrongHisto")) {
0271           MeanHistosSetRange.find(HistoName)->second->Fill((*channel).getBin());
0272           MeanHistosSetRange2D.find(HistoName)->second->Fill((*channel).getBin(), (*channel).getContents());
0273         }
0274         // fill the wheel summary histos if the SL has not passed the test
0275         if (abs((*channel).getContents()) < parameters.getUntrackedParameter<double>("meanMaxLimit"))
0276           wheelMeanHistos[(*hMean).first.first]->Fill(((*hMean).first.second) - 1, (*channel).getBin() - 1, 1);
0277         else
0278           wheelMeanHistos[(*hMean).first.first]->Fill(((*hMean).first.second) - 1, (*channel).getBin() - 1, 2);
0279         // fill the cms summary histo if the percentual of SL which have not passed the test
0280         // is more than a predefined treshold
0281         if (abs((*channel).getContents()) > parameters.getUntrackedParameter<double>("meanMaxLimit")) {
0282           cmsMeanHistos[make_pair((*hMean).first.first, (*hMean).first.second)]++;
0283           if (((*hMean).first.second < 13 &&
0284                double(cmsMeanHistos[make_pair((*hMean).first.first, (*hMean).first.second)]) / 11 >
0285                    double(percentual) / 100 &&
0286                MeanFilled[make_pair((*hMean).first.first, (*hMean).first.second)] == false) ||
0287               ((*hMean).first.first >= 13 &&
0288                double(cmsMeanHistos[make_pair((*hMean).first.first, (*hMean).first.second)]) / 2 >
0289                    double(percentual) / 100 &&
0290                MeanFilled[make_pair((*hMean).first.first, (*hMean).first.second)] == false)) {
0291             MeanFilled[make_pair((*hMean).first.first, (*hMean).first.second)] = true;
0292             wheelMeanHistos[3]->Fill(((*hMean).first.second) - 1, (*hMean).first.first);
0293           }
0294         }
0295       }
0296     }
0297   }
0298 
0299   // Sigma test
0300   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0301     string SigmaCriterionName = parameters.getUntrackedParameter<string>("sigmaTestName", "ResidualsSigmaInRange");
0302     for (map<pair<int, int>, MonitorElement*>::const_iterator hSigma = SigmaHistos.begin(); hSigma != SigmaHistos.end();
0303          hSigma++) {
0304       const QReport* theSigmaQReport = (*hSigma).second->getQReport(SigmaCriterionName);
0305       stringstream wheel;
0306       wheel << (*hSigma).first.first;
0307       stringstream sector;
0308       sector << (*hSigma).first.second;
0309       if (theSigmaQReport) {
0310         vector<dqm::me_util::Channel> badChannels = theSigmaQReport->getBadChannels();
0311         for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0312              channel++) {
0313           edm::LogError("resolution") << "Bad sigma: wh: " << wheel.str()
0314                                       << " st: " << stationFromBin((*channel).getBin()) << " sect: " << sector.str()
0315                                       << " sl: " << slFromBin((*channel).getBin())
0316                                       << " sigma (cm): " << (*channel).getContents();
0317           string HistoName = "W" + wheel.str() + "_Sec" + sector.str();
0318           SigmaHistosSetRange.find(HistoName)->second->Fill((*channel).getBin());
0319           SigmaHistosSetRange2D.find(HistoName)->second->Fill((*channel).getBin(), (*channel).getContents());
0320           // fill the wheel summary histos if the SL has not passed the test
0321           wheelSigmaHistos[(*hSigma).first.first]->Fill(((*hSigma).first.second) - 1, (*channel).getBin() - 1);
0322           // fill the cms summary histo if the percentual of SL which have not passed the test
0323           // is more than a predefined treshold
0324           cmsSigmaHistos[make_pair((*hSigma).first.first, (*hSigma).first.second)]++;
0325           if (((*hSigma).first.second < 13 &&
0326                double(cmsSigmaHistos[make_pair((*hSigma).first.first, (*hSigma).first.second)]) / 11 >
0327                    double(percentual) / 100 &&
0328                SigmaFilled[make_pair((*hSigma).first.first, (*hSigma).first.second)] == false) ||
0329               ((*hSigma).first.first >= 13 &&
0330                double(cmsSigmaHistos[make_pair((*hSigma).first.first, (*hSigma).first.second)]) / 2 >
0331                    double(percentual) / 100 &&
0332                SigmaFilled[make_pair((*hSigma).first.first, (*hSigma).first.second)] == false)) {
0333             SigmaFilled[make_pair((*hSigma).first.first, (*hSigma).first.second)] = true;
0334             wheelSigmaHistos[3]->Fill((*hSigma).first.second - 1, (*hSigma).first.first);
0335           }
0336         }
0337       }
0338     }
0339   }
0340 
0341   // Slope test
0342   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0343     string SlopeCriterionName = parameters.getUntrackedParameter<string>("slopeTestName", "ResidualsSlopeInRange");
0344     for (map<pair<int, int>, MonitorElement*>::const_iterator hSlope = SlopeHistos.begin(); hSlope != SlopeHistos.end();
0345          hSlope++) {
0346       const QReport* theSlopeQReport = (*hSlope).second->getQReport(SlopeCriterionName);
0347       stringstream wheel;
0348       wheel << (*hSlope).first.first;
0349       stringstream sector;
0350       sector << (*hSlope).first.second;
0351       if (theSlopeQReport) {
0352         vector<dqm::me_util::Channel> badChannels = theSlopeQReport->getBadChannels();
0353         for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0354              channel++) {
0355           edm::LogError("resolution") << "Bad slope: wh: " << wheel.str()
0356                                       << " st: " << stationFromBin((*channel).getBin()) << " sect: " << sector.str()
0357                                       << " sl: " << slFromBin((*channel).getBin())
0358                                       << " slope: " << (*channel).getContents();
0359           string HistoName = "W" + wheel.str() + "_Sec" + sector.str();
0360           SlopeHistosSetRange.find(HistoName)->second->Fill((*channel).getBin());
0361           SlopeHistosSetRange2D.find(HistoName)->second->Fill((*channel).getBin(), (*channel).getContents());
0362           // fill the wheel summary histos if the SL has not passed the test
0363           wheelSlopeHistos[(*hSlope).first.first]->Fill(((*hSlope).first.second) - 1, (*channel).getBin() - 1);
0364           // fill the cms summary histo if the percentual of SL which have not passed the test
0365           // is more than a predefined treshold
0366           cmsSlopeHistos[make_pair((*hSlope).first.first, (*hSlope).first.second)]++;
0367           if (((*hSlope).first.second < 13 &&
0368                double(cmsSlopeHistos[make_pair((*hSlope).first.first, (*hSlope).first.second)]) / 11 >
0369                    double(percentual) / 100 &&
0370                SlopeFilled[make_pair((*hSlope).first.first, (*hSlope).first.second)] == false) ||
0371               ((*hSlope).first.first >= 13 &&
0372                double(cmsSlopeHistos[make_pair((*hSlope).first.first, (*hSlope).first.second)]) / 2 >
0373                    double(percentual) / 100 &&
0374                SlopeFilled[make_pair((*hSlope).first.first, (*hSlope).first.second)] == false)) {
0375             SlopeFilled[make_pair((*hSlope).first.first, (*hSlope).first.second)] = true;
0376             wheelSlopeHistos[3]->Fill((*hSlope).first.second - 1, (*hSlope).first.first);
0377           }
0378         }
0379       }
0380     }
0381   }
0382 }
0383 
0384 void DTResolutionTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0385   edm::LogVerbatim("resolution") << "[DTResolutionTest] endjob called!";
0386 
0387   //FR: no idea if this output file needs to be written!! comment it for now
0388   //bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
0389   //if(outputMEsInRootFile){
0390   //    std::string outputFileName = parameters.getParameter<std::string>("OutputFileName");
0391   //    dbe->save(outputFileName,"DT/CalibrationSummary");
0392   //}
0393 }
0394 
0395 string DTResolutionTest::getMEName(const DTSuperLayerId& slID) {
0396   stringstream wheel;
0397   wheel << slID.wheel();
0398   stringstream station;
0399   station << slID.station();
0400   stringstream sector;
0401   sector << slID.sector();
0402   stringstream superLayer;
0403   superLayer << slID.superlayer();
0404 
0405   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
0406   string folderName;
0407 
0408   if (parameters.getUntrackedParameter<bool>("calibModule", false)) {
0409     folderName = folderRoot + "DT/DTCalibValidation/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" +
0410                  sector.str() + "/";
0411   } else {
0412     folderName = folderRoot + "DT/DTResolutionAnalysisTask/Wheel" + wheel.str() + "/Station" + station.str() +
0413                  "/Sector" + sector.str() + "/";
0414   }
0415 
0416   string histoTag = parameters.getUntrackedParameter<string>("histoTag", "hResDist");
0417 
0418   string histoname = folderName + histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +
0419                      "_SL" + superLayer.str();
0420 
0421   return histoname;
0422 }
0423 
0424 string DTResolutionTest::getMEName2D(const DTSuperLayerId& slID) {
0425   stringstream wheel;
0426   wheel << slID.wheel();
0427   stringstream station;
0428   station << slID.station();
0429   stringstream sector;
0430   sector << slID.sector();
0431   stringstream superLayer;
0432   superLayer << slID.superlayer();
0433 
0434   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
0435   string folderName;
0436 
0437   if (parameters.getUntrackedParameter<bool>("calibModule", false)) {
0438     folderName = folderRoot + "DT/DTCalibValidation/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" +
0439                  sector.str() + "/";
0440   } else {
0441     folderName = folderRoot + "DT/DTResolutionAnalysisTask/Wheel" + wheel.str() + "/Station" + station.str() +
0442                  "/Sector" + sector.str() + "/";
0443   }
0444 
0445   string histoTag2D = parameters.getUntrackedParameter<string>("histoTag2D", "hResDistVsDist");
0446 
0447   string histoname = folderName + histoTag2D + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +
0448                      "_SL" + superLayer.str();
0449 
0450   return histoname;
0451 }
0452 
0453 void DTResolutionTest::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& ch) {
0454   stringstream wheel;
0455   wheel << ch.wheel();
0456   stringstream sector;
0457   sector << ch.sector();
0458 
0459   string MeanHistoName = "MeanTest_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" + wheel.str() +
0460                          "_Sec" + sector.str();
0461   string SigmaHistoName = "SigmaTest_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0462                           wheel.str() + "_Sec" + sector.str();
0463   string SlopeHistoName = "SlopeTest_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0464                           wheel.str() + "_Sec" + sector.str();
0465 
0466   ibooker.setCurrentFolder("DT/Tests/DTResolution");
0467 
0468   // Book the histo for the mean value and set the axis labels
0469 
0470   MeanHistos[make_pair(ch.wheel(), ch.sector())] =
0471       ibooker.book1D(MeanHistoName.c_str(), MeanHistoName.c_str(), 11, 0, 11);
0472   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(1, "MB1_SL1", 1);
0473   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(2, "MB1_SL2", 1);
0474   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(3, "MB1_SL3", 1);
0475   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(4, "MB2_SL1", 1);
0476   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(5, "MB2_SL2", 1);
0477   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(6, "MB2_SL3", 1);
0478   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(7, "MB3_SL1", 1);
0479   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(8, "MB3_SL2", 1);
0480   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(9, "MB3_SL3", 1);
0481   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(10, "MB4_SL1", 1);
0482   (MeanHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(11, "MB4_SL3", 1);
0483 
0484   // Book the histo for the sigma value and set the axis labels
0485   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0486     SigmaHistos[make_pair(ch.wheel(), ch.sector())] =
0487         ibooker.book1D(SigmaHistoName.c_str(), SigmaHistoName.c_str(), 11, 0, 11);
0488     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(1, "MB1_SL1", 1);
0489     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(2, "MB1_SL2", 1);
0490     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(3, "MB1_SL3", 1);
0491     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(4, "MB2_SL1", 1);
0492     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(5, "MB2_SL2", 1);
0493     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(6, "MB2_SL3", 1);
0494     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(7, "MB3_SL1", 1);
0495     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(8, "MB3_SL2", 1);
0496     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(9, "MB3_SL3", 1);
0497     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(10, "MB4_SL1", 1);
0498     (SigmaHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(11, "MB4_SL3", 1);
0499   }
0500 
0501   // Book the histo for the slope value and set the axis labels
0502   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0503     SlopeHistos[make_pair(ch.wheel(), ch.sector())] =
0504         ibooker.book1D(SlopeHistoName.c_str(), SlopeHistoName.c_str(), 11, 0, 11);
0505     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(1, "MB1_SL1", 1);
0506     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(2, "MB1_SL2", 1);
0507     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(3, "MB1_SL3", 1);
0508     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(4, "MB2_SL1", 1);
0509     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(5, "MB2_SL2", 1);
0510     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(6, "MB2_SL3", 1);
0511     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(7, "MB3_SL1", 1);
0512     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(8, "MB3_SL2", 1);
0513     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(9, "MB3_SL3", 1);
0514     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(10, "MB4_SL1", 1);
0515     (SlopeHistos[make_pair(ch.wheel(), ch.sector())])->setBinLabel(11, "MB4_SL3", 1);
0516   }
0517 
0518   string HistoName = "W" + wheel.str() + "_Sec" + sector.str();
0519 
0520   if (parameters.getUntrackedParameter<bool>("meanWrongHisto")) {
0521     string MeanHistoNameSetRange = "MeanWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0522                                    wheel.str() + "_Sec" + sector.str() + "_SetRange";
0523 
0524     MeanHistosSetRange[HistoName] =
0525         ibooker.book1D(MeanHistoNameSetRange.c_str(), MeanHistoNameSetRange.c_str(), 11, 0.5, 11.5);
0526     string MeanHistoNameSetRange2D = "MeanWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0527                                      wheel.str() + "_Sec" + sector.str() + "_SetRange" + "_2D";
0528     MeanHistosSetRange2D[HistoName] = ibooker.book2D(
0529         MeanHistoNameSetRange2D.c_str(), MeanHistoNameSetRange2D.c_str(), 11, 0.5, 11.5, 100, -0.05, 0.05);
0530   }
0531 
0532   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0533     string SigmaHistoNameSetRange = "SigmaWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0534                                     wheel.str() + "_Sec" + sector.str() + "_SetRange";
0535 
0536     SigmaHistosSetRange[HistoName] =
0537         ibooker.book1D(SigmaHistoNameSetRange.c_str(), SigmaHistoNameSetRange.c_str(), 11, 0.5, 11.5);
0538     string SigmaHistoNameSetRange2D = "SigmaWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0539                                       wheel.str() + "_Sec" + sector.str() + "_SetRange" + "_2D";
0540 
0541     SigmaHistosSetRange2D[HistoName] =
0542         ibooker.book2D(SigmaHistoNameSetRange2D.c_str(), SigmaHistoNameSetRange2D.c_str(), 11, 0.5, 11.5, 500, 0, 0.5);
0543   }
0544 
0545   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0546     string SlopeHistoNameSetRange = "SlopeWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0547                                     wheel.str() + "_Sec" + sector.str() + "_SetRange";
0548 
0549     SlopeHistosSetRange[HistoName] =
0550         ibooker.book1D(SlopeHistoNameSetRange.c_str(), SlopeHistoNameSetRange.c_str(), 11, 0.5, 11.5);
0551     string SlopeHistoNameSetRange2D = "SlopeWrong_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" +
0552                                       wheel.str() + "_Sec" + sector.str() + "_SetRange" + "_2D";
0553 
0554     SlopeHistosSetRange2D[HistoName] = ibooker.book2D(
0555         SlopeHistoNameSetRange2D.c_str(), SlopeHistoNameSetRange2D.c_str(), 11, 0.5, 11.5, 200, -0.1, 0.1);
0556   }
0557 }
0558 
0559 void DTResolutionTest::bookHistos(DQMStore::IBooker& ibooker, int wh) {
0560   ibooker.setCurrentFolder("DT/CalibrationSummary");
0561 
0562   if (wheelMeanHistos.find(3) == wheelMeanHistos.end()) {
0563     string histoName =
0564         "MeanSummaryRes_testFailedByAtLeastBadSL_" + parameters.getUntrackedParameter<string>("STEP", "STEP3");
0565 
0566     wheelMeanHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 5, -2, 3);
0567     wheelMeanHistos[3]->setBinLabel(1, "Sector1", 1);
0568     wheelMeanHistos[3]->setBinLabel(1, "Sector1", 1);
0569     wheelMeanHistos[3]->setBinLabel(2, "Sector2", 1);
0570     wheelMeanHistos[3]->setBinLabel(3, "Sector3", 1);
0571     wheelMeanHistos[3]->setBinLabel(4, "Sector4", 1);
0572     wheelMeanHistos[3]->setBinLabel(5, "Sector5", 1);
0573     wheelMeanHistos[3]->setBinLabel(6, "Sector6", 1);
0574     wheelMeanHistos[3]->setBinLabel(7, "Sector7", 1);
0575     wheelMeanHistos[3]->setBinLabel(8, "Sector8", 1);
0576     wheelMeanHistos[3]->setBinLabel(9, "Sector9", 1);
0577     wheelMeanHistos[3]->setBinLabel(10, "Sector10", 1);
0578     wheelMeanHistos[3]->setBinLabel(11, "Sector11", 1);
0579     wheelMeanHistos[3]->setBinLabel(12, "Sector12", 1);
0580     wheelMeanHistos[3]->setBinLabel(13, "Sector13", 1);
0581     wheelMeanHistos[3]->setBinLabel(14, "Sector14", 1);
0582     wheelMeanHistos[3]->setBinLabel(1, "Wheel-2", 2);
0583     wheelMeanHistos[3]->setBinLabel(2, "Wheel-1", 2);
0584     wheelMeanHistos[3]->setBinLabel(3, "Wheel0", 2);
0585     wheelMeanHistos[3]->setBinLabel(4, "Wheel+1", 2);
0586     wheelMeanHistos[3]->setBinLabel(5, "Wheel+2", 2);
0587   }
0588 
0589   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0590     if (wheelSigmaHistos.find(3) == wheelSigmaHistos.end()) {
0591       string histoName =
0592           "SigmaSummaryRes_testFailedByAtLeastBadSL_" + parameters.getUntrackedParameter<string>("STEP", "STEP3");
0593 
0594       wheelSigmaHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 5, -2, 3);
0595       wheelSigmaHistos[3]->setBinLabel(1, "Sector1", 1);
0596       wheelSigmaHistos[3]->setBinLabel(1, "Sector1", 1);
0597       wheelSigmaHistos[3]->setBinLabel(2, "Sector2", 1);
0598       wheelSigmaHistos[3]->setBinLabel(3, "Sector3", 1);
0599       wheelSigmaHistos[3]->setBinLabel(4, "Sector4", 1);
0600       wheelSigmaHistos[3]->setBinLabel(5, "Sector5", 1);
0601       wheelSigmaHistos[3]->setBinLabel(6, "Sector6", 1);
0602       wheelSigmaHistos[3]->setBinLabel(7, "Sector7", 1);
0603       wheelSigmaHistos[3]->setBinLabel(8, "Sector8", 1);
0604       wheelSigmaHistos[3]->setBinLabel(9, "Sector9", 1);
0605       wheelSigmaHistos[3]->setBinLabel(10, "Sector10", 1);
0606       wheelSigmaHistos[3]->setBinLabel(11, "Sector11", 1);
0607       wheelSigmaHistos[3]->setBinLabel(12, "Sector12", 1);
0608       wheelSigmaHistos[3]->setBinLabel(13, "Sector13", 1);
0609       wheelSigmaHistos[3]->setBinLabel(14, "Sector14", 1);
0610       wheelSigmaHistos[3]->setBinLabel(1, "Wheel-2", 2);
0611       wheelSigmaHistos[3]->setBinLabel(2, "Wheel-1", 2);
0612       wheelSigmaHistos[3]->setBinLabel(3, "Wheel0", 2);
0613       wheelSigmaHistos[3]->setBinLabel(4, "Wheel+1", 2);
0614       wheelSigmaHistos[3]->setBinLabel(5, "Wheel+2", 2);
0615     }
0616   }
0617 
0618   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0619     if (wheelSlopeHistos.find(3) == wheelSlopeHistos.end()) {
0620       string histoName =
0621           "SlopeSummaryRes_testFailedByAtLeastBadSL_" + parameters.getUntrackedParameter<string>("STEP", "STEP3");
0622 
0623       wheelSlopeHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 5, -2, 3);
0624       wheelSlopeHistos[3]->setBinLabel(1, "Sector1", 1);
0625       wheelSlopeHistos[3]->setBinLabel(1, "Sector1", 1);
0626       wheelSlopeHistos[3]->setBinLabel(2, "Sector2", 1);
0627       wheelSlopeHistos[3]->setBinLabel(3, "Sector3", 1);
0628       wheelSlopeHistos[3]->setBinLabel(4, "Sector4", 1);
0629       wheelSlopeHistos[3]->setBinLabel(5, "Sector5", 1);
0630       wheelSlopeHistos[3]->setBinLabel(6, "Sector6", 1);
0631       wheelSlopeHistos[3]->setBinLabel(7, "Sector7", 1);
0632       wheelSlopeHistos[3]->setBinLabel(8, "Sector8", 1);
0633       wheelSlopeHistos[3]->setBinLabel(9, "Sector9", 1);
0634       wheelSlopeHistos[3]->setBinLabel(10, "Sector10", 1);
0635       wheelSlopeHistos[3]->setBinLabel(11, "Sector11", 1);
0636       wheelSlopeHistos[3]->setBinLabel(12, "Sector12", 1);
0637       wheelSlopeHistos[3]->setBinLabel(13, "Sector13", 1);
0638       wheelSlopeHistos[3]->setBinLabel(14, "Sector14", 1);
0639       wheelSlopeHistos[3]->setBinLabel(1, "Wheel-2", 2);
0640       wheelSlopeHistos[3]->setBinLabel(2, "Wheel-1", 2);
0641       wheelSlopeHistos[3]->setBinLabel(3, "Wheel0", 2);
0642       wheelSlopeHistos[3]->setBinLabel(4, "Wheel+1", 2);
0643       wheelSlopeHistos[3]->setBinLabel(5, "Wheel+2", 2);
0644     }
0645   }
0646 
0647   stringstream wheel;
0648   wheel << wh;
0649 
0650   if (wheelMeanHistos.find(wh) == wheelMeanHistos.end()) {
0651     string histoName =
0652         "MeanSummaryRes_testFailed_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") + "_W" + wheel.str();
0653 
0654     wheelMeanHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 11, 0, 11);
0655     wheelMeanHistos[wh]->setBinLabel(1, "Sector1", 1);
0656     wheelMeanHistos[wh]->setBinLabel(2, "Sector2", 1);
0657     wheelMeanHistos[wh]->setBinLabel(3, "Sector3", 1);
0658     wheelMeanHistos[wh]->setBinLabel(4, "Sector4", 1);
0659     wheelMeanHistos[wh]->setBinLabel(5, "Sector5", 1);
0660     wheelMeanHistos[wh]->setBinLabel(6, "Sector6", 1);
0661     wheelMeanHistos[wh]->setBinLabel(7, "Sector7", 1);
0662     wheelMeanHistos[wh]->setBinLabel(8, "Sector8", 1);
0663     wheelMeanHistos[wh]->setBinLabel(9, "Sector9", 1);
0664     wheelMeanHistos[wh]->setBinLabel(10, "Sector10", 1);
0665     wheelMeanHistos[wh]->setBinLabel(11, "Sector11", 1);
0666     wheelMeanHistos[wh]->setBinLabel(12, "Sector12", 1);
0667     wheelMeanHistos[wh]->setBinLabel(13, "Sector13", 1);
0668     wheelMeanHistos[wh]->setBinLabel(14, "Sector14", 1);
0669     wheelMeanHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0670     wheelMeanHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0671     wheelMeanHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0672     wheelMeanHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0673     wheelMeanHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0674     wheelMeanHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0675     wheelMeanHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0676     wheelMeanHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0677     wheelMeanHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0678     wheelMeanHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0679     wheelMeanHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0680   }
0681 
0682   if (parameters.getUntrackedParameter<bool>("sigmaTest")) {
0683     if (wheelSigmaHistos.find(wh) == wheelSigmaHistos.end()) {
0684       string histoName = "SigmaSummaryRes_testFailed_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") +
0685                          "_W" + wheel.str();
0686 
0687       wheelSigmaHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 11, 0, 11);
0688       wheelSigmaHistos[wh]->setBinLabel(1, "Sector1", 1);
0689       wheelSigmaHistos[wh]->setBinLabel(2, "Sector2", 1);
0690       wheelSigmaHistos[wh]->setBinLabel(3, "Sector3", 1);
0691       wheelSigmaHistos[wh]->setBinLabel(4, "Sector4", 1);
0692       wheelSigmaHistos[wh]->setBinLabel(5, "Sector5", 1);
0693       wheelSigmaHistos[wh]->setBinLabel(6, "Sector6", 1);
0694       wheelSigmaHistos[wh]->setBinLabel(7, "Sector7", 1);
0695       wheelSigmaHistos[wh]->setBinLabel(8, "Sector8", 1);
0696       wheelSigmaHistos[wh]->setBinLabel(9, "Sector9", 1);
0697       wheelSigmaHistos[wh]->setBinLabel(10, "Sector10", 1);
0698       wheelSigmaHistos[wh]->setBinLabel(11, "Sector11", 1);
0699       wheelSigmaHistos[wh]->setBinLabel(12, "Sector12", 1);
0700       wheelSigmaHistos[wh]->setBinLabel(13, "Sector13", 1);
0701       wheelSigmaHistos[wh]->setBinLabel(14, "Sector14", 1);
0702       wheelSigmaHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0703       wheelSigmaHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0704       wheelSigmaHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0705       wheelSigmaHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0706       wheelSigmaHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0707       wheelSigmaHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0708       wheelSigmaHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0709       wheelSigmaHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0710       wheelSigmaHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0711       wheelSigmaHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0712       wheelSigmaHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0713     }
0714   }
0715 
0716   if (parameters.getUntrackedParameter<bool>("slopeTest")) {
0717     if (wheelSlopeHistos.find(wh) == wheelSlopeHistos.end()) {
0718       string histoName = "SlopeSummaryRes_testFailed_" + parameters.getUntrackedParameter<string>("STEP", "STEP3") +
0719                          "_W" + wheel.str();
0720 
0721       wheelSlopeHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 11, 0, 11);
0722       wheelSlopeHistos[wh]->setBinLabel(1, "Sector1", 1);
0723       wheelSlopeHistos[wh]->setBinLabel(2, "Sector2", 1);
0724       wheelSlopeHistos[wh]->setBinLabel(3, "Sector3", 1);
0725       wheelSlopeHistos[wh]->setBinLabel(4, "Sector4", 1);
0726       wheelSlopeHistos[wh]->setBinLabel(5, "Sector5", 1);
0727       wheelSlopeHistos[wh]->setBinLabel(6, "Sector6", 1);
0728       wheelSlopeHistos[wh]->setBinLabel(7, "Sector7", 1);
0729       wheelSlopeHistos[wh]->setBinLabel(8, "Sector8", 1);
0730       wheelSlopeHistos[wh]->setBinLabel(9, "Sector9", 1);
0731       wheelSlopeHistos[wh]->setBinLabel(10, "Sector10", 1);
0732       wheelSlopeHistos[wh]->setBinLabel(11, "Sector11", 1);
0733       wheelSlopeHistos[wh]->setBinLabel(12, "Sector12", 1);
0734       wheelSlopeHistos[wh]->setBinLabel(13, "Sector13", 1);
0735       wheelSlopeHistos[wh]->setBinLabel(14, "Sector14", 1);
0736       wheelSlopeHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0737       wheelSlopeHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0738       wheelSlopeHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0739       wheelSlopeHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0740       wheelSlopeHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0741       wheelSlopeHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0742       wheelSlopeHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0743       wheelSlopeHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0744       wheelSlopeHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0745       wheelSlopeHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0746       wheelSlopeHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0747     }
0748   }
0749 }
0750 
0751 int DTResolutionTest::stationFromBin(int bin) const { return (int)(bin / 3.1) + 1; }
0752 
0753 int DTResolutionTest::slFromBin(int bin) const {
0754   int ret = bin % 3;
0755   if (ret == 0 || bin == 11)
0756     ret = 3;
0757 
0758   return ret;
0759 }