{ "cells": [ { "cell_type": "markdown", "id": "258b92c7", "metadata": {}, "source": [ "# line3Dfit\n", "Fitting of a TGraph2D with a 3D straight line\n", "\n", "run this macro by doing:\n", "\n", "```cpp\n", "root>.x line3Dfit.C+\n", "```\n", "\n", "\n", "\n", "\n", "**Author:** Lorenzo Moneta \n", "This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, May 19, 2026 at 08:25 PM." ] }, { "cell_type": "code", "execution_count": null, "id": "9c37dde9", "metadata": { "collapsed": false }, "outputs": [], "source": [ "%%cpp -d\n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "#include \n", "\n", "#include \n", "\n", "using namespace ROOT::Math;\n", "\n", "\n", "%%cpp -d\n", "bool first = true;" ] }, { "cell_type": "markdown", "id": "0571e4b2", "metadata": {}, "source": [ "function Object to be minimized" ] }, { "cell_type": "code", "execution_count": null, "id": "d1388458", "metadata": { "collapsed": false }, "outputs": [], "source": [ "struct SumDistance2 {\n", " // the TGraph is a data member of the object\n", " TGraph2D *fGraph;\n", "\n", " SumDistance2(TGraph2D *g) : fGraph(g) {}\n", "\n", " // calculate distance line-point\n", " double distance2(double x,double y,double z, const double *p) {\n", " // distance line point is D= | (xp-x0) cross ux |\n", " // where ux is direction of line and x0 is a point in the line (like t = 0)\n", " XYZVector xp(x,y,z);\n", " XYZVector x0(p[0], p[2], 0. );\n", " XYZVector x1(p[0] + p[1], p[2] + p[3], 1. );\n", " XYZVector u = (x1-x0).Unit();\n", " double d2 = ((xp-x0).Cross(u)).Mag2();\n", " return d2;\n", " }\n", "\n", " // implementation of the function to be minimized\n", " double operator() (const double *par) {\n", " assert(fGraph != nullptr);\n", " double * x = fGraph->GetX();\n", " double * y = fGraph->GetY();\n", " double * z = fGraph->GetZ();\n", " int npoints = fGraph->GetN();\n", " double sum = 0;\n", " for (int i = 0; i < npoints; ++i) {\n", " double d = distance2(x[i],y[i],z[i],par);\n", " sum += d;\n", " }\n", " if (first) {\n", " std::cout << \"Total Initial distance square = \" << sum << std::endl;\n", " }\n", " first = false;\n", " return sum;\n", " }\n", "\n", "};" ] }, { "cell_type": "markdown", "id": "59f468c9", "metadata": {}, "source": [ " define the parametric line equation\n", " " ] }, { "cell_type": "code", "execution_count": null, "id": "98ac044e", "metadata": { "collapsed": false }, "outputs": [], "source": [ "%%cpp -d\n", "void line(double t, const double *p, double &x, double &y, double &z) {\n", " // a parametric line is define from 6 parameters but 4 are independent\n", " // x0,y0,z0,z1,y1,z1 which are the coordinates of two points on the line\n", " // can choose z0 = 0 if line not parallel to x-y plane and z1 = 1;\n", " x = p[0] + p[1]*t;\n", " y = p[2] + p[3]*t;\n", " z = t;\n", "}" ] }, { "cell_type": "code", "execution_count": null, "id": "86e71891", "metadata": { "collapsed": false }, "outputs": [], "source": [ "gStyle->SetOptStat(0);\n", "gStyle->SetOptFit();" ] }, { "cell_type": "markdown", "id": "8767e40c", "metadata": {}, "source": [ "double e = 0.1;" ] }, { "cell_type": "code", "execution_count": null, "id": "48f99a81", "metadata": { "collapsed": false }, "outputs": [], "source": [ "int nd = 10000;" ] }, { "cell_type": "markdown", "id": "90125fac", "metadata": {}, "source": [ "double xmin = 0; double ymin = 0;\n", "double xmax = 10; double ymax = 10;" ] }, { "cell_type": "code", "execution_count": null, "id": "7746dc8b", "metadata": { "collapsed": false }, "outputs": [], "source": [ "TGraph2D * gr = new TGraph2D();" ] }, { "cell_type": "markdown", "id": "ee3909f7", "metadata": {}, "source": [ "Fill the 2D graph" ] }, { "cell_type": "code", "execution_count": null, "id": "586a9b22", "metadata": { "collapsed": false }, "outputs": [], "source": [ "double p0[4] = {10,20,1,2};" ] }, { "cell_type": "markdown", "id": "c794e107", "metadata": {}, "source": [ "generate graph with the 3d points" ] }, { "cell_type": "code", "execution_count": null, "id": "36d17d8a", "metadata": { "collapsed": false }, "outputs": [], "source": [ "for (int N=0; NUniform(0,10);\n", " line(t,p0,x,y,z);\n", " double err = 1;\n", " // do a Gaussian smearing around the points in all coordinates\n", " x += gRandom->Gaus(0,err);\n", " y += gRandom->Gaus(0,err);\n", " z += gRandom->Gaus(0,err);\n", " gr->SetPoint(N,x,y,z);\n", " //dt->SetPointError(N,0,0,err);\n", "}" ] }, { "cell_type": "markdown", "id": "52ad2c37", "metadata": {}, "source": [ "fit the graph now" ] }, { "cell_type": "code", "execution_count": null, "id": "7bbcb75d", "metadata": { "collapsed": false }, "outputs": [], "source": [ "ROOT::Fit::Fitter fitter;" ] }, { "cell_type": "markdown", "id": "71ecf9dd", "metadata": {}, "source": [ "make the functor objet" ] }, { "cell_type": "code", "execution_count": null, "id": "b5179f51", "metadata": { "collapsed": false }, "outputs": [], "source": [ "SumDistance2 sdist(gr);\n", "ROOT::Math::Functor fcn(sdist,4);" ] }, { "cell_type": "markdown", "id": "6fe6d416", "metadata": {}, "source": [ "set the function and the initial parameter values" ] }, { "cell_type": "code", "execution_count": null, "id": "fe914877", "metadata": { "collapsed": false }, "outputs": [], "source": [ "double pStart[4] = {1,1,1,1};\n", "fitter.SetFCN(fcn,pStart);" ] }, { "cell_type": "markdown", "id": "e2e7686c", "metadata": {}, "source": [ "set step sizes different than default ones (0.3 times parameter values)" ] }, { "cell_type": "code", "execution_count": null, "id": "1955793c", "metadata": { "collapsed": false }, "outputs": [], "source": [ "for (int i = 0; i < 4; ++i) fitter.Config().ParSettings(i).SetStepSize(0.01);\n", "\n", "bool ok = fitter.FitFCN();\n", "if (!ok) {\n", " Error(\"line3Dfit\",\"Line3D Fit failed\");\n", " return 1;\n", "}\n", "\n", "const ROOT::Fit::FitResult & result = fitter.Result();\n", "\n", "std::cout << \"Total final distance square \" << result.MinFcnValue() << std::endl;\n", "result.Print(std::cout);\n", "\n", "\n", "gr->Draw(\"p0\");" ] }, { "cell_type": "markdown", "id": "f91b0159", "metadata": {}, "source": [ "get fit parameters" ] }, { "cell_type": "code", "execution_count": null, "id": "a1e43f62", "metadata": { "collapsed": false }, "outputs": [], "source": [ "const double * parFit = result.GetParams();" ] }, { "cell_type": "markdown", "id": "1192eb31", "metadata": {}, "source": [ "draw the fitted line" ] }, { "cell_type": "code", "execution_count": null, "id": "2a7833a9", "metadata": { "collapsed": false }, "outputs": [], "source": [ "int n = 1000;\n", "double t0 = 0;\n", "double dt = 10;\n", "TPolyLine3D *l = new TPolyLine3D(n);\n", "for (int i = 0; i SetPoint(i,x,y,z);\n", "}\n", "l->SetLineColor(kRed);\n", "l->Draw(\"same\");" ] }, { "cell_type": "markdown", "id": "d81b9015", "metadata": {}, "source": [ "draw original line" ] }, { "cell_type": "code", "execution_count": null, "id": "51bb7d86", "metadata": { "collapsed": false }, "outputs": [], "source": [ "TPolyLine3D *l0 = new TPolyLine3D(n);\n", "for (int i = 0; i SetPoint(i,x,y,z);\n", "}\n", "l0->SetLineColor(kBlue);\n", "l0->Draw(\"same\");\n", "return 0;" ] }, { "cell_type": "markdown", "id": "c94bb9ab", "metadata": {}, "source": [ "Draw all canvases " ] }, { "cell_type": "code", "execution_count": null, "id": "62d7e147", "metadata": { "collapsed": false }, "outputs": [], "source": [ "gROOT->GetListOfCanvases()->Draw()" ] } ], "metadata": { "kernelspec": { "display_name": "ROOT C++", "language": "c++", "name": "root" }, "language_info": { "codemirror_mode": "text/x-c++src", "file_extension": ".C", "mimetype": " text/x-c++src", "name": "c++" } }, "nbformat": 4, "nbformat_minor": 5 }