From $ROOTSYS/tutorials/fit/fitLinear.C

#include "TGraphErrors.h"
#include "TF1.h"
#include "TRandom.h"
#include "TCanvas.h"
#include "TLegend.h"
#include "TMath.h"


void makePoints(Int_t n, Double_t *x, Double_t *y, Double_t *e, Int_t p);

void fitLinear()
{
   //Example of fitting with a linear function, using TLinearFitter
   //This example is for a TGraphErrors, but it can also be used
   //when fitting a histogram, a TGraph2D or a TMultiGraph
   //Author: Anna Kreshuk

   Int_t n = 40;
   Double_t *x = new Double_t[n];
   Double_t *y = new Double_t[n];
   Double_t *e = new Double_t[n];
   TCanvas *myc = new TCanvas("myc",
      "Fitting 3 TGraphErrors with linear functions");
   myc->SetFillColor(42);
   myc->SetGrid();

   //Generate points along a 3rd degree polynomial:
   makePoints(n, x, y, e, 3);
   TGraphErrors *gre3 = new TGraphErrors(n, x, y, 0, e);
   gre3->Draw("a*");
   //Fit the graph with the predefined "pol3" function
   gre3->Fit("pol3");
   //Access the fit resuts
   TF1 *f3 = gre3->GetFunction("pol3");
   f3->SetLineWidth(1);

   //Generate points along a sin(x)+sin(2x) function
   makePoints(n, x, y, e, 2);
   TGraphErrors *gre2=new TGraphErrors(n, x, y, 0, e);
   gre2->Draw("*same");
   gre2->SetMarkerColor(kBlue);
   gre2->SetLineColor(kBlue);
   //The fitting function can be predefined and passed to the Fit function
   //The "++" mean that the linear fitter should be used, and the following
   //formula is equivalent to "[0]*sin(x) + [1]*sin(2*x)"
   //A function, defined this way, is in no way different from any other TF1,
   //it can be evaluated, drawn, you can get its parameters, etc.
   //The fit result (parameter values, parameter errors, chisquare, etc) are
   //written into the fitting function.
   TF1 *f2 = new TF1("f2", "sin(x) ++ sin(2*x)", -2, 2);
   gre2->Fit(f2);
   f2 = gre2->GetFunction("f2");
   f2->SetLineColor(kBlue);
   f2->SetLineWidth(1);

   //Generate points along a -2+exp(-x) function
   makePoints(n, x, y, e, 4);
   TGraphErrors *gre4=new TGraphErrors(n, x, y, 0, e);
   gre4->Draw("*same");
   gre4->SetMarkerColor(kRed);
   gre4->SetLineColor(kRed);
   //If you don't want to define the function, you can just pass the string
   //with the the formula:
   gre4->Fit("1 ++ exp(-x)");
   //Access the fit results:
   TF1 *f4 = gre4->GetFunction("1 ++ exp(-x)");
   f4->SetName("f4");
   f4->SetLineColor(kRed);
   f4->SetLineWidth(1);

   TLegend *leg = new TLegend(0.3, 0.7, 0.65, 0.9);
   leg->AddEntry(gre3, " -7 + 2*x*x + x*x*x", "p");
   leg->AddEntry(gre2, "sin(x) + sin(2*x)", "p");
   leg->AddEntry(gre4, "-2 + exp(-x)", "p");
   leg->Draw();
   leg->SetFillColor(42);


}

void makePoints(Int_t n, Double_t *x, Double_t *y, Double_t *e, Int_t p)
{
  Int_t i;
  TRandom r;

  if (p==2) {
    for (i=0; i<n; i++) {
      x[i] = r.Uniform(-2, 2);
      y[i]=TMath::Sin(x[i]) + TMath::Sin(2*x[i]) + r.Gaus()*0.1;
      e[i] = 0.1;
    }
  }
  if (p==3) {
    for (i=0; i<n; i++) {
      x[i] = r.Uniform(-2, 2);
      y[i] = -7 + 2*x[i]*x[i] + x[i]*x[i]*x[i]+ r.Gaus()*0.1;
      e[i] = 0.1;
    }
  }
  if (p==4) {
    for (i=0; i<n; i++) {
      x[i] = r.Uniform(-2, 2);
      y[i]=-2 + TMath::Exp(-x[i]) + r.Gaus()*0.1;
      e[i] = 0.1;
    }
  }
}
 fitLinear.C:1
 fitLinear.C:2
 fitLinear.C:3
 fitLinear.C:4
 fitLinear.C:5
 fitLinear.C:6
 fitLinear.C:7
 fitLinear.C:8
 fitLinear.C:9
 fitLinear.C:10
 fitLinear.C:11
 fitLinear.C:12
 fitLinear.C:13
 fitLinear.C:14
 fitLinear.C:15
 fitLinear.C:16
 fitLinear.C:17
 fitLinear.C:18
 fitLinear.C:19
 fitLinear.C:20
 fitLinear.C:21
 fitLinear.C:22
 fitLinear.C:23
 fitLinear.C:24
 fitLinear.C:25
 fitLinear.C:26
 fitLinear.C:27
 fitLinear.C:28
 fitLinear.C:29
 fitLinear.C:30
 fitLinear.C:31
 fitLinear.C:32
 fitLinear.C:33
 fitLinear.C:34
 fitLinear.C:35
 fitLinear.C:36
 fitLinear.C:37
 fitLinear.C:38
 fitLinear.C:39
 fitLinear.C:40
 fitLinear.C:41
 fitLinear.C:42
 fitLinear.C:43
 fitLinear.C:44
 fitLinear.C:45
 fitLinear.C:46
 fitLinear.C:47
 fitLinear.C:48
 fitLinear.C:49
 fitLinear.C:50
 fitLinear.C:51
 fitLinear.C:52
 fitLinear.C:53
 fitLinear.C:54
 fitLinear.C:55
 fitLinear.C:56
 fitLinear.C:57
 fitLinear.C:58
 fitLinear.C:59
 fitLinear.C:60
 fitLinear.C:61
 fitLinear.C:62
 fitLinear.C:63
 fitLinear.C:64
 fitLinear.C:65
 fitLinear.C:66
 fitLinear.C:67
 fitLinear.C:68
 fitLinear.C:69
 fitLinear.C:70
 fitLinear.C:71
 fitLinear.C:72
 fitLinear.C:73
 fitLinear.C:74
 fitLinear.C:75
 fitLinear.C:76
 fitLinear.C:77
 fitLinear.C:78
 fitLinear.C:79
 fitLinear.C:80
 fitLinear.C:81
 fitLinear.C:82
 fitLinear.C:83
 fitLinear.C:84
 fitLinear.C:85
 fitLinear.C:86
 fitLinear.C:87
 fitLinear.C:88
 fitLinear.C:89
 fitLinear.C:90
 fitLinear.C:91
 fitLinear.C:92
 fitLinear.C:93
 fitLinear.C:94
 fitLinear.C:95
 fitLinear.C:96
 fitLinear.C:97
 fitLinear.C:98
 fitLinear.C:99
 fitLinear.C:100
 fitLinear.C:101
 fitLinear.C:102
 fitLinear.C:103
 fitLinear.C:104
 fitLinear.C:105
 fitLinear.C:106
 fitLinear.C:107
 fitLinear.C:108