Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SimpleFitting.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_r
3///
4/// Create an exponential fitting
5/// The idea is to create a set of numbers x,y with the function x^3 and some noise from ROOT,
6/// fit the function to get the exponent (which must be near 3) and plot the points with noise,
7/// the known function and the fitted function
8///
9/// \macro_code
10///
11/// \author Omar Zapata
12
13#include<TRInterface.h>
14#include<TRandom.h>
15
16TCanvas *SimpleFitting(){
17 TCanvas *c1 = new TCanvas("c1","Curve Fitting",700,500);
18 c1->SetGrid();
19
20 // draw a frame to define the range
21 TMultiGraph *mg = new TMultiGraph();
22
23 // create the first plot (points with gaussian noise)
24 const Int_t n = 24;
25 Double_t x1[n] ;
26 Double_t y1[n] ;
27 //Generate the points along a X^3 with noise
28 TRandom rg;
29 rg.SetSeed(520);
30 for (Int_t i = 0; i < n; i++) {
31 x1[i] = rg.Uniform(0, 1);
32 y1[i] = TMath::Power(x1[i], 3) + rg.Gaus() * 0.06;
33 }
34
35 TGraph *gr1 = new TGraph(n,x1,y1);
37 gr1->SetMarkerStyle(8);
38 gr1->SetMarkerSize(1);
39 mg->Add(gr1);
40
41 // create the second plot
42 TF1 *f_known=new TF1("f_known","pow(x,3)",0,1);
43 TGraph *gr2 = new TGraph(f_known);
44 gr2->SetMarkerColor(kRed);
45 gr2->SetMarkerStyle(8);
46 gr2->SetMarkerSize(1);
47 mg->Add(gr2);
48 //passing data to Rfot fitting
50 r["x"]<<TVectorD(n, x1);
51 r["y"]<<TVectorD(n, y1);
52 //creating a R data frame
53 r<<"ds<-data.frame(x=x,y=y)";
54 //fitting x and y to X^power using Nonlinear Least Squares
55 r<<"m <- nls(y ~ I(x^power),data = ds, start = list(power = 1),trace = T)";
56 //getting the exponent
57 Double_t power;
58 r["summary(m)$coefficients[1]"]>>power;
59
60 TF1 *f_fitted=new TF1("f_fitted","pow(x,[0])",0,1);
61 f_fitted->SetParameter(0,power);
62 //plotting the fitted function
63 TGraph *gr3 = new TGraph(f_fitted);
65 gr3->SetMarkerStyle(8);
66 gr3->SetMarkerSize(1);
67
68 mg->Add(gr3);
69 mg->Draw("ap");
70
71 //displaying basic results
72 TPaveText *pt = new TPaveText(0.1,0.6,0.5,0.9,"brNDC");
73 pt->SetFillColor(18);
74 pt->SetTextAlign(12);
75 pt->AddText("Fitting x^power ");
76 pt->AddText(" \"Blue\" Points with gaussian noise to be fitted");
77 pt->AddText(" \"Red\" Known function x^3");
78 TString fmsg;
79 fmsg.Form(" \"Green\" Fitted function with power=%.4lf",power);
80 pt->AddText(fmsg);
81 pt->Draw();
82 c1->Update();
83 return c1;
84}
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y1
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
ROOT R was implemented using the R Project library and the modules Rcpp and RInside
static TRInterface & Instance()
static method to get an TRInterface instance reference
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:42
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:660
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
void Draw(Option_t *chopt="") override
Draw this multigraph with its current attributes.
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
void Draw(Option_t *option="") override
Draw this pavetext with its current attributes.
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
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:275
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:682
Basic string class.
Definition TString.h:139
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2356
TPaveText * pt
return c1
Definition legend1.C:41
const Int_t n
Definition legend1.C:16
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721