CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
void setNumNeighbors (unsigned n)
 
void setTolerance (float percentage)
 
- Public Member Functions inherited from SimpleTest
virtual std::vector< DQMChannelgetBadChannels (void) const
 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 (void) const
 get algorithm name More...
 
std::string getMessage (void) const
 get message attached to test More...
 
std::string getName (void) const
 get name of quality test More...
 
int getStatus (void) const
 (class should be created by DQMStore class) More...
 
void setErrorProb (float prob)
 
void setWarningProb (float prob)
 set probability limit for warning and error (default: 90% and 50%) More...
 

Static Public Member Functions

static std::string getAlgoName (void)
 

Protected Member Functions

double getAverage (int bin, const TH1 *h) const
 
- Protected Member Functions inherited from SimpleTest
virtual void setMessage (void)
 set status & message after test has run More...
 
- Protected Member Functions inherited from QCriterion
void init (void)
 initialize values More...
 
 QCriterion (std::string qtname)
 
float runTest (const MonitorElement *me, QReport &qr, DQMNet::QValue &qv)
 
void setAlgoName (std::string name)
 set algorithm name More...
 
void setVerbose (int verbose)
 probability limits for warnings, errors More...
 
virtual ~QCriterion (void)
 

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...
 

Detailed Description

Check if any channels are noisy compared to neighboring ones.

Definition at line 303 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 306 of file QTest.h.

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

306  : SimpleTest(name,true)
307  {
308  rangeInitialized_ = false;
309  numNeighbors_ = 1;
311  }
static std::string getAlgoName(void)
Definition: QTest.h:312
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:77
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:137
bool rangeInitialized_
Definition: QTest.h:345
unsigned numNeighbors_
Definition: QTest.h:343

Member Function Documentation

static std::string NoisyChannel::getAlgoName ( void  )
inlinestatic
double NoisyChannel::getAverage ( int  bin,
const TH1 *  h 
) const
protected

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

do NOT use underflow bin

do NOT use overflow bin

use symmetric-to-bin bins to calculate average

check if need to consider bins on other side of spectrum (ie. if bins below 1 or above ncx)

average is sum over the # of bins used

Definition at line 843 of file QTest.cc.

References first, and i.

844 {
846  int first = 1;
848  int ncx = h->GetXaxis()->GetNbins();
849  double sum = 0; int bin_low, bin_hi;
850  for (unsigned i = 1; i <= numNeighbors_; ++i)
851  {
853  bin_low = bin-i; bin_hi = bin+i;
856  while (bin_low < first) // shift bin by +ncx
857  bin_low = ncx + bin_low;
858  while (bin_hi > ncx) // shift bin by -ncx
859  bin_hi = bin_hi - ncx;
860  sum += h->GetBinContent(bin_low) + h->GetBinContent(bin_hi);
861  }
863  return sum/(numNeighbors_ * 2);
864 }
int i
Definition: DBlmapReader.cc:9
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
bool first
Definition: L1TdeRCT.cc:75
unsigned numNeighbors_
Definition: QTest.h:343
float NoisyChannel::runTest ( const MonitorElement me)
virtual

Reimplemented from QCriterion.

Definition at line 747 of file QTest.cc.

References Abs(), PDRates::average, newFWLiteAna::bin, relmon_rootfiles_spy::contents, gather_cfg::cout, MonitorElement::DQM_KIND_TH1D, MonitorElement::DQM_KIND_TH1F, MonitorElement::DQM_KIND_TH1S, MonitorElement::DQM_KIND_TH2D, MonitorElement::DQM_KIND_TH2F, MonitorElement::DQM_KIND_TH2S, cmsPerfPublish::fail(), first, MonitorElement::getFullname(), MonitorElement::getRootObject(), MonitorElement::getTH1D(), MonitorElement::getTH1F(), MonitorElement::getTH1S(), MonitorElement::getTH2D(), MonitorElement::getTH2F(), MonitorElement::getTH2S(), h, MonitorElement::kind(), prof2calltree::last, and pileupCalc::nbins.

748 {
749  badChannels_.clear();
750  if (!me)
751  return -1;
752  if (!me->getRootObject())
753  return -1;
754  TH1* h=0;//initialize histogram pointer
755 
756  if (verbose_>1)
757  std::cout << "QTest:" << getAlgoName() << "::runTest called on "
758  << me-> getFullname() << "\n";
759 
760  int nbins=0;
761  //-- TH1F
763  {
764  nbins = me->getTH1F()->GetXaxis()->GetNbins();
765  h = me->getTH1F(); // access Test histo
766  }
767  //-- TH1S
768  else if (me->kind()==MonitorElement::DQM_KIND_TH1S)
769  {
770  nbins = me->getTH1S()->GetXaxis()->GetNbins();
771  h = me->getTH1S(); // access Test histo
772  }
773  //-- TH1D
774  else if (me->kind()==MonitorElement::DQM_KIND_TH1D)
775  {
776  nbins = me->getTH1D()->GetXaxis()->GetNbins();
777  h = me->getTH1D(); // access Test histo
778  }
779  //-- TH2
780  else if (me->kind()==MonitorElement::DQM_KIND_TH2F)
781  {
782  nbins = me->getTH2F()->GetXaxis()->GetNbins() *
783  me->getTH2F()->GetYaxis()->GetNbins();
784  h = me->getTH2F(); // access Test histo
785  }
786  //-- TH2
787  else if (me->kind()==MonitorElement::DQM_KIND_TH2S)
788  {
789  nbins = me->getTH2S()->GetXaxis()->GetNbins() *
790  me->getTH2S()->GetYaxis()->GetNbins();
791  h = me->getTH2S(); // access Test histo
792  }
793  //-- TH2
794  else if (me->kind()==MonitorElement::DQM_KIND_TH2D)
795  {
796  nbins = me->getTH2D()->GetXaxis()->GetNbins() *
797  me->getTH2D()->GetYaxis()->GetNbins();
798  h = me->getTH2D(); // access Test histo
799  }
800  else
801  {
802  if (verbose_>0)
803  std::cout << "QTest:NoisyChannel"
804  << " ME " << me->getFullname()
805  << " does not contain TH1F/TH1S/TH1D or TH2F/TH2S/TH2D, exiting\n";
806  return -1;
807  }
808 
809  //-- QUALITY TEST itself
810 
811  if ( !rangeInitialized_ || !h->GetXaxis() )
812  return 1; // all channels are accepted if tolerance has not been set
813 
814  // do NOT use underflow bin
815  int first = 1;
816  // do NOT use overflow bin
817  int last = nbins;
818  // bins outside Y-range
819  int fail = 0;
820  int bin;
821  for (bin = first; bin <= last; ++bin)
822  {
823  double contents = h->GetBinContent(bin);
824  double average = getAverage(bin, h);
825  bool failure = false;
826  if (average != 0)
827  failure = (((contents-average)/TMath::Abs(average)) > tolerance_);
828 
829  if (failure)
830  {
831  ++fail;
832  DQMChannel chan(bin, 0, 0, contents, h->GetBinError(bin));
833  badChannels_.push_back(chan);
834  }
835  }
836 
837  // return fraction of bins that passed test
838  return 1.*(nbins - fail)/nbins;
839 }
TH2S * getTH2S(void) const
TH1S * getTH1S(void) const
static std::string getAlgoName(void)
Definition: QTest.h:312
double getAverage(int bin, const TH1 *h) const
Definition: QTest.cc:843
int verbose_
Definition: QTest.h:117
TH1D * getTH1D(void) const
TH2D * getTH2D(void) const
T Abs(T a)
Definition: MathUtil.h:49
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Kind kind(void) const
Get the type of the monitor element.
bool first
Definition: L1TdeRCT.cc:75
const std::string getFullname(void) const
get full name of ME including Pathname
TObject * getRootObject(void) const
TH1F * getTH1F(void) const
std::vector< DQMChannel > badChannels_
Definition: QTest.h:159
int average
Definition: PDRates.py:137
float tolerance_
Definition: QTest.h:342
tuple cout
Definition: gather_cfg.py:121
TH2F * getTH2F(void) const
bool rangeInitialized_
Definition: QTest.h:345
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 322 of file QTest.h.

References n, and numNeighbors_.

Referenced by QTestConfigure::EnableNoisyChannelTest().

322 { if (n > 0) numNeighbors_ = n; }
unsigned numNeighbors_
Definition: QTest.h:343
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 328 of file QTest.h.

References rangeInitialized_, and tolerance_.

Referenced by QTestConfigure::EnableNoisyChannelTest().

329  {
330  if (percentage >=0)
331  {
332  tolerance_ = percentage;
333  rangeInitialized_ = true;
334  }
335  }
float tolerance_
Definition: QTest.h:342
bool rangeInitialized_
Definition: QTest.h:345

Member Data Documentation

unsigned NoisyChannel::numNeighbors_
protected

Definition at line 343 of file QTest.h.

Referenced by NoisyChannel(), and setNumNeighbors().

bool NoisyChannel::rangeInitialized_
protected

Definition at line 345 of file QTest.h.

Referenced by NoisyChannel(), and setTolerance().

float NoisyChannel::tolerance_
protected

Definition at line 342 of file QTest.h.

Referenced by setTolerance().