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 308 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 311 of file QTest.h.

References QCriterion::setAlgoName().

311  : SimpleTest(name,true)
312  {
313  rangeInitialized_ = false;
315  }
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:80
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:140
static std::string getAlgoName()
Definition: QTest.h:316
bool rangeInitialized_
Definition: QTest.h:328

Member Function Documentation

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

Definition at line 316 of file QTest.h.

References QCriterion::runTest().

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

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

loop over all channels

loop over all bins

Reimplemented from QCriterion.

Definition at line 738 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.

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

References Phase2TrackerMonitorDigi_cff::ymin.

Referenced by QTestConfigure::EnableDeadChannelTest().

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

Member Data Documentation

bool DeadChannel::rangeInitialized_
protected

Definition at line 328 of file QTest.h.

double DeadChannel::ymin_
protected

ymin - threshold

Definition at line 327 of file QTest.h.