Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author C. Battilana S. Marcellini - INFN Bologna
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah ncpp-um-my
0007  *
0008  */
0009 
0010 // This class header
0011 #include "DQM/DTMonitorClient/src/DTLocalTriggerLutTest.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 // Root
0023 #include "TF1.h"
0024 #include "TProfile.h"
0025 
0026 //C++ headers
0027 #include <iostream>
0028 #include <sstream>
0029 
0030 using namespace edm;
0031 using namespace std;
0032 
0033 DTLocalTriggerLutTest::DTLocalTriggerLutTest(const edm::ParameterSet& ps) {
0034   setConfig(ps, "DTLocalTriggerLut");
0035   baseFolderTM = "DT/03-LocalTrigger-TM/";
0036   thresholdPhiMean = ps.getUntrackedParameter<double>("thresholdPhiMean", 1.5);
0037   thresholdPhiRMS = ps.getUntrackedParameter<double>("thresholdPhiRMS", .5);
0038   thresholdPhibMean = ps.getUntrackedParameter<double>("thresholdPhibMean", 1.5);
0039   thresholdPhibRMS = ps.getUntrackedParameter<double>("thresholdPhibRMS", .8);
0040   doCorrStudy = ps.getUntrackedParameter<bool>("doCorrelationStudy", false);
0041 
0042   bookingdone = false;
0043 }
0044 
0045 DTLocalTriggerLutTest::~DTLocalTriggerLutTest() {}
0046 
0047 void DTLocalTriggerLutTest::Bookings(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0048   vector<string>::const_iterator iTr = trigSources.begin();
0049   vector<string>::const_iterator trEnd = trigSources.end();
0050   vector<string>::const_iterator iHw = hwSources.begin();
0051   vector<string>::const_iterator hwEnd = hwSources.end();
0052 
0053   //Booking
0054   if (parameters.getUntrackedParameter<bool>("staticBooking", true)) {
0055     for (; iTr != trEnd; ++iTr) {
0056       trigSource = (*iTr);
0057       for (; iHw != hwEnd; ++iHw) {
0058         hwSource = (*iHw);
0059         // Loop over the TriggerUnits
0060         for (int wh = -2; wh <= 2; ++wh) {
0061           bookWheelHistos(ibooker, wh, "PhiResidualMean");
0062           bookWheelHistos(ibooker, wh, "PhiResidualRMS");
0063           bookWheelHistos(ibooker, wh, "PhibResidualMean");
0064           bookWheelHistos(ibooker, wh, "PhibResidualRMS");
0065           if (doCorrStudy) {
0066             bookWheelHistos(ibooker, wh, "PhiTkvsTrigSlope");
0067             bookWheelHistos(ibooker, wh, "PhiTkvsTrigIntercept");
0068             bookWheelHistos(ibooker, wh, "PhiTkvsTrigCorr");
0069             bookWheelHistos(ibooker, wh, "PhibTkvsTrigSlope");
0070             bookWheelHistos(ibooker, wh, "PhibTkvsTrigIntercept");
0071             bookWheelHistos(ibooker, wh, "PhibTkvsTrigCorr");
0072           }
0073         }
0074       }
0075     }
0076   }
0077 
0078   // Summary test histo booking (only static)
0079   for (iTr = trigSources.begin(); iTr != trEnd; ++iTr) {
0080     trigSource = (*iTr);
0081     for (iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0082       hwSource = (*iHw);
0083       // Loop over the TriggerUnits
0084       for (int wh = -2; wh <= 2; ++wh) {
0085         bookWheelHistos(ibooker, wh, "PhiLutSummary");
0086         bookWheelHistos(ibooker, wh, "PhibLutSummary");
0087       }
0088 
0089       bookCmsHistos(ibooker, "PhiLutSummary");
0090       bookCmsHistos(ibooker, "PhibLutSummary");
0091     }
0092   }
0093 
0094   bookingdone = true;
0095 }
0096 
0097 void DTLocalTriggerLutTest::beginRun(const edm::Run& r, const edm::EventSetup& c) {
0098   DTLocalTriggerBaseTest::beginRun(r, c);
0099 }
0100 
0101 void DTLocalTriggerLutTest::runClientDiagnostic(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0102   if (!bookingdone)
0103     Bookings(ibooker, igetter);
0104 
0105   // Loop over Trig & Hw sources
0106   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0107     trigSource = (*iTr);
0108     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0109       hwSource = (*iHw);
0110       vector<const DTChamber*>::const_iterator chIt = muonGeom->chambers().begin();
0111       vector<const DTChamber*>::const_iterator chEnd = muonGeom->chambers().end();
0112       for (; chIt != chEnd; ++chIt) {
0113         DTChamberId chId((*chIt)->id());
0114         int wh = chId.wheel();
0115         int sect = chId.sector();
0116         int stat = chId.station();
0117 
0118         if (doCorrStudy) {
0119           // Perform Correlation Plots analysis (TM + segment Phi)
0120 
0121           TH2F* TrackPhitkvsPhitrig = getHisto<TH2F>(igetter.get(getMEName("PhitkvsPhitrig", "Segment", chId)));
0122 
0123           if (TrackPhitkvsPhitrig && TrackPhitkvsPhitrig->GetEntries() > 10) {
0124             // Fill client histos
0125             if (whME[wh].find(fullName("PhiTkvsTrigCorr")) == whME[wh].end()) {
0126               bookWheelHistos(ibooker, wh, "PhiTkvsTrigSlope");
0127               bookWheelHistos(ibooker, wh, "PhiTkvsTrigIntercept");
0128               bookWheelHistos(ibooker, wh, "PhiTkvsTrigCorr");
0129             }
0130 
0131             TProfile* PhitkvsPhitrigProf = TrackPhitkvsPhitrig->ProfileX();
0132             double phiInt = 0;
0133             double phiSlope = 0;
0134             double phiCorr = 0;
0135             try {
0136               TF1 ffPhi("mypol1", "pol1");
0137               PhitkvsPhitrigProf->Fit(&ffPhi, "CQO");
0138               phiInt = ffPhi.GetParameter(0);
0139               phiSlope = ffPhi.GetParameter(1);
0140               phiCorr = TrackPhitkvsPhitrig->GetCorrelationFactor();
0141             } catch (cms::Exception& iException) {
0142               edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhitkvsPhitrig for Wheel " << wh
0143                                         << " Sector " << sect << " Station " << stat;
0144             }
0145 
0146             std::map<std::string, MonitorElement*>& innerME = whME[wh];
0147             fillWhPlot(innerME.find(fullName("PhiTkvsTrigSlope"))->second, sect, stat, phiSlope - 1);
0148             fillWhPlot(innerME.find(fullName("PhiTkvsTrigIntercept"))->second, sect, stat, phiInt);
0149             fillWhPlot(innerME.find(fullName("PhiTkvsTrigCorr"))->second, sect, stat, phiCorr, false);
0150           }
0151 
0152           // Perform Correlation Plots analysis (TM + segment Phib)
0153           TH2F* TrackPhibtkvsPhibtrig = getHisto<TH2F>(igetter.get(getMEName("PhibtkvsPhibtrig", "Segment", chId)));
0154 
0155           if (stat != 3 && TrackPhibtkvsPhibtrig &&
0156               TrackPhibtkvsPhibtrig->GetEntries() > 10) {  // station 3 has no meaningful MB3 phi bending information
0157 
0158             // Fill client histos
0159             if (whME[wh].find(fullName("PhibTkvsTrigCorr")) == whME[wh].end()) {
0160               bookWheelHistos(ibooker, wh, "PhibTkvsTrigSlope");
0161               bookWheelHistos(ibooker, wh, "PhibTkvsTrigIntercept");
0162               bookWheelHistos(ibooker, wh, "PhibTkvsTrigCorr");
0163             }
0164 
0165             TProfile* PhibtkvsPhibtrigProf = TrackPhibtkvsPhibtrig->ProfileX();
0166             double phibInt = 0;
0167             double phibSlope = 0;
0168             double phibCorr = 0;
0169             try {
0170               TF1 ffPhib("ffPhib", "pol1");
0171               PhibtkvsPhibtrigProf->Fit(&ffPhib, "CQO");
0172               phibInt = ffPhib.GetParameter(0);
0173               phibSlope = ffPhib.GetParameter(1);
0174               phibCorr = TrackPhibtkvsPhibtrig->GetCorrelationFactor();
0175             } catch (cms::Exception& iException) {
0176               edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhibtkvsPhibtrig for Wheel " << wh
0177                                         << " Sector " << sect << " Station " << stat;
0178             }
0179 
0180             std::map<std::string, MonitorElement*>& innerME = whME[wh];
0181             fillWhPlot(innerME.find(fullName("PhibTkvsTrigSlope"))->second, sect, stat, phibSlope - 1);
0182             fillWhPlot(innerME.find(fullName("PhibTkvsTrigIntercept"))->second, sect, stat, phibInt);
0183             fillWhPlot(innerME.find(fullName("PhibTkvsTrigCorr"))->second, sect, stat, phibCorr, false);
0184           }
0185         }
0186 
0187         // Make Phi Residual Summary
0188 
0189         TH1F* PhiResidual = getHisto<TH1F>(igetter.get(getMEName("PhiResidual", "Segment", chId)));
0190         int phiSummary = 1;
0191 
0192         if (PhiResidual && PhiResidual->GetEffectiveEntries() > 10) {
0193           // Fill client histos
0194           if (whME[wh].find(fullName("PhiResidualMean")) == whME[wh].end()) {
0195             bookWheelHistos(ibooker, wh, "PhiResidualMean");
0196             bookWheelHistos(ibooker, wh, "PhiResidualRMS");
0197           }
0198 
0199           double peak = PhiResidual->GetBinCenter(PhiResidual->GetMaximumBin());
0200           double phiMean = 0;
0201           double phiRMS = 0;
0202           try {
0203             TF1 ffPhi("ffPhi", "gaus");
0204             PhiResidual->Fit(&ffPhi, "CQO", "", peak - 5, peak + 5);
0205             phiMean = ffPhi.GetParameter(1);
0206             phiRMS = ffPhi.GetParameter(2);
0207           } catch (cms::Exception& iException) {
0208             edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhiResidual for Wheel " << wh
0209                                       << " Sector " << sect << " Station " << stat;
0210           }
0211 
0212           std::map<std::string, MonitorElement*>& innerME = whME[wh];
0213           fillWhPlot(innerME.find(fullName("PhiResidualMean"))->second, sect, stat, phiMean);
0214           fillWhPlot(innerME.find(fullName("PhiResidualRMS"))->second, sect, stat, phiRMS);
0215 
0216           phiSummary = performLutTest(phiMean, phiRMS, thresholdPhiMean, thresholdPhiRMS);
0217         }
0218         fillWhPlot(whME[wh].find(fullName("PhiLutSummary"))->second, sect, stat, phiSummary);
0219 
0220         // Make Phib Residual Summary
0221         TH1F* PhibResidual = getHisto<TH1F>(igetter.get(getMEName("PhibResidual", "Segment", chId)));
0222         int phibSummary = stat == 3 ? 0 : 1;  // station 3 has no meaningful MB3 phi bending information
0223 
0224         if (stat != 3 && PhibResidual &&
0225             PhibResidual->GetEffectiveEntries() > 10) {  // station 3 has no meaningful MB3 phi bending information
0226 
0227           // Fill client histos
0228           if (whME[wh].find(fullName("PhibResidualMean")) == whME[wh].end()) {
0229             bookWheelHistos(ibooker, wh, "PhibResidualMean");
0230             bookWheelHistos(ibooker, wh, "PhibResidualRMS");
0231           }
0232 
0233           double peak = PhibResidual->GetBinCenter(PhibResidual->GetMaximumBin());
0234           double phibMean = 0;
0235           double phibRMS = 0;
0236           try {
0237             TF1 ffPhib("ffPhib", "gaus");
0238             PhibResidual->Fit(&ffPhib, "CQO", "", peak - 5, peak + 5);
0239             phibMean = ffPhib.GetParameter(1);
0240             phibRMS = ffPhib.GetParameter(2);
0241           } catch (cms::Exception& iException) {
0242             edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhibResidual for Wheel " << wh
0243                                       << " Sector " << sect << " Station " << stat;
0244           }
0245 
0246           std::map<std::string, MonitorElement*>& innerME = whME[wh];
0247           fillWhPlot(innerME.find(fullName("PhibResidualMean"))->second, sect, stat, phibMean);
0248           fillWhPlot(innerME.find(fullName("PhibResidualRMS"))->second, sect, stat, phibRMS);
0249 
0250           phibSummary = performLutTest(phibMean, phibRMS, thresholdPhibMean, thresholdPhibRMS);
0251         }
0252         fillWhPlot(whME[wh].find(fullName("PhibLutSummary"))->second, sect, stat, phibSummary);
0253       }
0254     }
0255   }
0256 
0257   // Barrel Summary Plots
0258   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0259     trigSource = (*iTr);
0260     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0261       hwSource = (*iHw);
0262       for (int wh = -2; wh <= 2; ++wh) {
0263         std::map<std::string, MonitorElement*>* innerME = &(whME[wh]);
0264 
0265         TH2F* phiWhSummary = getHisto<TH2F>(innerME->find(fullName("PhiLutSummary"))->second);
0266         TH2F* phibWhSummary = getHisto<TH2F>(innerME->find(fullName("PhibLutSummary"))->second);
0267         for (int sect = 1; sect <= 12; ++sect) {
0268           int phiErr = 0;
0269           int phibErr = 0;
0270           int phiNoData = 0;
0271           int phibNoData = 0;
0272           for (int stat = 1; stat <= 4; ++stat) {
0273             switch (static_cast<int>(phiWhSummary->GetBinContent(sect, stat))) {
0274               case 1:
0275                 phiNoData++;
0276                 [[fallthrough]];
0277               case 2:
0278               case 3:
0279                 phiErr++;
0280             }
0281             switch (static_cast<int>(phibWhSummary->GetBinContent(sect, stat))) {
0282               case 1:
0283                 phibNoData++;
0284                 [[fallthrough]];
0285               case 2:
0286               case 3:
0287                 phibErr++;
0288             }
0289           }
0290           if (phiNoData == 4)
0291             phiErr = 5;
0292           if (phibNoData == 3)
0293             phibErr = 5;  // MB3 has no phib information
0294           cmsME.find(fullName("PhiLutSummary"))->second->setBinContent(sect, wh + wheelArrayShift, phiErr);
0295           cmsME.find(fullName("PhibLutSummary"))->second->setBinContent(sect, wh + wheelArrayShift, phibErr);
0296         }
0297       }
0298     }
0299   }
0300 }
0301 
0302 int DTLocalTriggerLutTest::performLutTest(double mean, double RMS, double thresholdMean, double thresholdRMS) {
0303   bool meanErr = fabs(mean) > thresholdMean;
0304   bool rmsErr = RMS > thresholdRMS;
0305 
0306   return (meanErr || rmsErr) ? 2 + (meanErr != rmsErr) : 0;
0307 }
0308 
0309 void DTLocalTriggerLutTest::fillWhPlot(MonitorElement* plot, int sect, int stat, float value, bool lessIsBest) {
0310   if (sect > 12) {
0311     int scsect = sect == 13 ? 4 : 10;
0312     if ((fabs(value) > fabs(plot->getBinContent(scsect, stat))) == lessIsBest) {
0313       plot->setBinContent(scsect, stat, value);
0314     }
0315   } else {
0316     plot->setBinContent(sect, stat, value);
0317   }
0318 
0319   return;
0320 }