CMS 3D CMS Logo

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

List of all members.

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)

Static Public Member Functions

static std::string getAlgoName (void)

Protected Attributes

bool rangeInitialized_
double ymin_
 ymin - threshold

Detailed Description

test that histogram contents are above Ymin

Definition at line 275 of file QTest.h.


Constructor & Destructor Documentation

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

Definition at line 278 of file QTest.h.

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


Member Function Documentation

static std::string DeadChannel::getAlgoName ( void  ) [inline, static]
float DeadChannel::runTest ( const MonitorElement me) [virtual]

loop over all channels

loop over all bins

Reimplemented from QCriterion.

Definition at line 626 of file QTest.cc.

References newFWLiteAna::bin, cmsMakeMELists::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(), MonitorElement::kind(), prof2calltree::last, and NULL.

{
  badChannels_.clear();
  if (!me) 
    return -1;
  if (!me->getRootObject()) 
    return -1;
  TH1* h1=0;
  TH2* h2=0;//initialize histogram pointers

  if (verbose_>1) 
    std::cout << "QTest:" << getAlgoName() << "::runTest called on " 
              << me-> getFullname() << "\n";
  //TH1F
  if (me->kind()==MonitorElement::DQM_KIND_TH1F) 
  { 
    h1 = me->getTH1F(); //access Test histo
  } 
  //TH1S
  else if (me->kind()==MonitorElement::DQM_KIND_TH1S) 
  { 
    h1 = me->getTH1S(); //access Test histo
  } 
  //TH1D
  else if (me->kind()==MonitorElement::DQM_KIND_TH1D) 
  { 
    h1 = me->getTH1D(); //access Test histo
  } 
  //-- TH2F
  else if (me->kind()==MonitorElement::DQM_KIND_TH2F)
  { 
    h2  = me->getTH2F(); // access Test histo
  } 
  //-- TH2S
  else if (me->kind()==MonitorElement::DQM_KIND_TH2S)
  { 
    h2  = me->getTH2S(); // access Test histo
  } 
  //-- TH2D
  else if (me->kind()==MonitorElement::DQM_KIND_TH2D)
  { 
    h2  = me->getTH2D(); // access Test histo
  } 
  else 
  {
    if (verbose_>0) 
      std::cout << "QTest:DeadChannel"
         << " ME " << me->getFullname() 
         << " does not contain TH1F/TH1S/TH1D/TH2F/TH2S/TH2D, exiting\n"; 
    return -1;
  } 

  int fail = 0; // number of failed channels

  //--------- do the quality test for 1D histo ---------------//
  if (h1 != NULL)  
  {
    if (!rangeInitialized_ || !h1->GetXaxis() ) 
      return 1; // all bins are accepted if no initialization
    int ncx = h1->GetXaxis()->GetNbins();
    int first = 1;
    int last  = ncx;
    int bin;

    for (bin = first; bin <= last; ++bin)
    {
      double contents = h1->GetBinContent(bin);
      bool failure = false;
      failure = contents <= ymin_; // dead channel: equal to or less than ymin_
      if (failure)
      { 
        DQMChannel chan(bin, 0, 0, contents, h1->GetBinError(bin));
        badChannels_.push_back(chan);
        ++fail;
      }
    }
    //return fraction of alive channels
    return 1.*(ncx - fail)/ncx;
  }
  //----------------------------------------------------------//
 
  //--------- do the quality test for 2D -------------------//
  else if (h2 !=NULL )
  {
    int ncx = h2->GetXaxis()->GetNbins(); // get X bins
    int ncy = h2->GetYaxis()->GetNbins(); // get Y bins

    for (int cx = 1; cx <= ncx; ++cx)
    {
      for (int cy = 1; cy <= ncy; ++cy)
      {
        double contents = h2->GetBinContent(h2->GetBin(cx, cy));
        bool failure = false;
        failure = contents <= ymin_; // dead channel: equal to or less than ymin_
        if (failure)
        { 
          DQMChannel chan(cx, cy, 0, contents, h2->GetBinError(h2->GetBin(cx, cy)));
          badChannels_.push_back(chan);
          ++fail;
        }
      }
    }
    //return fraction of alive channels
    return 1.*(ncx*ncy - fail) / (ncx*ncy);
  }
  else 
  {
    if (verbose_>0) 
      std::cout << "QTest:DeadChannel"
                << " TH1/TH2F are NULL, exiting\n"; 
    return -1;
  }
}
void DeadChannel::setThreshold ( double  ymin) [inline]

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

Definition at line 287 of file QTest.h.

References rangeInitialized_, and ymin_.

Referenced by QTestConfigure::EnableDeadChannelTest().

  { 
    ymin_ = ymin;  
    rangeInitialized_ = true; 
  } 

Member Data Documentation

Definition at line 295 of file QTest.h.

Referenced by DeadChannel(), and setThreshold().

double DeadChannel::ymin_ [protected]

ymin - threshold

Definition at line 294 of file QTest.h.

Referenced by setThreshold().