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/DTEfficiencyTest.h"
0011 
0012 // Framework
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 
0015 // Geometry
0016 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0017 #include "Geometry/DTGeometry/interface/DTLayer.h"
0018 #include "Geometry/DTGeometry/interface/DTTopology.h"
0019 
0020 #include "DQMServices/Core/interface/DQMStore.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 
0023 #include <cstdio>
0024 #include <sstream>
0025 #include <cmath>
0026 
0027 using namespace edm;
0028 using namespace std;
0029 
0030 DTEfficiencyTest::DTEfficiencyTest(const edm::ParameterSet& ps)
0031     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0032   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest]: Constructor";
0033 
0034   parameters = ps;
0035 
0036   prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
0037 
0038   percentual = parameters.getUntrackedParameter<int>("BadSLpercentual", 10);
0039 }
0040 
0041 DTEfficiencyTest::~DTEfficiencyTest() {
0042   edm::LogVerbatim("efficiency") << "DTEfficiencyTest: analyzed " << nevents << " events";
0043 }
0044 
0045 void DTEfficiencyTest::beginRun(edm::Run const& run, edm::EventSetup const& context) {
0046   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest]: Begin run";
0047 
0048   nevents = 0;
0049 
0050   // Get the geometry
0051   muonGeom = &context.getData(muonGeomToken_);
0052 }
0053 
0054 void DTEfficiencyTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0055                                              DQMStore::IGetter& igetter,
0056                                              edm::LuminosityBlock const& lumiSeg,
0057                                              edm::EventSetup const& context) {
0058   for (map<int, MonitorElement*>::const_iterator histo = wheelHistos.begin(); histo != wheelHistos.end(); histo++) {
0059     (*histo).second->Reset();
0060   }
0061 
0062   for (map<int, MonitorElement*>::const_iterator histo = wheelUnassHistos.begin(); histo != wheelUnassHistos.end();
0063        histo++) {
0064     (*histo).second->Reset();
0065   }
0066 
0067   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest]: End of LS transition, performing the DQM client operation";
0068 
0069   // counts number of lumiSegs
0070   nLumiSegs = lumiSeg.id().luminosityBlock();
0071 
0072   // prescale factor
0073   if (nLumiSegs % prescaleFactor != 0)
0074     return;
0075 
0076   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest]: " << nLumiSegs << " updates";
0077 
0078   vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0079   vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0080 
0081   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest]: Efficiency tests results";
0082 
0083   map<DTLayerId, vector<double> > LayerBadCells;
0084   LayerBadCells.clear();
0085   map<DTLayerId, vector<double> > LayerUnassBadCells;
0086   LayerUnassBadCells.clear();
0087   map<DTSuperLayerId, vector<double> > SuperLayerBadCells;
0088   SuperLayerBadCells.clear();
0089   map<DTSuperLayerId, vector<double> > SuperLayerUnassBadCells;
0090   SuperLayerUnassBadCells.clear();
0091   map<pair<int, int>, int> cmsHistos;
0092   cmsHistos.clear();
0093   map<pair<int, int>, bool> filled;
0094   for (int i = -2; i < 3; i++) {
0095     for (int j = 1; j < 15; j++) {
0096       filled[make_pair(i, j)] = false;
0097     }
0098   }
0099   map<pair<int, int>, int> cmsUnassHistos;
0100   cmsUnassHistos.clear();
0101   map<pair<int, int>, bool> UnassFilled;
0102   for (int i = -2; i < 3; i++) {
0103     for (int j = 1; j < 15; j++) {
0104       UnassFilled[make_pair(i, j)] = false;
0105     }
0106   }
0107 
0108   // Loop over the chambers
0109   for (; ch_it != ch_end; ++ch_it) {
0110     DTChamberId chID = (*ch_it)->id();
0111     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
0112     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
0113 
0114     // Loop over the SuperLayers
0115     for (; sl_it != sl_end; ++sl_it) {
0116       DTSuperLayerId slID = (*sl_it)->id();
0117       vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin();
0118       vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
0119 
0120       // Loop over the layers
0121       for (; l_it != l_end; ++l_it) {
0122         DTLayerId lID = (*l_it)->id();
0123 
0124         stringstream wheel;
0125         wheel << chID.wheel();
0126         stringstream station;
0127         station << chID.station();
0128         stringstream sector;
0129         sector << chID.sector();
0130         stringstream superLayer;
0131         superLayer << slID.superlayer();
0132         stringstream layer;
0133         layer << lID.layer();
0134 
0135         string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" +
0136                            superLayer.str() + "_L" + layer.str();
0137 
0138         // Get the ME produced by EfficiencyTask Source
0139         MonitorElement* occupancy_histo = igetter.get(getMEName("hEffOccupancy", lID));
0140         MonitorElement* unassOccupancy_histo = igetter.get(getMEName("hEffUnassOccupancy", lID));
0141         MonitorElement* recSegmOccupancy_histo = igetter.get(getMEName("hRecSegmOccupancy", lID));
0142 
0143         // ME -> TH1F
0144         if (occupancy_histo && unassOccupancy_histo && recSegmOccupancy_histo) {
0145           TH1F* occupancy_histo_root = occupancy_histo->getTH1F();
0146           TH1F* unassOccupancy_histo_root = unassOccupancy_histo->getTH1F();
0147           TH1F* recSegmOccupancy_histo_root = recSegmOccupancy_histo->getTH1F();
0148 
0149           const int firstWire = muonGeom->layer(lID)->specificTopology().firstChannel();
0150           const int lastWire = muonGeom->layer(lID)->specificTopology().lastChannel();
0151 
0152           // Loop over the TH1F bin and fill the ME to be used for the Quality Test
0153           for (int bin = firstWire; bin <= lastWire; bin++) {
0154             if ((recSegmOccupancy_histo_root->GetBinContent(bin)) != 0) {
0155               if (EfficiencyHistos.find(lID) == EfficiencyHistos.end())
0156                 bookHistos(ibooker, lID, firstWire, lastWire);
0157               float efficiency =
0158                   occupancy_histo_root->GetBinContent(bin) / recSegmOccupancy_histo_root->GetBinContent(bin);
0159               float errorEff = sqrt(efficiency * (1 - efficiency) / recSegmOccupancy_histo_root->GetBinContent(bin));
0160               EfficiencyHistos.find(lID)->second->setBinContent(bin, efficiency);
0161               EfficiencyHistos.find(lID)->second->setBinError(bin, errorEff);
0162 
0163               if (UnassEfficiencyHistos.find(lID) == EfficiencyHistos.end())
0164                 bookHistos(ibooker, lID, firstWire, lastWire);
0165               float unassEfficiency =
0166                   unassOccupancy_histo_root->GetBinContent(bin) / recSegmOccupancy_histo_root->GetBinContent(bin);
0167               float errorUnassEff =
0168                   sqrt(unassEfficiency * (1 - unassEfficiency) / recSegmOccupancy_histo_root->GetBinContent(bin));
0169               UnassEfficiencyHistos.find(lID)->second->setBinContent(bin, unassEfficiency);
0170               UnassEfficiencyHistos.find(lID)->second->setBinError(bin, errorUnassEff);
0171             }
0172           }
0173         }
0174       }  // loop on layers
0175     }  // loop on superlayers
0176   }  //loop on chambers
0177 
0178   // Efficiency test
0179   //cout<<"[DTEfficiencyTest]: Efficiency Tests results"<<endl;
0180   string EfficiencyCriterionName = parameters.getUntrackedParameter<string>("EfficiencyTestName", "EfficiencyInRange");
0181   for (map<DTLayerId, MonitorElement*>::const_iterator hEff = EfficiencyHistos.begin(); hEff != EfficiencyHistos.end();
0182        hEff++) {
0183     const QReport* theEfficiencyQReport = (*hEff).second->getQReport(EfficiencyCriterionName);
0184     double counter = 0;
0185     if (theEfficiencyQReport) {
0186       vector<dqm::me_util::Channel> badChannels = theEfficiencyQReport->getBadChannels();
0187       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0188            channel++) {
0189         edm::LogError("efficiency") << "LayerID : " << getMEName("hEffOccupancy", (*hEff).first)
0190                                     << " Bad efficiency channels: " << (*channel).getBin()
0191                                     << "  Contents : " << (*channel).getContents();
0192         counter++;
0193       }
0194       LayerBadCells[(*hEff).first].push_back(counter);
0195       LayerBadCells[(*hEff).first].push_back(muonGeom->layer((*hEff).first)->specificTopology().channels());
0196       // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
0197       // edm::LogWarning ("efficiency") << "-------- "<<theEfficiencyQReport->getMessage()<<" ------- "<<theEfficiencyQReport->getStatus();
0198     }
0199   }
0200 
0201   // UnassEfficiency test
0202   //cout<<"[DTEfficiencyTest]: UnassEfficiency Tests results"<<endl;
0203   string UnassEfficiencyCriterionName =
0204       parameters.getUntrackedParameter<string>("UnassEfficiencyTestName", "UnassEfficiencyInRange");
0205   for (map<DTLayerId, MonitorElement*>::const_iterator hUnassEff = UnassEfficiencyHistos.begin();
0206        hUnassEff != UnassEfficiencyHistos.end();
0207        hUnassEff++) {
0208     const QReport* theUnassEfficiencyQReport = (*hUnassEff).second->getQReport(UnassEfficiencyCriterionName);
0209     double counter = 0;
0210     if (theUnassEfficiencyQReport) {
0211       vector<dqm::me_util::Channel> badChannels = theUnassEfficiencyQReport->getBadChannels();
0212       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0213            channel++) {
0214         edm::LogError("efficiency") << "Bad unassEfficiency channels: " << (*channel).getBin() << " "
0215                                     << (*channel).getContents();
0216         counter++;
0217       }
0218       LayerUnassBadCells[(*hUnassEff).first].push_back(counter);
0219       LayerUnassBadCells[(*hUnassEff).first].push_back(
0220           double(muonGeom->layer((*hUnassEff).first)->specificTopology().channels()));
0221       // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
0222       // edm::LogWarning ("efficiency") << theUnassEfficiencyQReport->getMessage()<<" ------- "<<theUnassEfficiencyQReport->getStatus();
0223     }
0224   }
0225 
0226   vector<const DTChamber*>::const_iterator ch2_it = muonGeom->chambers().begin();
0227   vector<const DTChamber*>::const_iterator ch2_end = muonGeom->chambers().end();
0228   for (; ch2_it != ch2_end; ++ch2_it) {
0229     vector<const DTSuperLayer*>::const_iterator sl2_it = (*ch2_it)->superLayers().begin();
0230     vector<const DTSuperLayer*>::const_iterator sl2_end = (*ch2_it)->superLayers().end();
0231     // Loop over the SLs
0232     for (; sl2_it != sl2_end; ++sl2_it) {
0233       DTSuperLayerId sl = (*sl2_it)->id();
0234       double superLayerBadC = 0;
0235       double superLayerTotC = 0;
0236       double superLayerUnassBadC = 0;
0237       double superLayerUnassTotC = 0;
0238       bool fill = false;
0239       vector<const DTLayer*>::const_iterator l2_it = (*sl2_it)->layers().begin();
0240       vector<const DTLayer*>::const_iterator l2_end = (*sl2_it)->layers().end();
0241       // Loop over the Ls
0242       for (; l2_it != l2_end; ++l2_it) {
0243         DTLayerId layerId = (*l2_it)->id();
0244         if (LayerBadCells.find(layerId) != LayerBadCells.end() &&
0245             LayerUnassBadCells.find(layerId) != LayerUnassBadCells.end()) {
0246           fill = true;
0247           superLayerBadC += LayerBadCells[layerId][0];
0248           superLayerTotC += LayerBadCells[layerId][1];
0249           superLayerUnassBadC += LayerUnassBadCells[layerId][0];
0250           superLayerUnassTotC += LayerUnassBadCells[layerId][1];
0251         }
0252       }
0253       if (fill) {
0254         SuperLayerBadCells[sl].push_back(superLayerBadC);
0255         SuperLayerBadCells[sl].push_back(superLayerTotC);
0256         SuperLayerUnassBadCells[sl].push_back(superLayerUnassBadC);
0257         SuperLayerUnassBadCells[sl].push_back(superLayerUnassTotC);
0258       }
0259     }
0260   }
0261 
0262   for (map<DTSuperLayerId, vector<double> >::const_iterator SLBCells = SuperLayerBadCells.begin();
0263        SLBCells != SuperLayerBadCells.end();
0264        SLBCells++) {
0265     if ((*SLBCells).second[0] / (*SLBCells).second[1] > double(percentual / 100)) {
0266       if (wheelHistos.find((*SLBCells).first.wheel()) == wheelHistos.end())
0267         bookHistos(ibooker, (*SLBCells).first.wheel());
0268       if (!((*SLBCells).first.station() == 4 && (*SLBCells).first.superlayer() == 3))
0269         wheelHistos[(*SLBCells).first.wheel()]->Fill(
0270             (*SLBCells).first.sector() - 1,
0271             ((*SLBCells).first.superlayer() - 1) + 3 * ((*SLBCells).first.station() - 1));
0272       else
0273         wheelHistos[(*SLBCells).first.wheel()]->Fill((*SLBCells).first.sector() - 1, 10);
0274       // fill the cms summary histo if the percentual of SL which have not passed the test
0275       // is more than a predefined treshold
0276       cmsHistos[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())]++;
0277       if (((*SLBCells).first.sector() < 13 &&
0278            double(cmsHistos[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())]) / 11 >
0279                double(percentual) / 100 &&
0280            filled[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())] == false) ||
0281           ((*SLBCells).first.sector() >= 13 &&
0282            double(cmsHistos[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())]) / 2 >
0283                double(percentual) / 100 &&
0284            filled[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())] == false)) {
0285         filled[make_pair((*SLBCells).first.wheel(), (*SLBCells).first.sector())] = true;
0286         wheelHistos[3]->Fill((*SLBCells).first.sector() - 1, (*SLBCells).first.wheel());
0287       }
0288     }
0289   }
0290 
0291   for (map<DTSuperLayerId, vector<double> >::const_iterator SLUBCells = SuperLayerUnassBadCells.begin();
0292        SLUBCells != SuperLayerUnassBadCells.end();
0293        SLUBCells++) {
0294     if ((*SLUBCells).second[0] / (*SLUBCells).second[1] > double(percentual / 100)) {
0295       if (wheelUnassHistos.find((*SLUBCells).first.wheel()) == wheelUnassHistos.end())
0296         bookHistos(ibooker, (*SLUBCells).first.wheel());
0297       if (!((*SLUBCells).first.station() == 4 && (*SLUBCells).first.superlayer() == 3))
0298         wheelUnassHistos[(*SLUBCells).first.wheel()]->Fill(
0299             (*SLUBCells).first.sector() - 1,
0300             ((*SLUBCells).first.superlayer() - 1) + 3 * ((*SLUBCells).first.station() - 1));
0301       else
0302         wheelUnassHistos[(*SLUBCells).first.wheel()]->Fill((*SLUBCells).first.sector() - 1, 10);
0303       // fill the cms summary histo if the percentual of SL which have not passed the test
0304       // is more than a predefined treshold
0305       cmsUnassHistos[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())]++;
0306       if (((*SLUBCells).first.sector() < 13 &&
0307            double(cmsUnassHistos[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())]) / 11 >
0308                double(percentual) / 100 &&
0309            UnassFilled[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())] == false) ||
0310           ((*SLUBCells).first.sector() >= 13 &&
0311            double(cmsUnassHistos[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())]) / 2 >
0312                double(percentual) / 100 &&
0313            UnassFilled[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())] == false)) {
0314         UnassFilled[make_pair((*SLUBCells).first.wheel(), (*SLUBCells).first.sector())] = true;
0315         wheelUnassHistos[3]->Fill((*SLUBCells).first.sector() - 1, (*SLUBCells).first.wheel());
0316       }
0317     }
0318   }
0319 }
0320 
0321 void DTEfficiencyTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0322   edm::LogVerbatim("efficiency") << "[DTEfficiencyTest] endjob called!";
0323 }
0324 
0325 string DTEfficiencyTest::getMEName(string histoTag, const DTLayerId& lID) {
0326   stringstream wheel;
0327   wheel << lID.superlayerId().wheel();
0328   stringstream station;
0329   station << lID.superlayerId().station();
0330   stringstream sector;
0331   sector << lID.superlayerId().sector();
0332   stringstream superLayer;
0333   superLayer << lID.superlayerId().superlayer();
0334   stringstream layer;
0335   layer << lID.layer();
0336 
0337   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
0338   string folderName = folderRoot + "DT/DTEfficiencyTask/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" +
0339                       sector.str() + "/SuperLayer" + superLayer.str() + "/";
0340 
0341   string histoname = folderName + histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +
0342                      "_SL" + superLayer.str() + "_L" + layer.str();
0343 
0344   return histoname;
0345 }
0346 
0347 void DTEfficiencyTest::bookHistos(DQMStore::IBooker& ibooker, const DTLayerId& lId, int firstWire, int lastWire) {
0348   stringstream wheel;
0349   wheel << lId.superlayerId().wheel();
0350   stringstream station;
0351   station << lId.superlayerId().station();
0352   stringstream sector;
0353   sector << lId.superlayerId().sector();
0354   stringstream superLayer;
0355   superLayer << lId.superlayerId().superlayer();
0356   stringstream layer;
0357   layer << lId.layer();
0358 
0359   string HistoName =
0360       "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str() + "_L" + layer.str();
0361   string EfficiencyHistoName = "Efficiency_" + HistoName;
0362   string UnassEfficiencyHistoName = "UnassEfficiency_" + HistoName;
0363 
0364   ibooker.setCurrentFolder("DT/Tests/DTEfficiency/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" +
0365                            sector.str());
0366 
0367   EfficiencyHistos[lId] = ibooker.book1D(EfficiencyHistoName.c_str(),
0368                                          EfficiencyHistoName.c_str(),
0369                                          lastWire - firstWire + 1,
0370                                          firstWire - 0.5,
0371                                          lastWire + 0.5);
0372   UnassEfficiencyHistos[lId] = ibooker.book1D(UnassEfficiencyHistoName.c_str(),
0373                                               UnassEfficiencyHistoName.c_str(),
0374                                               lastWire - firstWire + 1,
0375                                               firstWire - 0.5,
0376                                               lastWire + 0.5);
0377 }
0378 
0379 void DTEfficiencyTest::bookHistos(DQMStore::IBooker& ibooker, int wh) {
0380   ibooker.setCurrentFolder("DT/Tests/DTEfficiency/SummaryPlot");
0381 
0382   if (wheelHistos.find(3) == wheelHistos.end()) {
0383     string histoName = "ESummary_testFailedByAtLeastBadSL";
0384     wheelHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 5, -2, 2);
0385     wheelHistos[3]->setBinLabel(1, "Sector1", 1);
0386     wheelHistos[3]->setBinLabel(1, "Sector1", 1);
0387     wheelHistos[3]->setBinLabel(2, "Sector2", 1);
0388     wheelHistos[3]->setBinLabel(3, "Sector3", 1);
0389     wheelHistos[3]->setBinLabel(4, "Sector4", 1);
0390     wheelHistos[3]->setBinLabel(5, "Sector5", 1);
0391     wheelHistos[3]->setBinLabel(6, "Sector6", 1);
0392     wheelHistos[3]->setBinLabel(7, "Sector7", 1);
0393     wheelHistos[3]->setBinLabel(8, "Sector8", 1);
0394     wheelHistos[3]->setBinLabel(9, "Sector9", 1);
0395     wheelHistos[3]->setBinLabel(10, "Sector10", 1);
0396     wheelHistos[3]->setBinLabel(11, "Sector11", 1);
0397     wheelHistos[3]->setBinLabel(12, "Sector12", 1);
0398     wheelHistos[3]->setBinLabel(13, "Sector13", 1);
0399     wheelHistos[3]->setBinLabel(14, "Sector14", 1);
0400     wheelHistos[3]->setBinLabel(1, "Wheel-2", 2);
0401     wheelHistos[3]->setBinLabel(2, "Wheel-1", 2);
0402     wheelHistos[3]->setBinLabel(3, "Wheel0", 2);
0403     wheelHistos[3]->setBinLabel(4, "Wheel+1", 2);
0404     wheelHistos[3]->setBinLabel(5, "Wheel+2", 2);
0405   }
0406   if (wheelUnassHistos.find(3) == wheelUnassHistos.end()) {
0407     string histoName = "UESummary_testFailedByAtLeastBadSL";
0408     wheelUnassHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 5, -2, 2);
0409     wheelUnassHistos[3]->setBinLabel(1, "Sector1", 1);
0410     wheelUnassHistos[3]->setBinLabel(1, "Sector1", 1);
0411     wheelUnassHistos[3]->setBinLabel(2, "Sector2", 1);
0412     wheelUnassHistos[3]->setBinLabel(3, "Sector3", 1);
0413     wheelUnassHistos[3]->setBinLabel(4, "Sector4", 1);
0414     wheelUnassHistos[3]->setBinLabel(5, "Sector5", 1);
0415     wheelUnassHistos[3]->setBinLabel(6, "Sector6", 1);
0416     wheelUnassHistos[3]->setBinLabel(7, "Sector7", 1);
0417     wheelUnassHistos[3]->setBinLabel(8, "Sector8", 1);
0418     wheelUnassHistos[3]->setBinLabel(9, "Sector9", 1);
0419     wheelUnassHistos[3]->setBinLabel(10, "Sector10", 1);
0420     wheelUnassHistos[3]->setBinLabel(11, "Sector11", 1);
0421     wheelUnassHistos[3]->setBinLabel(12, "Sector12", 1);
0422     wheelUnassHistos[3]->setBinLabel(13, "Sector13", 1);
0423     wheelUnassHistos[3]->setBinLabel(14, "Sector14", 1);
0424     wheelUnassHistos[3]->setBinLabel(1, "Wheel-2", 2);
0425     wheelUnassHistos[3]->setBinLabel(2, "Wheel-1", 2);
0426     wheelUnassHistos[3]->setBinLabel(3, "Wheel0", 2);
0427     wheelUnassHistos[3]->setBinLabel(4, "Wheel+1", 2);
0428     wheelUnassHistos[3]->setBinLabel(5, "Wheel+2", 2);
0429   }
0430 
0431   stringstream wheel;
0432   wheel << wh;
0433 
0434   if (wheelHistos.find(wh) == wheelHistos.end()) {
0435     string histoName = "ESummary_testFailed_W" + wheel.str();
0436     wheelHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 11, 0, 11);
0437     wheelHistos[wh]->setBinLabel(1, "Sector1", 1);
0438     wheelHistos[wh]->setBinLabel(2, "Sector2", 1);
0439     wheelHistos[wh]->setBinLabel(3, "Sector3", 1);
0440     wheelHistos[wh]->setBinLabel(4, "Sector4", 1);
0441     wheelHistos[wh]->setBinLabel(5, "Sector5", 1);
0442     wheelHistos[wh]->setBinLabel(6, "Sector6", 1);
0443     wheelHistos[wh]->setBinLabel(7, "Sector7", 1);
0444     wheelHistos[wh]->setBinLabel(8, "Sector8", 1);
0445     wheelHistos[wh]->setBinLabel(9, "Sector9", 1);
0446     wheelHistos[wh]->setBinLabel(10, "Sector10", 1);
0447     wheelHistos[wh]->setBinLabel(11, "Sector11", 1);
0448     wheelHistos[wh]->setBinLabel(12, "Sector12", 1);
0449     wheelHistos[wh]->setBinLabel(13, "Sector13", 1);
0450     wheelHistos[wh]->setBinLabel(14, "Sector14", 1);
0451     wheelHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0452     wheelHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0453     wheelHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0454     wheelHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0455     wheelHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0456     wheelHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0457     wheelHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0458     wheelHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0459     wheelHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0460     wheelHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0461     wheelHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0462   }
0463   if (wheelUnassHistos.find(wh) == wheelUnassHistos.end()) {
0464     string histoName = "UESummary_testFailed_W" + wheel.str();
0465     wheelUnassHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 14, 0, 14, 11, 0, 11);
0466     wheelUnassHistos[wh]->setBinLabel(1, "Sector1", 1);
0467     wheelUnassHistos[wh]->setBinLabel(2, "Sector2", 1);
0468     wheelUnassHistos[wh]->setBinLabel(3, "Sector3", 1);
0469     wheelUnassHistos[wh]->setBinLabel(4, "Sector4", 1);
0470     wheelUnassHistos[wh]->setBinLabel(5, "Sector5", 1);
0471     wheelUnassHistos[wh]->setBinLabel(6, "Sector6", 1);
0472     wheelUnassHistos[wh]->setBinLabel(7, "Sector7", 1);
0473     wheelUnassHistos[wh]->setBinLabel(8, "Sector8", 1);
0474     wheelUnassHistos[wh]->setBinLabel(9, "Sector9", 1);
0475     wheelUnassHistos[wh]->setBinLabel(10, "Sector10", 1);
0476     wheelUnassHistos[wh]->setBinLabel(11, "Sector11", 1);
0477     wheelUnassHistos[wh]->setBinLabel(12, "Sector12", 1);
0478     wheelUnassHistos[wh]->setBinLabel(13, "Sector13", 1);
0479     wheelUnassHistos[wh]->setBinLabel(14, "Sector14", 1);
0480     wheelUnassHistos[wh]->setBinLabel(1, "MB1_SL1", 2);
0481     wheelUnassHistos[wh]->setBinLabel(2, "MB1_SL2", 2);
0482     wheelUnassHistos[wh]->setBinLabel(3, "MB1_SL3", 2);
0483     wheelUnassHistos[wh]->setBinLabel(4, "MB2_SL1", 2);
0484     wheelUnassHistos[wh]->setBinLabel(5, "MB2_SL2", 2);
0485     wheelUnassHistos[wh]->setBinLabel(6, "MB2_SL3", 2);
0486     wheelUnassHistos[wh]->setBinLabel(7, "MB3_SL1", 2);
0487     wheelUnassHistos[wh]->setBinLabel(8, "MB3_SL2", 2);
0488     wheelUnassHistos[wh]->setBinLabel(9, "MB3_SL3", 2);
0489     wheelUnassHistos[wh]->setBinLabel(10, "MB4_SL1", 2);
0490     wheelUnassHistos[wh]->setBinLabel(11, "MB4_SL3", 2);
0491   }
0492 }