Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fitLinear2.C File Reference

Detailed Description

View in nbviewer Open in SWAN Fit a 5d hyperplane by n points, using the linear fitter directly

This macro shows some features of the TLinearFitter class A 5-d hyperplane is fit, a constant term is assumed in the hyperplane equation (y = a0 + a1*x0 + a2*x1 + a3*x2 + a4*x3 + a5*x4)

par[0]=0.000069+-0.001011
par[1]=3.999934+-0.000164
par[2]=0.999835+-0.000172
par[3]=1.999892+-0.000178
par[4]=2.999967+-0.000185
par[5]=0.199823+-0.000174
chisquare=70.148012
More Points:
par[0]=0.000551+-0.000712
par[1]=3.999910+-0.000121
par[2]=0.999886+-0.000125
par[3]=2.000067+-0.000123
par[4]=2.999915+-0.000127
par[5]=0.199918+-0.000130
chisquare=145.050322490891915
Without Constant
par[0]=3.999913+-0.000121
par[1]=0.999890+-0.000125
par[2]=2.000057+-0.000123
par[3]=2.999919+-0.000127
par[4]=0.199918+-0.000130
chisquare=145.649621
Fixed Constant:
par[0]=0.000536+-0.000712
par[1]=4.000000+-1.000000
par[2]=0.999884+-0.000125
par[3]=2.000070+-0.000123
par[4]=2.999910+-0.000127
par[5]=0.199920+-0.000130
chisquare=145.602523231222818
#include "TLinearFitter.h"
#include "TF1.h"
#include "TRandom.h"
void fitLinear2()
{
Int_t n=100;
Int_t i;
TRandom randNum;
//The predefined "hypN" functions are the fastest to fit
lf->SetFormula("hyp5");
Double_t *x=new Double_t[n*10*5];
Double_t *y=new Double_t[n*10];
Double_t *e=new Double_t[n*10];
//Create the points and put them into the fitter
for (i=0; i<n; i++){
x[0 + i*5] = randNum.Uniform(-10, 10);
x[1 + i*5] = randNum.Uniform(-10, 10);
x[2 + i*5] = randNum.Uniform(-10, 10);
x[3 + i*5] = randNum.Uniform(-10, 10);
x[4 + i*5] = randNum.Uniform(-10, 10);
e[i] = 0.01;
y[i] = 4*x[0+i*5] + x[1+i*5] + 2*x[2+i*5] + 3*x[3+i*5] + 0.2*x[4+i*5] + randNum.Gaus()*e[i];
}
//To avoid copying the data into the fitter, the following function can be used:
lf->AssignData(n, 5, x, y, e);
//A different way to put the points into the fitter would be to use
//the AddPoint function for each point. This way the points are copied and stored
//inside the fitter
//Perform the fitting and look at the results
lf->Eval();
TVectorD params;
TVectorD errors;
lf->GetParameters(params);
lf->GetErrors(errors);
for (Int_t i=0; i<6; i++)
printf("par[%d]=%f+-%f\n", i, params(i), errors(i));
Double_t chisquare=lf->GetChisquare();
printf("chisquare=%f\n", chisquare);
//Now suppose you want to add some more points and see if the parameters will change
for (i=n; i<n*2; i++) {
x[0+i*5] = randNum.Uniform(-10, 10);
x[1+i*5] = randNum.Uniform(-10, 10);
x[2+i*5] = randNum.Uniform(-10, 10);
x[3+i*5] = randNum.Uniform(-10, 10);
x[4+i*5] = randNum.Uniform(-10, 10);
e[i] = 0.01;
y[i] = 4*x[0+i*5] + x[1+i*5] + 2*x[2+i*5] + 3*x[3+i*5] + 0.2*x[4+i*5] + randNum.Gaus()*e[i];
}
//Assign the data the same way as before
lf->AssignData(n*2, 5, x, y, e);
lf->Eval();
lf->GetParameters(params);
lf->GetErrors(errors);
printf("\nMore Points:\n");
for (Int_t i=0; i<6; i++)
printf("par[%d]=%f+-%f\n", i, params(i), errors(i));
chisquare=lf->GetChisquare();
printf("chisquare=%.15f\n", chisquare);
//Suppose, you are not satisfied with the result and want to try a different formula
//Without a constant:
//Since the AssignData() function was used, you don't have to add all points to the fitter again
lf->SetFormula("x0++x1++x2++x3++x4");
lf->Eval();
lf->GetParameters(params);
lf->GetErrors(errors);
printf("\nWithout Constant\n");
for (Int_t i=0; i<5; i++)
printf("par[%d]=%f+-%f\n", i, params(i), errors(i));
chisquare=lf->GetChisquare();
printf("chisquare=%f\n", chisquare);
//Now suppose that you want to fix the value of one of the parameters
//Let's fix the first parameter at 4:
lf->SetFormula("hyp5");
lf->FixParameter(1, 4);
lf->Eval();
lf->GetParameters(params);
lf->GetErrors(errors);
printf("\nFixed Constant:\n");
for (i=0; i<6; i++)
printf("par[%d]=%f+-%f\n", i, params(i), errors(i));
chisquare=lf->GetChisquare();
printf("chisquare=%.15f\n", chisquare);
//The fixed parameters can then be released by the ReleaseParameter method
delete lf;
}
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
The Linear Fitter - For fitting functions that are LINEAR IN PARAMETERS.
virtual Double_t GetChisquare()
Get the Chisquare.
virtual void GetErrors(TVectorD &vpar)
Returns parameter errors.
virtual Int_t Eval()
Perform the fit and evaluate the parameters Returns 0 if the fit is ok, 1 if there are errors.
virtual void GetParameters(TVectorD &vpar)
Returns parameter values.
virtual void FixParameter(Int_t ipar)
Fixes paramter #ipar at its current value.
virtual void SetFormula(const char *formula)
Additive parts should be separated by "++".
virtual void AssignData(Int_t npoints, Int_t xncols, Double_t *x, Double_t *y, Double_t *e=0)
This function is to use when you already have all the data in arrays and don't want to copy them into...
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition TRandom.cxx:274
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:672
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Author
Anna Kreshuk

Definition in file fitLinear2.C.