English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Uncertainty Modeler: URANIE::UncertModeler::TTester Class Reference
Uranie / Uncertainty Modeler v4.9.0
/* @license-end */
URANIE::UncertModeler::TTester Class Referenceabstract

Description of the class TTester. More...

#include <TTester.h>

Inheritance diagram for URANIE::UncertModeler::TTester:
Collaboration diagram for URANIE::UncertModeler::TTester:

Public Types

enum  ETestType { kKolmogorovSmirnov , kAndersonDarling , kCramerVonMises , kUnknown }
 
enum  ECase {
  kCase0 , kCase1 , kCase2 , kCase3 ,
  kCaseUnknown
}
 

Public Member Functions

Constructor and Destructor
 TTester (URANIE::DataServer::TDataServer *tds, ETestType nTestType, const char *varexpinput, Option_t *option="")
 Constructor by the type of Test.
 
virtual ~TTester ()
 Default constructor.
 
Setting and Getting attributes
Int_t getID ()
 
Double_t getScore ()
 Get the score of the last test.
 
Double_t getModifiedScore ()
 Get the modified score of the last test.
 
Double_t getPValue ()
 Get the pValue of the last test.
 
ECase getCase ()
 getCaseOfGoF
 
Computes the score
Double_t computeScore (Option_t *option="")
 Compute the score of the statistic.
 
virtual Double_t computeTheCurrentScore (Int_t i, Double_t dzi)=0
 The true computation with a vector of information and return the current value.
 
virtual void computeModifiedScore (URANIE::DataServer::TStochasticAttribute::ELawType tlaw)=0
 Compute the modified value of the statistic.
 
virtual const char * getTitleStatistic ()=0
 Return the name of the statistic.
 
void addUniformLaw (Double_t min, Double_t max, Color_t color=kBlack, TString stitle="")
 Add in the current graph the CDF of a uniform law defined by theses parameters $ \min, \max $.
 
void addNormalLaw (Double_t mu, Double_t sigma, Color_t color=kBlack, TString stitle="")
 Add in the current graph the CDF of a normal law defined by theses parameters $ {\cal N} ( \mu, \sigma ) $.
 
void addLogNormalLaw (Double_t mu, Double_t sigma, Color_t color=kBlack, TString stitle="")
 Adds in the current graph the CDF of a LogNormal law defined by theses parameters $ \mu, \sigma $ of the $ ln(X) \sim {\cal N} ( \mu, \sigma )$ law.
 
void addWeibullLaw (Double_t beta, Double_t m, Color_t color=kBlack, TString stitle="")
 Add in the current graph the CDF of a Weibull law defined by theses parameters $ {\cal W}( lambda, k) $.
 
void addGammaLaw (Double_t alpha, Double_t beta, Color_t color=kBlack, TString stitle="")
 Adds in the current graph the CDF of a Gamma law defined by theses parameters $ {\cal G}amma ( \alpha, \beta) $.
 
void addBetaLaw (Double_t alpha, Double_t beta, Double_t min, Double_t max, Color_t color=kBlack, TString stitle="")
 Add in the current graph the CDF of a Beta law defined by theses parameters $ \alpha, \beta, \min, \max $.
 
void treatOption (Option_t *option)
 Treats the options of the computeScore method.
 
Computes the PValues

from the book D'Agostino and Stephens "Goodness-of-Fit Techniques" Depends on the case :

  • Case0 : The parameters of the law are defined
  • Case3 : The parameters of the law are not defined; then they are estimated by the MLE ("Maximum LikeLihood Estimate"). and the type of the EDF test (KS, CvM ou AD)
void computePValue (URANIE::DataServer::TStochasticAttribute::ELawType tlaw)
 Compute the PValue of the test.
 
Double_t computePValue_KS_CaseAll (Double_t dz)
 Compute the PValue of the test KS in all cases (Case0 && Case3)
 
Double_t computePValue_AD_Case0 (Double_t dz)
 Compute the PValue of the test AD in case Case0 (all parameters of the PDF are given)
 
Double_t computePValue_AD_Case3 (Double_t dz, URANIE::DataServer::TStochasticAttribute::ELawType tlaw)
 Compute the PValue of the test AD in case Case3 (all parameters of the PDF are estimated by MLE)
 
Double_t computePValue_CvM_Case3 (Double_t dz, URANIE::DataServer::TStochasticAttribute::ELawType tlaw)
 Compute the PValue of the test CvM in case Case3 (all parameters of the PDF are estimated by MLE)
 
Log Printing
void setLog ()
 
void unsetLog ()
 
void changeLog ()
 
Bool_t getLog ()
 
virtual void printLog (Option_t *option="")
 

Public Attributes

ETestType _nTestType
 The type of test.
 
ECase _nCase
 
URANIE::DataServer::TDataServer * _tds
 Pointer to a TDS.
 
Bool_t _blog
 Boolean for edit the log.
 
Bool_t _bquiet
 Boolean for being quiet.
 

Protected Member Functions

Protected Methods
void findParametersLaw (TString soption, TVectorD &vec)
 Finds the parameters in the option string.
 
Int_t nTimesValueLTvalVector (Double_t dval, TVectorD vecValues)
 
TGraph * createGraph (Color_t color=kBlack)
 

Protected Attributes

Long64_t _nS
 the size of selected values
 
Double_t _pValue
 The pValue of the last test.
 
Double_t _dScore
 The score of the last test.
 
Double_t _dModifiedScore
 The Modified value score of the last test.
 
TString _sAttributes
 The string of attributes separated by the colon ":" character.
 
TList * _listOfLaws
 The list of law to test.
 
Bool_t _bDrawSame
 Boolean for drawing in the current pad.
 

Detailed Description

Description of the class TTester.

These class implements the tests of fit based on the Empirical Distribution Function ("EDF"). The EDF is a step function, calculated from the sample, wich estimates the population distribution.

See also
Goofness-of-fit techniques, R. D'Agostini and M. Stephens, Dekker

Member Enumeration Documentation

◆ ECase

Enumerator
kCase0 
kCase1 
kCase2 
kCase3 
kCaseUnknown 

◆ ETestType

Enumerator
kKolmogorovSmirnov 
kAndersonDarling 
kCramerVonMises 
kUnknown 

Constructor & Destructor Documentation

◆ TTester()

URANIE::UncertModeler::TTester::TTester ( URANIE::DataServer::TDataServer *  tds,
ETestType  nTestType,
const char *  varexpinput,
Option_t *  option = "" 
)

Constructor by the type of Test.

Referenced by ClassImp().

◆ ~TTester()

virtual URANIE::UncertModeler::TTester::~TTester ( )
virtual

Default constructor.

Referenced by ClassImp().

Member Function Documentation

◆ addBetaLaw()

void URANIE::UncertModeler::TTester::addBetaLaw ( Double_t  alpha,
Double_t  beta,
Double_t  min,
Double_t  max,
Color_t  color = kBlack,
TString  stitle = "" 
)

Add in the current graph the CDF of a Beta law defined by theses parameters $ \alpha, \beta, \min, \max $.

Parameters
alpha(Double_t) the alpha value of the law
beta(Double_t) the beta value of the law
min(Double_t) the minnimum value of the law
max(Double_t) the maximum value of the law
color(Color_t)[kBlack] the color of the line of the CDF's law
stitle(TString)[] the title of the entry in the legend box

Referenced by ClassImp().

◆ addGammaLaw()

void URANIE::UncertModeler::TTester::addGammaLaw ( Double_t  alpha,
Double_t  beta,
Color_t  color = kBlack,
TString  stitle = "" 
)

Adds in the current graph the CDF of a Gamma law defined by theses parameters $ {\cal G}amma ( \alpha, \beta) $.

where the CDF formula is : $ F(x) = inc_gamma( \alpha, x/scale )$

Parameters
alpha(Double_t) The scale parameter $ \alpha $ of the Gamma law
beta(Double_t) The shape parameter $ \beta $ of the the Gamma law
color(Color_t)[kBlack] The color of the line of the CDF's law
stitle(TString)[] The title of the entry in the legend box

Referenced by ClassImp().

◆ addLogNormalLaw()

void URANIE::UncertModeler::TTester::addLogNormalLaw ( Double_t  mu,
Double_t  sigma,
Color_t  color = kBlack,
TString  stitle = "" 
)

Adds in the current graph the CDF of a LogNormal law defined by theses parameters $ \mu, \sigma $ of the $ ln(X) \sim {\cal N} ( \mu, \sigma )$ law.

Parameters
mu(Double_t) The mean $ \mu $ of the $ ln(X)$ law
sigma(Double_t) The standard-deviation $ \sigma $ of the $ ln(X)$ law
color(Color_t)[kBlack] The color of the line of the CDF's law
stitle(TString)[] The title of the entry in the legend box

Referenced by ClassImp().

◆ addNormalLaw()

void URANIE::UncertModeler::TTester::addNormalLaw ( Double_t  mu,
Double_t  sigma,
Color_t  color = kBlack,
TString  stitle = "" 
)

Add in the current graph the CDF of a normal law defined by theses parameters $ {\cal N} ( \mu, \sigma ) $.

Parameters
mu(Double_t) the mean $ \mu $ of the law
sigma(Double_t) the standard-deviation $ \sigma $ of the law
color(Color_t)[kBlack] the color of the line of the CDF's law
stitle(TString)[] the title of the entry in the legend box

Referenced by ClassImp().

◆ addUniformLaw()

void URANIE::UncertModeler::TTester::addUniformLaw ( Double_t  min,
Double_t  max,
Color_t  color = kBlack,
TString  stitle = "" 
)

Add in the current graph the CDF of a uniform law defined by theses parameters $ \min, \max $.

Parameters
min(Double_t) the minnimum value of the law
max(Double_t) the maximum value of the law
color(Color_t)[kBlack] the color of the line of the CDF's law
stitle(TString)[] the title of the entry in the legend box

Referenced by ClassImp().

◆ addWeibullLaw()

void URANIE::UncertModeler::TTester::addWeibullLaw ( Double_t  beta,
Double_t  m,
Color_t  color = kBlack,
TString  stitle = "" 
)

Add in the current graph the CDF of a Weibull law defined by theses parameters $ {\cal W}( lambda, k) $.

where the CDF formula is : $ F(x) = 1 - \exp\left(-(x/\lambda)^k\right) \forall x \ge 0$

Parameters
lambda(Double_t) The scale parameter $ \beta $ of the Weibull law
k(Double_t) The shape parameter $ m $ of the the Weibull law
color(Color_t)[kBlack] The color of the line of the CDF's law
stitle(TString)[] The title of the entry in the legend box

Referenced by ClassImp().

◆ changeLog()

void URANIE::UncertModeler::TTester::changeLog ( )
inline

References _blog.

◆ computeModifiedScore()

virtual void URANIE::UncertModeler::TTester::computeModifiedScore ( URANIE::DataServer::TStochasticAttribute::ELawType  tlaw)
pure virtual

◆ computePValue()

void URANIE::UncertModeler::TTester::computePValue ( URANIE::DataServer::TStochasticAttribute::ELawType  tlaw)

Compute the PValue of the test.

Parameters
tlawType of Law to test

Referenced by ClassImp().

◆ computePValue_AD_Case0()

Double_t URANIE::UncertModeler::TTester::computePValue_AD_Case0 ( Double_t  dz)

Compute the PValue of the test AD in case Case0 (all parameters of the PDF are given)

Internal function. These function comes from the Table 4.3 Page 108 of the "Goodness-of-Fit Techniques" book Distribution of $ A^{2}$, Case 0 : The table gives $ q = I\!\!P [ A^{2} > z ] $. Here, we return the Pvalue = $ 1 -q$.

For the AD test, the AD value is not modified (

See also
Table 4.2 page 105 : $ A^{2} \forall n \ge 5$)
Parameters
dzThe $A^{2}$ value of the AD test (not modified)
Returns
the PValue

Referenced by ClassImp().

◆ computePValue_AD_Case3()

Double_t URANIE::UncertModeler::TTester::computePValue_AD_Case3 ( Double_t  dz,
URANIE::DataServer::TStochasticAttribute::ELawType  tlaw 
)

Compute the PValue of the test AD in case Case3 (all parameters of the PDF are estimated by MLE)

Parameters
dzThe $A^{2}$ value of the AD test (not modified)
tlawThe type of Law
Returns
the PValue

Referenced by ClassImp().

◆ computePValue_CvM_Case3()

Double_t URANIE::UncertModeler::TTester::computePValue_CvM_Case3 ( Double_t  dz,
URANIE::DataServer::TStochasticAttribute::ELawType  tlaw 
)

Compute the PValue of the test CvM in case Case3 (all parameters of the PDF are estimated by MLE)

Parameters
dzThe $W^{2}$ value of the CvM test (not modified)
tlawThe type of Law
Returns
the PValue

Referenced by ClassImp().

◆ computePValue_KS_CaseAll()

Double_t URANIE::UncertModeler::TTester::computePValue_KS_CaseAll ( Double_t  dz)

Compute the PValue of the test KS in all cases (Case0 && Case3)

Parameters
dzThe $ D^{\star} $ value of the KS test (modified score)
Returns
The PValue

Referenced by ClassImp().

◆ computeScore()

Double_t URANIE::UncertModeler::TTester::computeScore ( Option_t *  option = "")

Compute the score of the statistic.

A test to generalize the compute the score of the statistic.

Referenced by ClassImp().

◆ computeTheCurrentScore()

virtual Double_t URANIE::UncertModeler::TTester::computeTheCurrentScore ( Int_t  i,
Double_t  dzi 
)
pure virtual

The true computation with a vector of information and return the current value.

Parameters
i(Int_t) the current index
dzi(Double_t) the current value of Z = F(x)

Implemented in URANIE::UncertModeler::TTestAndersonDarling, URANIE::UncertModeler::TTestCramerVonMises, and URANIE::UncertModeler::TTestKolmogorovSmirnov.

Referenced by ClassImp().

◆ createGraph()

TGraph * URANIE::UncertModeler::TTester::createGraph ( Color_t  color = kBlack)
protected

Referenced by ClassImp().

◆ findParametersLaw()

void URANIE::UncertModeler::TTester::findParametersLaw ( TString  soption,
TVectorD &  vec 
)
protected

Finds the parameters in the option string.

Example of string parameters are

  • "normal(-125.025, 125)"
  • "normal(125.25, .235)"

Referenced by ClassImp().

◆ getCase()

ECase URANIE::UncertModeler::TTester::getCase ( )
inline

getCaseOfGoF

kCase0 : All the parameters of the law are specified kCase3 : All the parameters of the law are not specified. Then they are estimated by MLE.

The other cases are not treated in URANIE: kCase1 : Just the position parameters are specified, but not the scaled parameters kCase2 : Just the scaled parameters are specified, but not the position parameters

Returns
ECase the case of Goodness-Of-Fit

References _nCase.

◆ getID()

Int_t URANIE::UncertModeler::TTester::getID ( )
inline

Referenced by ClassImp().

◆ getLog()

Bool_t URANIE::UncertModeler::TTester::getLog ( )
inline

References _blog.

◆ getModifiedScore()

Double_t URANIE::UncertModeler::TTester::getModifiedScore ( )
inline

Get the modified score of the last test.

References _dModifiedScore.

◆ getPValue()

Double_t URANIE::UncertModeler::TTester::getPValue ( )
inline

Get the pValue of the last test.

References _pValue.

◆ getScore()

Double_t URANIE::UncertModeler::TTester::getScore ( )
inline

Get the score of the last test.

References _dScore.

◆ getTitleStatistic()

virtual const char * URANIE::UncertModeler::TTester::getTitleStatistic ( )
pure virtual

◆ nTimesValueLTvalVector()

Int_t URANIE::UncertModeler::TTester::nTimesValueLTvalVector ( Double_t  dval,
TVectorD  vecValues 
)
protected

Referenced by ClassImp().

◆ printLog()

virtual void URANIE::UncertModeler::TTester::printLog ( Option_t *  option = "")
virtual

Referenced by ClassImp().

◆ setLog()

void URANIE::UncertModeler::TTester::setLog ( )
inline

References _blog.

◆ treatOption()

void URANIE::UncertModeler::TTester::treatOption ( Option_t *  option)

Treats the options of the computeScore method.

The option of the computeScore method can be :

  • "normal"
  • "normal:normal(200,35)" to test the normal law with $ (\mu, \sigma) = (200, 35) $

Referenced by ClassImp().

◆ unsetLog()

void URANIE::UncertModeler::TTester::unsetLog ( )
inline

References _blog.

Member Data Documentation

◆ _bDrawSame

Bool_t URANIE::UncertModeler::TTester::_bDrawSame
protected

Boolean for drawing in the current pad.

Referenced by ClassImp().

◆ _blog

Bool_t URANIE::UncertModeler::TTester::_blog

Boolean for edit the log.

Referenced by changeLog(), ClassImp(), ClassImp(), ClassImp(), ClassImp(), getLog(), setLog(), and unsetLog().

◆ _bquiet

Bool_t URANIE::UncertModeler::TTester::_bquiet

Boolean for being quiet.

Referenced by ClassImp().

◆ _dModifiedScore

Double_t URANIE::UncertModeler::TTester::_dModifiedScore
protected

The Modified value score of the last test.

Referenced by ClassImp(), ClassImp(), ClassImp(), ClassImp(), and getModifiedScore().

◆ _dScore

Double_t URANIE::UncertModeler::TTester::_dScore
protected

The score of the last test.

Referenced by ClassImp(), ClassImp(), ClassImp(), ClassImp(), and getScore().

◆ _listOfLaws

TList* URANIE::UncertModeler::TTester::_listOfLaws
protected

The list of law to test.

Referenced by ClassImp().

◆ _nCase

ECase URANIE::UncertModeler::TTester::_nCase

◆ _nS

Long64_t URANIE::UncertModeler::TTester::_nS
protected

the size of selected values

Referenced by ClassImp(), ClassImp(), ClassImp(), and ClassImp().

◆ _nTestType

ETestType URANIE::UncertModeler::TTester::_nTestType

The type of test.

Referenced by ClassImp().

◆ _pValue

Double_t URANIE::UncertModeler::TTester::_pValue
protected

The pValue of the last test.

Referenced by ClassImp(), and getPValue().

◆ _sAttributes

TString URANIE::UncertModeler::TTester::_sAttributes
protected

The string of attributes separated by the colon ":" character.

Referenced by ClassImp().

◆ _tds

URANIE::DataServer::TDataServer* URANIE::UncertModeler::TTester::_tds

Pointer to a TDS.

Referenced by ClassImp().