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 G. Mila - INFN Torino
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah ncpp-um-my
0007  *
0008  */
0009 
0010 #include "DQM/DTMonitorClient/src/DTSegmentAnalysisTest.h"
0011 
0012 // Framework
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "FWCore/Framework/interface/ESHandle.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 
0019 // Geometry
0020 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0021 #include "Geometry/DTGeometry/interface/DTLayer.h"
0022 #include "Geometry/DTGeometry/interface/DTTopology.h"
0023 
0024 #include "DQMServices/Core/interface/DQMStore.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 
0027 #include <iostream>
0028 #include <cstdio>
0029 #include <string>
0030 #include <sstream>
0031 #include <cmath>
0032 
0033 using namespace edm;
0034 using namespace std;
0035 
0036 DTSegmentAnalysisTest::DTSegmentAnalysisTest(const ParameterSet& ps)
0037     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0038   LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest") << "[DTSegmentAnalysisTest]: Constructor";
0039   parameters = ps;
0040 
0041   // get the cfi parameters
0042   detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis", false);
0043   runOnline = parameters.getUntrackedParameter<bool>("runOnline", true);
0044   // top folder for the histograms in DQMStore
0045   topHistoFolder = ps.getUntrackedParameter<string>("topHistoFolder", "DT/02-Segments");
0046   nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000);
0047   maxPhiHit = ps.getUntrackedParameter<int>("maxPhiHit", 7);
0048   maxPhiZHit = ps.getUntrackedParameter<int>("maxPhiZHit", 11);
0049 
0050   nLSs = 0;
0051 
0052   bookingdone = false;
0053 }
0054 
0055 DTSegmentAnalysisTest::~DTSegmentAnalysisTest() {
0056   LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest") << "DTSegmentAnalysisTest: analyzed " << nLSs << " LS";
0057 }
0058 
0059 void DTSegmentAnalysisTest::beginRun(const Run& run, const EventSetup& context) {
0060   LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest") << "[DTSegmentAnalysisTest]: BeginRun";
0061   muonGeom = &context.getData(muonGeomToken_);
0062 }
0063 
0064 void DTSegmentAnalysisTest::dqmBeginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const&) {
0065   nLSs++;
0066   LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0067       << "DTSegmentAnalysisTest: analyzing LS" << lumiSeg.id().luminosityBlock() << " of " << nLSs << endl;
0068 }
0069 
0070 void DTSegmentAnalysisTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0071                                                   DQMStore::IGetter& igetter,
0072                                                   edm::LuminosityBlock const& lumiSeg,
0073                                                   edm::EventSetup const& context) {
0074   // book the histos
0075 
0076   if (!bookingdone)
0077     bookHistos(ibooker);
0078   bookingdone = true;
0079 
0080   // counts number of lumiSegs
0081   nLumiSegs = lumiSeg.id().luminosityBlock();
0082 
0083   if (runOnline) {
0084     LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0085         << "[DTSegmentAnalysisTest]: End of LS " << nLumiSegs
0086         << ". Client called in online mode , perform DQM client operation";
0087 
0088     performClientDiagnostic(igetter);
0089   }
0090 }
0091 
0092 void DTSegmentAnalysisTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0093   if (!runOnline) {
0094     LogTrace("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0095         << "[DTSegmentAnalysisTest]: endJob. Client called in offline mode , perform DQM client operation";
0096 
0097     performClientDiagnostic(igetter);
0098   }
0099 }
0100 
0101 void DTSegmentAnalysisTest::performClientDiagnostic(DQMStore::IGetter& igetter) {
0102   summaryHistos[3]->Reset();
0103   summaryHistos[4]->Reset();
0104   vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0105   vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0106 
0107   for (; ch_it != ch_end; ++ch_it) {
0108     DTChamberId chID = (*ch_it)->id();
0109 
0110     MonitorElement* hNHits = igetter.get(getMEName(chID, "h4DSegmNHits"));
0111     MonitorElement* hSegmOcc = igetter.get(getMEName(chID, "numberOfSegments"));
0112 
0113     if (hNHits && hSegmOcc) {
0114       TH1F* hNHits_root = hNHits->getTH1F();
0115       TH2F* hSegmOcc_root = hSegmOcc->getTH2F();
0116       TH2F* summary_histo_root = summaryHistos[3]->getTH2F();
0117 
0118       int sector = chID.sector();
0119       if (sector == 13)
0120         sector = 4;
0121       if (sector == 14)
0122         sector = 10;
0123 
0124       if ((chID.station() != 4 && hNHits_root->GetMaximumBin() < maxPhiZHit) ||
0125           (chID.station() == 4 && hNHits_root->GetMaximumBin() < maxPhiHit)) {
0126         summaryHistos[chID.wheel()]->setBinContent(sector, chID.station(), 1);
0127         if (summary_histo_root->GetBinContent(sector, chID.wheel() + 3) < 1)
0128           summaryHistos[3]->setBinContent(sector, chID.wheel() + 3, 1);
0129       } else
0130         summaryHistos[chID.wheel()]->setBinContent(sector, chID.station(), 0);
0131 
0132       if (detailedAnalysis) {
0133         if (chID.station() != 4)
0134           segmRecHitHistos[make_pair(chID.wheel(), chID.sector())]->Fill(chID.station(),
0135                                                                          abs(12 - hNHits_root->GetMaximumBin()));
0136         else
0137           segmRecHitHistos[make_pair(chID.wheel(), chID.sector())]->Fill(chID.station(),
0138                                                                          abs(8 - hNHits_root->GetMaximumBin()));
0139       }
0140 
0141       TH2F* summary2_histo_root = summaryHistos[3]->getTH2F();
0142       float weight = 0.001;
0143       if (hSegmOcc_root->GetBinContent(sector, chID.station()) == 0) {
0144         summaryHistos[chID.wheel()]->setBinContent(sector, chID.station(), 2);
0145         if (summary2_histo_root->GetBinContent(sector, chID.wheel() + 3) < 2)
0146           summaryHistos[3]->setBinContent(sector, chID.wheel() + 3, 2);
0147       } else {
0148         // Fill the percentage of segment occupancy
0149         weight = 1. / 4.;
0150         if ((sector == 4 || sector == 10) && chID.station() == 4)
0151           weight = 1. / 8.;
0152       }
0153       summaryHistos[4]->Fill(sector, chID.wheel(), weight);
0154     } else {
0155       LogVerbatim("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0156           << "[DTSegmentAnalysisTest]: histos not found!!";  // FIXME
0157     }
0158 
0159     if (detailedAnalysis) {  // switch on detailed analysis
0160 
0161       //test on chi2 segment quality
0162 
0163       MonitorElement* chi2_histo = igetter.get(getMEName(chID, "h4DChi2"));
0164       if (chi2_histo) {
0165         TH1F* chi2_histo_root = chi2_histo->getTH1F();
0166         double threshold = parameters.getUntrackedParameter<double>("chi2Threshold", 5);
0167         double maximum = chi2_histo_root->GetXaxis()->GetXmax();
0168         double minimum = chi2_histo_root->GetXaxis()->GetXmin();
0169         int nbins = chi2_histo_root->GetXaxis()->GetNbins();
0170         int thresholdBin = int(threshold / ((maximum - minimum) / nbins));
0171 
0172         double badSegments = 0;
0173         for (int bin = thresholdBin; bin <= nbins; bin++) {
0174           badSegments += chi2_histo_root->GetBinContent(bin);
0175         }
0176 
0177         if (chi2_histo_root->GetEntries() != 0) {
0178           double badSegmentsPercentual = badSegments / double(chi2_histo_root->GetEntries());
0179           chi2Histos[make_pair(chID.wheel(), chID.sector())]->Fill(chID.station(), badSegmentsPercentual);
0180         }
0181       } else {
0182         LogVerbatim("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0183             << "[DTSegmentAnalysisTest]: Histo: " << getMEName(chID, "h4DChi2") << " not found!" << endl;
0184       }
0185     }  // end of switch for detailed analysis
0186 
0187   }  //loop over all the chambers
0188 
0189   string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsSegment";
0190 
0191   MonitorElement* meProcEvts = igetter.get(nEvtsName);
0192 
0193   if (meProcEvts) {
0194     int nProcEvts = meProcEvts->getFloatValue();
0195     summaryHistos[4]->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
0196   } else {
0197     summaryHistos[4]->setEntries(nMinEvts + 1);
0198     LogVerbatim("DTDQM|DTMonitorClient|DTOccupancyTest")
0199         << "[DTOccupancyTest] ME: " << nEvtsName << " not found!" << endl;
0200   }
0201 
0202   if (detailedAnalysis) {
0203     string chi2CriterionName = parameters.getUntrackedParameter<string>("chi2TestName", "chi2InRange");
0204     for (map<pair<int, int>, MonitorElement*>::const_iterator histo = chi2Histos.begin(); histo != chi2Histos.end();
0205          histo++) {
0206       const QReport* theChi2QReport = (*histo).second->getQReport(chi2CriterionName);
0207       if (theChi2QReport) {
0208         vector<dqm::me_util::Channel> badChannels = theChi2QReport->getBadChannels();
0209         for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0210              channel++) {
0211           LogError("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0212               << "Wheel: " << (*histo).first.first << " Sector: " << (*histo).first.second
0213               << " Bad stations: " << (*channel).getBin() << "  Contents : " << (*channel).getContents();
0214         }
0215       }
0216     }
0217 
0218     string segmRecHitCriterionName =
0219         parameters.getUntrackedParameter<string>("segmRecHitTestName", "segmRecHitInRange");
0220     for (map<pair<int, int>, MonitorElement*>::const_iterator histo = segmRecHitHistos.begin();
0221          histo != segmRecHitHistos.end();
0222          histo++) {
0223       const QReport* theSegmRecHitQReport = (*histo).second->getQReport(segmRecHitCriterionName);
0224       if (theSegmRecHitQReport) {
0225         vector<dqm::me_util::Channel> badChannels = theSegmRecHitQReport->getBadChannels();
0226         for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0227              channel++) {
0228           LogError("DTDQM|DTMonitorClient|DTSegmentAnalysisTest")
0229               << "Wheel: " << (*histo).first.first << " Sector: " << (*histo).first.second
0230               << " Bad stations on recHit number: " << (*channel).getBin()
0231               << "  Contents : " << (*channel).getContents();
0232         }
0233       }
0234     }
0235 
0236   }  // end of detailedAnalysis
0237 }
0238 
0239 string DTSegmentAnalysisTest::getMEName(const DTChamberId& chID, string histoTag) {
0240   stringstream wheel;
0241   wheel << chID.wheel();
0242   stringstream station;
0243   station << chID.station();
0244   stringstream sector;
0245   sector << chID.sector();
0246 
0247   string folderName =
0248       topHistoFolder + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() + "/";
0249 
0250   string histoname = folderName + histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0251 
0252   if (histoTag == "numberOfSegments")
0253     histoname = topHistoFolder + "/Wheel" + wheel.str() + "/" + histoTag + +"_W" + wheel.str();
0254 
0255   return histoname;
0256 }
0257 
0258 void DTSegmentAnalysisTest::bookHistos(DQMStore::IBooker& ibooker) {
0259   for (int wh = -2; wh <= 2; wh++) {
0260     stringstream wheel;
0261     wheel << wh;
0262     string histoName = "segmentSummary_W" + wheel.str();
0263 
0264     ibooker.setCurrentFolder(topHistoFolder);
0265 
0266     summaryHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 12, 1, 13, 4, 1, 5);
0267     summaryHistos[wh]->setAxisTitle("Sector", 1);
0268     summaryHistos[wh]->setBinLabel(1, "MB1", 2);
0269     summaryHistos[wh]->setBinLabel(2, "MB2", 2);
0270     summaryHistos[wh]->setBinLabel(3, "MB3", 2);
0271     summaryHistos[wh]->setBinLabel(4, "MB4", 2);
0272 
0273     if (detailedAnalysis) {
0274       for (int sect = 1; sect <= 14; sect++) {
0275         stringstream sector;
0276         sector << sect;
0277         string chi2HistoName = "chi2BadSegmPercentual_W" + wheel.str() + "_Sec" + sector.str();
0278         ibooker.setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() + "/Tests");
0279         chi2Histos[make_pair(wh, sect)] = ibooker.book1D(chi2HistoName.c_str(), chi2HistoName.c_str(), 4, 1, 5);
0280         chi2Histos[make_pair(wh, sect)]->setBinLabel(1, "MB1");
0281         chi2Histos[make_pair(wh, sect)]->setBinLabel(2, "MB2");
0282         chi2Histos[make_pair(wh, sect)]->setBinLabel(3, "MB3");
0283         chi2Histos[make_pair(wh, sect)]->setBinLabel(4, "MB4");
0284 
0285         string segmHistoName = "residualsOnSegmRecHitNumber_W" + wheel.str() + "_Sec" + sector.str();
0286         segmRecHitHistos[make_pair(wh, sect)] = ibooker.book1D(segmHistoName.c_str(), segmHistoName.c_str(), 4, 1, 5);
0287         segmRecHitHistos[make_pair(wh, sect)]->setBinLabel(1, "MB1");
0288         segmRecHitHistos[make_pair(wh, sect)]->setBinLabel(2, "MB2");
0289         segmRecHitHistos[make_pair(wh, sect)]->setBinLabel(3, "MB3");
0290         segmRecHitHistos[make_pair(wh, sect)]->setBinLabel(4, "MB4");
0291       }
0292     }
0293   }
0294 
0295   string histoName = "segmentSummary";
0296 
0297   ibooker.setCurrentFolder(topHistoFolder);
0298 
0299   summaryHistos[3] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 12, 1, 13, 5, -2, 3);
0300   summaryHistos[3]->setAxisTitle("Sector", 1);
0301   summaryHistos[3]->setAxisTitle("Wheel", 2);
0302 
0303   summaryHistos[4] = ibooker.book2D("SegmentGlbSummary", histoName.c_str(), 12, 1, 13, 5, -2, 3);
0304   summaryHistos[4]->setAxisTitle("Sector", 1);
0305   summaryHistos[4]->setAxisTitle("Wheel", 2);
0306 }