Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author C. Battilana - CIEMAT
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
0007  *
0008  */
0009 
0010 // This class header
0011 #include "DQM/DTMonitorClient/src/DTTriggerEfficiencyTest.h"
0012 
0013 // Framework headers
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 
0018 // Geometry
0019 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0020 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0021 
0022 // Trigger
0023 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
0024 
0025 // Root
0026 #include "TF1.h"
0027 #include "TProfile.h"
0028 
0029 //C++ headers
0030 #include <iostream>
0031 #include <sstream>
0032 
0033 using namespace edm;
0034 using namespace std;
0035 
0036 DTTriggerEfficiencyTest::DTTriggerEfficiencyTest(const edm::ParameterSet& ps) {
0037   setConfig(ps, "DTTriggerEfficiency");
0038   baseFolderTM = "DT/03-LocalTrigger-TM/";
0039   detailedPlots = ps.getUntrackedParameter<bool>("detailedAnalysis", true);
0040 
0041   bookingdone = false;
0042 }
0043 
0044 DTTriggerEfficiencyTest::~DTTriggerEfficiencyTest() {}
0045 
0046 void DTTriggerEfficiencyTest::beginRun(const edm::Run& r, const edm::EventSetup& c) {
0047   DTLocalTriggerBaseTest::beginRun(r, c);
0048   trigGeomUtils = new DTTrigGeomUtils(muonGeom);
0049 }
0050 
0051 void DTTriggerEfficiencyTest::runClientDiagnostic(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0052   if (!bookingdone)
0053     Bookings(ibooker, igetter);
0054 
0055   // Loop over Trig & Hw sources
0056   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0057     trigSource = (*iTr);
0058     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0059       hwSource = (*iHw);
0060       // Loop over the TriggerUnits
0061       if (globalEffDistr.find(fullName("TrigEffPhi")) == globalEffDistr.end()) {
0062         bookHistos(ibooker, "TrigEffPhi", "");
0063         bookHistos(ibooker, "TrigEffCorrPhi", "");
0064       }
0065       for (int wh = -2; wh <= 2; ++wh) {
0066         TH2F* TrigEffDenum = getHisto<TH2F>(igetter.get(getMEName("TrigEffDenum", "Task", wh)));
0067         TH2F* TrigEffNum = getHisto<TH2F>(igetter.get(getMEName("TrigEffNum", "Task", wh)));
0068         TH2F* TrigEffCorrNum = getHisto<TH2F>(igetter.get(getMEName("TrigEffCorrNum", "Task", wh)));
0069 
0070         if (TrigEffDenum && TrigEffNum && TrigEffCorrNum && TrigEffDenum->GetEntries() > 1) {
0071           if (whME[wh].find(fullName("TrigEffPhi")) == whME[wh].end()) {
0072             bookWheelHistos(ibooker, wh, "TrigEffPhi", "");
0073             bookWheelHistos(ibooker, wh, "TrigEffCorrPhi", "");
0074           }
0075 
0076           MonitorElement* Eff1DAll_TrigEffPhi = (&globalEffDistr)->find(fullName("TrigEffPhi"))->second;
0077           MonitorElement* Eff1DAll_TrigEffCorrPhi = (&globalEffDistr)->find(fullName("TrigEffCorrPhi"))->second;
0078 
0079           MonitorElement* Eff1DWh_TrigEffPhi = (&(EffDistrPerWh[wh]))->find(fullName("TrigEffPhi"))->second;
0080           MonitorElement* Eff1DWh_TrigEffCorrPhi = (&(EffDistrPerWh[wh]))->find(fullName("TrigEffCorrPhi"))->second;
0081 
0082           MonitorElement* Eff2DWh_TrigEffPhi = (&(whME[wh]))->find(fullName("TrigEffPhi"))->second;
0083           MonitorElement* Eff2DWh_TrigEffCorrPhi = (&(whME[wh]))->find(fullName("TrigEffCorrPhi"))->second;
0084 
0085           makeEfficiencyME(TrigEffNum, TrigEffDenum, Eff2DWh_TrigEffPhi, Eff1DWh_TrigEffPhi, Eff1DAll_TrigEffPhi);
0086           makeEfficiencyME(
0087               TrigEffCorrNum, TrigEffDenum, Eff2DWh_TrigEffCorrPhi, Eff1DWh_TrigEffCorrPhi, Eff1DAll_TrigEffCorrPhi);
0088         }
0089 
0090         if (detailedPlots) {
0091           for (int stat = 1; stat <= 4; ++stat) {
0092             for (int sect = 1; sect <= 12; ++sect) {
0093               DTChamberId chId(wh, stat, sect);
0094               uint32_t indexCh = chId.rawId();
0095 
0096               // Perform Efficiency analysis (Phi+Segments 2D)
0097 
0098               TH2F* TrackPosvsAngle = getHisto<TH2F>(igetter.get(getMEName("TrackPosvsAngle", "Segment", chId)));
0099               TH2F* TrackPosvsAngleAnyQual =
0100                   getHisto<TH2F>(igetter.get(getMEName("TrackPosvsAngleAnyQual", "Segment", chId)));
0101               TH2F* TrackPosvsAngleCorr =
0102                   getHisto<TH2F>(igetter.get(getMEName("TrackPosvsAngleCorr", "Segment", chId)));
0103 
0104               if (TrackPosvsAngle && TrackPosvsAngleAnyQual && TrackPosvsAngleCorr &&
0105                   TrackPosvsAngle->GetEntries() > 1) {
0106                 if (chambME[indexCh].find(fullName("TrigEffAnglePhi")) == chambME[indexCh].end()) {
0107                   bookChambHistos(ibooker, chId, "TrigEffPosvsAnglePhi", "Segment");
0108                   bookChambHistos(ibooker, chId, "TrigEffPosvsAngleCorrPhi", "Segment");
0109                 }
0110 
0111                 std::map<std::string, MonitorElement*>* innerME = &(chambME[indexCh]);
0112                 makeEfficiencyME(
0113                     TrackPosvsAngleAnyQual, TrackPosvsAngle, innerME->find(fullName("TrigEffPosvsAnglePhi"))->second);
0114                 makeEfficiencyME(
0115                     TrackPosvsAngleCorr, TrackPosvsAngle, innerME->find(fullName("TrigEffPosvsAngleCorrPhi"))->second);
0116               }
0117             }
0118           }
0119         }
0120       }
0121     }
0122   }
0123 }
0124 
0125 void DTTriggerEfficiencyTest::makeEfficiencyME(TH2F* numerator,
0126                                                TH2F* denominator,
0127                                                MonitorElement* result2DWh,
0128                                                MonitorElement* result1DWh,
0129                                                MonitorElement* result1D) {
0130   TH2F* efficiency = result2DWh->getTH2F();
0131   efficiency->Divide(numerator, denominator, 1, 1, "");
0132 
0133   int nbinsx = efficiency->GetNbinsX();
0134   int nbinsy = efficiency->GetNbinsY();
0135   for (int binx = 1; binx <= nbinsx; ++binx) {
0136     for (int biny = 1; biny <= nbinsy; ++biny) {
0137       float error = 0;
0138       float bineff = efficiency->GetBinContent(binx, biny);
0139 
0140       result1DWh->Fill(bineff);
0141       result1D->Fill(bineff);
0142 
0143       if (denominator->GetBinContent(binx, biny)) {
0144         error = sqrt(bineff * (1 - bineff) / denominator->GetBinContent(binx, biny));
0145       } else {
0146         error = 1;
0147         efficiency->SetBinContent(binx, biny, 0.);
0148       }
0149 
0150       efficiency->SetBinError(binx, biny, error);
0151     }
0152   }
0153 }
0154 
0155 void DTTriggerEfficiencyTest::makeEfficiencyME(TH2F* numerator, TH2F* denominator, MonitorElement* result2DWh) {
0156   TH2F* efficiency = result2DWh->getTH2F();
0157   efficiency->Divide(numerator, denominator, 1, 1, "");
0158 
0159   int nbinsx = efficiency->GetNbinsX();
0160   int nbinsy = efficiency->GetNbinsY();
0161   for (int binx = 1; binx <= nbinsx; ++binx) {
0162     for (int biny = 1; biny <= nbinsy; ++biny) {
0163       float error = 0;
0164       float bineff = efficiency->GetBinContent(binx, biny);
0165 
0166       if (denominator->GetBinContent(binx, biny)) {
0167         error = sqrt(bineff * (1 - bineff) / denominator->GetBinContent(binx, biny));
0168       } else {
0169         error = 1;
0170         efficiency->SetBinContent(binx, biny, 0.);
0171       }
0172 
0173       efficiency->SetBinError(binx, biny, error);
0174     }
0175   }
0176 }
0177 
0178 string DTTriggerEfficiencyTest::getMEName(string histoTag, string folder, int wh) {
0179   stringstream wheel;
0180   wheel << wh;
0181 
0182   string folderName = topFolder() + folder + "/";
0183 
0184   string histoname = sourceFolder + folderName + fullName(histoTag) + "_W" + wheel.str();
0185 
0186   return histoname;
0187 }
0188 
0189 void DTTriggerEfficiencyTest::bookHistos(DQMStore::IBooker& ibooker, string hTag, string folder) {
0190   string basedir;
0191   basedir = topFolder();  //Book summary histo outside Task directory
0192 
0193   if (!folder.empty()) {
0194     basedir += folder + "/";
0195   }
0196 
0197   ibooker.setCurrentFolder(basedir);
0198 
0199   string fullTag = fullName(hTag);
0200   string hname = fullTag + "_All";
0201 
0202   globalEffDistr[fullTag] = ibooker.book1D(hname.c_str(), hname.c_str(), 51, 0., 1.02);
0203   globalEffDistr[fullTag]->setAxisTitle("Trig Eff", 1);
0204 }
0205 
0206 void DTTriggerEfficiencyTest::bookWheelHistos(DQMStore::IBooker& ibooker, int wheel, string hTag, string folder) {
0207   stringstream wh;
0208   wh << wheel;
0209   string basedir;
0210   if (hTag.find("Summary") != string::npos) {
0211     basedir = topFolder();  //Book summary histo outside wheel directories
0212   } else {
0213     basedir = topFolder() + "Wheel" + wh.str() + "/";
0214   }
0215   if (!folder.empty()) {
0216     basedir += folder + "/";
0217   }
0218 
0219   ibooker.setCurrentFolder(basedir);
0220 
0221   string fullTag = fullName(hTag);
0222   string hname = fullTag + "_W" + wh.str();
0223 
0224   string hnameAll = fullTag + "_All_W" + wh.str();
0225 
0226   LogTrace(category()) << "[" << testName << "Test]: booking " << basedir << hname;
0227 
0228   (EffDistrPerWh[wheel])[fullTag] = ibooker.book1D(hnameAll.c_str(), hnameAll.c_str(), 51, 0., 1.02);
0229 
0230   if (hTag.find("Phi") != string::npos || hTag.find("Summary") != string::npos) {
0231     MonitorElement* me = ibooker.book2D(hname.c_str(), hname.c_str(), 12, 1, 13, 4, 1, 5);
0232 
0233     //     setLabelPh(me);
0234     me->setBinLabel(1, "MB1", 2);
0235     me->setBinLabel(2, "MB2", 2);
0236     me->setBinLabel(3, "MB3", 2);
0237     me->setBinLabel(4, "MB4", 2);
0238     me->setAxisTitle("Sector", 1);
0239 
0240     whME[wheel][fullTag] = me;
0241     return;
0242   }
0243 
0244   if (hTag.find("Theta") != string::npos) {
0245     MonitorElement* me = ibooker.book2D(hname.c_str(), hname.c_str(), 12, 1, 13, 3, 1, 4);
0246 
0247     //     setLabelTh(me);
0248     me->setBinLabel(1, "MB1", 2);
0249     me->setBinLabel(2, "MB2", 2);
0250     me->setBinLabel(3, "MB3", 2);
0251     me->setAxisTitle("Sector", 1);
0252 
0253     whME[wheel][fullTag] = me;
0254     return;
0255   }
0256 }
0257 
0258 void DTTriggerEfficiencyTest::bookChambHistos(DQMStore::IBooker& ibooker,
0259                                               DTChamberId chambId,
0260                                               string htype,
0261                                               string folder) {
0262   stringstream wheel;
0263   wheel << chambId.wheel();
0264   stringstream station;
0265   station << chambId.station();
0266   stringstream sector;
0267   sector << chambId.sector();
0268 
0269   string fullType = fullName(htype);
0270   string HistoName = fullType + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0271 
0272   ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() +
0273                            "/" + folder + "/");
0274 
0275   LogTrace(category()) << "[" << testName << "Test]: booking " + topFolder() + "Wheel" << wheel.str() << "/Sector"
0276                        << sector.str() << "/Station" << station.str() << "/" + folder + "/" << HistoName;
0277 
0278   uint32_t indexChId = chambId.rawId();
0279   float min, max;
0280   int nbins;
0281   trigGeomUtils->phiRange(chambId, min, max, nbins, 20);
0282   if (htype.find("TrigEffPosvsAnglePhi") == 0) {
0283     chambME[indexChId][fullType] = ibooker.book2D(
0284         HistoName.c_str(), "Trigger efficiency (any qual.) position vs angle (Phi)", 12, -30., 30., nbins, min, max);
0285     return;
0286   }
0287   if (htype.find("TrigEffPosvsAngleCorrPhi") == 0) {
0288     chambME[indexChId][fullType] = ibooker.book2D(
0289         HistoName.c_str(), "Trigger efficiency (correlated) pos vs angle (Phi)", 12, -30., 30., nbins, min, max);
0290     return;
0291   }
0292 }
0293 
0294 void DTTriggerEfficiencyTest::Bookings(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0295   vector<string>::const_iterator iTr = trigSources.begin();
0296   vector<string>::const_iterator trEnd = trigSources.end();
0297   vector<string>::const_iterator iHw = hwSources.begin();
0298   vector<string>::const_iterator hwEnd = hwSources.end();
0299 
0300   //Booking
0301   if (parameters.getUntrackedParameter<bool>("staticBooking", true)) {
0302     for (; iTr != trEnd; ++iTr) {
0303       trigSource = (*iTr);
0304       for (; iHw != hwEnd; ++iHw) {
0305         hwSource = (*iHw);
0306         // Loop over the TriggerUnits
0307 
0308         bookHistos(ibooker, "TrigEffPhi", "");
0309         bookHistos(ibooker, "TrigEffCorrPhi", "");
0310         for (int wh = -2; wh <= 2; ++wh) {
0311           if (detailedPlots) {
0312             for (int sect = 1; sect <= 12; ++sect) {
0313               for (int stat = 1; stat <= 4; ++stat) {
0314                 DTChamberId chId(wh, stat, sect);
0315 
0316                 bookChambHistos(ibooker, chId, "TrigEffPosvsAnglePhi", "Segment");
0317                 bookChambHistos(ibooker, chId, "TrigEffPosvsAngleCorrPhi", "Segment");
0318               }
0319             }
0320           }
0321 
0322           bookWheelHistos(ibooker, wh, "TrigEffPhi", "");
0323           bookWheelHistos(ibooker, wh, "TrigEffCorrPhi", "");
0324         }
0325       }
0326     }
0327   }
0328   bookingdone = true;
0329 }