line3Dfit.C

Go to the documentation of this file.
00001 //Fitting of a TGraph2D with a 3D straight line
00002 //
00003 // run this macro by doing: 
00004 // 
00005 // root>.x line3Dfit.C+
00006 //
00007 
00008 //Author: L. Moneta
00009 //
00010 
00011    
00012 #include <TMath.h>
00013 #include <TGraph2D.h>
00014 #include <TRandom2.h>
00015 #include <TStyle.h>
00016 #include <TCanvas.h>
00017 #include <TF2.h>
00018 #include <TH1.h>
00019 #include <TVirtualFitter.h>
00020 #include <TPolyLine3D.h>
00021 #include <Math/Vector3D.h>
00022 
00023 using namespace ROOT::Math;
00024 
00025 
00026 // define the parameteric line equation 
00027 void line(double t, double *p, double &x, double &y, double &z) { 
00028    // a parameteric line is define from 6 parameters but 4 are independent
00029    // x0,y0,z0,z1,y1,z1 which are the coordinates of two points on the line
00030    // can choose z0 = 0 if line not parallel to x-y plane and z1 = 1; 
00031    x = p[0] + p[1]*t; 
00032    y = p[2] + p[3]*t;
00033    z = t; 
00034 } 
00035 
00036 // calculate distance line-point 
00037 double distance2(double x,double y,double z, double *p) { 
00038    // distance line point is D= | (xp-x0) cross  ux | 
00039    // where ux is direction of line and x0 is a point in the line (like t = 0) 
00040    XYZVector xp(x,y,z); 
00041    XYZVector x0(p[0], p[2], 0. ); 
00042    XYZVector x1(p[0] + p[1], p[2] + p[3], 1. ); 
00043    XYZVector u = (x1-x0).Unit(); 
00044    double d2 = ((xp-x0).Cross(u)) .Mag2(); 
00045    return d2; 
00046 }
00047 bool first = true; 
00048 
00049 
00050 // function to be minimized 
00051 void SumDistance2(int &, double *, double & sum, double * par, int ) { 
00052    // the TGraph must be a global variable
00053    TGraph2D * gr = dynamic_cast<TGraph2D*>( (TVirtualFitter::GetFitter())->GetObjectFit() );
00054    assert(gr != 0);
00055    double * x = gr->GetX();
00056    double * y = gr->GetY();
00057    double * z = gr->GetZ();
00058    int npoints = gr->GetN();
00059    sum = 0;
00060    for (int i  = 0; i < npoints; ++i) { 
00061       double d = distance2(x[i],y[i],z[i],par); 
00062       sum += d;
00063 #ifdef DEBUG
00064       if (first) std::cout << "point " << i << "\t" 
00065                            << x[i] << "\t" 
00066                            << y[i] << "\t" 
00067                            << z[i] << "\t" 
00068                            << std::sqrt(d) << std::endl; 
00069 #endif
00070    }
00071    if (first) 
00072       std::cout << "Total sum2 = " << sum << std::endl;
00073    first = false;
00074 }
00075 
00076 void line3Dfit()
00077 {
00078    gStyle->SetOptStat(0);
00079    gStyle->SetOptFit();
00080 
00081 
00082    //double e = 0.1;
00083    Int_t nd = 10000;
00084 
00085 
00086 //    double xmin = 0; double ymin = 0;
00087 //    double xmax = 10; double ymax = 10;
00088 
00089    TGraph2D * gr = new TGraph2D();
00090 
00091    // Fill the 2D graph
00092    double p0[4] = {10,20,1,2};
00093 
00094    // generate graph with the 3d points
00095    for (Int_t N=0; N<nd; N++) {
00096       double x,y,z = 0;
00097       // Generate a random number 
00098       double t = gRandom->Uniform(0,10);
00099       line(t,p0,x,y,z);
00100       double err = 1;
00101     // do a gaussian smearing around the points in all coordinates
00102       x += gRandom->Gaus(0,err);  
00103       y += gRandom->Gaus(0,err);  
00104       z += gRandom->Gaus(0,err);  
00105       gr->SetPoint(N,x,y,z);
00106       //dt->SetPointError(N,0,0,err);
00107    }
00108    // fit the graph now 
00109    
00110    TVirtualFitter *min = TVirtualFitter::Fitter(0,4);
00111    min->SetObjectFit(gr);
00112    min->SetFCN(SumDistance2);
00113    
00114   
00115    Double_t arglist[10];
00116    arglist[0] = 3;
00117    min->ExecuteCommand("SET PRINT",arglist,1);
00118   
00119    double pStart[4] = {1,1,1,1};
00120    min->SetParameter(0,"x0",pStart[0],0.01,0,0);
00121    min->SetParameter(1,"Ax",pStart[1],0.01,0,0);
00122    min->SetParameter(2,"y0",pStart[2],0.01,0,0);
00123    min->SetParameter(3,"Ay",pStart[3],0.01,0,0);
00124     
00125    arglist[0] = 1000; // number of function calls 
00126    arglist[1] = 0.001; // tolerance 
00127    min->ExecuteCommand("MIGRAD",arglist,2);
00128 
00129   //if (minos) min->ExecuteCommand("MINOS",arglist,0);
00130    int nvpar,nparx; 
00131    double amin,edm, errdef;
00132    min->GetStats(amin,edm,errdef,nvpar,nparx);
00133    min->PrintResults(1,amin);
00134    gr->Draw("p0");
00135 
00136    // get fit parameters
00137    double parFit[4];
00138    for (int i = 0; i <4; ++i) 
00139       parFit[i] = min->GetParameter(i);
00140    
00141    // draw the fitted line
00142    int n = 1000;
00143    double t0 = 0;
00144    double dt = 10;
00145    TPolyLine3D *l = new TPolyLine3D(n);
00146    for (int i = 0; i <n;++i) {
00147       double t = t0+ dt*i/n;
00148       double x,y,z;
00149       line(t,parFit,x,y,z);
00150       l->SetPoint(i,x,y,z);
00151    }
00152    l->SetLineColor(kRed);
00153    l->Draw("same");
00154    
00155    // draw original line
00156    TPolyLine3D *l0 = new TPolyLine3D(n);
00157    for (int i = 0; i <n;++i) {
00158       double t = t0+ dt*i/n;
00159       double x,y,z;
00160       line(t,p0,x,y,z);
00161       l0->SetPoint(i,x,y,z);
00162    }
00163    l0->SetLineColor(kBlue);
00164    l0->Draw("same");
00165    
00166 }
00167 
00168 int main() { 
00169    line3Dfit();
00170 }

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