CMS 3D CMS Logo

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

Check the sigma of each bin against the rest of the chamber by a factor of tolerance/. More...

#include <QTest.h>

Inheritance diagram for ContentSigma:
SimpleTest QCriterion

Public Member Functions

 ContentSigma (const std::string &name)
 
float runTest (const MonitorElement *me) override
 
void setDead (bool dead)
 
void setNoisy (bool noisy)
 
void setNumNeighborsX (unsigned ncx)
 
void setNumNeighborsY (unsigned ncy)
 
void setNumXblocks (unsigned ncx)
 
void setNumYblocks (unsigned ncy)
 
void setToleranceDead (float factorDead)
 
void setToleranceNoisy (float factorNoisy)
 
void setXMax (unsigned xMax)
 
void setXMin (unsigned xMin)
 
void setYMax (unsigned yMax)
 
void setYMin (unsigned yMin)
 
- 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
 get test status (see Core/interface/DQMDefinitions.h) 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 Member Functions

double getNeighborSigma (double average, unsigned groupx, unsigned groupy, unsigned Xblocks, unsigned Yblocks, unsigned neighborsX, unsigned neighborsY, const TH1 *h) const
 
double getNeighborSum (unsigned groupx, unsigned groupy, unsigned Xblocks, unsigned Yblocks, unsigned neighborsX, unsigned neighborsY, const TH1 *h) const
 for each bin get sum of the surrounding neighbors More...
 
- 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
 

Protected Attributes

bool dead_
 
bool noisy_
 
unsigned numNeighborsX_
 
unsigned numNeighborsY_
 
unsigned numXblocks_
 
unsigned numYblocks_
 
bool rangeInitialized_
 
float toleranceDead_
 
float toleranceNoisy_
 
unsigned xMax_
 
unsigned xMin_
 
unsigned yMax_
 
unsigned yMin_
 
- 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

- Public Types inherited from QCriterion
typedef dqm::legacy::MonitorElement MonitorElement
 (class should be created by DQMStore class) More...
 

Detailed Description

Check the sigma of each bin against the rest of the chamber by a factor of tolerance/.

Definition at line 301 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 303 of file QTest.h.

References QCriterion::setAlgoName().

303  : SimpleTest(name, true) {
304  rangeInitialized_ = false;
305  numXblocks_ = 1;
306  numYblocks_ = 1;
307  numNeighborsX_ = 1;
308  numNeighborsY_ = 1;
310  }
unsigned numYblocks_
Definition: QTest.h:385
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:95
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:158
unsigned numXblocks_
Definition: QTest.h:384
bool rangeInitialized_
Definition: QTest.h:390
unsigned numNeighborsY_
Definition: QTest.h:388
unsigned numNeighborsX_
Definition: QTest.h:386
static std::string getAlgoName()
Definition: QTest.h:311

Member Function Documentation

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

Definition at line 311 of file QTest.h.

References hlt_dqm_clientPB-live_cfg::me, and QCriterion::runTest().

Referenced by QTestConfigure::EnableContentSigmaTest(), and QTestConfigure::enableTests().

311 { return "ContentSigma"; }
double ContentSigma::getNeighborSigma ( double  average,
unsigned  groupx,
unsigned  groupy,
unsigned  Xblocks,
unsigned  Yblocks,
unsigned  neighborsX,
unsigned  neighborsY,
const TH1 *  h 
) const
protected

Definition at line 779 of file QTest.cc.

References hlt_dqm_clientPB-live_cfg::nbinsX, hlt_dqm_clientPB-live_cfg::nbinsY, multiplicitycorr_cfi::xMax, photonAnalyzer_cfi::xMin, multiplicitycorr_cfi::yMax, and photonAnalyzer_cfi::yMin.

786  {
787  double sigma = 0;
788  unsigned nbinsX = h->GetXaxis()->GetNbins();
789  unsigned nbinsY = h->GetYaxis()->GetNbins();
790  unsigned xMin = 1;
791  unsigned yMin = 1;
792  unsigned xMax = nbinsX;
793  unsigned yMax = nbinsY;
794  unsigned Xbinnum = 1;
795  unsigned Ybinnum = 1;
796  double block_sum;
797 
798  if (xMin_ != 0 && xMax_ != 0) {
799  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) {
800  nbinsX = xMax_ - xMin_ + 1;
801  xMax = xMax_;
802  xMin = xMin_;
803  }
804  }
805  if (yMin_ != 0 && yMax_ != 0) {
806  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
807  nbinsY = yMax_ - yMin_ + 1;
808  yMax = yMax_;
809  yMin = yMin_;
810  }
811  }
812  if (Xblocks == 999) {
813  Xblocks = nbinsX;
814  }
815  if (Yblocks == 999) {
816  Yblocks = nbinsY;
817  }
818 
819  Xbinnum = nbinsX / Xblocks;
820  Ybinnum = nbinsY / Yblocks;
821 
822  unsigned xLow, xHi, yLow, yHi;
823  for (unsigned x_block_count = 0; x_block_count <= 2 * neighborsX; ++x_block_count) {
824  for (unsigned y_block_count = 0; y_block_count <= 2 * neighborsY; ++y_block_count) {
825  //Sum over blocks. Need to find standard deviation of average of blocksums. Set up low and hi values similar to sum but for blocks now.
826  if (groupx + x_block_count > neighborsX) {
827  xLow = (groupx + x_block_count - neighborsX) * Xbinnum + xMin;
828  if (xLow < xMin) {
829  xLow = xMin;
830  }
831  } else {
832  xLow = xMin;
833  }
834  xHi = xLow + Xbinnum;
835  if (xHi > xMax) {
836  xHi = xMax;
837  }
838  if (groupy + y_block_count > neighborsY) {
839  yLow = (groupy + y_block_count - neighborsY) * Ybinnum + yMin;
840  if (yLow < yMin) {
841  yLow = yMin;
842  }
843  } else {
844  yLow = yMin;
845  }
846  yHi = yLow + Ybinnum;
847  if (yHi > yMax) {
848  yHi = yMax;
849  }
850  block_sum = 0;
851  for (unsigned x_block_bin = xLow; x_block_bin <= xHi; ++x_block_bin) {
852  for (unsigned y_block_bin = yLow; y_block_bin <= yHi; ++y_block_bin) {
853  block_sum += h->GetBinContent(x_block_bin, y_block_bin);
854  }
855  }
856  sigma += (average - block_sum) * (average - block_sum); //will sqrt and divide by sqrt(N-1) outside of function
857  }
858  }
859  return sigma;
860 }
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
unsigned yMax_
Definition: QTest.h:394
unsigned yMin_
Definition: QTest.h:393
unsigned xMin_
Definition: QTest.h:391
unsigned xMax_
Definition: QTest.h:392
double ContentSigma::getNeighborSum ( unsigned  groupx,
unsigned  groupy,
unsigned  Xblocks,
unsigned  Yblocks,
unsigned  neighborsX,
unsigned  neighborsY,
const TH1 *  h 
) const
protected

for each bin get sum of the surrounding neighbors

Definition at line 699 of file QTest.cc.

References hlt_dqm_clientPB-live_cfg::nbinsX, hlt_dqm_clientPB-live_cfg::nbinsY, multiplicitycorr_cfi::xMax, photonAnalyzer_cfi::xMin, multiplicitycorr_cfi::yMax, and photonAnalyzer_cfi::yMin.

705  {
706  double sum = 0;
707  unsigned nbinsX = h->GetXaxis()->GetNbins();
708  unsigned nbinsY = h->GetYaxis()->GetNbins();
709  unsigned xMin = 1;
710  unsigned yMin = 1;
711  unsigned xMax = nbinsX;
712  unsigned yMax = nbinsY;
713  unsigned Xbinnum = 1;
714  unsigned Ybinnum = 1;
715 
716  //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
717  // check that user's parameters are completely in agreement with histogram
718  // for instance, if inputted xMax is out of range xMin will automatically be ignored
719  if (xMin_ != 0 && xMax_ != 0) {
720  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) {
721  nbinsX = xMax_ - xMin_ + 1;
722  xMax = xMax_; // do NOT use overflow bin
723  xMin = xMin_; // do NOT use underflow bin
724  }
725  }
726  if (yMin_ != 0 && yMax_ != 0) {
727  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
728  nbinsY = yMax_ - yMin_ + 1;
729  yMax = yMax_;
730  yMin = yMin_;
731  }
732  }
733 
734  if (Xblocks == 999) { //Check to see if blocks should be ignored
735  Xblocks = nbinsX;
736  }
737  if (Yblocks == 999) {
738  Yblocks = nbinsY;
739  }
740 
741  Xbinnum = nbinsX / Xblocks;
742  Ybinnum = nbinsY / Yblocks;
743 
744  unsigned xLow, xHi, yLow, yHi; //Define the neighbor blocks edges to be summed
745  if (groupx > neighborsX) {
746  xLow = (groupx - neighborsX) * Xbinnum + xMin;
747  if (xLow < xMin) {
748  xLow = xMin; //If the neigbor block would go outside the histogram edge, set it the edge
749  }
750  } else {
751  xLow = xMin;
752  }
753  xHi = (groupx + 1 + neighborsX) * Xbinnum + xMin;
754  if (xHi > xMax) {
755  xHi = xMax;
756  }
757  if (groupy > neighborsY) {
758  yLow = (groupy - neighborsY) * Ybinnum + yMin;
759  if (yLow < yMin) {
760  yLow = yMin;
761  }
762  } else {
763  yLow = yMin;
764  }
765  yHi = (groupy + 1 + neighborsY) * Ybinnum + yMin;
766  if (yHi > yMax) {
767  yHi = yMax;
768  }
769 
770  for (unsigned xbin = xLow; xbin <= xHi; ++xbin) { //now sum over all the bins
771  for (unsigned ybin = yLow; ybin <= yHi; ++ybin) {
772  sum += h->GetBinContent(xbin, ybin);
773  }
774  }
775  return sum;
776 }
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
unsigned yMax_
Definition: QTest.h:394
unsigned yMin_
Definition: QTest.h:393
unsigned xMin_
Definition: QTest.h:391
unsigned xMax_
Definition: QTest.h:392
float ContentSigma::runTest ( const MonitorElement me)
overridevirtual

Reimplemented from QCriterion.

Definition at line 492 of file QTest.cc.

References funct::abs(), gather_cfg::cout, dqm::impl::MonitorElement::getFullname(), dqm::legacy::MonitorElement::getRootObject(), dqm::legacy::MonitorElement::getTH1D(), dqm::legacy::MonitorElement::getTH1F(), dqm::legacy::MonitorElement::getTH1S(), dqm::legacy::MonitorElement::getTH2D(), dqm::legacy::MonitorElement::getTH2F(), dqm::legacy::MonitorElement::getTH2S(), h, dqm::impl::MonitorElement::kind(), hlt_dqm_clientPB-live_cfg::nbinsX, hlt_dqm_clientPB-live_cfg::nbinsY, mathSSE::sqrt(), MonitorElementData::TH1D, MonitorElementData::TH1F, MonitorElementData::TH1S, MonitorElementData::TH2D, MonitorElementData::TH2F, MonitorElementData::TH2S, multiplicitycorr_cfi::xMax, photonAnalyzer_cfi::xMin, multiplicitycorr_cfi::yMax, and photonAnalyzer_cfi::yMin.

492  {
493  badChannels_.clear();
494  if (!me)
495  return -1;
496  if (!me->getRootObject())
497  return -1;
498  TH1* h = nullptr; //initialize histogram pointer
499 
500  if (verbose_ > 1)
501  std::cout << "QTest:" << getAlgoName() << "::runTest called on " << me->getFullname() << "\n";
502 
503  unsigned nbinsX;
504  unsigned nbinsY;
505 
506  //-- TH1F
507  if (me->kind() == MonitorElement::Kind::TH1F) {
508  nbinsX = me->getTH1F()->GetXaxis()->GetNbins();
509  nbinsY = me->getTH1F()->GetYaxis()->GetNbins();
510  h = me->getTH1F(); // access Test histo
511  }
512  //-- TH1S
513  else if (me->kind() == MonitorElement::Kind::TH1S) {
514  nbinsX = me->getTH1S()->GetXaxis()->GetNbins();
515  nbinsY = me->getTH1S()->GetYaxis()->GetNbins();
516  h = me->getTH1S(); // access Test histo
517  }
518  //-- TH1D
519  else if (me->kind() == MonitorElement::Kind::TH1D) {
520  nbinsX = me->getTH1D()->GetXaxis()->GetNbins();
521  nbinsY = me->getTH1D()->GetYaxis()->GetNbins();
522  h = me->getTH1D(); // access Test histo
523  }
524  //-- TH2
525  else if (me->kind() == MonitorElement::Kind::TH2F) {
526  nbinsX = me->getTH2F()->GetXaxis()->GetNbins();
527  nbinsY = me->getTH2F()->GetYaxis()->GetNbins();
528  h = me->getTH2F(); // access Test histo
529  }
530  //-- TH2
531  else if (me->kind() == MonitorElement::Kind::TH2S) {
532  nbinsX = me->getTH2S()->GetXaxis()->GetNbins();
533  nbinsY = me->getTH2S()->GetYaxis()->GetNbins();
534  h = me->getTH2S(); // access Test histo
535  }
536  //-- TH2
537  else if (me->kind() == MonitorElement::Kind::TH2D) {
538  nbinsX = me->getTH2D()->GetXaxis()->GetNbins();
539  nbinsY = me->getTH2D()->GetYaxis()->GetNbins();
540  h = me->getTH2D(); // access Test histo
541  } else {
542  if (verbose_ > 0)
543  std::cout << "QTest:ContentSigma"
544  << " ME " << me->getFullname() << " does not contain TH1F/TH1S/TH1D or TH2F/TH2S/TH2D, exiting\n";
545  return -1;
546  }
547 
548  //-- QUALITY TEST itself
549 
550  if (!rangeInitialized_ || !h->GetXaxis())
551  return 1; // all channels are accepted if tolerance has not been set
552 
553  int fail = 0; // initialize bin failure count
554  unsigned xMin = 1; //initialize minimums and maximums with expected values
555  unsigned yMin = 1;
556  unsigned xMax = nbinsX;
557  unsigned yMax = nbinsY;
558  unsigned XBlocks = numXblocks_; //Initialize xml inputs blocks and neighbors
559  unsigned YBlocks = numYblocks_;
560  unsigned neighborsX = numNeighborsX_;
561  unsigned neighborsY = numNeighborsY_;
562  unsigned Xbinnum = 1;
563  unsigned Ybinnum = 1;
564  unsigned XWidth = 0;
565  unsigned YWidth = 0;
566 
567  if (neighborsX == 999) {
568  neighborsX = 0;
569  }
570  if (neighborsY == 999) {
571  neighborsY = 0;
572  }
573 
574  //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
575  // check that user's parameters are completely in agreement with histogram
576  // for instance, if inputted xMax is out of range xMin will automatically be ignored
577  if (xMin_ != 0 && xMax_ != 0) {
578  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) { // rescale area of histogram being analyzed
579  nbinsX = xMax_ - xMin_ + 1;
580  xMax = xMax_; // do NOT use overflow bin
581  xMin = xMin_; // do NOT use underflow bin
582  }
583  }
584  //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
585  if (yMin_ != 0 && yMax_ != 0) {
586  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
587  nbinsY = yMax_ - yMin_ + 1;
588  yMax = yMax_;
589  yMin = yMin_;
590  }
591  }
592 
593  if (neighborsX * 2 >= nbinsX) { //make sure neighbor check does not overlap with bin under consideration
594  if (nbinsX % 2 == 0) {
595  neighborsX = nbinsX / 2 - 1; //set neighbors for no overlap
596  } else {
597  neighborsX = (nbinsX - 1) / 2;
598  }
599  }
600 
601  if (neighborsY * 2 >= nbinsY) {
602  if (nbinsY % 2 == 0) {
603  neighborsY = nbinsY / 2 - 1;
604  } else {
605  neighborsY = (nbinsY - 1) / 2;
606  }
607  }
608 
609  if (XBlocks == 999) { //Setting 999 prevents blocks and does quality tests by bins only
610  XBlocks = nbinsX;
611  }
612  if (YBlocks == 999) {
613  YBlocks = nbinsY;
614  }
615 
616  Xbinnum = nbinsX / XBlocks;
617  Ybinnum = nbinsY / YBlocks;
618  for (unsigned groupx = 0; groupx < XBlocks; ++groupx) { //Test over all the blocks
619  for (unsigned groupy = 0; groupy < YBlocks; ++groupy) {
620  double blocksum = 0;
621  for (unsigned binx = 0; binx < Xbinnum; ++binx) { //Sum the contents of the block in question
622  for (unsigned biny = 0; biny < Ybinnum; ++biny) {
623  if (groupx * Xbinnum + xMin + binx <= xMax && groupy * Ybinnum + yMin + biny <= yMax) {
624  blocksum += abs(h->GetBinContent(groupx * Xbinnum + xMin + binx, groupy * Ybinnum + yMin + biny));
625  }
626  }
627  }
628 
629  double sum = getNeighborSum(groupx, groupy, XBlocks, YBlocks, neighborsX, neighborsY, h);
630  sum -= blocksum; //remove center block to test
631 
632  if (neighborsX > groupx) { //Find correct average at the edges
633  XWidth = neighborsX + groupx + 1;
634  } else if (neighborsX > (XBlocks - (groupx + 1))) {
635  XWidth = (XBlocks - groupx) + neighborsX;
636  } else {
637  XWidth = 2 * neighborsX + 1;
638  }
639  if (neighborsY > groupy) {
640  YWidth = neighborsY + groupy + 1;
641  } else if (neighborsY > (YBlocks - (groupy + 1))) {
642  YWidth = (YBlocks - groupy) + neighborsY;
643  } else {
644  YWidth = 2 * neighborsY + 1;
645  }
646 
647  double average = sum / (XWidth * YWidth - 1);
648  double sigma = getNeighborSigma(average, groupx, groupy, XBlocks, YBlocks, neighborsX, neighborsY, h);
649  //get rid of block being tested just like we did with the average
650  sigma -= (average - blocksum) * (average - blocksum);
651  double sigma_2 = sqrt(sigma) / sqrt(XWidth * YWidth - 2); //N-1 where N=XWidth*YWidth - 1
652  double sigma_real = sigma_2 / (XWidth * YWidth - 1);
653  //double avg_uncrt = average*sqrt(sum)/sum;//Obsolete now(Chad Freer)
654  double avg_uncrt = 3 * sigma_real;
655 
656  double probNoisy = ROOT::Math::poisson_cdf_c(blocksum - 1, average + avg_uncrt);
657  double probDead = ROOT::Math::poisson_cdf(blocksum, average - avg_uncrt);
658  double thresholdNoisy = ROOT::Math::normal_cdf_c(toleranceNoisy_);
659  double thresholdDead = ROOT::Math::normal_cdf(-1 * toleranceDead_);
660 
661  int failureNoisy = 0;
662  int failureDead = 0;
663 
664  if (average != 0) {
665  if (noisy_ && dead_) {
666  if (blocksum > average) {
667  failureNoisy = probNoisy < thresholdNoisy;
668  } else {
669  failureDead = probDead < thresholdDead;
670  }
671  } else if (noisy_) {
672  if (blocksum > average) {
673  failureNoisy = probNoisy < thresholdNoisy;
674  }
675  } else if (dead_) {
676  if (blocksum < average) {
677  failureDead = probDead < thresholdDead;
678  }
679  } else {
680  std::cout << "No test type selected!\n";
681  }
682  //Following lines useful for debugging using verbose (Chad Freer)
683  //string histName = h->GetName();
684  //if (histName == "emtfTrackBX") {
685  // std::printf("Chad says: %i XBlocks, %i XBlocks, %f Blocksum, %f Average", XBlocks,YBlocks,blocksum,average);}
686  }
687 
688  if (failureNoisy || failureDead) {
689  ++fail;
690  //DQMChannel chan(groupx*Xbinnum+xMin+binx, 0, 0, blocksum, h->GetBinError(groupx*Xbinnum+xMin+binx));
691  //badChannels_.push_back(chan);
692  }
693  }
694  }
695  return 1. * ((XBlocks * YBlocks) - fail) / (XBlocks * YBlocks);
696 }
virtual TH2D * getTH2D() const
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
unsigned numYblocks_
Definition: QTest.h:385
virtual TH2F * getTH2F() const
int verbose_
Definition: QTest.h:138
virtual TH1F * getTH1F() const
unsigned yMax_
Definition: QTest.h:394
Kind kind() const
Get the type of the monitor element.
virtual TH1S * getTH1S() const
const std::string getFullname() const
get full name of ME including Pathname
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned numXblocks_
Definition: QTest.h:384
float toleranceNoisy_
Definition: QTest.h:382
virtual TH2S * getTH2S() const
bool rangeInitialized_
Definition: QTest.h:390
unsigned numNeighborsY_
Definition: QTest.h:388
bool dead_
Definition: QTest.h:381
float toleranceDead_
Definition: QTest.h:383
unsigned numNeighborsX_
Definition: QTest.h:386
double getNeighborSigma(double average, unsigned groupx, unsigned groupy, unsigned Xblocks, unsigned Yblocks, unsigned neighborsX, unsigned neighborsY, const TH1 *h) const
Definition: QTest.cc:779
unsigned yMin_
Definition: QTest.h:393
virtual TH1D * getTH1D() const
std::vector< DQMChannel > badChannels_
Definition: QTest.h:173
unsigned xMin_
Definition: QTest.h:391
double getNeighborSum(unsigned groupx, unsigned groupy, unsigned Xblocks, unsigned Yblocks, unsigned neighborsX, unsigned neighborsY, const TH1 *h) const
for each bin get sum of the surrounding neighbors
Definition: QTest.cc:699
TObject * getRootObject() const override
static std::string getAlgoName()
Definition: QTest.h:311
unsigned xMax_
Definition: QTest.h:392
bool noisy_
Definition: QTest.h:380
void ContentSigma::setDead ( bool  dead)
inline

Definition at line 354 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

354 { dead_ = dead; }
bool dead_
Definition: QTest.h:381
void ContentSigma::setNoisy ( bool  noisy)
inline

Definition at line 353 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

353 { noisy_ = noisy; }
bool noisy_
Definition: QTest.h:380
void ContentSigma::setNumNeighborsX ( unsigned  ncx)
inline

Definition at line 329 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

329  {
330  if (ncx > 0)
331  numNeighborsX_ = ncx;
332  }
unsigned numNeighborsX_
Definition: QTest.h:386
void ContentSigma::setNumNeighborsY ( unsigned  ncy)
inline

Definition at line 333 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

333  {
334  if (ncy > 0)
335  numNeighborsY_ = ncy;
336  }
unsigned numNeighborsY_
Definition: QTest.h:388
void ContentSigma::setNumXblocks ( unsigned  ncx)
inline

set # of neighboring channels for calculating average to be used for comparison with channel under consideration; use 1 for considering bin+1 and bin-1 (default), use 2 for considering bin+1,bin-1, bin+2,bin-2, etc; Will use rollover when bin+i or bin-i is beyond histogram limits (e.g. for histogram with N bins, bin N+1 corresponds to bin 1, and bin -1 corresponds to bin N)

Definition at line 321 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

321  {
322  if (ncx > 0)
323  numXblocks_ = ncx;
324  }
unsigned numXblocks_
Definition: QTest.h:384
void ContentSigma::setNumYblocks ( unsigned  ncy)
inline

Definition at line 325 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

325  {
326  if (ncy > 0)
327  numYblocks_ = ncy;
328  }
unsigned numYblocks_
Definition: QTest.h:385
void ContentSigma::setToleranceDead ( float  factorDead)
inline

Definition at line 347 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

347  {
348  if (factorDead >= 0) {
349  toleranceDead_ = factorDead;
350  rangeInitialized_ = true;
351  }
352  }
bool rangeInitialized_
Definition: QTest.h:390
float toleranceDead_
Definition: QTest.h:383
void ContentSigma::setToleranceNoisy ( float  factorNoisy)
inline

set factor tolerance for considering a channel noisy or dead; eg. if tolerance = 1, channel will be noisy if (content - 1 x sigma) > chamber_avg or channel will be dead if (content - 1 x sigma) < chamber_avg

Definition at line 341 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

341  {
342  if (factorNoisy >= 0) {
343  toleranceNoisy_ = factorNoisy;
344  rangeInitialized_ = true;
345  }
346  }
float toleranceNoisy_
Definition: QTest.h:382
bool rangeInitialized_
Definition: QTest.h:390
void ContentSigma::setXMax ( unsigned  xMax)
inline

Definition at line 357 of file QTest.h.

References multiplicitycorr_cfi::xMax.

Referenced by QTestConfigure::EnableContentSigmaTest().

void ContentSigma::setXMin ( unsigned  xMin)
inline

Definition at line 356 of file QTest.h.

References photonAnalyzer_cfi::xMin.

Referenced by QTestConfigure::EnableContentSigmaTest().

356 { xMin_ = xMin; }
unsigned xMin_
Definition: QTest.h:391
void ContentSigma::setYMax ( unsigned  yMax)
inline

Definition at line 359 of file QTest.h.

References h, and multiplicitycorr_cfi::yMax.

Referenced by QTestConfigure::EnableContentSigmaTest().

void ContentSigma::setYMin ( unsigned  yMin)
inline

Definition at line 358 of file QTest.h.

References photonAnalyzer_cfi::yMin.

Referenced by QTestConfigure::EnableContentSigmaTest().

358 { yMin_ = yMin; }
unsigned yMin_
Definition: QTest.h:393

Member Data Documentation

bool ContentSigma::dead_
protected

Definition at line 381 of file QTest.h.

bool ContentSigma::noisy_
protected

Definition at line 380 of file QTest.h.

unsigned ContentSigma::numNeighborsX_
protected

Definition at line 386 of file QTest.h.

unsigned ContentSigma::numNeighborsY_
protected

Definition at line 388 of file QTest.h.

unsigned ContentSigma::numXblocks_
protected

Definition at line 384 of file QTest.h.

unsigned ContentSigma::numYblocks_
protected

Definition at line 385 of file QTest.h.

bool ContentSigma::rangeInitialized_
protected

Definition at line 390 of file QTest.h.

float ContentSigma::toleranceDead_
protected

Definition at line 383 of file QTest.h.

float ContentSigma::toleranceNoisy_
protected

Definition at line 382 of file QTest.h.

unsigned ContentSigma::xMax_
protected

Definition at line 392 of file QTest.h.

unsigned ContentSigma::xMin_
protected

Definition at line 391 of file QTest.h.

unsigned ContentSigma::yMax_
protected

Definition at line 394 of file QTest.h.

unsigned ContentSigma::yMin_
protected

Definition at line 393 of file QTest.h.