fitLinearRobust.C

Go to the documentation of this file.
00001 #include "TRandom.h"
00002 #include "TGraphErrors.h"
00003 #include "TF1.h"
00004 #include "TCanvas.h"
00005 #include "TLegend.h"
00006 
00007 void fitLinearRobust()
00008 {
00009 //This tutorial shows how the least trimmed squares regression,
00010 //included in the TLinearFitter class, can be used for fitting
00011 //in cases when the data contains outliers.
00012 //Here the fitting is done via the TGraph::Fit function with option "rob":
00013 //If you want to use the linear fitter directly for computing
00014 //the robust fitting coefficients, just use the TLinearFitter::EvalRobust
00015 //function instead of TLinearFitter::Eval
00016 //Author: Anna Kreshuk
00017    
00018    //First generate a dataset, where 20% of points are spoiled by large
00019    //errors
00020    Int_t npoints = 250;
00021    Int_t fraction = Int_t(0.8*npoints);
00022    Double_t *x = new Double_t[npoints];
00023    Double_t *y = new Double_t[npoints];
00024    Double_t *e = new Double_t[npoints];
00025    TRandom r;
00026    Int_t i;
00027    for (i=0; i<fraction; i++){
00028       //the good part of the sample
00029       x[i]=r.Uniform(-2, 2);
00030       e[i]=1;
00031       y[i]=1 + 2*x[i] + 3*x[i]*x[i] + 4*x[i]*x[i]*x[i] + e[i]*r.Gaus();
00032    }
00033    for (i=fraction; i<npoints; i++){
00034       //the bad part of the sample
00035       x[i]=r.Uniform(-1, 1);
00036       e[i]=1;
00037       y[i] = 1 + 2*x[i] + 3*x[i]*x[i] + 4*x[i]*x[i]*x[i] + r.Landau(10, 5);
00038    }
00039 
00040    TGraphErrors *grr = new TGraphErrors(npoints, x, y, 0, e);
00041    grr->SetMinimum(-30);
00042    grr->SetMaximum(80);
00043    TF1 *ffit1 = new TF1("ffit1", "pol3", -5, 5);
00044    TF1 *ffit2 = new TF1("ffit2", "pol3", -5, 5);
00045    ffit1->SetLineColor(kBlue);
00046    ffit2->SetLineColor(kRed);
00047    TCanvas *myc = new TCanvas("myc", "Linear and robust linear fitting");
00048    myc->SetFillColor(42);
00049    myc->SetGrid();
00050    grr->Draw("ap");
00051    //first, let's try to see the result sof ordinary least-squares fit:
00052    printf("Ordinary least squares:\n");
00053    grr->Fit(ffit1);
00054    //the fitted function doesn't really follow the pattern of the data
00055    //and the coefficients are far from the real ones
00056 
00057    printf("Resistant Least trimmed squares fit:\n");
00058    //Now let's try the resistant regression
00059    //The option "rob=0.75" means that we want to use robust fitting and
00060    //we know that at least 75% of data is good points (at least 50% of points
00061    //should be good to use this algorithm). If you don't specify any number
00062    //and just use "rob" for the option, default value of (npoints+nparameters+1)/2
00063    //will be taken
00064    grr->Fit(ffit2, "+rob=0.75");
00065    //
00066    TLegend *leg = new TLegend(0.6, 0.8, 0.89, 0.89);
00067    leg->AddEntry(ffit1, "Ordinary least squares", "l");
00068    leg->AddEntry(ffit2, "LTS regression", "l");
00069    leg->SetFillColor(42);
00070    leg->Draw();
00071 
00072    delete [] x;
00073    delete [] y;
00074    delete [] e;
00075 
00076 }

Generated on Tue Jul 5 15:44:07 2011 for ROOT_528-00b_version by  doxygen 1.5.1