Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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