CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes
NoisyChannel Class Reference

Check if any channels are noisy compared to neighboring ones. More...

#include <QTest.h>

Inheritance diagram for NoisyChannel:
SimpleTest QCriterion

Public Member Functions

 NoisyChannel (const std::string &name)
 
float runTest (const MonitorElement *me) override
 
void setNumNeighbors (unsigned n)
 
void setTolerance (float percentage)
 
- Public Member Functions inherited from SimpleTest
std::vector< DQMChannelgetBadChannels () const override
 get vector of channels that failed test (not always relevant!) More...
 
void setMinimumEntries (unsigned n)
 set minimum # of entries needed More...
 
 SimpleTest (const std::string &name, bool keepBadChannels=false)
 
- Public Member Functions inherited from QCriterion
std::string algoName () const
 get algorithm name More...
 
std::string getMessage () const
 get message attached to test More...
 
std::string getName () const
 get name of quality test More...
 
int getStatus () const
 get test status More...
 
void init ()
 initialize values More...
 
 QCriterion (std::string qtname)
 
float runTest (const MonitorElement *me, QReport &qr, DQMNet::QValue &qv)
 
void setErrorProb (float prob)
 
void setWarningProb (float prob)
 set probability limit for warning and error (default: 90% and 50%) More...
 
virtual ~QCriterion ()=default
 

Static Public Member Functions

static std::string getAlgoName ()
 

Protected Member Functions

double getAverage (int bin, const TH1 *h) const
 
double getAverage2D (int binX, int binY, const TH2 *h) const
 
- Protected Member Functions inherited from SimpleTest
void setMessage () override
 set status & message after test has run More...
 
- Protected Member Functions inherited from QCriterion
void setAlgoName (std::string name)
 set algorithm name More...
 
void setVerbose (int verbose)
 probability limits for warnings, errors More...
 

Protected Attributes

unsigned numNeighbors_
 
bool rangeInitialized_
 
float tolerance_
 
- Protected Attributes inherited from SimpleTest
std::vector< DQMChannelbadChannels_
 
bool keepBadChannels_
 
unsigned minEntries_
 
- Protected Attributes inherited from QCriterion
std::string algoName_
 name of quality test More...
 
float errorProb_
 
std::string message_
 quality test status More...
 
float prob_
 name of algorithm More...
 
std::string qtname_
 
int status_
 
int verbose_
 
float warningProb_
 message attached to test More...
 

Additional Inherited Members

- Public Types inherited from QCriterion
typedef dqm::legacy::MonitorElement MonitorElement
 (class should be created by DQMStore class) More...
 
- Static Public Attributes inherited from QCriterion
static const float ERROR_PROB_THRESHOLD = 0.50
 
static const float WARNING_PROB_THRESHOLD = 0.90
 default "probability" values for setting warnings & errors when running tests More...
 

Detailed Description

Check if any channels are noisy compared to neighboring ones.

Definition at line 220 of file QTest.h.

Constructor & Destructor Documentation

◆ NoisyChannel()

NoisyChannel::NoisyChannel ( const std::string &  name)
inline

Definition at line 222 of file QTest.h.

222  : SimpleTest(name, true) {
223  rangeInitialized_ = false;
224  numNeighbors_ = 1;
226  }

References getAlgoName(), numNeighbors_, rangeInitialized_, and QCriterion::setAlgoName().

Member Function Documentation

◆ getAlgoName()

static std::string NoisyChannel::getAlgoName ( )
inlinestatic

Definition at line 227 of file QTest.h.

227 { return "NoisyChannel"; }

Referenced by QualityTester::makeQCriterion(), and NoisyChannel().

◆ getAverage()

double NoisyChannel::getAverage ( int  bin,
const TH1 *  h 
) const
protected

get average for bin under consideration (see description of method setNumNeighbors)

Definition at line 375 of file QTest.cc.

375  {
376  int first = 1; // Do NOT use underflow bin
377  int ncx = h->GetXaxis()->GetNbins(); // Do NOT use overflow bin
378  double sum = 0;
379  int bin_start, bin_end;
380  int add_right = 0;
381  int add_left = 0;
382 
383  bin_start = bin - numNeighbors_; // First bin in integral
384  bin_end = bin + numNeighbors_; // Last bin in integral
385 
386  if (bin_start < first) { // If neighbors take you outside of histogram range shift integral right
387  add_right = first - bin_start; // How much to shift remembering we are not using underflow
388  bin_start = first; // Remember to reset the starting bin
389  bin_end += add_right;
390  if (bin_end > ncx)
391  bin_end = ncx; // If the test would be larger than histogram just sum histogram without overflow
392  }
393 
394  if (bin_end > ncx) { // Now we make sure doesn't run off right edge of histogram
395  add_left = bin_end - ncx;
396  bin_end = ncx;
397  bin_start -= add_left;
398  if (bin_start < first)
399  bin_start = first; // If the test would be larger than histogram just sum histogram without underflow
400  }
401 
402  sum += h->Integral(bin_start, bin_end);
403  sum -= h->GetBinContent(bin);
404 
405  int dimension = 2 * numNeighbors_ + 1;
406  if (dimension > h->GetNbinsX())
407  dimension = h->GetNbinsX();
408 
409  return sum / (dimension - 1);
410 }

References newFWLiteAna::bin, pat::helper::ParametrizationHelper::dimension(), and dqmdumpme::first.

◆ getAverage2D()

double NoisyChannel::getAverage2D ( int  binX,
int  binY,
const TH2 *  h 
) const
protected

Do NOT use underflow or overflow bins

Definition at line 412 of file QTest.cc.

412  {
414  int firstX = 1;
415  int firstY = 1;
416  double sum = 0;
417  int ncx = h2->GetXaxis()->GetNbins();
418  int ncy = h2->GetYaxis()->GetNbins();
419 
420  int neighborsX, neighborsY; // Convert unsigned input to int so we can use comparators
421  neighborsX = numNeighbors_;
422  neighborsY = numNeighbors_;
423  int bin_startX, bin_endX;
424  int add_rightX = 0; // Start shifts at 0
425  int add_leftX = 0;
426  int bin_startY, bin_endY;
427  int add_topY = 0;
428  int add_downY = 0;
429 
430  bin_startX = binX - neighborsX; // First bin in X
431  bin_endX = binX + neighborsX; // Last bin in X
432 
433  if (bin_startX < firstX) { // If neighbors take you outside of histogram range shift integral right
434  add_rightX = firstX - bin_startX; // How much to shift remembering we are no using underflow
435  bin_startX = firstX; // Remember to reset the starting bin
436  bin_endX += add_rightX;
437  if (bin_endX > ncx)
438  bin_endX = ncx;
439  }
440 
441  if (bin_endX > ncx) { // Now we make sure doesn't run off right edge of histogram
442  add_leftX = bin_endX - ncx;
443  bin_endX = ncx;
444  bin_startX -= add_leftX;
445  if (bin_startX < firstX)
446  bin_startX = firstX; // If the test would be larger than histogram just sum histogram without underflow
447  }
448 
449  bin_startY = binY - neighborsY; // First bin in Y
450  bin_endY = binY + neighborsY; // Last bin in Y
451 
452  if (bin_startY < firstY) { // If neighbors take you outside of histogram range shift integral up
453  add_topY = firstY - bin_startY; // How much to shift remembering we are no using underflow
454  bin_startY = firstY; // Remember to reset the starting bin
455  bin_endY += add_topY;
456  if (bin_endY > ncy)
457  bin_endY = ncy;
458  }
459 
460  if (bin_endY > ncy) { // Now we make sure doesn't run off top edge of histogram
461  add_downY = bin_endY - ncy;
462  bin_endY = ncy;
463  bin_startY -= add_downY;
464  if (bin_startY < firstY)
465  bin_startY = firstY; // If the test would be larger than histogram just sum histogram without underflow
466  }
467 
468  sum += h2->Integral(bin_startX, bin_endX, bin_startY, bin_endY);
469  sum -= h2->GetBinContent(binX, binY);
470 
471  int dimensionX = 2 * neighborsX + 1;
472  int dimensionY = 2 * neighborsY + 1;
473 
474  if (dimensionX > h2->GetNbinsX())
475  dimensionX = h2->GetNbinsX();
476  if (dimensionY > h2->GetNbinsY())
477  dimensionY = h2->GetNbinsY();
478 
479  return sum / (dimensionX * dimensionY - 1); // Average is sum over the # of bins used
480 
481 } // End getAverage2D

◆ runTest()

float NoisyChannel::runTest ( const MonitorElement me)
overridevirtual

Reimplemented from QCriterion.

Definition at line 263 of file QTest.cc.

263  {
264  badChannels_.clear();
265  if (!me)
266  return -1;
267  if (!me->getRootObject())
268  return -1;
269  TH1* h = nullptr; //initialize histogram pointer
270  TH2* h2 = nullptr; //initialize histogram pointer
271 
272  if (verbose_ > 1)
273  std::cout << "QTest:" << getAlgoName() << "::runTest called on " << me->getFullname() << "\n";
274 
275  int nbins = 0;
276  int nbinsX = 0, nbinsY = 0;
277  //-- TH1F
278  if (me->kind() == MonitorElement::Kind::TH1F) {
279  nbins = me->getTH1F()->GetXaxis()->GetNbins();
280  h = me->getTH1F(); // access Test histo
281  }
282  //-- TH1S
283  else if (me->kind() == MonitorElement::Kind::TH1S) {
284  nbins = me->getTH1S()->GetXaxis()->GetNbins();
285  h = me->getTH1S(); // access Test histo
286  }
287  //-- TH1D
288  else if (me->kind() == MonitorElement::Kind::TH1D) {
289  nbins = me->getTH1D()->GetXaxis()->GetNbins();
290  h = me->getTH1D(); // access Test histo
291  }
292  //-- TH2
293  else if (me->kind() == MonitorElement::Kind::TH2F) {
294  nbinsX = me->getTH2F()->GetXaxis()->GetNbins();
295  nbinsY = me->getTH2F()->GetYaxis()->GetNbins();
296  h2 = me->getTH2F(); // access Test histo
297  }
298  //-- TH2
299  else if (me->kind() == MonitorElement::Kind::TH2S) {
300  nbinsX = me->getTH2S()->GetXaxis()->GetNbins();
301  nbinsY = me->getTH2S()->GetYaxis()->GetNbins();
302  h2 = me->getTH2S(); // access Test histo
303  }
304  //-- TH2
305  else if (me->kind() == MonitorElement::Kind::TH2D) {
306  nbinsX = me->getTH2F()->GetXaxis()->GetNbins();
307  nbinsY = me->getTH2F()->GetYaxis()->GetNbins();
308  h2 = me->getTH2D(); // access Test histo
309  } else {
310  if (verbose_ > 0)
311  std::cout << "QTest:NoisyChannel"
312  << " ME " << me->getFullname() << " does not contain TH1F/TH1S/TH1D or TH2F/TH2S/TH2D, exiting\n";
313  return -1;
314  }
315 
316  //-- QUALITY TEST itself
317 
318  // do NOT use underflow bin
319  int first = 1;
320  // do NOT use overflow bin
321  int last = nbins;
322  int lastX = nbinsX, lastY = nbinsY;
323  // bins outside Y-range
324  int fail = 0;
325  int bin;
326  int binX, binY;
327  if (h != nullptr) {
328  if (!rangeInitialized_ || !h->GetXaxis()) {
329  return 1; // all channels are accepted if tolerance has not been set
330  }
331  for (bin = first; bin <= last; ++bin) {
332  double contents = h->GetBinContent(bin);
333  double average = getAverage(bin, h);
334  bool failure = false;
335  if (average != 0)
336  failure = (((contents - average) / std::abs(average)) > tolerance_);
337 
338  if (failure) {
339  ++fail;
340  DQMChannel chan(bin, 0, 0, contents, h->GetBinError(bin));
341  badChannels_.push_back(chan);
342  }
343  }
344 
345  // return fraction of bins that passed test
346  return 1. * (nbins - fail) / nbins;
347  } else if (h2 != nullptr) {
348  for (binY = first; binY <= lastY; ++binY) {
349  for (binX = first; binX <= lastX; ++binX) {
350  double contents = h2->GetBinContent(binX, binY);
351  double average = getAverage2D(binX, binY, h2);
352  bool failure = false;
353  if (average != 0)
354  failure = (((contents - average) / std::abs(average)) > tolerance_);
355  if (failure) {
356  ++fail;
357  DQMChannel chan(binX, 0, 0, contents, h2->GetBinError(binX));
358  badChannels_.push_back(chan);
359  }
360  } //end x loop
361  } //end y loop
362  // return fraction of bins that passed test
363  return 1. * ((nbinsX * nbinsY) - fail) / (nbinsX * nbinsY);
364  } //end nullptr conditional
365  else {
366  if (verbose_ > 0)
367  std::cout << "QTest:NoisyChannel"
368  << " TH1/TH2F are NULL, exiting\n";
369  return -1;
370  }
371 }

References funct::abs(), PDRates::average, newFWLiteAna::bin, officialStyle::chan, relmon_rootfiles_spy::contents, gather_cfg::cout, dqmdumpme::first, dqmdumpme::last, hlt_dqm_clientPB-live_cfg::me, LaserClient_cfi::nbins, hlt_dqm_clientPB-live_cfg::nbinsX, hlt_dqm_clientPB-live_cfg::nbinsY, MonitorElementData::TH1D, MonitorElementData::TH1F, MonitorElementData::TH1S, MonitorElementData::TH2D, MonitorElementData::TH2F, and MonitorElementData::TH2S.

◆ setNumNeighbors()

void NoisyChannel::setNumNeighbors ( unsigned  n)
inline

set # of neighboring channels for calculating average to be used for comparison with channel under consideration; use 1 for considering bin+1 and bin-1 (default), use 2 for considering bin+1,bin-1, bin+2,bin-2, etc; Will use rollover when bin+i or bin-i is beyond histogram limits (e.g. for histogram with N bins, bin N+1 corresponds to bin 1, and bin -1 corresponds to bin N)

Definition at line 237 of file QTest.h.

237  {
238  if (n > 0)
239  numNeighbors_ = n;
240  }

References dqmiodumpmetadata::n, and numNeighbors_.

◆ setTolerance()

void NoisyChannel::setTolerance ( float  percentage)
inline

set (percentage) tolerance for considering a channel noisy; eg. if tolerance = 20%, a channel will be noisy if (contents-average)/|average| > 20%; average is calculated from neighboring channels (also see method setNumNeighbors)

Definition at line 246 of file QTest.h.

246  {
247  if (percentage >= 0) {
248  tolerance_ = percentage;
249  rangeInitialized_ = true;
250  }
251  }

References rangeInitialized_, and tolerance_.

Member Data Documentation

◆ numNeighbors_

unsigned NoisyChannel::numNeighbors_
protected

Definition at line 260 of file QTest.h.

Referenced by NoisyChannel(), and setNumNeighbors().

◆ rangeInitialized_

bool NoisyChannel::rangeInitialized_
protected

Definition at line 262 of file QTest.h.

Referenced by NoisyChannel(), and setTolerance().

◆ tolerance_

float NoisyChannel::tolerance_
protected

Definition at line 259 of file QTest.h.

Referenced by setTolerance().

pat::helper::ParametrizationHelper::dimension
uint32_t dimension(pat::CandKinResolution::Parametrization parametrization)
Returns the number of free parameters in a parametrization (3 or 4)
Definition: ParametrizationHelper.h:12
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
NoisyChannel::getAlgoName
static std::string getAlgoName()
Definition: QTest.h:227
MonitorElementData::Kind::TH1S
MonitorElementData::Kind::TH1F
gather_cfg.cout
cout
Definition: gather_cfg.py:144
hlt_dqm_clientPB-live_cfg.nbinsX
nbinsX
Definition: hlt_dqm_clientPB-live_cfg.py:60
NoisyChannel::tolerance_
float tolerance_
Definition: QTest.h:259
dqmdumpme.first
first
Definition: dqmdumpme.py:55
MonitorElementData::Kind::TH2D
MonitorElementData::Kind::TH2F
PDRates.average
average
Definition: PDRates.py:139
relmon_rootfiles_spy.contents
contents
Definition: relmon_rootfiles_spy.py:129
SimpleTest::SimpleTest
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:124
dqmdumpme.last
last
Definition: dqmdumpme.py:56
NoisyChannel::getAverage2D
double getAverage2D(int binX, int binY, const TH2 *h) const
Definition: QTest.cc:412
h
QCriterion::verbose_
int verbose_
Definition: QTest.h:110
QCriterion::setAlgoName
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:96
LaserClient_cfi.nbins
nbins
Definition: LaserClient_cfi.py:51
hlt_dqm_clientPB-live_cfg.nbinsY
nbinsY
Definition: hlt_dqm_clientPB-live_cfg.py:64
NoisyChannel::numNeighbors_
unsigned numNeighbors_
Definition: QTest.h:260
average
Definition: average.py:1
MonitorElementData::Kind::TH1D
SimpleTest::badChannels_
std::vector< DQMChannel > badChannels_
Definition: QTest.h:139
newFWLiteAna.bin
bin
Definition: newFWLiteAna.py:161
MonitorElementData::Kind::TH2S
NoisyChannel::getAverage
double getAverage(int bin, const TH1 *h) const
Definition: QTest.cc:375
MonitorElementData::QReport::DQMChannel
Definition: MonitorElementCollection.h:64
officialStyle.chan
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi....
Definition: officialStyle.py:106
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
NoisyChannel::rangeInitialized_
bool rangeInitialized_
Definition: QTest.h:262
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
hlt_dqm_clientPB-live_cfg.me
me
Definition: hlt_dqm_clientPB-live_cfg.py:56