CMS 3D CMS Logo

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

test that histogram contents are above Ymin More...

#include <QTest.h>

Inheritance diagram for DeadChannel:
SimpleTest QCriterion

Public Member Functions

 DeadChannel (const std::string &name)
 
float runTest (const MonitorElement *me) override
 
void setThreshold (double ymin)
 set Ymin (inclusive) threshold for "dead" channel (default: 0) More...
 
- 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
 (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 ()
 

Protected Attributes

bool rangeInitialized_
 
double ymin_
 ymin - threshold More...
 
- 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

- Protected Member Functions inherited from SimpleTest
void setMessage () override
 set status & message after test has run More...
 
- Protected Member Functions inherited from QCriterion
void init ()
 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 ()=default
 

Detailed Description

test that histogram contents are above Ymin

Definition at line 307 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 310 of file QTest.h.

References QCriterion::setAlgoName().

310  : SimpleTest(name,true)
311  {
312  rangeInitialized_ = false;
314  }
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:79
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:139
static std::string getAlgoName()
Definition: QTest.h:315
bool rangeInitialized_
Definition: QTest.h:327

Member Function Documentation

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

Definition at line 315 of file QTest.h.

References QCriterion::runTest().

Referenced by QTestConfigure::EnableDeadChannelTest(), QTestConfigure::enableTests(), and QTestParameterNames::QTestParameterNames().

315 { return "DeadChannel"; }
float DeadChannel::runTest ( const MonitorElement me)
overridevirtual

loop over all channels

loop over all bins

Reimplemented from QCriterion.

Definition at line 737 of file QTest.cc.

References stringResolutionProvider_cfi::bin, officialStyle::chan, 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(), plotBeamSpotDB::first, MonitorElement::getFullname(), MonitorElement::getRootObject(), MonitorElement::getTH1D(), MonitorElement::getTH1F(), MonitorElement::getTH1S(), MonitorElement::getTH2D(), MonitorElement::getTH2F(), MonitorElement::getTH2S(), MonitorElement::kind(), and plotBeamSpotDB::last.

738 {
739  badChannels_.clear();
740  if (!me)
741  return -1;
742  if (!me->getRootObject())
743  return -1;
744  TH1* h1=nullptr;
745  TH2* h2=nullptr;//initialize histogram pointers
746 
747  if (verbose_>1)
748  std::cout << "QTest:" << getAlgoName() << "::runTest called on "
749  << me-> getFullname() << "\n";
750  //TH1F
752  {
753  h1 = me->getTH1F(); //access Test histo
754  }
755  //TH1S
756  else if (me->kind()==MonitorElement::DQM_KIND_TH1S)
757  {
758  h1 = me->getTH1S(); //access Test histo
759  }
760  //TH1D
761  else if (me->kind()==MonitorElement::DQM_KIND_TH1D)
762  {
763  h1 = me->getTH1D(); //access Test histo
764  }
765  //-- TH2F
766  else if (me->kind()==MonitorElement::DQM_KIND_TH2F)
767  {
768  h2 = me->getTH2F(); // access Test histo
769  }
770  //-- TH2S
771  else if (me->kind()==MonitorElement::DQM_KIND_TH2S)
772  {
773  h2 = me->getTH2S(); // access Test histo
774  }
775  //-- TH2D
776  else if (me->kind()==MonitorElement::DQM_KIND_TH2D)
777  {
778  h2 = me->getTH2D(); // access Test histo
779  }
780  else
781  {
782  if (verbose_>0)
783  std::cout << "QTest:DeadChannel"
784  << " ME " << me->getFullname()
785  << " does not contain TH1F/TH1S/TH1D/TH2F/TH2S/TH2D, exiting\n";
786  return -1;
787  }
788 
789  int fail = 0; // number of failed channels
790 
791  //--------- do the quality test for 1D histo ---------------//
792  if (h1 != nullptr)
793  {
794  if (!rangeInitialized_ || !h1->GetXaxis() )
795  return 1; // all bins are accepted if no initialization
796  int ncx = h1->GetXaxis()->GetNbins();
797  int first = 1;
798  int last = ncx;
799  int bin;
800 
802  for (bin = first; bin <= last; ++bin)
803  {
804  double contents = h1->GetBinContent(bin);
805  bool failure = false;
806  failure = contents <= ymin_; // dead channel: equal to or less than ymin_
807  if (failure)
808  {
809  DQMChannel chan(bin, 0, 0, contents, h1->GetBinError(bin));
810  badChannels_.push_back(chan);
811  ++fail;
812  }
813  }
814  //return fraction of alive channels
815  return 1.*(ncx - fail)/ncx;
816  }
817  //----------------------------------------------------------//
818 
819  //--------- do the quality test for 2D -------------------//
820  else if (h2 !=nullptr )
821  {
822  int ncx = h2->GetXaxis()->GetNbins(); // get X bins
823  int ncy = h2->GetYaxis()->GetNbins(); // get Y bins
824 
826  for (int cx = 1; cx <= ncx; ++cx)
827  {
828  for (int cy = 1; cy <= ncy; ++cy)
829  {
830  double contents = h2->GetBinContent(h2->GetBin(cx, cy));
831  bool failure = false;
832  failure = contents <= ymin_; // dead channel: equal to or less than ymin_
833  if (failure)
834  {
835  DQMChannel chan(cx, cy, 0, contents, h2->GetBinError(h2->GetBin(cx, cy)));
836  badChannels_.push_back(chan);
837  ++fail;
838  }
839  }
840  }
841  //return fraction of alive channels
842  return 1.*(ncx*ncy - fail) / (ncx*ncy);
843  }
844  else
845  {
846  if (verbose_>0)
847  std::cout << "QTest:DeadChannel"
848  << " TH1/TH2F are NULL, exiting\n";
849  return -1;
850  }
851 }
double ymin_
ymin - threshold
Definition: QTest.h:326
TH1F * getTH1F() const
TH2S * getTH2S() const
int verbose_
Definition: QTest.h:119
TObject * getRootObject() const
TH2D * getTH2D() const
TH2F * getTH2F() const
bin
set the eta bin as selection string.
const std::string getFullname() const
get full name of ME including Pathname
std::vector< DQMChannel > badChannels_
Definition: QTest.h:161
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
static std::string getAlgoName()
Definition: QTest.h:315
TH1D * getTH1D() const
TH1S * getTH1S() const
def fail(errstr="")
bool rangeInitialized_
Definition: QTest.h:327
Kind kind() const
Get the type of the monitor element.
void DeadChannel::setThreshold ( double  ymin)
inline

set Ymin (inclusive) threshold for "dead" channel (default: 0)

Definition at line 319 of file QTest.h.

References Phase2TrackerMonitorDigi_cff::ymin.

Referenced by QTestConfigure::EnableDeadChannelTest().

320  {
321  ymin_ = ymin;
322  rangeInitialized_ = true;
323  }
double ymin_
ymin - threshold
Definition: QTest.h:326
bool rangeInitialized_
Definition: QTest.h:327

Member Data Documentation

bool DeadChannel::rangeInitialized_
protected

Definition at line 327 of file QTest.h.

double DeadChannel::ymin_
protected

ymin - threshold

Definition at line 326 of file QTest.h.