Logo ROOT  
Reference Guide
line3Dfit.C File Reference

Detailed Description

View in nbviewer Open in SWAN Fitting of a TGraph2D with a 3D straight line

run this macro by doing:

root>.x line3Dfit.C+
Total Initial distance square = 8.65172e+07
Total final distance square 19988
****************************************
Minimizer is Minuit / Migrad
MinFCN = 19988
NDf = 0
Edm = 6.31385e-10
NCalls = 257
Par_0 = 10.3552 +/- 0.39454
Par_1 = 19.9602 +/- 0.0689483
Par_2 = 0.999198 +/- 0.0444211
Par_3 = 2.00379 +/- 0.00780946
(int) 0
#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>
#include <cassert>
using namespace ROOT::Math;
// define the parametric line equation
void line(double t, const double *p, double &x, double &y, double &z) {
// a parametric 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;
}
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) {}
// 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;
}
// 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;
}
if (first) {
std::cout << "Total Initial distance square = " << sum << std::endl;
}
first = false;
return sum;
}
};
Int_t line3Dfit()
{
//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
// make the functor objet
SumDistance2 sdist(gr);
ROOT::Math::Functor fcn(sdist,4);
// 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;
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
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->Draw("same");
return 0;
}
int main() {
return line3Dfit();
}
#define d(i)
Definition: RSha256.hxx:102
#define g(i)
Definition: RSha256.hxx:105
static const double x1[5]
int Int_t
Definition: RtypesCore.h:43
@ kRed
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
void Error(const char *location, const char *msgfmt,...)
#define N
TGLVector3 Cross(const TGLVector3 &v1, const TGLVector3 &v2)
Definition: TGLUtil.h:322
TRObject operator()(const T1 &t1) const
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
R__EXTERN TStyle * gStyle
Definition: TStyle.h:410
const ParameterSettings & ParSettings(unsigned int i) const
get the parameter settings for the i-th parameter (const method)
Definition: FitConfig.h:75
class containg the result of the fit and all the related information (fitted parameter values,...
Definition: FitResult.h:47
const double * GetParams() const
parameter values (return const pointer)
Definition: FitResult.h:176
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
Definition: FitResult.cxx:428
double MinFcnValue() const
Return value of the objective function (chi2 or likelihood) used in the fit.
Definition: FitResult.h:120
Fitter class, entry point for performing all type of fits.
Definition: Fitter.h:77
bool FitFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Fit using the a generic FCN function as a C++ callable object implementing double () (const double *)...
Definition: Fitter.h:610
bool SetFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Set a generic FCN function as a C++ callable object implementing double () (const double *) Note that...
Definition: Fitter.h:615
const FitResult & Result() const
get fit result
Definition: Fitter.h:384
const FitConfig & Config() const
access to the fit configuration (const method)
Definition: Fitter.h:412
void SetStepSize(double err)
set the step size
Documentation for class Functor class.
Definition: Functor.h:398
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:41
Double_t * GetY() const
Definition: TGraph2D.h:121
Double_t * GetX() const
Definition: TGraph2D.h:120
Int_t GetN() const
Definition: TGraph2D.h:119
Double_t * GetZ() const
Definition: TGraph2D.h:122
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2269
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:760
A 3-dimensional polyline.
Definition: TPolyLine3D.h:32
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
virtual void Draw(Option_t *option="")
Draw this 3-D polyline with its current attributes.
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:263
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1590
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1542
TLine * line
SVector< T, D > Unit(const SVector< T, D > &rhs)
Unit.
Definition: Functions.h:381
T Mag2(const SVector< T, D > &rhs)
Vector magnitude square Template to compute .
Definition: Functions.h:229
int main(int argc, char **argv)
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
TGraphErrors * gr
Definition: legend1.C:25
DisplacementVector3D< Cartesian3D< double >, DefaultCoordinateSystemTag > XYZVector
3D Vector based on the cartesian coordinates x,y,z in double precision
Definition: Vector3Dfwd.h:34
Definition: first.py:1
auto * l
Definition: textangle.C:4
static long int sum(long int i)
Definition: Factory.cxx:2275
Author
Lorenzo Moneta

Definition in file line3Dfit.C.