AliRoot Core  3dc7879 (3dc7879)
AliTMinuitToolkit Class Reference

The AliTMinuitToolkit serves as an easy to use interface for the TMinuit. More...

#include <AliTMinuitToolkit.h>

Inheritance diagram for AliTMinuitToolkit:

Public Types

enum  EVerboseFlags {
  kPrintMinuit =0x00001, kSysInfo =0x00002, kPrintAll =0x00004, kStreamFcn =0x00010,
  kStreamFcnPoint =0x00020, kFCNInfo =0x000400
}
 
enum  EFitFlags { kCovarFailed =0x0001 }
 

Public Member Functions

 AliTMinuitToolkit (const char *streamerName=nullptr, const char *mode="recreate")
 
virtual ~AliTMinuitToolkit ()
 
void SetStreamer (const char *streamerName, const char *mode="recreate")
 Set stream for writing info into TTree. see examples(STEER/STEERBase/TTreeStream.cxx) More...
 
TTreeSRedirectorGetStreamer ()
 
void ClearData ()
 
void Fit (Option_t *option="")
 Internal function that calls the fitter. More...
 
void FitHistogram (const TH1 *his, Option_t *option="")
 Fit a one dimensional histogram. More...
 
void FitGraph (const TGraph *gr, Option_t *option="")
 Fit a TGraph object. More...
 
Long64_t FillFitter (TTree *inputTree, TString values, TString variables, TString selection, Int_t firstEntry, Int_t nEntries, Bool_t doReset=kTRUE)
 Fill input matrices of the fitter. More...
 
TString GetFitFunctionAsAlias (Option_t *option="", TTree *tree=nullptr)
 Format alias string for the for tree alias or for fit title. More...
 
void Bootstrap (ULong_t nIter, const char *reportName, Option_t *option=nullptr)
 Bootstrap minimization. More...
 
void TwoFoldCrossValidation (UInt_t nIter, const char *reportName, Option_t *option=nullptr)
 
void MISAC (Int_t nFitPoints, UInt_t nIter, const char *reportName, Option_t *option=nullptr)
 Fitting with outlier removal inspired by RANSAC - see https://en.wikipedia.org/w/index.php?title=Random_sample_consensus&oldid=767657742. More...
 
void SetFitFunction (TF1 *formula, Bool_t doReset)
 Set formula of the fitted function and reset parameter. See doReset. More...
 
void SetInitialParam (const TMatrixD *initialParam)
 Set initial parameters matrix. More...
 
TMatrixD * GetInitialParam ()
 
void SetLogLikelihoodFunction (TF1 *logLikelihood)
 
void SetFitAlgorithm (Char_t *const name)
 
void SetMaxCalls (Int_t calls)
 
void SetTolerance (Double_t tol)
 
void SetPoints (TMatrixD *points)
 
void SetValues (TMatrixD *values)
 
void EnableRobust (Bool_t b)
 
const TMatrixD * GetPoints () const
 
const TMatrixD * GetValues () const
 
TF1 * GetFormula () const
 
const TF1 * GetLogLikelihoodFunction () const
 
const TVectorD * GetParameters () const
 
const TVectorD * GetRMSEstimator () const
 
const TMatrixD * GetCovarianceMatrix () const
 
const TMatrixD * GetMISACEstimators () const
 
Double_t GetChisquare ()
 
Int_t GetStatus ()
 
Bool_t IsHuberCost () const
 
TObjArrayGetListOfVariables ()
 
void SetVerbose (Int_t verbose)
 
Bool_t IsVerbose (Int_t flag)
 

Static Public Member Functions

static void FitterFCN (int &nParam, double *info, double &chi2, double *gin, int iflag)
 Internal function used in the TMinuit minimization. More...
 
static AliTMinuitToolkitFit (TH1 *his, const char *fitterName, Option_t *fitOption="", Option_t *hisOption="", Option_t *funOption="", Double_t xMin=0, Double_t xMax=0)
 Make histogram fit. More...
 
static AliTMinuitToolkitFit (TGraph *graph, const char *fitterName, Option_t *fitOption="", Option_t *graphOption="", Option_t *funOption="", Double_t xMin=0, Double_t xMax=0)
 Make graph fit. More...
 
static AliTMinuitToolkitGetPredefinedFitter (const std::string &fitterName)
 
static void SetPredefinedFitter (const std::string &fitterName, AliTMinuitToolkit *fitter)
 
static void RegisterDefaultFitters ()
 
static AliTMinuitToolkitRegisterPlaneFitter (Int_t nPlanes, Int_t logLikeType)
 Register default fitters. More...
 
static Double_t GaussCachyLogLike (const Double_t *x, const Double_t *p)
 
static Double_t HuberLogLike (const Double_t *x, const Double_t *p)
 
static Double_t PseudoHuberLogLike (const Double_t *x, const Double_t *p)
 
static Double_t GausKurtosisSkewness (const Double_t *x, const Double_t *p)
 Gaus convoluted with Erf to emulate kurtosis and skewness. More...
 
static Double_t RndmGaus (Double_t mean=0, Double_t sigma=1)
 Random gaus generator as static function (to use in root TFormula, TTreeFormula) More...
 
static Double_t RndmLandau (Double_t mean=0, Double_t sigma=1)
 Random gaus generator as static function (to use in root TFormula, TTreeFormula) More...
 
static void SetFunctionDrawOption (TF1 *fun, Option_t *option)
 

Private Member Functions

 AliTMinuitToolkit (const AliTMinuitToolkit &)
 
AliTMinuitToolkitoperator= (const AliTMinuitToolkit &)
 

Private Attributes

TTreeSRedirectorfStreamer
 
Int_t fVerbose
 
TF1 * fFormula
 
TF1 * fLogLikelihoodFunction
 
TString fFitAlgorithm
 
TMatrixD * fPoints
 
TMatrixD * fValues
 
TObjArrayfVarNames
 
TObjArrayfValueNames
 
TArrayI fPointIndex
 
TMatrixD * fInitialParam
 
TVectorD * fParam
 
TVectorD * fRMSEstimator
 
TMatrixD * fMISACEstimators
 
TMatrixD * fCovar
 
Int_t fStatus
 
Double_t fChi2
 
Int_t fMaxCalls
 
Double_t fPrecision
 
Bool_t fIsHuberCost
 
Int_t fFCNCounter
 

Static Private Attributes

static std::map< std::string, AliTMinuitToolkit * > fPredefinedFitters
 

Detailed Description

The AliTMinuitToolkit serves as an easy to use interface for the TMinuit.

Usage:

  • 1.) Setting the fit formula:
    • Set fit function:
      void SetFitFunction(TFormula * formula, Bool_t doReset)".
  • 2.) Adding the data points
    • Direct specification of the points is possible via
      void SetPoints(TMatrixD * points)
      • Each point corresponds to one row in the matrix.
    • The fitter is then started with the command
      void Fit()
    • In order to fit a histogram, graphs in once use:
      AliTMinuitToolkit::FitHistogram(const TH1 * his, Option_t *option)
      AliTMinuitToolkit::FitGraph(const TGraph * gr, Option_t *option)
      AliTMinuitToolkit::Fit(TH1* his, const char *fitterName, Option_t* option, Option_t* goption,Option_t* foption, Double_t xMin, Double_t xMax)
      AliTMinuitToolkit::Fit(TGraph *graph, const char *fitterName, Option_t* option, Option_t* goption, Option_t* foption, Double_t xMin, Double_t xMax)
    • Filling input points using TTree queries (example)
      tool1D->FillFitter(inputTree,"test1D+noise:1/sqrt(12.+0)","testx0", "", 0,fitEntries);
  • 3.) Accessing the fit results
    • The N parameters of the formula are stored in a N-dimensional vector which is returned by
    • The covariance (NxN) matrix of the fit
  • 4.) Log likelihood minimization
    • see example
      $ALICE_ROOT/../src/STAT/test/AliTMinuitToolkitTest.C:Test1D()

      Even a few outliers can lead to wrong results of a least-squares fitting procedure. In this case the use of robust(resistant) methods can be helpful, but a stronger dependence on starting values or convergence to local minimum can occur.

    • Cost (log likelihood) functions
      1. By default chi2 minimization is invoked.
      2. Optionally more robust log likelihood- Huber cost function can be specified instead of chi^2. For small deviations the function behaves like x^2 and for larger deviations like |x| - the so called Huber estimator:
        h(x) = x^2 , for x < 2.5*sigma
        h(x) = 2*(2.5*sigma)*x - (2.5*sigma)^2 , for x > 2.5*sigma
      3. User defined log likelihood function can be specified
        TF1 * f1Cost = new TF1("1","abs(x)<10?-log(0.8*exp(-x**2)+0.2/(1+x**2)):-log(0.2/(1+x**2))",-20,20); // 80 % gaus + 20% cachy
        tool1D->SetLogLikelihoodFunction(f1Cost);
        tool1D->Fit();
  • 5.) Additional fitting strategy (options)
  • 6.) Work in progress

Example usage in the $AliRoot_SRC/STAT/test/

  • TODO - save example images in the doxygen/html directory (STEERING part ?) - currently saved by hand - to check with OFFLINE
  • $AliRoot_SRC/STAT/test/AliTMinuitToolkitTest.C+
  • $AliRoot_SRC/STAT/test/AliTMinuitToolkitTestLinear.C+
    • compare fits example with different likelihood/resp fits strategies for data with outliers
      AliTMinuiToolkit.TestLinearFitExample.png
    • compare performance fits example with different likelihood/resp fits strategies for data with outliers
      AliTMinuiToolkit.TestLinearFitStatExample.png

Definition at line 104 of file AliTMinuitToolkit.h.

Member Enumeration Documentation

Enumerator
kCovarFailed 

Definition at line 115 of file AliTMinuitToolkit.h.

Enumerator
kPrintMinuit 
kSysInfo 
kPrintAll 
kStreamFcn 
kStreamFcnPoint 
kFCNInfo 

Definition at line 106 of file AliTMinuitToolkit.h.

Constructor & Destructor Documentation

AliTMinuitToolkit::AliTMinuitToolkit ( const char *  streamerName = nullptr,
const char *  mode = "recreate" 
)

Definition at line 43 of file AliTMinuitToolkit.cxx.

Referenced by MISAC(), RegisterPlaneFitter(), and SetPredefinedFitter().

AliTMinuitToolkit::~AliTMinuitToolkit ( )
virtual

Definition at line 78 of file AliTMinuitToolkit.cxx.

AliTMinuitToolkit::AliTMinuitToolkit ( const AliTMinuitToolkit )
private

Member Function Documentation

void AliTMinuitToolkit::Bootstrap ( ULong_t  nIter,
const char *  reportName,
Option_t *  option = nullptr 
)

Bootstrap minimization.

  • fitting parameters done several times on modified data sample (random samples with replacement)
    • to emulate different data population
    • reduce problem with local minimum in the M-estimator estimated parameters - robust mean/median
    • obtained RMS used to estimate error (take into account also model error)
  • TODO - parameters robust estimator setting
    Parameters
    nIter- number of iteration
    reportName- report name - in case results streamed
    option- fit option

Definition at line 582 of file AliTMinuitToolkit.cxx.

Referenced by Fit(), GetStreamer(), and Test1D().

void AliTMinuitToolkit::ClearData ( )

Definition at line 107 of file AliTMinuitToolkit.cxx.

Referenced by FillFitter(), FitGraph(), FitHistogram(), and GetStreamer().

void AliTMinuitToolkit::EnableRobust ( Bool_t  b)
inline

Definition at line 144 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkit_TestHistogram(), and Test1D().

Long64_t AliTMinuitToolkit::FillFitter ( TTree *  inputTree,
TString  values,
TString  variables,
TString  selection,
Int_t  firstEntry,
Int_t  nEntries,
Bool_t  doReset = kTRUE 
)

Fill input matrices of the fitter.

  • TODO - add better working example (using code/trees from the test directory)
    Parameters
    inputTree- pointer to input tree
    values- values description (parameter:<error>)
    variables- variables array (var0:<var1>)
    selection- tree selection
    firstEntry- first entry to query
    nEntries- number of entries to query
    doReset- switch reset previous/or append
    Returns
    Example: TPC occupancy 6 dimensional hyperplane fit
    hyp6R->FillFitter(treeQA, "NoThreshold.fElements:0.1*NoThreshold.fElements", "1:dRNorm:Phi2Norm:QNorm:gxNorm:gyNorm", "region&&occFitCut&&rndm<0.2", 0, 10000000);
    hyp6R->Fit();

Definition at line 415 of file AliTMinuitToolkit.cxx.

Referenced by GetStreamer(), and Test1D().

void AliTMinuitToolkit::Fit ( Option_t *  option = "")

Internal function that calls the fitter.

Parameters
option
  • bootstrap[0-9]* -
  • twofold[0-9]* -
  • misac([0-9]*,[0-9]*) -

Definition at line 275 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), AliTMinuitToolkitTestLinear(), Bootstrap(), FitGraph(), FitHistogram(), GetStreamer(), IsVerbose(), MISAC(), Test1D(), and TwoFoldCrossValidation().

AliTMinuitToolkit * AliTMinuitToolkit::Fit ( TH1 *  his,
const char *  fitterName,
Option_t *  fitOption = "",
Option_t *  hisOption = "",
Option_t *  funOption = "",
Double_t  xMin = 0,
Double_t  xMax = 0 
)
static

Make histogram fit.

  • TODO- add example code snippet - using the code from test directory
    Parameters
    his- pointer to histogram
    fitterName- fitter name ()
    fitOption- fit Option
    hisOption- graphic option for histogram
    funOption- graphic option for function drawing
    xMin- fit range min
    xMax- fit range max
    Returns
    - pointer to the fitter

Definition at line 943 of file AliTMinuitToolkit.cxx.

AliTMinuitToolkit * AliTMinuitToolkit::Fit ( TGraph *  graph,
const char *  fitterName,
Option_t *  fitOption = "",
Option_t *  graphOption = "",
Option_t *  funOption = "",
Double_t  xMin = 0,
Double_t  xMax = 0 
)
static

Make graph fit.

  • TODO- add example code snippet - using the code from test directory
    Parameters
    graph- pointer to graph
    fitterName- fitter name ()
    fitOption- fitter setting option
    graphOption- graphic option
    funOption- graphic option for function
    xMin- fit range min
    xMax- fit range max
    Returns
    - pointer to the fitter

Definition at line 979 of file AliTMinuitToolkit.cxx.

void AliTMinuitToolkit::FitGraph ( const TGraph *  gr,
Option_t *  option = "" 
)

Fit a TGraph object.

Parameters
gr- const TGraph object
option- fit option e.g (default, bootstrap20,misac .... )
  • In alternative usage more option and graph range can be specified
    Fit(TGraph* graph, const char *fitterName, Option_t* fitOption = "", Option_t* graphOption = "", Option_t* funOption="",Double_t xMin = 0, Double_t xMax = 0)

Definition at line 169 of file AliTMinuitToolkit.cxx.

Referenced by Fit(), and GetStreamer().

void AliTMinuitToolkit::FitHistogram ( const TH1 *  his,
Option_t *  option = "" 
)

Fit a one dimensional histogram.

Parameters
his- const TH1 object
option- fit option e.g (default, bootstrap20,misac ... )
  • In alternative usage more option and histogram range can be specified
    Fit(TH1* his, const char *fitterName, Option_t* fitOption = "", Option_t* hisOption = "", Option_t *funOption="", Double_t xMin = 0, Double_t xMax = 0)

Definition at line 150 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), Fit(), and GetStreamer().

void AliTMinuitToolkit::FitterFCN ( int &  nParam,
double *  info,
double &  chi2,
double *  gin,
int  iflag 
)
static

Internal function used in the TMinuit minimization.

  • see documentation for FCN parameters g.g in ROOT- TMinuit::Eval
Parameters
nPar- not used in our interface
grad- gradient of cost function
chi2(fVal))- cost function to minimize (chi2/log-likelihood/)
gin(par)- function parameters
iflag(flag)- flag used by root (1.) print/ 2.) derivative/ 3. final fit)
  • in our interface also possible to specify stream data (kStreamFcnPoint,kStreamFcn) typedef void(* FCNFunc_t) (Int_t &nPar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
  • Original documentation in TMinuit::Eval documentation in root
    • Input parameters:
      • nPar: number of currently variable parameters
      • par: array of (constant and variable) parameters
      • flag: Indicates what is to be calculated (see example below) ()
      • grad: array of gradients
    • Output parameters:
      • fVal: The calculated function value.
      • grad: The (optional) vector of first derivatives).

Definition at line 205 of file AliTMinuitToolkit.cxx.

Referenced by Bootstrap(), Fit(), MISAC(), SetValues(), and TwoFoldCrossValidation().

Double_t AliTMinuitToolkit::GausKurtosisSkewness ( const Double_t *  x,
const Double_t *  p 
)
static

Gaus convoluted with Erf to emulate kurtosis and skewness.

Parameters
x
p
  • p[0],p[1],p[2] - as for gauss
  • p[3] - Kurtosis - tails - proportional to 4 moment
  • p[4] - Asymmetry - proportional to 3 moment
Returns

Definition at line 562 of file AliTMinuitToolkit.cxx.

Referenced by SetPredefinedFitter().

Double_t AliTMinuitToolkit::GaussCachyLogLike ( const Double_t *  x,
const Double_t *  p 
)
static

Log likelihood for the gaus+cauchy - used for Log likelihood minimization of distribution with tails

Parameters
x-
p- p[0] - gaus fraction
  • p[1] - cauchy 0.5*FWHM
Returns
log likelihood for the mixture of gaus and cauchy distribution representation as in https://en.wikipedia.org/wiki/Cauchy_distribution

Definition at line 524 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and SetPredefinedFitter().

Double_t AliTMinuitToolkit::GetChisquare ( )
inline

Definition at line 155 of file AliTMinuitToolkit.h.

const TMatrixD* AliTMinuitToolkit::GetCovarianceMatrix ( ) const
inline

Definition at line 152 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear().

TString AliTMinuitToolkit::GetFitFunctionAsAlias ( Option_t *  option = "",
TTree *  tree = nullptr 
)

Format alias string for the for tree alias or for fit title.

Parameters
option- use latex in case of title request
Returns
- string description

Definition at line 465 of file AliTMinuitToolkit.cxx.

Referenced by GetStreamer(), and Test1D().

TF1* AliTMinuitToolkit::GetFormula ( ) const
inline

Definition at line 148 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkit_TestHistogram(), Fit(), and FitterFCN().

TMatrixD* AliTMinuitToolkit::GetInitialParam ( )
inline

Definition at line 136 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear().

TObjArray* AliTMinuitToolkit::GetListOfVariables ( )
inline

Definition at line 158 of file AliTMinuitToolkit.h.

Referenced by GetFitFunctionAsAlias().

const TF1* AliTMinuitToolkit::GetLogLikelihoodFunction ( ) const
inline

Definition at line 149 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear().

const TMatrixD* AliTMinuitToolkit::GetMISACEstimators ( ) const
inline

Definition at line 153 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear().

const TVectorD* AliTMinuitToolkit::GetParameters ( ) const
inline
const TMatrixD* AliTMinuitToolkit::GetPoints ( ) const
inline

Definition at line 146 of file AliTMinuitToolkit.h.

Referenced by FitterFCN().

static AliTMinuitToolkit* AliTMinuitToolkit::GetPredefinedFitter ( const std::string &  fitterName)
inlinestatic

Definition at line 164 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear(), and Fit().

const TVectorD* AliTMinuitToolkit::GetRMSEstimator ( ) const
inline

Definition at line 151 of file AliTMinuitToolkit.h.

Referenced by AliTMinuitToolkitTestLinear(), and GetFitFunctionAsAlias().

Int_t AliTMinuitToolkit::GetStatus ( )
inline

Definition at line 156 of file AliTMinuitToolkit.h.

TTreeSRedirector* AliTMinuitToolkit::GetStreamer ( )
inline

Definition at line 122 of file AliTMinuitToolkit.h.

const TMatrixD* AliTMinuitToolkit::GetValues ( ) const
inline

Definition at line 147 of file AliTMinuitToolkit.h.

Referenced by FitterFCN().

Double_t AliTMinuitToolkit::HuberLogLike ( const Double_t *  x,
const Double_t *  p 
)
static

Huber loss is a loss function used in robust regression - representation as in https://en.wikipedia.org/w/index.php?title=Huber_loss&diff=809252384&oldid=798150659

Parameters
x-
p- p[0] - gaus region boundary
Returns
- loss function (used as a log likelihood)

Definition at line 536 of file AliTMinuitToolkit.cxx.

Referenced by SetPredefinedFitter().

Bool_t AliTMinuitToolkit::IsHuberCost ( ) const
inline

Definition at line 157 of file AliTMinuitToolkit.h.

Referenced by FitterFCN().

Bool_t AliTMinuitToolkit::IsVerbose ( Int_t  flag)
inline

Definition at line 160 of file AliTMinuitToolkit.h.

Referenced by FitterFCN().

void AliTMinuitToolkit::MISAC ( Int_t  nFitPoints,
UInt_t  nIter,
const char *  reportName,
Option_t *  option = nullptr 
)

Fitting with outlier removal inspired by RANSAC - see https://en.wikipedia.org/w/index.php?title=Random_sample_consensus&oldid=767657742.

  • TODO: Make robust n-dimensional estimator
    Parameters
    nFitPoints- number of points used in one fit
    nIter- number of iterations
    reportName- report name (!=NULL- specified stream to save intermediate data)
    option- fit option
  • Algorithm:
    • 0.) Generate permutations of data
    • 1.) Calculate fit parameters for small subset of data
    • 2.) Calculate normalized distance between each fit vector
    • 3.) Reject outliers and calculate robust mean
    • 4.) Extract combined statistic
    • 5.) Dump results into tree in verbose mode

TODO- check failure of the fit

Definition at line 720 of file AliTMinuitToolkit.cxx.

Referenced by Fit(), and GetStreamer().

AliTMinuitToolkit& AliTMinuitToolkit::operator= ( const AliTMinuitToolkit )
private

Referenced by SetPredefinedFitter().

Double_t AliTMinuitToolkit::PseudoHuberLogLike ( const Double_t *  x,
const Double_t *  p 
)
static

Pseudo Huber loss is a loss function used in robust regression - representation as in https://en.wikipedia.org/w/index.php?title=Huber_loss&diff=809252384&oldid=798150659

Parameters
x-
p- p[0] - gaus region boundary
Returns
- loss function (used as a log likelihood)

Definition at line 547 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), MISAC(), RegisterPlaneFitter(), and SetPredefinedFitter().

static void AliTMinuitToolkit::RegisterDefaultFitters ( )
static
AliTMinuitToolkit * AliTMinuitToolkit::RegisterPlaneFitter ( Int_t  nPlanes,
Int_t  logLikeType 
)
static

Register default fitters.

  • Registered function can be used in standard minimization using sting identifier, e.g:
    • AliTMinuitToolkit::Fit(TH1* his, const char *fitterName, Option_t* option, Option_t* goption,Option_t* foption, Double_t xMin, Double_t xMax)
    • AliTMinuitToolkit::Fit(TH1* his, const char *fitterName, Option_t* option, Option_t* goption,Option_t* foption, Double_t xMin, Double_t xMax)
  • Registered fitter can be accessed - eventually modified (e.g to set initial values)
    • * Default fitters:
      * formula: pol<N>, gauss
      * log-likelihood: chi2, likeGausCachy, huber norm
      TF1 *likeGausCachy = new TF1("likeGausCachy", AliTMinuitToolkit::GaussCachyLogLike,-10,10,2);
      TF1 *likePseudoHuber = new TF1("likePseudoHuber", AliTMinuitToolkit::PseudoHuberLogLike,-10,10,1);
      likeGausCachy->SetParameters(0.8,1);
      likePseudoHuber->SetParameter(0,3);
      likeGausCachy->GetXaxis()->SetTitle("#Delta");
      likeGausCachy->GetYaxis()->SetTitle("-logLikelihood");
      likeGausCachy->SetLineColor(2);
      likePseudoHuber->GetXaxis()->SetTitle("#Delta");
      likePseudoHuber->GetYaxis()->SetTitle("-logLikelihood");
      likePseudoHuber->SetLineColor(4);
      //
      for (Int_t iPol=0; iPol<10; iPol++){ //register polynomial fitters
      TMatrixD initPar(iPol+1,4);
      for (Int_t iPar=0; iPar<iPol+1; iPar++) initPar(iPar,1)=1;
      //
      TF1 *fpol = new TF1(TString::Format("fpol%d",iPol).Data(),TString::Format("pol%d",iPol).Data());
      fitter1D->SetVerbose(0x1); fitter1D->SetFitFunction(fpol,kTRUE);
      fitter1D->SetInitialParam(&initPar);
      fitter1D->SetName(TString::Format("pol%d",iPol).Data());
      AliTMinuitToolkit::SetPredefinedFitter(fitter1D->GetName(),fitter1D);
      // gaus log likelihood cost function
      AliTMinuitToolkit * fitter1DR = new AliTMinuitToolkit();
      fitter1DR->SetVerbose(0x1); fitter1DR->SetFitFunction(fpol,kTRUE);
      fitter1DR->SetLogLikelihoodFunction(likeGausCachy);
      fitter1DR->SetInitialParam(&initPar);
      fitter1DR->SetName(TString::Format("pol%dR",iPol).Data());
      AliTMinuitToolkit::SetPredefinedFitter(fitter1DR->GetName(),fitter1DR);
      // huber cost function
      AliTMinuitToolkit * fitter1DH = new AliTMinuitToolkit();
      fitter1DH->SetVerbose(0x1); fitter1DH->SetFitFunction(fpol,kTRUE);
      fitter1DH->SetInitialParam(&initPar);
      //fitter1DH->EnableRobust(true);
      fitter1DH->SetLogLikelihoodFunction(likePseudoHuber);
      fitter1DH->SetName(TString::Format("pol%dH",iPol).Data());
      AliTMinuitToolkit::SetPredefinedFitter(fitter1DH->GetName(),fitter1DH);
      }
      //
      TF1 *fGaus = new TF1("fgaus","gaus");
      TMatrixD initPar(3,4);
      initPar(0,0)=0; initPar(0,1)=1;
      initPar(1,0)=0; initPar(1,1)=1;
      initPar(2,0)=1; initPar(2,1)=1;
      fitterGR->SetVerbose(0x1); fitterGR->SetFitFunction(fGaus,kTRUE);
      fitterGR->SetLogLikelihoodFunction(likeGausCachy);
      fitterGR->SetInitialParam(&initPar);
      //
      //
      // TFormulaPrimitive::AddFormula(new TFormulaPrimitive("GausKS","GausKS",(TFormulaPrimitive::GenFuncG)AliTMinuitToolkit::GausKurtosisSkewness,5));
      // hack number of arguments in the primitive (protected member) - //TODO - fix ROOT
      // TFormulaPrimitive * p = TFormulaPrimitive::FindFormula("GausKS");
      // TDataMember *m = (TDataMember*)p->IsA()->GetListOfDataMembers()->FindObject("fNArguments");
      // *((Int_t*)(((char*)p)+m->GetOffset()))=1;
      }
      \breief Register plane fitters (As in the TLinearFitter)
      Global registered fitters can be obtained e.g :
      Parameters
      nPlanes- number of planes
      logLikeType- type of the log likelihood
      • 0) chi2, 1.) gaus+cauchy, 2.) pseudo-huber
      Returns
      return registered fitter

Definition at line 902 of file AliTMinuitToolkit.cxx.

Referenced by SetPredefinedFitter().

Double_t AliTMinuitToolkit::RndmGaus ( Double_t  mean = 0,
Double_t  sigma = 1 
)
static

Random gaus generator as static function (to use in root TFormula, TTreeFormula)

Parameters
mean- mean Gaus
sigma- sigma gaus
Returns

Definition at line 505 of file AliTMinuitToolkit.cxx.

Referenced by SetPredefinedFitter().

Double_t AliTMinuitToolkit::RndmLandau ( Double_t  mean = 0,
Double_t  sigma = 1 
)
static

Random gaus generator as static function (to use in root TFormula, TTreeFormula)

Parameters
mean
sigma
Returns

Definition at line 513 of file AliTMinuitToolkit.cxx.

Referenced by SetPredefinedFitter().

void AliTMinuitToolkit::SetFitAlgorithm ( Char_t *const  name)
inline

Definition at line 138 of file AliTMinuitToolkit.h.

void AliTMinuitToolkit::SetFitFunction ( TF1 *  formula,
Bool_t  doReset 
)

Set formula of the fitted function and reset parameter. See doReset.

Parameters
formula- formula of the fitted function
doReset- if true reset parameters related to the fit function

Definition at line 120 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), GetStreamer(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and Test1D().

void AliTMinuitToolkit::SetFunctionDrawOption ( TF1 *  fun,
Option_t *  option 
)
static

Parse draw option for function

Parameters
fun- pointer to function to change
option- draw option for function

Definition at line 1014 of file AliTMinuitToolkit.cxx.

Referenced by Fit(), and SetPredefinedFitter().

void AliTMinuitToolkit::SetInitialParam ( const TMatrixD *  paramLimits)

Set initial parameters matrix.

Parameters
paramLimits- Initialization matrix
  • number of lines = number of parameters
  • 4 columns
    • [0] - initial value
    • [1] - initial error
    • [2],[3] - optional parameters to specify minimum and maximum of the parameter (see TMinuit)

Definition at line 140 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkit_TestHistogram(), GetStreamer(), MISAC(), RegisterGaussFitters(), RegisterPlaneFitter(), and Test1D().

void AliTMinuitToolkit::SetLogLikelihoodFunction ( TF1 *  logLikelihood)
inline
void AliTMinuitToolkit::SetMaxCalls ( Int_t  calls)
inline

Definition at line 139 of file AliTMinuitToolkit.h.

void AliTMinuitToolkit::SetPoints ( TMatrixD *  points)
inline

Definition at line 141 of file AliTMinuitToolkit.h.

static void AliTMinuitToolkit::SetPredefinedFitter ( const std::string &  fitterName,
AliTMinuitToolkit fitter 
)
inlinestatic
void AliTMinuitToolkit::SetStreamer ( const char *  streamerName,
const char *  mode = "recreate" 
)

Set stream for writing info into TTree. see examples(STEER/STEERBase/TTreeStream.cxx)

Parameters
streamerName- Name of .root file
mode- option for .root file open
  • recreate
  • update
  • ....

Definition at line 99 of file AliTMinuitToolkit.cxx.

Referenced by AliTMinuitToolkitTestLinear().

void AliTMinuitToolkit::SetTolerance ( Double_t  tol)
inline

Definition at line 140 of file AliTMinuitToolkit.h.

void AliTMinuitToolkit::SetValues ( TMatrixD *  values)
inline

Definition at line 142 of file AliTMinuitToolkit.h.

void AliTMinuitToolkit::SetVerbose ( Int_t  verbose)
inline

Definition at line 159 of file AliTMinuitToolkit.h.

Referenced by MISAC(), RegisterGaussFitters(), and Test1D().

void AliTMinuitToolkit::TwoFoldCrossValidation ( UInt_t  nIter,
const char *  reportName,
Option_t *  option = nullptr 
)

Two fold cross validation - https://en.wikipedia.org/wiki/Cross-validation_(statistics)

  • Sample randomly divided to the independent training set and complementary test set (Repeated random sub-sampling validation)
  • log likelihood cost function calculated for training fit set
  • TODO - Used for fit model validation. Report - Post processing using streamer not available for the moment
  • TODO - Add: The results are then averaged over the splits as described in the wiki section - (Repeated random sub-sampling validation)
    Parameters
    nIter- number of iterations
    reportName- report name (identifier in the streamer)
    option- fit option

Definition at line 637 of file AliTMinuitToolkit.cxx.

Referenced by Fit(), GetStreamer(), and Test1D().

Member Data Documentation

Double_t AliTMinuitToolkit::fChi2
private

Definition at line 197 of file AliTMinuitToolkit.h.

Referenced by Bootstrap(), Fit(), GetChisquare(), and MISAC().

TMatrixD* AliTMinuitToolkit::fCovar
private
Int_t AliTMinuitToolkit::fFCNCounter
private

Definition at line 201 of file AliTMinuitToolkit.h.

Referenced by Bootstrap(), FitterFCN(), and TwoFoldCrossValidation().

TString AliTMinuitToolkit::fFitAlgorithm
private

Definition at line 185 of file AliTMinuitToolkit.h.

Referenced by Fit(), and SetFitAlgorithm().

TF1* AliTMinuitToolkit::fFormula
private
TMatrixD* AliTMinuitToolkit::fInitialParam
private

Definition at line 191 of file AliTMinuitToolkit.h.

Referenced by Fit(), GetInitialParam(), SetInitialParam(), and ~AliTMinuitToolkit().

Bool_t AliTMinuitToolkit::fIsHuberCost
private

Definition at line 200 of file AliTMinuitToolkit.h.

Referenced by EnableRobust(), and IsHuberCost().

TF1* AliTMinuitToolkit::fLogLikelihoodFunction
private
Int_t AliTMinuitToolkit::fMaxCalls
private

Definition at line 198 of file AliTMinuitToolkit.h.

Referenced by Fit(), and SetMaxCalls().

TMatrixD* AliTMinuitToolkit::fMISACEstimators
private

Definition at line 194 of file AliTMinuitToolkit.h.

Referenced by GetMISACEstimators(), MISAC(), and ~AliTMinuitToolkit().

TVectorD* AliTMinuitToolkit::fParam
private
TArrayI AliTMinuitToolkit::fPointIndex
private

Definition at line 190 of file AliTMinuitToolkit.h.

Referenced by Bootstrap(), FitterFCN(), MISAC(), and TwoFoldCrossValidation().

TMatrixD* AliTMinuitToolkit::fPoints
private
Double_t AliTMinuitToolkit::fPrecision
private

Definition at line 199 of file AliTMinuitToolkit.h.

Referenced by Fit(), and SetTolerance().

std::map< std::string, AliTMinuitToolkit * > AliTMinuitToolkit::fPredefinedFitters
staticprivate

Definition at line 202 of file AliTMinuitToolkit.h.

Referenced by GetPredefinedFitter(), and SetPredefinedFitter().

TVectorD* AliTMinuitToolkit::fRMSEstimator
private
Int_t AliTMinuitToolkit::fStatus
private

Definition at line 196 of file AliTMinuitToolkit.h.

Referenced by Bootstrap(), Fit(), GetStatus(), MISAC(), and TwoFoldCrossValidation().

TTreeSRedirector* AliTMinuitToolkit::fStreamer
private
TObjArray* AliTMinuitToolkit::fValueNames
private

Definition at line 189 of file AliTMinuitToolkit.h.

Referenced by FillFitter().

TMatrixD* AliTMinuitToolkit::fValues
private
TObjArray* AliTMinuitToolkit::fVarNames
private

Definition at line 188 of file AliTMinuitToolkit.h.

Referenced by FillFitter(), GetFitFunctionAsAlias(), and GetListOfVariables().

Int_t AliTMinuitToolkit::fVerbose
private

Definition at line 182 of file AliTMinuitToolkit.h.

Referenced by Fit(), IsVerbose(), and SetVerbose().


The documentation for this class was generated from the following files: