Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fit1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_fit
3/// \notebook
4/// Simple fitting example (1-d histogram with an interpreted function)
5///
6/// \macro_image
7/// \macro_output
8/// \macro_code
9///
10/// \author Rene Brun
11
12#include "TCanvas.h"
13#include "TFrame.h"
14#include "TBenchmark.h"
15#include "TString.h"
16#include "TF1.h"
17#include "TH1.h"
18#include "TFile.h"
19#include "TROOT.h"
20#include "TError.h"
21#include "TInterpreter.h"
22#include "TSystem.h"
23#include "TPaveText.h"
24
25void fit1() {
26 TCanvas *c1 = new TCanvas("c1_fit1","The Fit Canvas",200,10,700,500);
27 c1->SetGridx();
28 c1->SetGridy();
29 c1->GetFrame()->SetFillColor(21);
30 c1->GetFrame()->SetBorderMode(-1);
31 c1->GetFrame()->SetBorderSize(5);
32
33 gBenchmark->Start("fit1");
34 //
35 // We connect the ROOT file generated in a previous tutorial
36 // (see <a href="fillrandom.C.nbconvert.ipynb">Filling histograms with random numbers from a function</a>)
37 //
38 TString dir = gROOT->GetTutorialDir();
39 dir.Append("/fit/");
40 TFile *file = TFile::Open("fillrandom.root");
41 if (!file) {
42 gROOT->ProcessLine(Form(".x %s../hist/fillrandom.C(0)",dir.Data()));
43 file = TFile::Open("fillrandom.root");
44 if (!file) return;
45 }
46
47 //
48 // The function "ls()" lists the directory contents of this file
49 //
50 file->ls();
51
52 //
53 // Get object "sqroot" from the file. Undefined objects are searched
54 // for using gROOT->FindObject("xxx"), e.g.:
55 // TF1 *sqroot = (TF1*) gROOT.FindObject("sqroot")
56 //
57 TF1 * sqroot = 0;
58 file->GetObject("sqroot",sqroot);
59 if (!sqroot){
60 Error("fit1.C","Cannot find object sqroot of type TF1\n");
61 return;
62 }
63 sqroot->Print();
64
65 //
66 // Now get and fit histogram h1f with the function sqroot
67 //
68 TH1F* h1f = 0;
69 file->GetObject("h1f",h1f);
70 if (!h1f){
71 Error("fit1.C","Cannot find object h1f of type TH1F\n");
72 return;
73 }
74 h1f->SetFillColor(45);
75 h1f->Fit("sqroot");
76
77 // We now annotate the picture by creating a PaveText object
78 // and displaying the list of commands in this macro
79 //
80 TPaveText * fitlabel = new TPaveText(0.6,0.4,0.9,0.75,"NDC");
81 fitlabel->SetTextAlign(12);
82 fitlabel->SetFillColor(42);
83 fitlabel->ReadFile(Form("%sfit1_C.txt",dir.Data()));
84 fitlabel->Draw();
85 c1->Update();
86 gBenchmark->Show("fit1");
87}
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:404
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:2887
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:4025
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:3893
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