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 M. Pelliccioni - INFN Torino
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah ncpp-um-my
0007  *
0008  */
0009 
0010 #include "DQM/DTMonitorClient/src/DTChamberEfficiencyClient.h"
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0017 
0018 #include <cstdio>
0019 #include <sstream>
0020 #include <cmath>
0021 
0022 using namespace edm;
0023 using namespace std;
0024 
0025 //two words about conventions: "All" histograms are those made for all segments
0026 //while "Qual" histograms are those for segments with at least 12 hits
0027 
0028 DTChamberEfficiencyClient::DTChamberEfficiencyClient(const ParameterSet& pSet)
0029     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0030   LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyClient") << "DTChamberEfficiencyClient: Constructor called";
0031 
0032   prescaleFactor = pSet.getUntrackedParameter<int>("diagnosticPrescale", 1);
0033 }
0034 
0035 DTChamberEfficiencyClient::~DTChamberEfficiencyClient() {
0036   LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyClient") << "DTChamberEfficiencyClient: Destructor called";
0037 }
0038 
0039 void DTChamberEfficiencyClient::beginRun(const edm::Run& run, const edm::EventSetup& setup) {
0040   // Get the DT Geometry
0041   muonGeom = &setup.getData(muonGeomToken_);
0042 }
0043 
0044 void DTChamberEfficiencyClient::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0045                                                       DQMStore::IGetter& igetter,
0046                                                       edm::LuminosityBlock const& lumiSeg,
0047                                                       edm::EventSetup const& setup) {
0048   LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyClient") << "DTChamberEfficiencyClient: endluminosityBlock";
0049 }
0050 
0051 void DTChamberEfficiencyClient::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0052   LogVerbatim("DTDQM|DTMonitorClient|DTChamberEfficiencyClient") << "DTChamberEfficiencyClient: endRun";
0053 
0054   bookHistos(ibooker);
0055 
0056   // reset the global summary
0057   globalEffSummary->Reset();
0058 
0059   //Loop over the wheels
0060   for (int wheel = -2; wheel <= 2; wheel++) {
0061     stringstream wheel_str;
0062     wheel_str << wheel;
0063 
0064     // Get the ME produced by EfficiencyTask Source
0065     // All means no selection on segments, Qual means segments with at least 12 hits
0066 
0067     MonitorElement* MECountAll = igetter.get("DT/05-ChamberEff/Task/hCountSectVsChamb_All_W" + wheel_str.str());
0068     MonitorElement* MECountQual = igetter.get("DT/05-ChamberEff/Task/hCountSectVsChamb_Qual_W" + wheel_str.str());
0069     MonitorElement* MEExtrap = igetter.get("DT/05-ChamberEff/Task/hExtrapSectVsChamb_W" + wheel_str.str());
0070 
0071     //get the TH2F
0072     if (!MECountAll || !(MECountAll->getTH2F())) {
0073       edm::LogWarning("DTChamberEfficiencyClient") << "ME not available" << std::endl;
0074       return;
0075     }
0076 
0077     TH2F* hCountAll = MECountAll->getTH2F();
0078     TH2F* hCountQual = MECountQual->getTH2F();
0079     TH2F* hExtrap = MEExtrap->getTH2F();
0080 
0081     const int nBinX = summaryHistos[wheel + 2][0]->getNbinsX();
0082     const int nBinY = summaryHistos[wheel + 2][0]->getNbinsY();
0083 
0084     for (int j = 1; j <= nBinX; j++) {
0085       for (int k = 1; k <= nBinY; k++) {
0086         summaryHistos[wheel + 2][0]->setBinContent(j, k, 0.);
0087         summaryHistos[wheel + 2][1]->setBinContent(j, k, 0.);
0088 
0089         const float numerAll = hCountAll->GetBinContent(j, k);
0090         const float numerQual = hCountQual->GetBinContent(j, k);
0091         const float denom = hExtrap->GetBinContent(j, k);
0092 
0093         if (denom != 0.) {
0094           const float effAll = numerAll / denom;
0095           const float eff_error_All = sqrt((effAll + effAll * effAll) / denom);
0096 
0097           const float effQual = numerQual / denom;
0098           const float eff_error_Qual = sqrt((effQual + effQual * effQual) / denom);
0099 
0100           summaryHistos[wheel + 2][0]->setBinContent(j, k, effAll);
0101           summaryHistos[wheel + 2][0]->setBinError(j, k, eff_error_All);
0102 
0103           summaryHistos[wheel + 2][1]->setBinContent(j, k, effQual);
0104           summaryHistos[wheel + 2][1]->setBinError(j, k, eff_error_Qual);
0105 
0106           // Fill 1D eff distributions
0107           globalEffDistr->Fill(effAll);
0108           EffDistrPerWh[wheel + 2]->Fill(effAll);
0109         }
0110       }
0111     }
0112   }
0113 
0114   // fill the global eff. summary
0115   // problems at a granularity smaller than the chamber are ignored
0116   for (int wheel = -2; wheel <= 2; wheel++) {  // loop over wheels
0117     // retrieve the chamber efficiency summary
0118     MonitorElement* segmentWheelSummary = summaryHistos[wheel + 2][0];
0119     if (segmentWheelSummary != nullptr) {
0120       for (int sector = 1; sector <= 12; sector++) {  // loop over sectors
0121         float nFailingChambers = 0.;
0122 
0123         double meaneff = 0.;
0124         double errorsum = 0.;
0125 
0126         for (int station = 1; station != 5; ++station) {  // loop over stations
0127 
0128           const double tmpefficiency = segmentWheelSummary->getBinContent(sector, station);
0129           const double tmpvariance = pow(segmentWheelSummary->getBinError(sector, station), 2);
0130 
0131           if (tmpefficiency < 0.2 || tmpvariance == 0) {
0132             nFailingChambers++;
0133             continue;
0134           }
0135 
0136           meaneff += tmpefficiency / tmpvariance;
0137           errorsum += 1. / tmpvariance;
0138 
0139           LogTrace("DTDQM|DTMonitorClient|DTChamberEfficiencyClient")
0140               << "Wheel: " << wheel << " Stat: " << station << " Sect: " << sector << " status: " << meaneff / errorsum
0141               << endl;
0142         }
0143 
0144         if (sector == 4 || sector == 10) {
0145           int whichSector = (sector == 4) ? 13 : 14;
0146 
0147           const double tmpefficiency = segmentWheelSummary->getBinContent(whichSector, 4);
0148           const double tmpvariance = pow(segmentWheelSummary->getBinError(whichSector, 4), 2);
0149 
0150           if (tmpefficiency > 0.2 && tmpvariance != 0) {
0151             meaneff += tmpefficiency / tmpvariance;
0152             errorsum += 1. / tmpvariance;
0153           } else
0154             nFailingChambers++;
0155         }
0156 
0157         double eff_result = 0;
0158         if (errorsum != 0)
0159           eff_result = meaneff / errorsum;
0160 
0161         if (nFailingChambers != 0) {
0162           if (sector != 4 && sector != 10)
0163             eff_result = eff_result * (4. - nFailingChambers) / 4.;
0164           else
0165             eff_result = eff_result * (5. - nFailingChambers) / 5.;
0166         }
0167 
0168         if (eff_result > 0.7)
0169           globalEffSummary->Fill(sector, wheel, 1.);
0170         else if (eff_result < 0.7 && eff_result > 0.5)
0171           globalEffSummary->Fill(sector, wheel, 0.6);
0172         else if (eff_result < 0.5 && eff_result > 0.3)
0173           globalEffSummary->Fill(sector, wheel, 0.4);
0174         else if (eff_result < 0.3 && eff_result > 0.)
0175           globalEffSummary->Fill(sector, wheel, 0.15);
0176         else
0177           globalEffSummary->Fill(sector, wheel, 0.01);
0178       }
0179     }
0180   }
0181   return;
0182 }
0183 
0184 void DTChamberEfficiencyClient::bookHistos(DQMStore::IBooker& ibooker) {
0185   ibooker.setCurrentFolder("DT/05-ChamberEff");
0186 
0187   globalEffSummary = ibooker.book2D("EfficiencyGlbSummary", "Efficiency Summary", 12, 1, 13, 5, -2, 3);
0188   globalEffSummary->setAxisTitle("sector", 1);
0189   globalEffSummary->setAxisTitle("wheel", 2);
0190 
0191   globalEffDistr = ibooker.book1D("TotalEfficiency", "Total efficiency", 51, 0., 1.02);
0192   globalEffDistr->setAxisTitle("Eff", 1);
0193 
0194   for (int wh = -2; wh <= 2; wh++) {
0195     stringstream wheel;
0196     wheel << wh;
0197     string histoNameAll = "EfficiencyMap_All_W" + wheel.str();
0198     string histoTitleAll = "Efficiency map for all segments for wheel " + wheel.str();
0199 
0200     string histoNameQual = "EfficiencyMap_Qual_W" + wheel.str();
0201     string histoTitleQual = "Efficiency map for quality segments for wheel " + wheel.str();
0202 
0203     string histoNameEff = "Efficiency_W" + wheel.str();
0204     string histoTitleEff = "Segment efficiency, wheel " + wheel.str();
0205 
0206     ibooker.setCurrentFolder("DT/05-ChamberEff");
0207 
0208     summaryHistos[wh + 2][0] = ibooker.book2D(histoNameAll.c_str(), histoTitleAll.c_str(), 14, 1., 15., 4, 1., 5.);
0209     summaryHistos[wh + 2][0]->setAxisTitle("Sector", 1);
0210     summaryHistos[wh + 2][0]->setBinLabel(1, "MB1", 2);
0211     summaryHistos[wh + 2][0]->setBinLabel(2, "MB2", 2);
0212     summaryHistos[wh + 2][0]->setBinLabel(3, "MB3", 2);
0213     summaryHistos[wh + 2][0]->setBinLabel(4, "MB4", 2);
0214 
0215     EffDistrPerWh[wh + 2] = ibooker.book1D(histoNameEff.c_str(), histoTitleEff.c_str(), 51, 0., 1.02);
0216     EffDistrPerWh[wh + 2]->setAxisTitle("Eff", 1);
0217 
0218     ibooker.setCurrentFolder("DT/05-ChamberEff/HighQual");
0219 
0220     summaryHistos[wh + 2][1] = ibooker.book2D(histoNameQual.c_str(), histoTitleQual.c_str(), 14, 1., 15., 4, 1., 5.);
0221     summaryHistos[wh + 2][1]->setAxisTitle("Sector", 1);
0222     summaryHistos[wh + 2][1]->setBinLabel(1, "MB1", 2);
0223     summaryHistos[wh + 2][1]->setBinLabel(2, "MB2", 2);
0224     summaryHistos[wh + 2][1]->setBinLabel(3, "MB3", 2);
0225     summaryHistos[wh + 2][1]->setBinLabel(4, "MB4", 2);
0226   }
0227 
0228   return;
0229 }