CMS 3D CMS Logo

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

#include <QTest.h>

Inheritance diagram for CompareLastFilledBin:
SimpleTest QCriterion

Public Member Functions

 CompareLastFilledBin (const std::string &name)
 
float runTest (const MonitorElement *me)
 
void setAverage (float average)
 
void setMax (float max)
 
void setMin (float min)
 
 ~CompareLastFilledBin ()
 
- 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

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)
 

Private Attributes

float _average
 
float _max
 
float _min
 

Additional Inherited Members

- 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

Definition at line 676 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 680 of file QTest.h.

References QCriterion::setAlgoName().

680  : SimpleTest(name,true){
681  this->_min = 0.0;
682  this->_max = 1.0;
683  this->_average = 0.0;
685  };
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:689
CompareLastFilledBin::~CompareLastFilledBin ( )
inline

Definition at line 687 of file QTest.h.

687 {};

Member Function Documentation

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

Definition at line 689 of file QTest.h.

References QCriterion::runTest().

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

689 { return "CompareLastFilledBin"; }
float CompareLastFilledBin::runTest ( const MonitorElement me)
virtual

Reimplemented from QCriterion.

Definition at line 1496 of file QTest.cc.

References gather_cfg::cout, MonitorElement::DQM_KIND_TH1F, MonitorElement::DQM_KIND_TH2F, MonitorElement::getRootObject(), MonitorElement::getTH1F(), MonitorElement::getTH2F(), MonitorElement::kind(), and NULL.

1496  {
1497  if (!me)
1498  return -1;
1499  if (!me->getRootObject())
1500  return -1;
1501  TH1* h1=0;
1502  TH2* h2=0;
1503  if (verbose_>1){
1504  std::cout << "QTest:" << getAlgoName() << "::runTest called on "
1505  << me-> getFullname() << "\n";
1506  std::cout << "\tMin = " << _min << "; Max = " << _max << "\n";
1507  }
1509  {
1510  h1 = me->getTH1F(); // access Test histo
1511  }
1512  else if (me->kind()==MonitorElement::DQM_KIND_TH2F)
1513  {
1514  h2 = me->getTH2F(); // access Test histo
1515  }
1516  else
1517  {
1518  if (verbose_>0)
1519  std::cout << "QTest:ContentsWithinExpected"
1520  << " ME does not contain TH1F or TH2F, exiting\n";
1521  return -1;
1522  }
1523  int lastBinX = 0;
1524  int lastBinY = 0;
1525  float lastBinVal;
1526 
1527  //--------- do the quality test for 1D histo ---------------//
1528  if (h1 != NULL)
1529  {
1530  lastBinX = h1->FindLastBinAbove(_average,1);
1531  lastBinVal = h1->GetBinContent(lastBinX);
1532  if (h1->GetEntries() == 0 || lastBinVal < 0) return 1;
1533  }
1534  else if (h2 != NULL)
1535  {
1536 
1537  lastBinX = h2->FindLastBinAbove(_average,1);
1538  lastBinY = h2->FindLastBinAbove(_average,2);
1539  if ( h2->GetEntries() == 0 || lastBinX < 0 || lastBinY < 0 ) return 1;
1540  lastBinVal = h2->GetBinContent(h2->GetBin(lastBinX,lastBinY));
1541  } else {
1542  if (verbose_ > 0) std::cout << "QTest:"<< getAlgoName() << " Histogram does not exist" << std::endl;
1543  return 1;
1544  }
1545  if (verbose_ > 0) std::cout << "Min and Max values " << _min << " " << _max << " Av value " << _average << " lastBinX " << lastBinX<< " lastBinY " << lastBinY << " lastBinVal " << lastBinVal << std::endl;
1546  if (lastBinVal > _min && lastBinVal <= _max)
1547  return 1;
1548  else
1549  return 0;
1550 }
int verbose_
Definition: QTest.h:118
#define NULL
Definition: scimark2.h:8
Kind kind(void) const
Get the type of the monitor element.
TObject * getRootObject(void) const
TH1F * getTH1F(void) const
TH2F * getTH2F(void) const
static std::string getAlgoName(void)
Definition: QTest.h:689
void CompareLastFilledBin::setAverage ( float  average)
inline

Definition at line 692 of file QTest.h.

Referenced by QTestConfigure::EnableCompareLastFilledBinTest().

692 {_average = average;};
void CompareLastFilledBin::setMax ( float  max)
inline
void CompareLastFilledBin::setMessage ( void  )
inlineprotectedvirtual

set status & message after test has run

Reimplemented from SimpleTest.

Definition at line 698 of file QTest.h.

References QCriterion::algoName_, python.rootplot.argparse::message, QCriterion::message_, QCriterion::prob_, and QCriterion::qtname_.

698  {
699  std::ostringstream message;
700  message << "Test " << qtname_ << " (" << algoName_
701  << "): Last Bin filled with desired value = " << prob_;
702  message_ = message.str();
703  }
std::string algoName_
name of quality test
Definition: QTest.h:112
std::string qtname_
Definition: QTest.h:111
float prob_
name of algorithm
Definition: QTest.h:113
std::string message_
quality test status
Definition: QTest.h:115
void CompareLastFilledBin::setMin ( float  min)
inline

Definition at line 693 of file QTest.h.

References min().

Referenced by QTestConfigure::EnableCompareLastFilledBinTest().

693 {_min = min;};
T min(T a, T b)
Definition: MathUtil.h:58

Member Data Documentation

float CompareLastFilledBin::_average
private

Definition at line 707 of file QTest.h.

float CompareLastFilledBin::_max
private

Definition at line 706 of file QTest.h.

Referenced by generateEDF.LumiInfoCont::max().

float CompareLastFilledBin::_min
private

Definition at line 706 of file QTest.h.

Referenced by generateEDF.LumiInfoCont::min().