fit2dHist.C: Example to fit two histograms at the same time via TVirtualFitter | Fitting tutorials | fitExclude.C: Illustrate how to fit excluding points in a given range |
//Generate points distributed with some errors around a circle //Fit a circle through the points and draw //To run the script, do, eg // root > .x fitCircle.C (10000 points by default) // root > .x fitCircle.C(100); (with only 100 points // root > .x fitCircle.C(100000); with ACLIC // //Author: Rene Brun #include "TCanvas.h" #include "TRandom3.h" #include "TGraph.h" #include "TMath.h" #include "TArc.h" #include "TVirtualFitter.h" TGraph *gr; //____________________________________________________________________ void myfcn(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t) { //minimisation function computing the sum of squares of residuals Int_t np = gr->GetN(); f = 0; Double_t *x = gr->GetX(); Double_t *y = gr->GetY(); for (Int_t i=0;i<np;i++) { Double_t u = x[i] - par[0]; Double_t v = y[i] - par[1]; Double_t dr = par[2] - TMath::Sqrt(u*u+v*v); f += dr*dr; } } //____________________________________________________________________ void fitCircle(Int_t n=10000) { //generates n points around a circle and fit them TCanvas *c1 = new TCanvas("c1","c1",600,600); c1->SetGrid(); gr = new TGraph(n); if (n> 999) gr->SetMarkerStyle(1); else gr->SetMarkerStyle(3); TRandom3 r; Double_t x,y; for (Int_t i=0;i<n;i++) { r.Circle(x,y,r.Gaus(4,0.3)); gr->SetPoint(i,x,y); } c1->DrawFrame(-5,-5,5,5); gr->Draw("p"); //Fit a circle to the graph points TVirtualFitter::SetDefaultFitter("Minuit"); //default is Minuit TVirtualFitter *fitter = TVirtualFitter::Fitter(0, 3); fitter->SetFCN(myfcn); fitter->SetParameter(0, "x0", 0, 0.1, 0,0); fitter->SetParameter(1, "y0", 0, 0.1, 0,0); fitter->SetParameter(2, "R", 1, 0.1, 0,0); Double_t arglist[1] = {0}; fitter->ExecuteCommand("MIGRAD", arglist, 0); //Draw the circle on top of the points TArc *arc = new TArc(fitter->GetParameter(0), fitter->GetParameter(1),fitter->GetParameter(2)); arc->SetLineColor(kRed); arc->SetLineWidth(4); arc->Draw(); } fitCircle.C:1 fitCircle.C:2 fitCircle.C:3 fitCircle.C:4 fitCircle.C:5 fitCircle.C:6 fitCircle.C:7 fitCircle.C:8 fitCircle.C:9 fitCircle.C:10 fitCircle.C:11 fitCircle.C:12 fitCircle.C:13 fitCircle.C:14 fitCircle.C:15 fitCircle.C:16 fitCircle.C:17 fitCircle.C:18 fitCircle.C:19 fitCircle.C:20 fitCircle.C:21 fitCircle.C:22 fitCircle.C:23 fitCircle.C:24 fitCircle.C:25 fitCircle.C:26 fitCircle.C:27 fitCircle.C:28 fitCircle.C:29 fitCircle.C:30 fitCircle.C:31 fitCircle.C:32 fitCircle.C:33 fitCircle.C:34 fitCircle.C:35 fitCircle.C:36 fitCircle.C:37 fitCircle.C:38 fitCircle.C:39 fitCircle.C:40 fitCircle.C:41 fitCircle.C:42 fitCircle.C:43 fitCircle.C:44 fitCircle.C:45 fitCircle.C:46 fitCircle.C:47 fitCircle.C:48 fitCircle.C:49 fitCircle.C:50 fitCircle.C:51 fitCircle.C:52 fitCircle.C:53 fitCircle.C:54 fitCircle.C:55 fitCircle.C:56 fitCircle.C:57 fitCircle.C:58 fitCircle.C:59 fitCircle.C:60 fitCircle.C:61 fitCircle.C:62 fitCircle.C:63 fitCircle.C:64 fitCircle.C:65 fitCircle.C:66 fitCircle.C:67 fitCircle.C:68 fitCircle.C:69 fitCircle.C:70 |
|