Loading [MathJax]/extensions/tex2jax.js
Logo ROOT   6.12/07
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
fit1.C File Reference

Detailed Description

View in nbviewer Open in SWAN Simple fitting example (1-d histogram with an interpreted function)

pict1_fit1.C.png
Processing /mnt/build/workspace/root-makedoc-v612/rootspi/rdoc/src/v6-12-00-patches/tutorials/fit/fit1.C...
TFile** fillrandom.root
TFile* fillrandom.root
KEY: TFormula form1;1 abs(sin(x)/x)
KEY: TF1 sqroot;1 x*gaus(0) + [3]*form1
KEY: TH1F h1f;1 Test random numbers
Formula based function: sqroot
sqroot : x*gaus(0) + [3]*form1 Ndim= 1, Npar= 4, Number= 0
Formula expression:
x*[p0]*exp(-0.5*((x-[p1])/[p2])*((x-[p1])/[p2]))+[p3]*(abs(sin(x)/x))
FCN=208.946 FROM MIGRAD STATUS=CONVERGED 150 CALLS 151 TOTAL
EDM=1.0681e-06 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 3.28169e+01 5.39754e-01 3.06400e-03 -2.20763e-03
2 p1 4.00921e+00 1.71703e-02 9.79184e-05 -3.06960e-02
3 p2 9.84098e-01 1.30567e-02 6.17324e-05 -4.04514e-02
4 p3 6.46663e+01 1.34272e+00 9.06115e-03 -8.26048e-04
fit1 : Real Time = 0.21 seconds Cpu Time = 0.21 seconds
#include "TCanvas.h"
#include "TFrame.h"
#include "TBenchmark.h"
#include "TString.h"
#include "TF1.h"
#include "TH1.h"
#include "TFile.h"
#include "TROOT.h"
#include "TError.h"
#include "TInterpreter.h"
#include "TSystem.h"
#include "TPaveText.h"
void fit1() {
TCanvas *c1 = new TCanvas("c1_fit1","The Fit Canvas",200,10,700,500);
c1->SetGridx();
c1->SetGridy();
c1->GetFrame()->SetFillColor(21);
c1->GetFrame()->SetBorderMode(-1);
gBenchmark->Start("fit1");
//
// We connect the ROOT file generated in a previous tutorial
// (see <a href="fillrandom.C.nbconvert.ipynb">Filling histograms with random numbers from a function</a>)
//
TString dir = gROOT->GetTutorialDir();
dir.Append("/fit/");
TFile *file = TFile::Open("fillrandom.root");
if (!file) {
gROOT->ProcessLine(Form(".x %s../hist/fillrandom.C",dir.Data()));
file = TFile::Open("fillrandom.root");
if (!file) return;
}
//
// The function "ls()" lists the directory contents of this file
//
file->ls();
//
// Get object "sqroot" from the file. Undefined objects are searched
// for using gROOT->FindObject("xxx"), e.g.:
// TF1 *sqroot = (TF1*) gROOT.FindObject("sqroot")
//
TF1 * sqroot = 0;
file->GetObject("sqroot",sqroot);
if (!sqroot){
Error("fit1.C","Cannot find object sqroot of type TF1\n");
return;
}
sqroot->Print();
//
// Now get and fit histogram h1f with the function sqroot
//
TH1F* h1f = 0;
file->GetObject("h1f",h1f);
if (!h1f){
Error("fit1.C","Cannot find object h1f of type TH1F\n");
return;
}
h1f->SetFillColor(45);
h1f->Fit("sqroot");
// We now annotate the picture by creating a PaveText object
// and displaying the list of commands in this macro
//
TPaveText * fitlabel = new TPaveText(0.6,0.4,0.9,0.75,"NDC");
fitlabel->SetTextAlign(12);
fitlabel->SetFillColor(42);
fitlabel->ReadFile(Form("%sfit1_C.txt",dir.Data()));
fitlabel->Draw();
c1->Update();
gBenchmark->Show("fit1");
}
Author
Rene Brun

Definition in file fit1.C.