Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:54:10

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/DTNoiseAnalysisTest.h"
0011 
0012 // Framework
0013 #include "FWCore/ServiceRegistry/interface/Service.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 
0017 // Geometry
0018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0019 #include "Geometry/DTGeometry/interface/DTLayer.h"
0020 #include "Geometry/DTGeometry/interface/DTTopology.h"
0021 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0022 
0023 #include "DQMServices/Core/interface/DQMStore.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025 
0026 #include <iostream>
0027 #include <sstream>
0028 
0029 using namespace edm;
0030 using namespace std;
0031 
0032 DTNoiseAnalysisTest::DTNoiseAnalysisTest(const edm::ParameterSet& ps)
0033     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0034   LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "[DTNoiseAnalysisTest]: Constructor";
0035 
0036   // get the cfi parameters
0037   noisyCellDef = ps.getUntrackedParameter<int>("noisyCellDef", 500);
0038 
0039   // switch on/off the summaries for the Synchronous noise
0040   doSynchNoise = ps.getUntrackedParameter<bool>("doSynchNoise", false);
0041   detailedAnalysis = ps.getUntrackedParameter<bool>("detailedAnalysis", false);
0042   maxSynchNoiseRate = ps.getUntrackedParameter<double>("maxSynchNoiseRate", 0.001);
0043   nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000);
0044 
0045   nevents = 0;
0046 
0047   bookingdone = false;
0048 }
0049 
0050 DTNoiseAnalysisTest::~DTNoiseAnalysisTest() {
0051   LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "DTNoiseAnalysisTest: analyzed " << nevents << " events";
0052 }
0053 
0054 void DTNoiseAnalysisTest::beginRun(Run const& run, EventSetup const& context) {
0055   // Get the geometry
0056   muonGeom = &context.getData(muonGeomToken_);
0057 }
0058 
0059 void DTNoiseAnalysisTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0060                                                 DQMStore::IGetter& igetter,
0061                                                 edm::LuminosityBlock const& lumiSeg,
0062                                                 edm::EventSetup const& context) {
0063   float chRate;
0064 
0065   if (!bookingdone) {
0066     // book the histos
0067     bookHistos(ibooker);
0068   }
0069   bookingdone = true;
0070 
0071   LogVerbatim("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
0072       << "[DTNoiseAnalysisTest]: End of LS transition, performing the DQM client operation";
0073 
0074   // Reset the summary plots
0075   for (map<int, MonitorElement*>::iterator plot = noiseHistos.begin(); plot != noiseHistos.end(); ++plot) {
0076     (*plot).second->Reset();
0077   }
0078 
0079   for (map<int, MonitorElement*>::iterator plot = noisyCellHistos.begin(); plot != noisyCellHistos.end(); ++plot) {
0080     (*plot).second->Reset();
0081   }
0082 
0083   summaryNoiseHisto->Reset();
0084 
0085   vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
0086   vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
0087 
0088   LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "[DTNoiseAnalysisTest]: Fill the summary histos";
0089 
0090   for (; ch_it != ch_end; ++ch_it) {  // loop over chambers
0091     DTChamberId chID = (*ch_it)->id();
0092 
0093     MonitorElement* histo = igetter.get(getMEName(chID));
0094 
0095     if (histo) {  // check the pointer
0096 
0097       TH2F* histo_root = histo->getTH2F();
0098 
0099       for (int sl = 1; sl != 4; ++sl) {  // loop over SLs
0100         // skip theta SL in MB4 chambers
0101         if (chID.station() == 4 && sl == 2)
0102           continue;
0103 
0104         int binYlow = ((sl - 1) * 4) + 1;
0105 
0106         for (int layer = 1; layer <= 4; ++layer) {  // loop over layers
0107 
0108           // Get the layer ID
0109           DTLayerId layID(chID, sl, layer);
0110 
0111           int nWires = muonGeom->layer(layID)->specificTopology().channels();
0112           int firstWire = muonGeom->layer(layID)->specificTopology().firstChannel();
0113 
0114           int binY = binYlow + (layer - 1);
0115 
0116           for (int wire = firstWire; wire != (nWires + firstWire); wire++) {  // loop over wires
0117 
0118             double noise = histo_root->GetBinContent(wire, binY);
0119             // fill the histos
0120             noiseHistos[chID.wheel()]->Fill(noise);
0121             noiseHistos[3]->Fill(noise);
0122             int sector = chID.sector();
0123             if (noise > noisyCellDef) {
0124               if (sector == 13) {
0125                 sector = 4;
0126               } else if (sector == 14) {
0127                 sector = 10;
0128               }
0129               noisyCellHistos[chID.wheel()]->Fill(sector, chID.station());
0130               summaryNoiseHisto->Fill(sector, chID.wheel());
0131             }
0132           }
0133         }
0134       }
0135     }
0136   }
0137 
0138   if (detailedAnalysis) {
0139     threshChannelsHisto->Reset();
0140     TH1F* histo = noiseHistos[3]->getTH1F();
0141     for (int step = 0; step != 15; step++) {
0142       int threshBin = step + 1;
0143       int minBin = 26 + step * 5;
0144       int nNoisyCh = histo->Integral(minBin, 101);
0145       threshChannelsHisto->setBinContent(threshBin, nNoisyCh);
0146     }
0147   }
0148 
0149   // build the summary of synch noise
0150 
0151   if (doSynchNoise) {
0152     LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
0153         << "[DTNoiseAnalysisTest]: fill summaries for synch noise" << endl;
0154     summarySynchNoiseHisto->Reset();
0155     glbSummarySynchNoiseHisto->Reset();
0156     for (int wheel = -2; wheel != 3; ++wheel) {
0157       // Get the histo produced by DTDigiTask
0158 
0159       MonitorElement* histoNoiseSynch = igetter.get(getSynchNoiseMEName(wheel));
0160       if (histoNoiseSynch != nullptr) {
0161         for (int sect = 1; sect != 13; ++sect) {  // loop over sectors
0162           TH2F* histo = histoNoiseSynch->getTH2F();
0163           float maxSectRate = 0;
0164           for (int sta = 1; sta != 5; ++sta) {
0165             if (nevents > 0)
0166               chRate = histo->GetBinContent(sect, sta) / (float)nevents;
0167             else
0168               chRate = -1.0;
0169             // in case nevents 0 e.g. counting not done
0170             LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
0171                 << "   Wheel: " << wheel << " sect: " << sect << " station: " << sta << " rate is: " << chRate << endl;
0172             if (chRate > maxSectRate)
0173               maxSectRate = chRate;
0174           }
0175           summarySynchNoiseHisto->Fill(sect, wheel, maxSectRate > maxSynchNoiseRate ? 1 : 0);
0176           float glbBinValue = 1 - 0.15 * maxSectRate / maxSynchNoiseRate;
0177           glbSummarySynchNoiseHisto->Fill(sect, wheel, glbBinValue > 0 ? glbBinValue : 0);
0178         }
0179       } else {
0180         LogWarning("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
0181             << "   Histo: " << getSynchNoiseMEName(wheel) << " not found!" << endl;
0182       }
0183     }
0184   }
0185 
0186   string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsNoise";
0187 
0188   MonitorElement* meProcEvts = igetter.get(nEvtsName);
0189 
0190   if (meProcEvts) {
0191     int nProcEvts = meProcEvts->getFloatValue();
0192     glbSummarySynchNoiseHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
0193     summarySynchNoiseHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
0194   } else {
0195     glbSummarySynchNoiseHisto->setEntries(nMinEvts + 1);
0196     summarySynchNoiseHisto->setEntries(nMinEvts + 1);
0197     LogVerbatim("DTDQM|DTMonitorClient|DTnoiseAnalysisTest")
0198         << "[DTNoiseAnalysisTest] ME: " << nEvtsName << " not found!" << endl;
0199   }
0200 }
0201 
0202 string DTNoiseAnalysisTest::getMEName(const DTChamberId& chID) {
0203   stringstream wheel;
0204   wheel << chID.wheel();
0205   stringstream station;
0206   station << chID.station();
0207   stringstream sector;
0208   sector << chID.sector();
0209 
0210   string folderName = "DT/05-Noise/Wheel" + wheel.str() + "/Sector" + sector.str() + "/";
0211 
0212   string histoname =
0213       folderName + string("NoiseRate") + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0214 
0215   return histoname;
0216 }
0217 
0218 void DTNoiseAnalysisTest::bookHistos(DQMStore::IBooker& ibooker) {
0219   ibooker.setCurrentFolder("DT/05-Noise");
0220   string histoName;
0221 
0222   for (int wh = -2; wh <= 2; wh++) {
0223     stringstream wheel;
0224     wheel << wh;
0225     histoName = "NoiseRateSummary_W" + wheel.str();
0226 
0227     noiseHistos[wh] = ibooker.book1D(histoName.c_str(), histoName.c_str(), 100, 0, 2000);
0228     noiseHistos[wh]->setAxisTitle("rate (Hz)", 1);
0229     noiseHistos[wh]->setAxisTitle("entries", 2);
0230   }
0231   histoName = "NoiseRateSummary";
0232 
0233   noiseHistos[3] = ibooker.book1D(histoName.c_str(), histoName.c_str(), 100, 0, 2000);
0234   noiseHistos[3]->setAxisTitle("rate (Hz)", 1);
0235   noiseHistos[3]->setAxisTitle("entries", 2);
0236 
0237   for (int wh = -2; wh <= 2; wh++) {
0238     stringstream wheel;
0239     wheel << wh;
0240     histoName = "NoiseSummary_W" + wheel.str();
0241 
0242     noisyCellHistos[wh] = ibooker.book2D(histoName.c_str(), "# of noisy channels", 12, 1, 13, 4, 1, 5);
0243     noisyCellHistos[wh]->setBinLabel(1, "MB1", 2);
0244     noisyCellHistos[wh]->setBinLabel(2, "MB2", 2);
0245     noisyCellHistos[wh]->setBinLabel(3, "MB3", 2);
0246     noisyCellHistos[wh]->setBinLabel(4, "MB4", 2);
0247     noisyCellHistos[wh]->setAxisTitle("Sector", 1);
0248   }
0249 
0250   histoName = "NoiseSummary";
0251 
0252   summaryNoiseHisto = ibooker.book2D(histoName.c_str(), "# of noisy channels", 12, 1, 13, 5, -2, 3);
0253   summaryNoiseHisto->setAxisTitle("Sector", 1);
0254   summaryNoiseHisto->setAxisTitle("Wheel", 2);
0255 
0256   if (detailedAnalysis) {
0257     histoName = "NoisyChannels";
0258 
0259     threshChannelsHisto = ibooker.book1D(histoName.c_str(), "# of noisy channels vs threshold", 15, 500, 2000);
0260     threshChannelsHisto->setAxisTitle("threshold", 1);
0261     threshChannelsHisto->setAxisTitle("# noisy channels", 2);
0262   }
0263 
0264   if (doSynchNoise) {
0265     ibooker.setCurrentFolder("DT/05-Noise/SynchNoise/");
0266     histoName = "SynchNoiseSummary";
0267 
0268     summarySynchNoiseHisto = ibooker.book2D(histoName.c_str(), "Summary Synch. Noise", 12, 1, 13, 5, -2, 3);
0269     summarySynchNoiseHisto->setAxisTitle("Sector", 1);
0270     summarySynchNoiseHisto->setAxisTitle("Wheel", 2);
0271     histoName = "SynchNoiseGlbSummary";
0272 
0273     glbSummarySynchNoiseHisto = ibooker.book2D(histoName.c_str(), "Summary Synch. Noise", 12, 1, 13, 5, -2, 3);
0274     glbSummarySynchNoiseHisto->setAxisTitle("Sector", 1);
0275     glbSummarySynchNoiseHisto->setAxisTitle("Wheel", 2);
0276   }
0277 }
0278 
0279 string DTNoiseAnalysisTest::getSynchNoiseMEName(int wheelId) const {
0280   stringstream wheel;
0281   wheel << wheelId;
0282   string folderName = "DT/05-Noise/SynchNoise/";
0283   string histoname = folderName + string("SyncNoiseEvents") + "_W" + wheel.str();
0284 
0285   return histoname;
0286 }
0287 
0288 void DTNoiseAnalysisTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {}