fitLinear.C: Example of fitting with a linear function, using TLinearFitter | Fitting tutorials | fitLinearRobust.C: This tutorial shows how the least trimmed squares regression, |
#include "TLinearFitter.h" #include "TF1.h" #include "TRandom.h" void fitLinear2() { //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) //Author: Anna Kreshuk Int_t n=100; Int_t i; TRandom rand; TLinearFitter *lf=new TLinearFitter(5); //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] = rand.Uniform(-10, 10); x[1 + i*5] = rand.Uniform(-10, 10); x[2 + i*5] = rand.Uniform(-10, 10); x[3 + i*5] = rand.Uniform(-10, 10); x[4 + i*5] = rand.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] + rand.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] = rand.Uniform(-10, 10); x[1+i*5] = rand.Uniform(-10, 10); x[2+i*5] = rand.Uniform(-10, 10); x[3+i*5] = rand.Uniform(-10, 10); x[4+i*5] = rand.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] + rand.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; } fitLinear2.C:1 fitLinear2.C:2 fitLinear2.C:3 fitLinear2.C:4 fitLinear2.C:5 fitLinear2.C:6 fitLinear2.C:7 fitLinear2.C:8 fitLinear2.C:9 fitLinear2.C:10 fitLinear2.C:11 fitLinear2.C:12 fitLinear2.C:13 fitLinear2.C:14 fitLinear2.C:15 fitLinear2.C:16 fitLinear2.C:17 fitLinear2.C:18 fitLinear2.C:19 fitLinear2.C:20 fitLinear2.C:21 fitLinear2.C:22 fitLinear2.C:23 fitLinear2.C:24 fitLinear2.C:25 fitLinear2.C:26 fitLinear2.C:27 fitLinear2.C:28 fitLinear2.C:29 fitLinear2.C:30 fitLinear2.C:31 fitLinear2.C:32 fitLinear2.C:33 fitLinear2.C:34 fitLinear2.C:35 fitLinear2.C:36 fitLinear2.C:37 fitLinear2.C:38 fitLinear2.C:39 fitLinear2.C:40 fitLinear2.C:41 fitLinear2.C:42 fitLinear2.C:43 fitLinear2.C:44 fitLinear2.C:45 fitLinear2.C:46 fitLinear2.C:47 fitLinear2.C:48 fitLinear2.C:49 fitLinear2.C:50 fitLinear2.C:51 fitLinear2.C:52 fitLinear2.C:53 fitLinear2.C:54 fitLinear2.C:55 fitLinear2.C:56 fitLinear2.C:57 fitLinear2.C:58 fitLinear2.C:59 fitLinear2.C:60 fitLinear2.C:61 fitLinear2.C:62 fitLinear2.C:63 fitLinear2.C:64 fitLinear2.C:65 fitLinear2.C:66 fitLinear2.C:67 fitLinear2.C:68 fitLinear2.C:69 fitLinear2.C:70 fitLinear2.C:71 fitLinear2.C:72 fitLinear2.C:73 fitLinear2.C:74 fitLinear2.C:75 fitLinear2.C:76 fitLinear2.C:77 fitLinear2.C:78 fitLinear2.C:79 fitLinear2.C:80 fitLinear2.C:81 fitLinear2.C:82 fitLinear2.C:83 fitLinear2.C:84 fitLinear2.C:85 fitLinear2.C:86 fitLinear2.C:87 fitLinear2.C:88 fitLinear2.C:89 fitLinear2.C:90 fitLinear2.C:91 fitLinear2.C:92 fitLinear2.C:93 fitLinear2.C:94 fitLinear2.C:95 fitLinear2.C:96 fitLinear2.C:97 fitLinear2.C:98 fitLinear2.C:99 fitLinear2.C:100 fitLinear2.C:101 fitLinear2.C:102 fitLinear2.C:103 fitLinear2.C:104 fitLinear2.C:105 fitLinear2.C:106 fitLinear2.C:107 fitLinear2.C:108 fitLinear2.C:109 fitLinear2.C:110 |
|