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)

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=198.935 FROM MIGRAD STATUS=CONVERGED 148 CALLS 149 TOTAL
EDM=2.98567e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 3.31658e+01 5.45703e-01 3.00376e-03 -1.11540e-03
2 p1 4.00667e+00 1.65304e-02 9.48491e-05 -3.06425e-02
3 p2 9.84663e-01 1.28238e-02 6.05976e-05 -3.04244e-02
4 p3 6.34464e+01 1.33233e+00 8.77483e-03 -3.96109e-04
fit1 : Real Time = 0.24 seconds Cpu Time = 0.24 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);
c1->GetFrame()->SetBorderSize(5);
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");
}
R__EXTERN TBenchmark * gBenchmark
Definition TBenchmark.h:59
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:187
#define gROOT
Definition TROOT.h:406
char * Form(const char *fmt,...)
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void Start(const char *name)
Starts Benchmark with the specified name.
virtual void Show(const char *name)
Stops Benchmark name and Prints results.
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:213
virtual void Print(Option_t *option="") const
Print TNamed name and title.
Definition TF1.cxx:2881
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:3997
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition TH1.cxx:3892
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
TString & Append(const char *cs)
Definition TString.h:564
return c1
Definition legend1.C:41
Definition file.py:1
Definition fit1.py:1
Author
Rene Brun

Definition in file fit1.C.