Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305
/*
 *  See header file for a description of this class.
 *
 *  \author G. Mila - INFN Torino
 *
 *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
 *
 */

#include "DQM/DTMonitorClient/src/DTNoiseAnalysisTest.h"

// Framework
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"

// Geometry
#include "Geometry/DTGeometry/interface/DTGeometry.h"
#include "Geometry/DTGeometry/interface/DTLayer.h"
#include "Geometry/DTGeometry/interface/DTTopology.h"
#include "DataFormats/MuonDetId/interface/DTLayerId.h"

#include "DQMServices/Core/interface/DQMStore.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include <iostream>
#include <sstream>

using namespace edm;
using namespace std;

DTNoiseAnalysisTest::DTNoiseAnalysisTest(const edm::ParameterSet& ps)
    : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
  LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "[DTNoiseAnalysisTest]: Constructor";

  // get the cfi parameters
  noisyCellDef = ps.getUntrackedParameter<int>("noisyCellDef",
                                               500);  //value set in DQM/DTMonitorClient/python/dtNoiseAnalysis_cfi.py
  isCosmics = ps.getUntrackedParameter<bool>("isCosmics", false);

  // switch on/off the summaries for the Synchronous noise
  doSynchNoise = ps.getUntrackedParameter<bool>("doSynchNoise", false);
  detailedAnalysis = ps.getUntrackedParameter<bool>("detailedAnalysis", false);
  maxSynchNoiseRate = ps.getUntrackedParameter<double>("maxSynchNoiseRate", 0.001);
  noiseSafetyFactor = ps.getUntrackedParameter<double>("noiseSafetyFactor", 5.);  //for collisions
  nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000);

  nevents = 0;

  bookingdone = false;
}

DTNoiseAnalysisTest::~DTNoiseAnalysisTest() {
  LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "DTNoiseAnalysisTest: analyzed " << nevents << " events";
}

void DTNoiseAnalysisTest::beginRun(Run const& run, EventSetup const& context) {
  // Get the geometry
  muonGeom = &context.getData(muonGeomToken_);
}

void DTNoiseAnalysisTest::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
                                                DQMStore::IGetter& igetter,
                                                edm::LuminosityBlock const& lumiSeg,
                                                edm::EventSetup const& context) {
  float chRate;

  if (!bookingdone) {
    // book the histos
    bookHistos(ibooker);
  }
  bookingdone = true;

  LogVerbatim("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
      << "[DTNoiseAnalysisTest]: End of LS transition, performing the DQM client operation";

  // Reset the summary plots
  for (map<int, MonitorElement*>::iterator plot = noiseHistos.begin(); plot != noiseHistos.end(); ++plot) {
    (*plot).second->Reset();
  }

  for (map<int, MonitorElement*>::iterator plot = noisyCellHistos.begin(); plot != noisyCellHistos.end(); ++plot) {
    (*plot).second->Reset();
  }

  summaryNoiseHisto->Reset();

  vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
  vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();

  LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest") << "[DTNoiseAnalysisTest]: Fill the summary histos";

  for (; ch_it != ch_end; ++ch_it) {  // loop over chambers
    DTChamberId chID = (*ch_it)->id();

    MonitorElement* histo = igetter.get(getMEName(chID));

    if (histo) {  // check the pointer

      TH2D* histo_root = histo->getTH2D();

      for (int sl = 1; sl != 4; ++sl) {  // loop over SLs
        // skip theta SL in MB4 chambers
        if (chID.station() == 4 && sl == 2)
          continue;

        int binYlow = ((sl - 1) * 4) + 1;

        for (int layer = 1; layer <= 4; ++layer) {  // loop over layers

          // Get the layer ID
          DTLayerId layID(chID, sl, layer);

          int nWires = muonGeom->layer(layID)->specificTopology().channels();
          int firstWire = muonGeom->layer(layID)->specificTopology().firstChannel();

          int binY = binYlow + (layer - 1);

          for (int wire = firstWire; wire != (nWires + firstWire); wire++) {  // loop over wires

            double noise = histo_root->GetBinContent(wire, binY);
            // fill the histos
            noiseHistos[chID.wheel()]->Fill(noise);
            noiseHistos[3]->Fill(noise);
            int sector = chID.sector();

            if (!isCosmics) {  // for collisions
              float k = kW_MB[2 - abs(chID.wheel())][chID.station() - 1];
              if (chID.station() == 4) {  //special geometry cases for MB4
                if (sector == 9 || sector == 10 || sector == 11)
                  k = 0.05;
                else if (sector == 4 && chID.wheel() == 0)
                  k = 0.25;
              }
              noisyCellDef =
                  cellW * instLumi * k * lenghtSL_MB[(sl % 2)][chID.station() - 1];  // background expected per chamber
              noisyCellDef = noisyCellDef * noiseSafetyFactor;
            }  //else value read from DQM/DTMonitorClient/python/dtNoiseAnalysis_cfi.py

            if (noise > noisyCellDef) {
              if (sector == 13) {
                sector = 4;
              } else if (sector == 14) {
                sector = 10;
              }
              noisyCellHistos[chID.wheel()]->Fill(sector, chID.station());
              summaryNoiseHisto->Fill(sector, chID.wheel());
            }
          }
        }
      }
    }
  }

  if (detailedAnalysis) {
    threshChannelsHisto->Reset();
    TH1F* histo = noiseHistos[3]->getTH1F();
    for (int step = 0; step != 15; step++) {
      int threshBin = step + 1;
      int minBin = 26 + step * 5;
      int nNoisyCh = histo->Integral(minBin, 101);
      threshChannelsHisto->setBinContent(threshBin, nNoisyCh);
    }
  }

  // build the summary of synch noise

  if (doSynchNoise) {
    LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
        << "[DTNoiseAnalysisTest]: fill summaries for synch noise" << endl;
    summarySynchNoiseHisto->Reset();
    glbSummarySynchNoiseHisto->Reset();
    for (int wheel = -2; wheel != 3; ++wheel) {
      // Get the histo produced by DTDigiTask

      MonitorElement* histoNoiseSynch = igetter.get(getSynchNoiseMEName(wheel));
      if (histoNoiseSynch != nullptr) {
        for (int sect = 1; sect != 13; ++sect) {  // loop over sectors
          TH2F* histo = histoNoiseSynch->getTH2F();
          float maxSectRate = 0;
          for (int sta = 1; sta != 5; ++sta) {
            if (nevents > 0)
              chRate = histo->GetBinContent(sect, sta) / (float)nevents;
            else
              chRate = -1.0;
            // in case nevents 0 e.g. counting not done
            LogTrace("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
                << "   Wheel: " << wheel << " sect: " << sect << " station: " << sta << " rate is: " << chRate << endl;
            if (chRate > maxSectRate)
              maxSectRate = chRate;
          }
          summarySynchNoiseHisto->Fill(sect, wheel, maxSectRate > maxSynchNoiseRate ? 1 : 0);
          float glbBinValue = 1 - 0.15 * maxSectRate / maxSynchNoiseRate;
          glbSummarySynchNoiseHisto->Fill(sect, wheel, glbBinValue > 0 ? glbBinValue : 0);
        }
      } else {
        LogWarning("DTDQM|DTMonitorClient|DTNoiseAnalysisTest")
            << "   Histo: " << getSynchNoiseMEName(wheel) << " not found!" << endl;
      }
    }
  }

  string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsNoise";

  MonitorElement* meProcEvts = igetter.get(nEvtsName);

  if (meProcEvts) {
    int nProcEvts = meProcEvts->getFloatValue();
    glbSummarySynchNoiseHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
    summarySynchNoiseHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
  } else {
    glbSummarySynchNoiseHisto->setEntries(nMinEvts + 1);
    summarySynchNoiseHisto->setEntries(nMinEvts + 1);
    LogVerbatim("DTDQM|DTMonitorClient|DTnoiseAnalysisTest")
        << "[DTNoiseAnalysisTest] ME: " << nEvtsName << " not found!" << endl;
  }
}

string DTNoiseAnalysisTest::getMEName(const DTChamberId& chID) {
  stringstream wheel;
  wheel << chID.wheel();
  stringstream station;
  station << chID.station();
  stringstream sector;
  sector << chID.sector();

  string folderName = "DT/05-Noise/Wheel" + wheel.str() + "/Sector" + sector.str() + "/";

  string histoname =
      folderName + string("NoiseRate") + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();

  return histoname;
}

void DTNoiseAnalysisTest::bookHistos(DQMStore::IBooker& ibooker) {
  ibooker.setCurrentFolder("DT/05-Noise");
  string histoName;

  for (int wh = -2; wh <= 2; wh++) {
    stringstream wheel;
    wheel << wh;
    histoName = "NoiseRateSummary_W" + wheel.str();

    noiseHistos[wh] = ibooker.book1D(histoName.c_str(), histoName.c_str(), 500, 0, 10000);
    noiseHistos[wh]->setAxisTitle("rate (Hz)", 1);
    noiseHistos[wh]->setAxisTitle("entries", 2);
  }
  histoName = "NoiseRateSummary";

  noiseHistos[3] = ibooker.book1D(histoName.c_str(), histoName.c_str(), 500, 0, 10000);
  noiseHistos[3]->setAxisTitle("rate (Hz)", 1);
  noiseHistos[3]->setAxisTitle("entries", 2);

  for (int wh = -2; wh <= 2; wh++) {
    stringstream wheel;
    wheel << wh;
    histoName = "NoiseSummary_W" + wheel.str();

    noisyCellHistos[wh] = ibooker.book2D(histoName.c_str(), "# of noisy channels", 12, 1, 13, 4, 1, 5);
    noisyCellHistos[wh]->setBinLabel(1, "MB1", 2);
    noisyCellHistos[wh]->setBinLabel(2, "MB2", 2);
    noisyCellHistos[wh]->setBinLabel(3, "MB3", 2);
    noisyCellHistos[wh]->setBinLabel(4, "MB4", 2);
    noisyCellHistos[wh]->setAxisTitle("Sector", 1);
  }

  histoName = "NoiseSummary";

  summaryNoiseHisto = ibooker.book2D(histoName.c_str(), "# of noisy channels", 12, 1, 13, 5, -2, 3);
  summaryNoiseHisto->setAxisTitle("Sector", 1);
  summaryNoiseHisto->setAxisTitle("Wheel", 2);

  if (detailedAnalysis) {
    histoName = "NoisyChannels";

    threshChannelsHisto = ibooker.book1D(histoName.c_str(), "# of noisy channels vs threshold", 15, 500, 2000);
    threshChannelsHisto->setAxisTitle("threshold", 1);
    threshChannelsHisto->setAxisTitle("# noisy channels", 2);
  }

  if (doSynchNoise) {
    ibooker.setCurrentFolder("DT/05-Noise/SynchNoise/");
    histoName = "SynchNoiseSummary";

    summarySynchNoiseHisto = ibooker.book2D(histoName.c_str(), "Summary Synch. Noise", 12, 1, 13, 5, -2, 3);
    summarySynchNoiseHisto->setAxisTitle("Sector", 1);
    summarySynchNoiseHisto->setAxisTitle("Wheel", 2);
    histoName = "SynchNoiseGlbSummary";

    glbSummarySynchNoiseHisto = ibooker.book2D(histoName.c_str(), "Summary Synch. Noise", 12, 1, 13, 5, -2, 3);
    glbSummarySynchNoiseHisto->setAxisTitle("Sector", 1);
    glbSummarySynchNoiseHisto->setAxisTitle("Wheel", 2);
  }
}

string DTNoiseAnalysisTest::getSynchNoiseMEName(int wheelId) const {
  stringstream wheel;
  wheel << wheelId;
  string folderName = "DT/05-Noise/SynchNoise/";
  string histoname = folderName + string("SyncNoiseEvents") + "_W" + wheel.str();

  return histoname;
}

void DTNoiseAnalysisTest::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {}