Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:27

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author G. Mila - INFN Torino
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
0007  *
0008  */
0009 
0010 #include "DQM/DTMonitorClient/src/DTResolutionAnalysisTest.h"
0011 
0012 // Framework
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 
0017 // Geometry
0018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0019 
0020 #include "FWCore/ServiceRegistry/interface/Service.h"
0021 #include "DQMServices/Core/interface/DQMStore.h"
0022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0023 
0024 #include <string>
0025 #include <sstream>
0026 #include <cmath>
0027 
0028 using namespace edm;
0029 using namespace std;
0030 
0031 DTResolutionAnalysisTest::DTResolutionAnalysisTest(const ParameterSet& ps)
0032     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0033   LogTrace("DTDQM|DTMonitorClient|DTResolutionAnalysisTest") << "[DTResolutionAnalysisTest]: Constructor";
0034 
0035   prescaleFactor = ps.getUntrackedParameter<int>("diagnosticPrescale", 1);
0036   // permitted test range
0037   maxGoodMeanValue = ps.getUntrackedParameter<double>("maxGoodMeanValue", 0.02);
0038   minBadMeanValue = ps.getUntrackedParameter<double>("minBadMeanValue", 0.04);
0039   maxGoodSigmaValue = ps.getUntrackedParameter<double>("maxGoodSigmaValue", 0.08);
0040   minBadSigmaValue = ps.getUntrackedParameter<double>("minBadSigmaValue", 0.16);
0041   // top folder for the histograms in DQMStore
0042   topHistoFolder = ps.getUntrackedParameter<string>("topHistoFolder", "DT/02-Segments");
0043 
0044   doCalibAnalysis = ps.getUntrackedParameter<bool>("doCalibAnalysis", false);
0045 }
0046 
0047 DTResolutionAnalysisTest::~DTResolutionAnalysisTest() {
0048   LogTrace("DTDQM|DTMonitorClient|DTResolutionAnalysisTest")
0049       << "DTResolutionAnalysisTest: analyzed " << nevents << " events";
0050 }
0051 
0052 void DTResolutionAnalysisTest::beginRun(const Run& run, const EventSetup& context) {
0053   LogTrace("DTDQM|DTMonitorClient|DTResolutionAnalysisTest") << "[DTResolutionAnalysisTest]: BeginRun";
0054 
0055   nevents = 0;
0056 
0057   // Get the geometry
0058   muonGeom = &context.getData(muonGeomToken_);
0059 }
0060 
0061 void DTResolutionAnalysisTest::bookHistos(DQMStore::IBooker& ibooker) {
0062   // global residual summary
0063 
0064   ibooker.setCurrentFolder(topHistoFolder);
0065   globalResSummary =
0066       ibooker.book2D("ResidualsGlbSummary", "# of SLs with good mean and good sigma of residuals", 12, 1, 13, 5, -2, 3);
0067 
0068   // book summaries for mean and sigma
0069   ibooker.setCurrentFolder(topHistoFolder + "/00-MeanRes");
0070   meanDistr[-2] = ibooker.book1D("MeanDistr", "Mean value of the residuals all (cm)", 100, -0.1, 0.1);
0071   meanDistr[-1] = ibooker.book1D("MeanDistr_Phi", "Mean value of the residuals #phi SL (cm)", 100, -0.1, 0.1);
0072   meanDistr[0] =
0073       ibooker.book1D("MeanDistr_ThetaWh0", "Mean values of the residuals #theta SL Wh 0 (cm)", 100, -0.1, 0.1);
0074   meanDistr[1] =
0075       ibooker.book1D("MeanDistr_ThetaWh1", "Mean value of the residuals #theta SL Wh +/-1 (cm)", 100, -0.1, 0.1);
0076   meanDistr[2] =
0077       ibooker.book1D("MeanDistr_ThetaWh2", "Mean value of the residuals #theta SL Wh +/-2 (cm)", 100, -0.1, 0.1);
0078 
0079   for (int wh = -2; wh <= 2; wh++) {
0080     stringstream wheel;
0081     wheel << wh;
0082     for (int st = 1; st <= 4; st++) {
0083       stringstream station;
0084       station << st;
0085       for (int sl = 1; sl <= 3; sl++) {
0086         if (sl == 2 && st == 4)
0087           continue;
0088         ibooker.setCurrentFolder(topHistoFolder + "/00-MeanRes/Wheel" + wheel.str());
0089         string histoLabel = "MeanDistr";
0090         string histoName = histoLabel + "_W" + wheel.str() + "_MB" + station.str() + "_Phi";
0091         string histoType = histoLabel + "_Phi";
0092         if (sl == 2) {
0093           histoName = histoLabel + "_W" + wheel.str() + "_MB" + station.str() + "_Theta";
0094           histoType = histoLabel + "_Theta";
0095         }
0096         (wheelRingHistos[wh][st])[histoType] = ibooker.book1D(histoName, histoLabel, 100, -0.1, 0.1);
0097 
0098         ibooker.setCurrentFolder(topHistoFolder + "/01-SigmaRes/Wheel" + wheel.str());
0099         histoLabel = "SigmaDistr";
0100         histoName = histoLabel + "_W" + wheel.str() + "_MB" + station.str() + "_Phi";
0101         histoType = histoLabel + "_Phi";
0102         if (sl == 2) {
0103           histoName = histoLabel + "_W" + wheel.str() + "_MB" + station.str() + "_Theta";
0104           histoType = histoLabel + "_Theta";
0105         }
0106         (wheelRingHistos[wh][st])[histoType] = ibooker.book1D(histoName, histoLabel, 50, 0.0, 0.2);
0107       }
0108     }
0109   }
0110 
0111   ibooker.setCurrentFolder(topHistoFolder + "/00-MeanRes");
0112   string histoTitle = "# of SLs with good mean of residuals";
0113   wheelMeanHistos[3] = ibooker.book2D("MeanResGlbSummary", histoTitle.c_str(), 12, 1, 13, 5, -2, 3);
0114   wheelMeanHistos[3]->setAxisTitle("Sector", 1);
0115   wheelMeanHistos[3]->setAxisTitle("Wheel", 2);
0116 
0117   ibooker.setCurrentFolder(topHistoFolder + "/01-SigmaRes");
0118   sigmaDistr[-2] = ibooker.book1D("SigmaDistr", "Sigma value of the residuals all (cm)", 50, 0.0, 0.2);
0119   sigmaDistr[-1] = ibooker.book1D("SigmaDistr_Phi", "Sigma value of the residuals #phi SL (cm)", 50, 0.0, 0.2);
0120   sigmaDistr[0] =
0121       ibooker.book1D("SigmaDistr_ThetaWh0", "Sigma value of the residuals #theta SL Wh 0 (cm)", 50, 0.0, 0.2);
0122   sigmaDistr[1] =
0123       ibooker.book1D("SigmaDistr_ThetaWh1", "Sigma value of the residuals #theta SL Wh +/-1 (cm)", 50, 0.0, 0.2);
0124   sigmaDistr[2] =
0125       ibooker.book1D("SigmaDistr_ThetaWh2", "Sigma value of the residuals #theta SL Wh +/-2 (cm)", 50, 0.0, 0.2);
0126 
0127   histoTitle = "# of SLs with good sigma of residuals";
0128   wheelSigmaHistos[3] = ibooker.book2D("SigmaResGlbSummary", histoTitle.c_str(), 12, 1, 13, 5, -2, 3);
0129   wheelSigmaHistos[3]->setAxisTitle("Sector", 1);
0130   wheelSigmaHistos[3]->setAxisTitle("Wheel", 2);
0131 
0132   // loop over all the CMS wheels, sectors & book the summary histos
0133   for (int wheel = -2; wheel <= 2; wheel++) {
0134     bookHistos(ibooker, wheel);
0135     for (int sector = 1; sector <= 12; sector++) {
0136       bookHistos(ibooker, wheel, sector);
0137     }
0138   }
0139 }
0140 
0141 void DTResolutionAnalysisTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0142   if (!igetter.dirExists(topHistoFolder)) {
0143     LogTrace("DTDQM|DTMonitorClient|DTResolutionAnalysisTest")
0144         << "[DTResolutionAnalysisTest]: Base folder " << topHistoFolder << " does not exist. Skipping client operation."
0145         << endl;
0146     return;
0147   }
0148 
0149   bookHistos(ibooker);  // histos booked only if top histo folder exist
0150   // as Standard/AlcaReco Harvest is performed in the same step
0151 
0152   LogTrace("DTDQM|DTMonitorClient|DTResolutionAnalysisTest")
0153       << "[DTResolutionAnalysisTest]: End of Run transition, performing the DQM client operation" << endl;
0154 
0155   // reset the ME with fixed scale
0156   resetMEs();
0157 
0158   for (vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0159        ch_it != muonGeom->chambers().end();
0160        ++ch_it) {  // loop over the chambers
0161 
0162     DTChamberId chID = (*ch_it)->id();
0163 
0164     // Fill the test histos
0165     for (vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
0166          sl_it != (*ch_it)->superLayers().end();
0167          ++sl_it) {  // loop over SLs
0168 
0169       DTSuperLayerId slID = (*sl_it)->id();
0170       MonitorElement* res_histo = igetter.get(getMEName(slID));
0171 
0172       if (res_histo) {  // Gaussian Fit
0173         float statMean = res_histo->getMean(1);
0174         float statSigma = res_histo->getRMS(1);
0175         double mean = -1;
0176         double sigma = -1;
0177         TH1F* histo_root = res_histo->getTH1F();
0178 
0179         // fill the summaries
0180         int entry = (chID.station() - 1) * 3;
0181         int binSect = slID.sector();
0182         if (slID.sector() == 13)
0183           binSect = 4;
0184         else if (slID.sector() == 14)
0185           binSect = 10;
0186         int binSL = entry + slID.superLayer();
0187         if (chID.station() == 4 && slID.superLayer() == 3)
0188           binSL--;
0189         if ((slID.sector() == 13 || slID.sector() == 14) && slID.superLayer() == 1)
0190           binSL = 12;
0191         if ((slID.sector() == 13 || slID.sector() == 14) && slID.superLayer() == 3)
0192           binSL = 13;
0193 
0194         if (histo_root->GetEntries() > 20) {
0195           TF1* gfit = new TF1("Gaussian", "gaus", (statMean - (2 * statSigma)), (statMean + (2 * statSigma)));
0196           try {
0197             histo_root->Fit(gfit, "Q0 SERIAL", "", -0.1, 0.1);
0198           } catch (cms::Exception& iException) {
0199             LogWarning("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
0200                 << "[DTResolutionAnalysisTask]: Exception when fitting SL : " << slID;
0201             // FIXME: the SL is set as OK in the summary
0202             double weight = 1 / 11.;
0203             if ((binSect == 4 || binSect == 10) && slID.station() == 4)
0204               weight = 1 / 22.;
0205             globalResSummary->Fill(binSect, slID.wheel(), weight);
0206             continue;
0207           }
0208 
0209           if (gfit) {
0210             // get the mean and the sigma of the distribution
0211             mean = gfit->GetParameter(1);
0212             sigma = gfit->GetParameter(2);
0213 
0214             // fill the distributions
0215             meanDistr[-2]->Fill(mean);
0216             sigmaDistr[-2]->Fill(sigma);
0217             if (slID.superlayer() == 2) {
0218               meanDistr[abs(slID.wheel())]->Fill(mean);
0219               sigmaDistr[abs(slID.wheel())]->Fill(sigma);
0220             } else {
0221               meanDistr[-1]->Fill(mean);
0222               sigmaDistr[-1]->Fill(sigma);
0223             }
0224 
0225             string histoType = "MeanDistr_Phi";
0226             if (slID.superlayer() == 2)
0227               histoType = "MeanDistr_Theta";
0228             (wheelRingHistos[slID.wheel()][slID.station()])[histoType]->Fill(mean);
0229 
0230             histoType = "SigmaDistr_Phi";
0231             if (slID.superlayer() == 2)
0232               histoType = "SigmaDistr_Theta";
0233             (wheelRingHistos[slID.wheel()][slID.station()])[histoType]->Fill(sigma);
0234 
0235             // sector summaries
0236             MeanHistos[make_pair(slID.wheel(), binSect)]->setBinContent(binSL, mean);
0237             SigmaHistos[make_pair(slID.wheel(), binSect)]->setBinContent(binSL, sigma);
0238 
0239             if ((slID.sector() == 13 || slID.sector() == 14) && binSL == 12)
0240               binSL = 10;
0241             if ((slID.sector() == 13 || slID.sector() == 14) && binSL == 13)
0242               binSL = 11;
0243 
0244             if ((slID.sector() == 13 || slID.sector() == 14)) {
0245               double MeanVal = wheelMeanHistos[slID.wheel()]->getBinContent(binSect, binSL);
0246               double MeanBinVal = (MeanVal > 0. && MeanVal < meanInRange(mean)) ? MeanVal : meanInRange(mean);
0247               wheelMeanHistos[slID.wheel()]->setBinContent(binSect, binSL, MeanBinVal);
0248 
0249               double SigmaVal = wheelSigmaHistos[slID.wheel()]->getBinContent(binSect, binSL);
0250               double SigmaBinVal = (SigmaVal > 0. && SigmaVal < sigmaInRange(sigma)) ? SigmaVal : sigmaInRange(sigma);
0251               wheelSigmaHistos[slID.wheel()]->setBinContent(binSect, binSL, SigmaBinVal);
0252 
0253             } else {
0254               wheelMeanHistos[slID.wheel()]->setBinContent(binSect, binSL, meanInRange(mean));
0255               wheelSigmaHistos[slID.wheel()]->setBinContent(binSect, binSL, sigmaInRange(sigma));
0256             }
0257 
0258             // set the weight
0259             double weight = 1 / 11.;
0260             if ((binSect == 4 || binSect == 10) && slID.station() == 4)
0261               weight = 1 / 22.;
0262 
0263             // test the values of mean and sigma
0264             if ((meanInRange(mean) > 0.85) && (sigmaInRange(sigma) > 0.85)) {  // sigma and mean ok
0265               globalResSummary->Fill(binSect, slID.wheel(), weight);
0266               wheelMeanHistos[3]->Fill(binSect, slID.wheel(), weight);
0267               wheelSigmaHistos[3]->Fill(binSect, slID.wheel(), weight);
0268             } else {
0269               if ((meanInRange(mean) < 0.85) && (sigmaInRange(sigma) > 0.85)) {  // only sigma ok
0270                 wheelSigmaHistos[3]->Fill(binSect, slID.wheel(), weight);
0271               }
0272               if ((meanInRange(mean) > 0.85) && (sigmaInRange(sigma) < 0.85)) {  // only mean ok
0273                 wheelMeanHistos[3]->Fill(binSect, slID.wheel(), weight);
0274               }
0275             }
0276           }
0277           delete gfit;
0278         } else {
0279           LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
0280               << "[DTResolutionAnalysisTask] Fit of " << slID << " not performed because # entries < 20 ";
0281           // FIXME: the SL is set as OK in the summary
0282           double weight = 1 / 11.;
0283           if ((binSect == 4 || binSect == 10) && slID.station() == 4)
0284             weight = 1 / 22.;
0285           globalResSummary->Fill(binSect, slID.wheel(), weight);
0286           wheelMeanHistos[3]->Fill(binSect, slID.wheel(), weight);
0287           wheelSigmaHistos[3]->Fill(binSect, slID.wheel(), weight);
0288           wheelMeanHistos[slID.wheel()]->setBinContent(binSect, binSL, 1.);
0289           wheelSigmaHistos[slID.wheel()]->setBinContent(binSect, binSL, 1.);
0290         }
0291       } else {
0292         LogWarning("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
0293             << "[DTResolutionAnalysisTask] Histo: " << getMEName(slID) << " not found" << endl;
0294       }
0295     }  // loop on SLs
0296   }  // Loop on Stations
0297 }
0298 
0299 void DTResolutionAnalysisTest::bookHistos(DQMStore::IBooker& ibooker, int wh) {
0300   stringstream wheel;
0301   wheel << wh;
0302 
0303   ibooker.setCurrentFolder(topHistoFolder + "/00-MeanRes");
0304   string histoName = "MeanSummaryRes_W" + wheel.str();
0305   string histoTitle = "# of SLs with wrong mean of residuals (Wheel " + wheel.str() + ")";
0306 
0307   wheelMeanHistos[wh] = ibooker.book2D(histoName.c_str(), histoTitle.c_str(), 12, 1, 13, 11, 1, 12);
0308   wheelMeanHistos[wh]->setAxisTitle("Sector", 1);
0309   wheelMeanHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0310   wheelMeanHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0311   wheelMeanHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0312   wheelMeanHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0313   wheelMeanHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0314   wheelMeanHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0315   wheelMeanHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0316   wheelMeanHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0317   wheelMeanHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0318   wheelMeanHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0319   wheelMeanHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0320 
0321   ibooker.setCurrentFolder(topHistoFolder + "/01-SigmaRes");
0322   histoName = "SigmaSummaryRes_W" + wheel.str();
0323   histoTitle = "# of SLs with wrong sigma of residuals (Wheel " + wheel.str() + ")";
0324 
0325   wheelSigmaHistos[wh] = ibooker.book2D(histoName.c_str(), histoTitle.c_str(), 12, 1, 13, 11, 1, 12);
0326   wheelSigmaHistos[wh]->setAxisTitle("Sector", 1);
0327   wheelSigmaHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0328   wheelSigmaHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0329   wheelSigmaHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0330   wheelSigmaHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0331   wheelSigmaHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0332   wheelSigmaHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0333   wheelSigmaHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0334   wheelSigmaHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0335   wheelSigmaHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0336   wheelSigmaHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0337   wheelSigmaHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0338 }
0339 
0340 void DTResolutionAnalysisTest::bookHistos(DQMStore::IBooker& ibooker, int wh, int sect) {
0341   stringstream wheel;
0342   wheel << wh;
0343   stringstream sector;
0344   sector << sect;
0345 
0346   string MeanHistoName = "MeanTest_W" + wheel.str() + "_Sec" + sector.str();
0347   string SigmaHistoName = "SigmaTest_W" + wheel.str() + "_Sec" + sector.str();
0348 
0349   string folder = topHistoFolder + "/Wheel" + wheel.str() + "/Sector" + sector.str();
0350   ibooker.setCurrentFolder(folder);
0351 
0352   if (sect != 4 && sect != 10) {
0353     MeanHistos[make_pair(wh, sect)] =
0354 
0355         ibooker.book1D(MeanHistoName.c_str(), "Mean (from gaussian fit) of the residuals distribution", 11, 1, 12);
0356   } else {
0357     MeanHistos[make_pair(wh, sect)] =
0358         ibooker.book1D(MeanHistoName.c_str(), "Mean (from gaussian fit) of the residuals distribution", 13, 1, 14);
0359   }
0360   (MeanHistos[make_pair(wh, sect)])->setBinLabel(1, "MB1_SL1", 1);
0361   (MeanHistos[make_pair(wh, sect)])->setBinLabel(2, "MB1_SL2", 1);
0362   (MeanHistos[make_pair(wh, sect)])->setBinLabel(3, "MB1_SL3", 1);
0363   (MeanHistos[make_pair(wh, sect)])->setBinLabel(4, "MB2_SL1", 1);
0364   (MeanHistos[make_pair(wh, sect)])->setBinLabel(5, "MB2_SL2", 1);
0365   (MeanHistos[make_pair(wh, sect)])->setBinLabel(6, "MB2_SL3", 1);
0366   (MeanHistos[make_pair(wh, sect)])->setBinLabel(7, "MB3_SL1", 1);
0367   (MeanHistos[make_pair(wh, sect)])->setBinLabel(8, "MB3_SL2", 1);
0368   (MeanHistos[make_pair(wh, sect)])->setBinLabel(9, "MB3_SL3", 1);
0369   (MeanHistos[make_pair(wh, sect)])->setBinLabel(10, "MB4_SL1", 1);
0370   (MeanHistos[make_pair(wh, sect)])->setBinLabel(11, "MB4_SL3", 1);
0371   if (sect == 4) {
0372     (MeanHistos[make_pair(wh, sect)])->setBinLabel(12, "MB4S13_SL1", 1);
0373     (MeanHistos[make_pair(wh, sect)])->setBinLabel(13, "MB4S13_SL3", 1);
0374   }
0375   if (sect == 10) {
0376     (MeanHistos[make_pair(wh, sect)])->setBinLabel(12, "MB4S14_SL1", 1);
0377     (MeanHistos[make_pair(wh, sect)])->setBinLabel(13, "MB4S14_SL3", 1);
0378   }
0379 
0380   if (sect != 4 && sect != 10) {
0381     SigmaHistos[make_pair(wh, sect)] =
0382         ibooker.book1D(SigmaHistoName.c_str(), "Sigma (from gaussian fit) of the residuals distribution", 11, 1, 12);
0383   } else {
0384     SigmaHistos[make_pair(wh, sect)] =
0385         ibooker.book1D(SigmaHistoName.c_str(), "Sigma (from gaussian fit) of the residuals distribution", 13, 1, 14);
0386   }
0387   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(1, "MB1_SL1", 1);
0388   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(2, "MB1_SL2", 1);
0389   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(3, "MB1_SL3", 1);
0390   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(4, "MB2_SL1", 1);
0391   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(5, "MB2_SL2", 1);
0392   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(6, "MB2_SL3", 1);
0393   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(7, "MB3_SL1", 1);
0394   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(8, "MB3_SL2", 1);
0395   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(9, "MB3_SL3", 1);
0396   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(10, "MB4_SL1", 1);
0397   (SigmaHistos[make_pair(wh, sect)])->setBinLabel(11, "MB4_SL3", 1);
0398   if (sect == 4) {
0399     (SigmaHistos[make_pair(wh, sect)])->setBinLabel(12, "MB4S13_SL1", 1);
0400     (SigmaHistos[make_pair(wh, sect)])->setBinLabel(13, "MB4S13_SL3", 1);
0401   }
0402   if (sect == 10) {
0403     (SigmaHistos[make_pair(wh, sect)])->setBinLabel(12, "MB4S14_SL1", 1);
0404     (SigmaHistos[make_pair(wh, sect)])->setBinLabel(13, "MB4S14_SL3", 1);
0405   }
0406 }
0407 
0408 string DTResolutionAnalysisTest::getMEName(const DTSuperLayerId& slID) {
0409   stringstream wheel;
0410   wheel << slID.wheel();
0411   stringstream station;
0412   station << slID.station();
0413   stringstream sector;
0414   sector << slID.sector();
0415   stringstream superLayer;
0416   superLayer << slID.superlayer();
0417 
0418   string folderName =
0419       topHistoFolder + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() + "/";
0420 
0421   if (doCalibAnalysis)
0422     folderName =
0423         "DT/DTCalibValidation/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" + sector.str() + "/";
0424 
0425   string histoname = folderName + "hResDist" + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +
0426                      "_SL" + superLayer.str();
0427 
0428   if (doCalibAnalysis)
0429     histoname = folderName + "hResDist_STEP3" + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +
0430                 "_SL" + superLayer.str();
0431 
0432   return histoname;
0433 }
0434 
0435 int DTResolutionAnalysisTest::stationFromBin(int bin) const { return (int)(bin / 3.1) + 1; }
0436 
0437 int DTResolutionAnalysisTest::slFromBin(int bin) const {
0438   int ret = bin % 3;
0439   if (ret == 0 || bin == 11)
0440     ret = 3;
0441 
0442   return ret;
0443 }
0444 
0445 double DTResolutionAnalysisTest::meanInRange(double mean) const {
0446   double value(0.);
0447   if (fabs(mean) <= maxGoodMeanValue) {
0448     value = 1.;
0449   } else if (fabs(mean) > maxGoodMeanValue && fabs(mean) < minBadMeanValue) {
0450     value = 0.9;
0451   } else if (fabs(mean) >= minBadMeanValue) {
0452     value = 0.1;
0453   }
0454   return value;
0455 }
0456 
0457 double DTResolutionAnalysisTest::sigmaInRange(double sigma) const {
0458   double value(0.);
0459   if (sigma <= maxGoodSigmaValue) {
0460     value = 1.;
0461   } else if (sigma > maxGoodSigmaValue && sigma < minBadSigmaValue) {
0462     value = 0.9;
0463   } else if (sigma >= minBadSigmaValue) {
0464     value = 0.1;
0465   }
0466   return value;
0467 }
0468 
0469 void DTResolutionAnalysisTest::resetMEs() {
0470   globalResSummary->Reset();
0471   // Reset the summary histo
0472   for (map<int, MonitorElement*>::const_iterator histo = wheelMeanHistos.begin(); histo != wheelMeanHistos.end();
0473        histo++) {
0474     (*histo).second->Reset();
0475   }
0476   for (map<int, MonitorElement*>::const_iterator histo = wheelSigmaHistos.begin(); histo != wheelSigmaHistos.end();
0477        histo++) {
0478     (*histo).second->Reset();
0479   }
0480 
0481   for (int indx = -2; indx != 3; ++indx) {
0482     meanDistr[indx]->Reset();
0483     sigmaDistr[indx]->Reset();
0484   }
0485 }