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
 (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 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...
 

Detailed Description

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

Definition at line 382 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 385 of file QTest.h.

References QCriterion::setAlgoName().

385  : SimpleTest(name,true)
386  {
387  rangeInitialized_ = false;
388  numXblocks_ = 1;
389  numYblocks_ = 1;
390  numNeighborsX_ = 1;
391  numNeighborsY_ = 1;
393  }
unsigned numYblocks_
Definition: QTest.h:458
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:80
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:140
unsigned numXblocks_
Definition: QTest.h:457
bool rangeInitialized_
Definition: QTest.h:463
unsigned numNeighborsY_
Definition: QTest.h:461
unsigned numNeighborsX_
Definition: QTest.h:459
static std::string getAlgoName()
Definition: QTest.h:394

Member Function Documentation

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

Definition at line 394 of file QTest.h.

References QCriterion::runTest().

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

394 { 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 1379 of file QTest.cc.

References anotherprimaryvertexanalyzer_cfi::xMax, anotherprimaryvertexanalyzer_cfi::xMin, CMSBoostedTauSeedingParameters_cfi::yMax, and CMSBoostedTauSeedingParameters_cfi::yMin.

1379  {
1380  double sigma = 0;
1381  unsigned nbinsX = h->GetXaxis()->GetNbins();
1382  unsigned nbinsY = h->GetYaxis()->GetNbins();
1383  unsigned xMin = 1;
1384  unsigned yMin = 1;
1385  unsigned xMax = nbinsX;
1386  unsigned yMax = nbinsY;
1387  unsigned Xbinnum = 1;
1388  unsigned Ybinnum = 1;
1389  double block_sum;
1390 
1391  if (xMin_ != 0 && xMax_ != 0) {
1392  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) {
1393  nbinsX = xMax_ - xMin_ + 1;
1394  xMax = xMax_;
1395  xMin = xMin_;
1396  }
1397  }
1398  if (yMin_ != 0 && yMax_ != 0) {
1399  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
1400  nbinsY = yMax_ - yMin_ + 1;
1401  yMax = yMax_;
1402  yMin = yMin_;
1403  }
1404  }
1405  if (Xblocks==999){
1406  Xblocks=nbinsX;
1407  }
1408  if (Yblocks==999){
1409  Yblocks=nbinsY;
1410  }
1411 
1412  Xbinnum = nbinsX/Xblocks;
1413  Ybinnum = nbinsY/Yblocks;
1414 
1415  unsigned xLow,xHi,yLow,yHi;
1416  for (unsigned x_block_count=0; x_block_count<=2*neighborsX; ++x_block_count){
1417  for (unsigned y_block_count=0; y_block_count<=2*neighborsY; ++y_block_count){
1418  //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.
1419  if (groupx + x_block_count > neighborsX){
1420  xLow=(groupx + x_block_count - neighborsX) * Xbinnum + xMin;
1421  if (xLow < xMin){
1422  xLow = xMin;
1423  }
1424  }else{
1425  xLow = xMin;
1426  }
1427  xHi = xLow + Xbinnum;
1428  if (xHi > xMax){
1429  xHi = xMax;
1430  }
1431  if (groupy + y_block_count > neighborsY){
1432  yLow=(groupy + y_block_count - neighborsY) * Ybinnum + yMin;
1433  if (yLow < yMin){
1434  yLow = yMin;
1435  }
1436  }else{
1437  yLow = yMin;
1438  }
1439  yHi = yLow + Ybinnum;
1440  if (yHi > yMax){
1441  yHi = yMax;
1442  }
1443  block_sum = 0;
1444  for (unsigned x_block_bin=xLow; x_block_bin<=xHi; ++x_block_bin){
1445  for (unsigned y_block_bin=yLow; y_block_bin<=yHi; ++y_block_bin){
1446  block_sum += h->GetBinContent(x_block_bin,y_block_bin);
1447  }
1448  }
1449  sigma += (average-block_sum)*(average-block_sum); //will sqrt and divide by sqrt(N-1) outside of function
1450  }
1451  }
1452  return sigma;
1453 }
unsigned yMax_
Definition: QTest.h:467
unsigned yMin_
Definition: QTest.h:466
unsigned xMin_
Definition: QTest.h:464
unsigned xMax_
Definition: QTest.h:465
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 1306 of file QTest.cc.

References anotherprimaryvertexanalyzer_cfi::xMax, anotherprimaryvertexanalyzer_cfi::xMin, CMSBoostedTauSeedingParameters_cfi::yMax, and CMSBoostedTauSeedingParameters_cfi::yMin.

1306  {
1307  double sum = 0;
1308  unsigned nbinsX = h->GetXaxis()->GetNbins();
1309  unsigned nbinsY = h->GetYaxis()->GetNbins();
1310  unsigned xMin = 1;
1311  unsigned yMin = 1;
1312  unsigned xMax = nbinsX;
1313  unsigned yMax = nbinsY;
1314  unsigned Xbinnum = 1;
1315  unsigned Ybinnum = 1;
1316 
1317  if (xMin_ != 0 && xMax_ != 0) { //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
1318  // check that user's parameters are completely in agreement with histogram
1319  // for instance, if inputted xMax is out of range xMin will automatically be ignored
1320  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) {
1321  nbinsX = xMax_ - xMin_ + 1;
1322  xMax = xMax_; // do NOT use overflow bin
1323  xMin = xMin_; // do NOT use underflow bin
1324  }
1325  }
1326  if (yMin_ != 0 && yMax_ != 0) {
1327  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
1328  nbinsY = yMax_ - yMin_ + 1;
1329  yMax = yMax_;
1330  yMin = yMin_;
1331  }
1332  }
1333 
1334  if (Xblocks==999){ //Check to see if blocks should be ignored
1335  Xblocks=nbinsX;
1336  }
1337  if (Yblocks==999){
1338  Yblocks=nbinsY;
1339  }
1340 
1341  Xbinnum = nbinsX/Xblocks;
1342  Ybinnum = nbinsY/Yblocks;
1343 
1344  unsigned xLow,xHi,yLow,yHi; //Define the neighbor blocks edges to be summed
1345  if (groupx>neighborsX){
1346  xLow=(groupx-neighborsX)*Xbinnum+xMin;
1347  if (xLow<xMin){
1348  xLow=xMin; //If the neigbor block would go outside the histogram edge, set it the edge
1349  }
1350  }else{
1351  xLow=xMin;
1352  }
1353  xHi=(groupx+1+neighborsX)*Xbinnum+xMin;
1354  if (xHi>xMax){
1355  xHi=xMax;
1356  }
1357  if (groupy>neighborsY){
1358  yLow=(groupy-neighborsY)*Ybinnum+yMin;
1359  if (yLow<yMin){
1360  yLow=yMin;
1361  }
1362  }else{
1363  yLow=yMin;
1364  }
1365  yHi=(groupy+1+neighborsY)*Ybinnum+yMin;
1366  if (yHi>yMax){
1367  yHi=yMax;
1368  }
1369 
1370  for (unsigned xbin=xLow; xbin<=xHi; ++xbin){ //now sum over all the bins
1371  for (unsigned ybin=yLow; ybin<=yHi; ++ybin){
1372  sum += h->GetBinContent(xbin,ybin);
1373  }
1374  }
1375  return sum;
1376 }
unsigned yMax_
Definition: QTest.h:467
unsigned yMin_
Definition: QTest.h:466
unsigned xMin_
Definition: QTest.h:464
unsigned xMax_
Definition: QTest.h:465
float ContentSigma::runTest ( const MonitorElement me)
overridevirtual

Reimplemented from QCriterion.

Definition at line 1093 of file QTest.cc.

References funct::abs(), 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(), MonitorElement::getFullname(), MonitorElement::getRootObject(), MonitorElement::getTH1D(), MonitorElement::getTH1F(), MonitorElement::getTH1S(), MonitorElement::getTH2D(), MonitorElement::getTH2F(), MonitorElement::getTH2S(), MonitorElement::kind(), mathSSE::sqrt(), anotherprimaryvertexanalyzer_cfi::xMax, anotherprimaryvertexanalyzer_cfi::xMin, CMSBoostedTauSeedingParameters_cfi::yMax, and CMSBoostedTauSeedingParameters_cfi::yMin.

1094 {
1095  badChannels_.clear();
1096  if (!me)
1097  return -1;
1098  if (!me->getRootObject())
1099  return -1;
1100  TH1* h=nullptr;//initialize histogram pointer
1101 
1102  if (verbose_>1)
1103  std::cout << "QTest:" << getAlgoName() << "::runTest called on "
1104  << me-> getFullname() << "\n";
1105 
1106  unsigned nbinsX;
1107  unsigned nbinsY;
1108 
1109  //-- TH1F
1111  {
1112  nbinsX = me->getTH1F()->GetXaxis()->GetNbins();
1113  nbinsY = me->getTH1F()->GetYaxis()->GetNbins();
1114  h = me->getTH1F(); // access Test histo
1115  }
1116  //-- TH1S
1117  else if (me->kind()==MonitorElement::DQM_KIND_TH1S)
1118  {
1119  nbinsX = me->getTH1S()->GetXaxis()->GetNbins();
1120  nbinsY = me->getTH1S()->GetYaxis()->GetNbins();
1121  h = me->getTH1S(); // access Test histo
1122  }
1123  //-- TH1D
1124  else if (me->kind()==MonitorElement::DQM_KIND_TH1D)
1125  {
1126  nbinsX = me->getTH1D()->GetXaxis()->GetNbins();
1127  nbinsY = me->getTH1D()->GetYaxis()->GetNbins();
1128  h = me->getTH1D(); // access Test histo
1129  }
1130  //-- TH2
1131  else if (me->kind()==MonitorElement::DQM_KIND_TH2F)
1132  {
1133  nbinsX = me->getTH2F()->GetXaxis()->GetNbins();
1134  nbinsY = me->getTH2F()->GetYaxis()->GetNbins();
1135  h = me->getTH2F(); // access Test histo
1136  }
1137  //-- TH2
1138  else if (me->kind()==MonitorElement::DQM_KIND_TH2S)
1139  {
1140  nbinsX = me->getTH2S()->GetXaxis()->GetNbins();
1141  nbinsY = me->getTH2S()->GetYaxis()->GetNbins();
1142  h = me->getTH2S(); // access Test histo
1143  }
1144  //-- TH2
1145  else if (me->kind()==MonitorElement::DQM_KIND_TH2D)
1146  {
1147  nbinsX = me->getTH2D()->GetXaxis()->GetNbins();
1148  nbinsY = me->getTH2D()->GetYaxis()->GetNbins();
1149  h = me->getTH2D(); // access Test histo
1150  }
1151  else
1152  {
1153  if (verbose_>0)
1154  std::cout << "QTest:ContentSigma"
1155  << " ME " << me->getFullname()
1156  << " does not contain TH1F/TH1S/TH1D or TH2F/TH2S/TH2D, exiting\n";
1157  return -1;
1158  }
1159 
1160  //-- QUALITY TEST itself
1161 
1162  if ( !rangeInitialized_ || !h->GetXaxis() )
1163  return 1; // all channels are accepted if tolerance has not been set
1164 
1165  int fail = 0; // initialize bin failure count
1166  unsigned xMin = 1; //initialize minimums and maximums with expected values
1167  unsigned yMin =1;
1168  unsigned xMax = nbinsX;
1169  unsigned yMax = nbinsY;
1170  unsigned XBlocks = numXblocks_; //Initialize xml inputs blocks and neighbors
1171  unsigned YBlocks = numYblocks_;
1172  unsigned neighborsX = numNeighborsX_;
1173  unsigned neighborsY = numNeighborsY_;
1174  unsigned Xbinnum = 1;
1175  unsigned Ybinnum =1;
1176  unsigned XWidth = 0;
1177  unsigned YWidth = 0;
1178 
1179  if (neighborsX==999){
1180  neighborsX=0;
1181  }
1182  if (neighborsY==999){
1183  neighborsY=0;
1184  }
1185 
1186  if (xMin_ != 0 && xMax_ != 0) { //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
1187  // check that user's parameters are completely in agreement with histogram
1188  // for instance, if inputted xMax is out of range xMin will automatically be ignored
1189  if ((xMax_ <= nbinsX) && (xMin_ <= xMax_)) { // rescale area of histogram being analyzed
1190  nbinsX = xMax_ - xMin_ + 1;
1191  xMax = xMax_; // do NOT use overflow bin
1192  xMin = xMin_; // do NOT use underflow bin
1193  }
1194  }
1195  if (yMin_ != 0 && yMax_ != 0) { //give users option for automatic mininum and maximum selection by inputting 0 to any of the parameters
1196  if ((yMax_ <= nbinsY) && (yMin_ <= yMax_)) {
1197  nbinsY = yMax_ - yMin_ + 1;
1198  yMax = yMax_;
1199  yMin = yMin_;
1200  }
1201  }
1202 
1203  if (neighborsX*2 >= nbinsX) { //make sure neighbor check does not overlap with bin under consideration
1204  if (nbinsX%2 == 0) {
1205  neighborsX = nbinsX/2 - 1; //set neighbors for no overlap
1206  }
1207  else { neighborsX = (nbinsX - 1)/2; }
1208  }
1209 
1210  if (neighborsY*2 >= nbinsY) {
1211  if (nbinsY%2 == 0) {
1212  neighborsY = nbinsY/2 - 1;
1213  }
1214  else { neighborsY = (nbinsY - 1)/2; }
1215  }
1216 
1217  if (XBlocks==999){ //Setting 999 prevents blocks and does quality tests by bins only
1218  XBlocks=nbinsX;
1219  }
1220  if (YBlocks==999){
1221  YBlocks=nbinsY;
1222  }
1223 
1224  Xbinnum = nbinsX/XBlocks;
1225  Ybinnum = nbinsY/YBlocks;
1226  for (unsigned groupx=0; groupx<XBlocks; ++groupx){ //Test over all the blocks
1227  for (unsigned groupy=0; groupy<YBlocks; ++groupy){
1228 
1229  double blocksum = 0;
1230  for (unsigned binx=0; binx<Xbinnum; ++binx){ //Sum the contents of the block in question
1231  for (unsigned biny=0; biny<Ybinnum; ++biny){
1232  if (groupx*Xbinnum+xMin+binx<=xMax && groupy*Ybinnum+yMin+biny<=yMax){
1233  blocksum += abs(h->GetBinContent(groupx*Xbinnum+xMin+binx,groupy*Ybinnum+yMin+biny));
1234  }
1235  }
1236  }
1237 
1238  double sum = getNeighborSum(groupx, groupy, XBlocks, YBlocks, neighborsX, neighborsY, h);
1239  sum -= blocksum; //remove center block to test
1240 
1241  if (neighborsX>groupx){ //Find correct average at the edges
1242  XWidth = neighborsX + groupx + 1;
1243  }else if (neighborsX>(XBlocks-(groupx+1))){
1244  XWidth = (XBlocks-groupx)+neighborsX;
1245  }else {
1246  XWidth = 2*neighborsX+1;
1247  }
1248  if (neighborsY>groupy){
1249  YWidth = neighborsY + groupy + 1;
1250  }else if (neighborsY>(YBlocks-(groupy+1))){
1251  YWidth = (YBlocks-groupy)+neighborsY;
1252  }else {
1253  YWidth = 2*neighborsY+1;
1254  }
1255 
1256  double average = sum/(XWidth*YWidth - 1);
1257  double sigma = getNeighborSigma(average, groupx, groupy, XBlocks, YBlocks, neighborsX, neighborsY, h);
1258  sigma -= (average-blocksum) * (average-blocksum); //get rid of block being tested just like we did with the average
1259  double sigma_2 = sqrt(sigma)/sqrt(XWidth*YWidth - 2); //N-1 where N=XWidth*YWidth - 1
1260  double sigma_real = sigma_2/(XWidth*YWidth - 1);
1261  //double avg_uncrt = average*sqrt(sum)/sum;//Obsolete now(Chad Freer)
1262  double avg_uncrt = 3*sigma_real;
1263 
1264  double probNoisy = ROOT::Math::poisson_cdf_c(blocksum - 1, average + avg_uncrt);
1265  double probDead = ROOT::Math::poisson_cdf(blocksum, average - avg_uncrt);
1266  double thresholdNoisy = ROOT::Math::normal_cdf_c(toleranceNoisy_);
1267  double thresholdDead = ROOT::Math::normal_cdf(-1 * toleranceDead_);
1268 
1269  int failureNoisy = 0;
1270  int failureDead = 0;
1271 
1272  if (average != 0){
1273  if (noisy_ && dead_) {
1274  if (blocksum > average) {
1275  failureNoisy = probNoisy < thresholdNoisy;
1276  }else {
1277  failureDead = probDead < thresholdDead;
1278  }
1279  }else if (noisy_) {
1280  if (blocksum > average) {
1281  failureNoisy = probNoisy < thresholdNoisy;
1282  }
1283  }else if (dead_) {
1284  if (blocksum < average) {
1285  failureDead = probDead < thresholdDead;
1286  }
1287  }else { std::cout<<"No test type selected!\n"; }
1288  //Following lines useful for debugging using verbose (Chad Freer)
1289  //string histName = h->GetName();
1290  //if (histName == "emtfTrackBX") {
1291  // std::printf("Chad says: %i XBlocks, %i XBlocks, %f Blocksum, %f Average", XBlocks,YBlocks,blocksum,average);}
1292  }
1293 
1294  if (failureNoisy || failureDead) {
1295  ++fail;
1296  //DQMChannel chan(groupx*Xbinnum+xMin+binx, 0, 0, blocksum, h->GetBinError(groupx*Xbinnum+xMin+binx));
1297  //badChannels_.push_back(chan);
1298  }
1299  }
1300  }
1301  return 1.*((XBlocks*YBlocks)-fail)/(XBlocks*YBlocks);
1302 }
unsigned numYblocks_
Definition: QTest.h:458
TH1F * getTH1F() const
TH2S * getTH2S() const
int verbose_
Definition: QTest.h:120
TObject * getRootObject() const
unsigned yMax_
Definition: QTest.h:467
T sqrt(T t)
Definition: SSEVec.h:18
TH2D * getTH2D() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned numXblocks_
Definition: QTest.h:457
float toleranceNoisy_
Definition: QTest.h:455
TH2F * getTH2F() const
bool rangeInitialized_
Definition: QTest.h:463
unsigned numNeighborsY_
Definition: QTest.h:461
const std::string getFullname() const
get full name of ME including Pathname
bool dead_
Definition: QTest.h:454
float toleranceDead_
Definition: QTest.h:456
unsigned numNeighborsX_
Definition: QTest.h:459
double getNeighborSigma(double average, unsigned groupx, unsigned groupy, unsigned Xblocks, unsigned Yblocks, unsigned neighborsX, unsigned neighborsY, const TH1 *h) const
Definition: QTest.cc:1379
unsigned yMin_
Definition: QTest.h:466
std::vector< DQMChannel > badChannels_
Definition: QTest.h:162
unsigned xMin_
Definition: QTest.h:464
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:1306
TH1D * getTH1D() const
TH1S * getTH1S() const
static std::string getAlgoName()
Definition: QTest.h:394
def fail(errstr="")
unsigned xMax_
Definition: QTest.h:465
bool noisy_
Definition: QTest.h:454
Kind kind() const
Get the type of the monitor element.
void ContentSigma::setDead ( bool  dead)
inline

Definition at line 431 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

431  {
432  dead_ = dead;
433  }
bool dead_
Definition: QTest.h:454
void ContentSigma::setNoisy ( bool  noisy)
inline

Definition at line 428 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

428  {
429  noisy_ = noisy;
430  }
bool noisy_
Definition: QTest.h:454
void ContentSigma::setNumNeighborsX ( unsigned  ncx)
inline

Definition at line 406 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

406 { if (ncx > 0) numNeighborsX_ = ncx; }
unsigned numNeighborsX_
Definition: QTest.h:459
void ContentSigma::setNumNeighborsY ( unsigned  ncy)
inline

Definition at line 407 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

407 { if (ncy > 0) numNeighborsY_ = ncy; }
unsigned numNeighborsY_
Definition: QTest.h:461
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 404 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

404 { if (ncx > 0) numXblocks_ = ncx; }
unsigned numXblocks_
Definition: QTest.h:457
void ContentSigma::setNumYblocks ( unsigned  ncy)
inline

Definition at line 405 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

405 { if (ncy > 0) numYblocks_ = ncy; }
unsigned numYblocks_
Definition: QTest.h:458
void ContentSigma::setToleranceDead ( float  factorDead)
inline

Definition at line 420 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

421  {
422  if (factorDead >=0)
423  {
424  toleranceDead_ = factorDead;
425  rangeInitialized_ = true;
426  }
427  }
bool rangeInitialized_
Definition: QTest.h:463
float toleranceDead_
Definition: QTest.h:456
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 412 of file QTest.h.

Referenced by QTestConfigure::EnableContentSigmaTest().

413  {
414  if (factorNoisy >=0)
415  {
416  toleranceNoisy_ = factorNoisy;
417  rangeInitialized_ = true;
418  }
419  }
float toleranceNoisy_
Definition: QTest.h:455
bool rangeInitialized_
Definition: QTest.h:463
void ContentSigma::setXMax ( unsigned  xMax)
inline
void ContentSigma::setXMin ( unsigned  xMin)
inline
void ContentSigma::setYMax ( unsigned  yMax)
inline
void ContentSigma::setYMin ( unsigned  yMin)
inline

Member Data Documentation

bool ContentSigma::dead_
protected

Definition at line 454 of file QTest.h.

bool ContentSigma::noisy_
protected

Definition at line 454 of file QTest.h.

unsigned ContentSigma::numNeighborsX_
protected

Definition at line 459 of file QTest.h.

unsigned ContentSigma::numNeighborsY_
protected

Definition at line 461 of file QTest.h.

unsigned ContentSigma::numXblocks_
protected

Definition at line 457 of file QTest.h.

unsigned ContentSigma::numYblocks_
protected

Definition at line 458 of file QTest.h.

bool ContentSigma::rangeInitialized_
protected

Definition at line 463 of file QTest.h.

float ContentSigma::toleranceDead_
protected

Definition at line 456 of file QTest.h.

float ContentSigma::toleranceNoisy_
protected

Definition at line 455 of file QTest.h.

unsigned ContentSigma::xMax_
protected

Definition at line 465 of file QTest.h.

unsigned ContentSigma::xMin_
protected

Definition at line 464 of file QTest.h.

unsigned ContentSigma::yMax_
protected

Definition at line 467 of file QTest.h.

unsigned ContentSigma::yMin_
protected

Definition at line 466 of file QTest.h.