Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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/DTChamberEfficiencyTest.h"
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 
0013 // Framework
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 
0016 // Geometry
0017 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0018 
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 
0021 #include <cstdio>
0022 #include <sstream>
0023 #include <cmath>
0024 
0025 using namespace edm;
0026 using namespace std;
0027 
0028 DTChamberEfficiencyTest::DTChamberEfficiencyTest(const edm::ParameterSet& ps)
0029     : muonGeomToken_(esConsumes<edm::Transition::EndLuminosityBlock>()) {
0030   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest]: Constructor";
0031 
0032   parameters = ps;
0033 
0034   prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
0035 
0036   nevents = 0;
0037 
0038   bookingdone = false;
0039 }
0040 
0041 DTChamberEfficiencyTest::~DTChamberEfficiencyTest() {
0042   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0043       << "DTChamberEfficiencyTest: analyzed " << nevents << " events";
0044 }
0045 
0046 void DTChamberEfficiencyTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0047                                                     DQMStore::IGetter& igetter,
0048                                                     edm::LuminosityBlock const& lumiSeg,
0049                                                     edm::EventSetup const& context) {
0050   if (!bookingdone) {
0051     // Get the DT Geometry
0052     muonGeom = &context.getData(muonGeomToken_);
0053 
0054     // Loop over all the chambers
0055     vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0056     vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0057     for (; ch_it != ch_end; ++ch_it) {
0058       // histo booking
0059       bookHistos(ibooker, (*ch_it)->id());
0060     }
0061 
0062     //summary histo booking
0063     bookHistos(ibooker);
0064   }
0065 
0066   bookingdone = true;
0067 
0068   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0069       << "[DTChamberEfficiencyTest]: End of LS transition, performing the DQM client operation";
0070 
0071   // counts number of lumiSegs
0072   nLumiSegs = lumiSeg.id().luminosityBlock();
0073 
0074   // prescale factor
0075   if (nLumiSegs % prescaleFactor != 0)
0076     return;
0077 
0078   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0079       << "[DTChamberEfficiencyTest]: " << nLumiSegs << " updates";
0080 
0081   vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0082   vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0083 
0084   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0085       << "[DTChamberEfficiencyTest]: ChamberEfficiency tests results";
0086 
0087   // Loop over the chambers
0088   for (; ch_it != ch_end; ++ch_it) {
0089     DTChamberId chID = (*ch_it)->id();
0090 
0091     stringstream wheel;
0092     wheel << chID.wheel();
0093     stringstream station;
0094     station << chID.station();
0095     stringstream sector;
0096     sector << chID.sector();
0097 
0098     string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0099 
0100     // Get the ME produced by EfficiencyTask Source
0101 
0102     MonitorElement* GoodSegDen_histo = igetter.get(getMEName("hEffGoodSegVsPosDen", chID));
0103     MonitorElement* GoodCloseSegNum_histo = igetter.get(getMEName("hEffGoodCloseSegVsPosNum", chID));
0104 
0105     // ME -> TH1F
0106     if (GoodSegDen_histo && GoodCloseSegNum_histo) {
0107       TH2F* GoodSegDen_histo_root = GoodSegDen_histo->getTH2F();
0108       TH2F* GoodCloseSegNum_histo_root = GoodCloseSegNum_histo->getTH2F();
0109 
0110       int lastBinX = (*GoodSegDen_histo_root).GetNbinsX();
0111       TH1D* proxN = GoodCloseSegNum_histo_root->ProjectionX();
0112       TH1D* proxD = GoodSegDen_histo_root->ProjectionX();
0113 
0114       int lastBinY = (*GoodSegDen_histo_root).GetNbinsY();
0115       TH1D* proyN = GoodCloseSegNum_histo_root->ProjectionY();
0116       TH1D* proyD = GoodSegDen_histo_root->ProjectionY();
0117 
0118       for (int xBin = 1; xBin <= lastBinX; xBin++) {
0119         if (proxD->GetBinContent(xBin) != 0) {
0120           float Xefficiency = proxN->GetBinContent(xBin) / proxD->GetBinContent(xBin);
0121           xEfficiencyHistos.find(HistoName)->second->setBinContent(xBin, Xefficiency);
0122         }
0123 
0124         for (int yBin = 1; yBin <= lastBinY; yBin++) {
0125           if (GoodSegDen_histo_root->GetBinContent(xBin, yBin) != 0) {
0126             float XvsYefficiency = GoodCloseSegNum_histo_root->GetBinContent(xBin, yBin) /
0127                                    GoodSegDen_histo_root->GetBinContent(xBin, yBin);
0128             xVSyEffHistos.find(HistoName)->second->setBinContent(xBin, yBin, XvsYefficiency);
0129           }
0130         }
0131       }
0132 
0133       for (int yBin = 1; yBin <= lastBinY; yBin++) {
0134         if (proyD->GetBinContent(yBin) != 0) {
0135           float Yefficiency = proyN->GetBinContent(yBin) / proyD->GetBinContent(yBin);
0136           yEfficiencyHistos.find(HistoName)->second->setBinContent(yBin, Yefficiency);
0137         }
0138       }
0139     }
0140   }  // loop on chambers
0141 
0142   // ChamberEfficiency test on X axis
0143   string XEfficiencyCriterionName =
0144       parameters.getUntrackedParameter<string>("XEfficiencyTestName", "ChEfficiencyInRangeX");
0145   for (map<string, MonitorElement*>::const_iterator hXEff = xEfficiencyHistos.begin(); hXEff != xEfficiencyHistos.end();
0146        hXEff++) {
0147     const QReport* theXEfficiencyQReport = (*hXEff).second->getQReport(XEfficiencyCriterionName);
0148     if (theXEfficiencyQReport) {
0149       vector<dqm::me_util::Channel> badChannels = theXEfficiencyQReport->getBadChannels();
0150       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0151            channel++) {
0152         edm::LogError("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0153             << "Chamber : " << (*hXEff).first << " Bad XChamberEfficiency channels: " << (*channel).getBin()
0154             << "  Contents : " << (*channel).getContents();
0155       }
0156     }
0157   }
0158 
0159   // ChamberEfficiency test on Y axis
0160   string YEfficiencyCriterionName =
0161       parameters.getUntrackedParameter<string>("YEfficiencyTestName", "ChEfficiencyInRangeY");
0162   for (map<string, MonitorElement*>::const_iterator hYEff = yEfficiencyHistos.begin(); hYEff != yEfficiencyHistos.end();
0163        hYEff++) {
0164     const QReport* theYEfficiencyQReport = (*hYEff).second->getQReport(YEfficiencyCriterionName);
0165     if (theYEfficiencyQReport) {
0166       vector<dqm::me_util::Channel> badChannels = theYEfficiencyQReport->getBadChannels();
0167       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0168            channel++) {
0169         edm::LogError("DTDQM|DTMonitorClient|DTChamberEfficiencyTest")
0170             << "Chamber : " << (*hYEff).first << " Bad YChamberEfficiency channels: " << (*channel).getBin()
0171             << "  Contents : " << (*channel).getContents();
0172       }
0173     }
0174   }
0175 
0176   //Fill the report summary histos
0177   for (int wh = -2; wh <= 2; wh++) {
0178     for (int sec = 1; sec <= 12; sec++) {
0179       for (int st = 1; st <= 4; st++) {
0180         summaryHistos[wh]->Fill(sec, st, 1);
0181       }
0182     }
0183   }
0184 }
0185 
0186 void DTChamberEfficiencyTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0187   edm::LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest] endjob called!";
0188 }
0189 
0190 string DTChamberEfficiencyTest::getMEName(string histoTag, const DTChamberId& chID) {
0191   stringstream wheel;
0192   wheel << chID.wheel();
0193   stringstream station;
0194   station << chID.station();
0195   stringstream sector;
0196   sector << chID.sector();
0197 
0198   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
0199   string folderName = folderRoot + "DT/01-DTChamberEfficiency/Task/Wheel" + wheel.str() + "/Sector" + sector.str() +
0200                       "/Station" + station.str() + "/";
0201 
0202   string histoname = folderName + histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0203 
0204   return histoname;
0205 }
0206 
0207 void DTChamberEfficiencyTest::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& chId) {
0208   stringstream wheel;
0209   wheel << chId.wheel();
0210   stringstream station;
0211   station << chId.station();
0212   stringstream sector;
0213   sector << chId.sector();
0214 
0215   string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0216   string xEfficiencyHistoName = "xEfficiency_" + HistoName;
0217   string yEfficiencyHistoName = "yEfficiency_" + HistoName;
0218   string xVSyEffHistoName = "xVSyEff_" + HistoName;
0219 
0220   ibooker.setCurrentFolder("DT/01-DTChamberEfficiency/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
0221                            station.str());
0222 
0223   xEfficiencyHistos[HistoName] =
0224       ibooker.book1D(xEfficiencyHistoName.c_str(), xEfficiencyHistoName.c_str(), 25, -250., 250.);
0225   yEfficiencyHistos[HistoName] =
0226       ibooker.book1D(yEfficiencyHistoName.c_str(), yEfficiencyHistoName.c_str(), 25, -250., 250.);
0227   xVSyEffHistos[HistoName] =
0228       ibooker.book2D(xVSyEffHistoName.c_str(), xVSyEffHistoName.c_str(), 25, -250., 250., 25, -250., 250.);
0229 }
0230 
0231 void DTChamberEfficiencyTest::bookHistos(DQMStore::IBooker& ibooker) {
0232   for (int wh = -2; wh <= 2; wh++) {
0233     stringstream wheel;
0234     wheel << wh;
0235     string histoName = "chEfficiencySummary_W" + wheel.str();
0236 
0237     ibooker.setCurrentFolder("DT/01-DTChamberEfficiency");
0238     summaryHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 12, 1, 13, 4, 1, 5);
0239     summaryHistos[wh]->setAxisTitle("Sector", 1);
0240     summaryHistos[wh]->setBinLabel(1, "MB1", 2);
0241     summaryHistos[wh]->setBinLabel(2, "MB2", 2);
0242     summaryHistos[wh]->setBinLabel(3, "MB3", 2);
0243     summaryHistos[wh]->setBinLabel(4, "MB4", 2);
0244   }
0245 }