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)
 
void setThreshold (double ymin)
 set Ymin (inclusive) threshold for "dead" channel (default: 0) More...
 
- 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 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
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)
 

Detailed Description

test that histogram contents are above Ymin

Definition at line 306 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 309 of file QTest.h.

References QCriterion::setAlgoName().

309  : SimpleTest(name,true)
310  {
311  rangeInitialized_ = false;
313  }
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:78
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:138
static std::string getAlgoName(void)
Definition: QTest.h:314
bool rangeInitialized_
Definition: QTest.h:326

Member Function Documentation

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

Definition at line 314 of file QTest.h.

References QCriterion::runTest().

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

314 { return "DeadChannel"; }
float DeadChannel::runTest ( const MonitorElement me)
virtual

loop over all channels

loop over all bins

Reimplemented from QCriterion.

Definition at line 727 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(), plotBeamSpotDB::last, and NULL.

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

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

Definition at line 318 of file QTest.h.

References Phase2TrackerMonitorDigi_cff::ymin.

Referenced by QTestConfigure::EnableDeadChannelTest().

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

Member Data Documentation

bool DeadChannel::rangeInitialized_
protected

Definition at line 326 of file QTest.h.

double DeadChannel::ymin_
protected

ymin - threshold

Definition at line 325 of file QTest.h.