From $ROOTSYS/tutorials/fit/line3Dfit.C

//Fitting of a TGraph2D with a 3D straight line
//
// run this macro by doing: 
// 
// root>.x line3Dfit.C+
//

//Author: L. Moneta
//

   
#include <TMath.h>
#include <TGraph2D.h>
#include <TRandom2.h>
#include <TStyle.h>
#include <TCanvas.h>
#include <TF2.h>
#include <TH1.h>
#include <Math/Functor.h>
#include <TPolyLine3D.h>
#include <Math/Vector3D.h>
#include <Fit/Fitter.h>

using namespace ROOT::Math;


// define the parameteric line equation 
void line(double t, const double *p, double &x, double &y, double &z) {
   // a parameteric line is define from 6 parameters but 4 are independent
   // x0,y0,z0,z1,y1,z1 which are the coordinates of two points on the line
   // can choose z0 = 0 if line not parallel to x-y plane and z1 = 1; 
   x = p[0] + p[1]*t; 
   y = p[2] + p[3]*t;
   z = t; 
} 

// calculate distance line-point 
double distance2(double x,double y,double z, const double *p) {
   // distance line point is D= | (xp-x0) cross  ux | 
   // where ux is direction of line and x0 is a point in the line (like t = 0) 
   XYZVector xp(x,y,z); 
   XYZVector x0(p[0], p[2], 0. ); 
   XYZVector x1(p[0] + p[1], p[2] + p[3], 1. ); 
   XYZVector u = (x1-x0).Unit(); 
   double d2 = ((xp-x0).Cross(u)) .Mag2(); 
   return d2; 
}
bool first = true; 


// function Object to be minimized
struct SumDistance2 {
   // the TGraph is a data member of the object
   TGraph2D * fGraph;
   
   SumDistance2(TGraph2D * g) : fGraph(g) {}
   
   // implementation of the function to be minimized
   double operator() (const double * par) {
      
      assert(fGraph    != 0);
      double * x = fGraph->GetX();
      double * y = fGraph->GetY();
      double * z = fGraph->GetZ();
      int npoints = fGraph->GetN();
      double sum = 0;
      for (int i  = 0; i < npoints; ++i) {
         double d = distance2(x[i],y[i],z[i],par);
         sum += d;
#ifdef DEBUG
         if (first) std::cout << "point " << i << "\t"
            << x[i] << "\t"
            << y[i] << "\t"
            << z[i] << "\t"
            << std::sqrt(d) << std::endl;
#endif
      }
      if (first)
         std::cout << "Total Initial distance square = " << sum << std::endl;
      first = false;
      return sum;
   }
   
};

Int_t line3Dfit()
{
   gStyle->SetOptStat(0);
   gStyle->SetOptFit();


   //double e = 0.1;
   Int_t nd = 10000;


//    double xmin = 0; double ymin = 0;
//    double xmax = 10; double ymax = 10;

   TGraph2D * gr = new TGraph2D();

   // Fill the 2D graph
   double p0[4] = {10,20,1,2};

   // generate graph with the 3d points
   for (Int_t N=0; N<nd; N++) {
      double x,y,z = 0;
      // Generate a random number 
      double t = gRandom->Uniform(0,10);
      line(t,p0,x,y,z);
      double err = 1;
    // do a gaussian smearing around the points in all coordinates
      x += gRandom->Gaus(0,err);  
      y += gRandom->Gaus(0,err);  
      z += gRandom->Gaus(0,err);  
      gr->SetPoint(N,x,y,z);
      //dt->SetPointError(N,0,0,err);
   }
   // fit the graph now 
   
   ROOT::Fit::Fitter  fitter;
   
   
   // make the functor objet
   SumDistance2 sdist(gr);
#ifdef __CINT__
   ROOT::Math::Functor fcn(&sdist,4,"SumDistance2");
#else
   ROOT::Math::Functor fcn(sdist,4);
#endif
   // set the function and the initial parameter values
   double pStart[4] = {1,1,1,1};
   fitter.SetFCN(fcn,pStart);
   // set step sizes different than default ones (0.3 times parameter values)
   for (int i = 0; i < 4; ++i) fitter.Config().ParSettings(i).SetStepSize(0.01);
  
   bool ok = fitter.FitFCN();
   if (!ok) {
      Error("line3Dfit","Line3D Fit failed");
      return 1;
   }
   
   const ROOT::Fit::FitResult & result = fitter.Result();
   
   std::cout << "Total final distance square " << result.MinFcnValue() << std::endl;
   result.Print(std::cout);
   

   gr->Draw("p0");

   // get fit parameters
   const double * parFit = result.GetParams();

   
   // draw the fitted line
   int n = 1000;
   double t0 = 0;
   double dt = 10;
   TPolyLine3D *l = new TPolyLine3D(n);
   for (int i = 0; i <n;++i) {
      double t = t0+ dt*i/n;
      double x,y,z;
      line(t,parFit,x,y,z);
      l->SetPoint(i,x,y,z);
   }
   l->SetLineColor(kRed);
   l->Draw("same");
   
   // draw original line
   TPolyLine3D *l0 = new TPolyLine3D(n);
   for (int i = 0; i <n;++i) {
      double t = t0+ dt*i/n;
      double x,y,z;
      line(t,p0,x,y,z);
      l0->SetPoint(i,x,y,z);
   }
   l0->SetLineColor(kBlue);
   l0->Draw("same");
   return 0;
}

int main() { 
   return line3Dfit();
}
 line3Dfit.C:1
 line3Dfit.C:2
 line3Dfit.C:3
 line3Dfit.C:4
 line3Dfit.C:5
 line3Dfit.C:6
 line3Dfit.C:7
 line3Dfit.C:8
 line3Dfit.C:9
 line3Dfit.C:10
 line3Dfit.C:11
 line3Dfit.C:12
 line3Dfit.C:13
 line3Dfit.C:14
 line3Dfit.C:15
 line3Dfit.C:16
 line3Dfit.C:17
 line3Dfit.C:18
 line3Dfit.C:19
 line3Dfit.C:20
 line3Dfit.C:21
 line3Dfit.C:22
 line3Dfit.C:23
 line3Dfit.C:24
 line3Dfit.C:25
 line3Dfit.C:26
 line3Dfit.C:27
 line3Dfit.C:28
 line3Dfit.C:29
 line3Dfit.C:30
 line3Dfit.C:31
 line3Dfit.C:32
 line3Dfit.C:33
 line3Dfit.C:34
 line3Dfit.C:35
 line3Dfit.C:36
 line3Dfit.C:37
 line3Dfit.C:38
 line3Dfit.C:39
 line3Dfit.C:40
 line3Dfit.C:41
 line3Dfit.C:42
 line3Dfit.C:43
 line3Dfit.C:44
 line3Dfit.C:45
 line3Dfit.C:46
 line3Dfit.C:47
 line3Dfit.C:48
 line3Dfit.C:49
 line3Dfit.C:50
 line3Dfit.C:51
 line3Dfit.C:52
 line3Dfit.C:53
 line3Dfit.C:54
 line3Dfit.C:55
 line3Dfit.C:56
 line3Dfit.C:57
 line3Dfit.C:58
 line3Dfit.C:59
 line3Dfit.C:60
 line3Dfit.C:61
 line3Dfit.C:62
 line3Dfit.C:63
 line3Dfit.C:64
 line3Dfit.C:65
 line3Dfit.C:66
 line3Dfit.C:67
 line3Dfit.C:68
 line3Dfit.C:69
 line3Dfit.C:70
 line3Dfit.C:71
 line3Dfit.C:72
 line3Dfit.C:73
 line3Dfit.C:74
 line3Dfit.C:75
 line3Dfit.C:76
 line3Dfit.C:77
 line3Dfit.C:78
 line3Dfit.C:79
 line3Dfit.C:80
 line3Dfit.C:81
 line3Dfit.C:82
 line3Dfit.C:83
 line3Dfit.C:84
 line3Dfit.C:85
 line3Dfit.C:86
 line3Dfit.C:87
 line3Dfit.C:88
 line3Dfit.C:89
 line3Dfit.C:90
 line3Dfit.C:91
 line3Dfit.C:92
 line3Dfit.C:93
 line3Dfit.C:94
 line3Dfit.C:95
 line3Dfit.C:96
 line3Dfit.C:97
 line3Dfit.C:98
 line3Dfit.C:99
 line3Dfit.C:100
 line3Dfit.C:101
 line3Dfit.C:102
 line3Dfit.C:103
 line3Dfit.C:104
 line3Dfit.C:105
 line3Dfit.C:106
 line3Dfit.C:107
 line3Dfit.C:108
 line3Dfit.C:109
 line3Dfit.C:110
 line3Dfit.C:111
 line3Dfit.C:112
 line3Dfit.C:113
 line3Dfit.C:114
 line3Dfit.C:115
 line3Dfit.C:116
 line3Dfit.C:117
 line3Dfit.C:118
 line3Dfit.C:119
 line3Dfit.C:120
 line3Dfit.C:121
 line3Dfit.C:122
 line3Dfit.C:123
 line3Dfit.C:124
 line3Dfit.C:125
 line3Dfit.C:126
 line3Dfit.C:127
 line3Dfit.C:128
 line3Dfit.C:129
 line3Dfit.C:130
 line3Dfit.C:131
 line3Dfit.C:132
 line3Dfit.C:133
 line3Dfit.C:134
 line3Dfit.C:135
 line3Dfit.C:136
 line3Dfit.C:137
 line3Dfit.C:138
 line3Dfit.C:139
 line3Dfit.C:140
 line3Dfit.C:141
 line3Dfit.C:142
 line3Dfit.C:143
 line3Dfit.C:144
 line3Dfit.C:145
 line3Dfit.C:146
 line3Dfit.C:147
 line3Dfit.C:148
 line3Dfit.C:149
 line3Dfit.C:150
 line3Dfit.C:151
 line3Dfit.C:152
 line3Dfit.C:153
 line3Dfit.C:154
 line3Dfit.C:155
 line3Dfit.C:156
 line3Dfit.C:157
 line3Dfit.C:158
 line3Dfit.C:159
 line3Dfit.C:160
 line3Dfit.C:161
 line3Dfit.C:162
 line3Dfit.C:163
 line3Dfit.C:164
 line3Dfit.C:165
 line3Dfit.C:166
 line3Dfit.C:167
 line3Dfit.C:168
 line3Dfit.C:169
 line3Dfit.C:170
 line3Dfit.C:171
 line3Dfit.C:172
 line3Dfit.C:173
 line3Dfit.C:174
 line3Dfit.C:175
 line3Dfit.C:176
 line3Dfit.C:177
 line3Dfit.C:178
 line3Dfit.C:179
 line3Dfit.C:180
 line3Dfit.C:181
 line3Dfit.C:182
 line3Dfit.C:183
 line3Dfit.C:184