CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes
Comp2RefChi2 Class Reference

#include <QTest.h>

Inheritance diagram for Comp2RefChi2:
SimpleTest QCriterion

Public Member Functions

 Comp2RefChi2 (const std::string &name)
 
float runTest (const MonitorElement *me)
 
- 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)
 

Protected Attributes

double chi2_
 
int Ndof_
 
- 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 183 of file QTest.h.

Constructor & Destructor Documentation

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

Definition at line 186 of file QTest.h.

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

186  :SimpleTest(name)
187  {
189  }
void setAlgoName(std::string name)
set algorithm name
Definition: QTest.h:77
SimpleTest(const std::string &name, bool keepBadChannels=false)
Definition: QTest.h:137
static std::string getAlgoName(void)
Definition: QTest.h:190

Member Function Documentation

static std::string Comp2RefChi2::getAlgoName ( void  )
inlinestatic
float Comp2RefChi2::runTest ( const MonitorElement me)
virtual

Reimplemented from QCriterion.

Definition at line 163 of file QTest.cc.

References HLT_25ns14e33_v1_cff::constraint, gather_cfg::cout, MonitorElement::DQM_KIND_TH1D, MonitorElement::DQM_KIND_TH1F, MonitorElement::DQM_KIND_TH1S, MonitorElement::DQM_KIND_TPROFILE, MonitorElement::getRefRootObject(), MonitorElement::getRefTH1D(), MonitorElement::getRefTH1F(), MonitorElement::getRefTH1S(), MonitorElement::getRefTProfile(), MonitorElement::getRootObject(), MonitorElement::getTH1D(), MonitorElement::getTH1F(), MonitorElement::getTH1S(), MonitorElement::getTProfile(), h, i, MonitorElement::kind(), ndof, and groupFilesInBlocks::temp.

164 {
165  if (!me)
166  return -1;
167  if (!me->getRootObject() || !me->getRefRootObject())
168  return -1;
169  TH1* h=0;
170  TH1* ref_=0;
171 
172  if (verbose_>1)
173  std::cout << "QTest:" << getAlgoName() << "::runTest called on "
174  << me-> getFullname() << "\n";
175  //-- TH1F
177  {
178  h = me->getTH1F(); // access Test histo
179  ref_ = me->getRefTH1F(); //access Ref histo
180  }
181  //-- TH1S
182  else if (me->kind()==MonitorElement::DQM_KIND_TH1S)
183  {
184  h = me->getTH1S(); // access Test histo
185  ref_ = me->getRefTH1S(); //access Ref histo
186  }
187  //-- TH1D
188  else if (me->kind()==MonitorElement::DQM_KIND_TH1D)
189  {
190  h = me->getTH1D(); // access Test histo
191  ref_ = me->getRefTH1D(); //access Ref histo
192  }
193  //-- TProfile
194  else if (me->kind()==MonitorElement::DQM_KIND_TPROFILE)
195  {
196  h = me->getTProfile(); // access Test histo
197  ref_ = me->getRefTProfile(); //access Ref histo
198  }
199  else
200  {
201  if (verbose_>0)
202  std::cout << "QTest::Comp2RefChi2"
203  << " ME does not contain TH1F/TH1S/TH1D/TProfile, exiting\n";
204  return -1;
205  }
206 
207  //-- isInvalid ? - Check consistency in number of channels
208  int ncx1 = h->GetXaxis()->GetNbins();
209  int ncx2 = ref_->GetXaxis()->GetNbins();
210  if ( ncx1 != ncx2)
211  {
212  if (verbose_>0)
213  std::cout << "QTest:Comp2RefChi2"
214  << " different number of channels! ("
215  << ncx1 << ", " << ncx2 << "), exiting\n";
216  return -1;
217  }
218 
219  //-- QUALITY TEST itself
220  //reset Results
221  Ndof_ = 0; chi2_ = -1.; ncx1 = ncx2 = -1;
222 
223  int i, i_start, i_end;
224  double chi2 = 0.; int ndof = 0; int constraint = 0;
225 
226  i_start = 1;
227  i_end = ncx1;
228  // if (fXaxis.TestBit(TAxis::kAxisRange)) {
229  i_start = h->GetXaxis()->GetFirst();
230  i_end = h->GetXaxis()->GetLast();
231  // }
232  ndof = i_end-i_start+1-constraint;
233 
234  //Compute the normalisation factor
235  double sum1=0, sum2=0;
236  for (i=i_start; i<=i_end; i++)
237  {
238  sum1 += h->GetBinContent(i);
239  sum2 += ref_->GetBinContent(i);
240  }
241 
242  //check that the histograms are not empty
243  if (sum1 == 0)
244  {
245  if (verbose_>0)
246  std::cout << "QTest:Comp2RefChi2"
247  << " Test Histogram " << h->GetName()
248  << " is empty, exiting\n";
249  return -1;
250  }
251  if (sum2 == 0)
252  {
253  if (verbose_>0)
254  std::cout << "QTest:Comp2RefChi2"
255  << " Ref Histogram " << ref_->GetName()
256  << " is empty, exiting\n";
257  return -1;
258  }
259 
260  double bin1, bin2, err1, err2, temp;
261  for (i=i_start; i<=i_end; i++)
262  {
263  bin1 = h->GetBinContent(i)/sum1;
264  bin2 = ref_->GetBinContent(i)/sum2;
265  if (bin1 ==0 && bin2==0)
266  {
267  --ndof; //no data means one less degree of freedom
268  }
269  else
270  {
271  temp = bin1-bin2;
272  err1=h->GetBinError(i); err2=ref_->GetBinError(i);
273  if (err1 == 0 && err2 == 0)
274  {
275  if (verbose_>0)
276  std::cout << "QTest:Comp2RefChi2"
277  << " bins with non-zero content and zero error, exiting\n";
278  return -1;
279  }
280  err1*=err1 ; err2*=err2;
281  err1/=sum1*sum1 ; err2/=sum2*sum2;
282  chi2 +=temp*temp/(err1+err2);
283  }
284  }
285  chi2_ = chi2; Ndof_ = ndof;
286  return TMath::Prob(0.5*chi2, int(0.5*ndof));
287 }
TH1F * getRefTH1F(void) const
TH1S * getTH1S(void) const
int i
Definition: DBlmapReader.cc:9
int verbose_
Definition: QTest.h:117
TH1D * getTH1D(void) const
int Ndof_
Definition: QTest.h:206
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Kind kind(void) const
Get the type of the monitor element.
TObject * getRootObject(void) const
TH1F * getTH1F(void) const
TProfile * getRefTProfile(void) const
static std::string getAlgoName(void)
Definition: QTest.h:190
TProfile * getTProfile(void) const
tuple cout
Definition: gather_cfg.py:121
double chi2_
Definition: QTest.h:206
TObject * getRefRootObject(void) const
TH1S * getRefTH1S(void) const
TH1D * getRefTH1D(void) const
void Comp2RefChi2::setMessage ( void  )
inlineprotectedvirtual

set status & message after test has run

Reimplemented from SimpleTest.

Definition at line 195 of file QTest.h.

References chi2_, QCriterion::errorProb_, python.rootplot.argparse::message, QCriterion::message_, SimpleTest::minEntries_, Ndof_, and QCriterion::warningProb_.

196  {
197  std::ostringstream message;
198  message << "chi2/Ndof = " << chi2_ << "/" << Ndof_
199  << ", minimum needed statistics = " << minEntries_
200  << " warning threshold = " << this->warningProb_
201  << " error threshold = " << this->errorProb_;
202  message_ = message.str();
203  }
int Ndof_
Definition: QTest.h:206
float errorProb_
Definition: QTest.h:115
std::string message_
quality test status
Definition: QTest.h:114
unsigned minEntries_
Definition: QTest.h:158
double chi2_
Definition: QTest.h:206
float warningProb_
message attached to test
Definition: QTest.h:115

Member Data Documentation

double Comp2RefChi2::chi2_
protected

Definition at line 206 of file QTest.h.

Referenced by setMessage().

int Comp2RefChi2::Ndof_
protected

Definition at line 206 of file QTest.h.

Referenced by setMessage().