Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fit1.C File Reference

Detailed Description

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

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))
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 198.935
NDf = 190
Edm = 1.49283e-07
NCalls = 149
p0 = 33.1658 +/- 0.545703
p1 = 4.00667 +/- 0.0165304
p2 = 0.984663 +/- 0.0128238
p3 = 63.4464 +/- 1.33233
#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);
c1->GetFrame()->SetBorderSize(5);
// (for more details, see
// <a href="hist001_TH1_fillrandom_userfunc.C.nbconvert.ipynb">filling histograms with random numbers from a function</a>)
TFormula *form1 = new TFormula("form1", "abs(sin(x)/x)");
TF1 *sqroot = new TF1("sqroot", "x*gaus(0) + [3]*form1", 0.0, 10.0);
sqroot->SetLineColor(4);
sqroot->SetLineWidth(6);
// Set parameters to the functions "gaus" and "form1".
sqroot->SetParameters(10.0, 4.0, 1.0, 20.0);
sqroot->Print();
TH1D *h1d = new TH1D("h1d", "Test random numbers", 200, 0.0, 10.0);
h1d->FillRandom("sqroot", 10000);
//
// Now fit histogram h1d with the function sqroot
//
h1d->SetFillColor(45);
h1d->Fit("sqroot");
// We now annotate the picture by creating a PaveText object
// and displaying the list of commands in this macro
//
TString dir = gROOT->GetTutorialDir();
dir.Append("/math/fit/");
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();
}
#define gROOT
Definition TROOT.h:406
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
The Formula class.
Definition TFormula.h:89
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:693
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
TString & Append(const char *cs)
Definition TString.h:572
return c1
Definition legend1.C:41
Definition fit1.py:1
Author
Rene Brun

Definition in file fit1.C.