Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
myfit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_fit
3/// \notebook -js
4/// Get in memory an histogram from a root file and fit a user defined function.
5/// Note that a user defined function must always be defined
6/// as in this example:
7/// - first parameter: array of variables (in this example only 1-dimension)
8/// - second parameter: array of parameters
9/// Note also that in case of user defined functions, one must set
10/// an initial value for each parameter.
11///
12/// \macro_image
13/// \macro_output
14/// \macro_code
15///
16/// \author Rene Brun
17
18#include <TCanvas.h>
19#include <TF1.h>
20#include <TFile.h>
21#include <TH1F.h>
22#include <TInterpreter.h>
23#include <TROOT.h>
24
25#include <cmath>
26
27double fitf(double *x, double *par)
28{
29 double arg = 0;
30 if (par[2] != 0) arg = (x[0] - par[1])/par[2];
31
32 double fitval = par[0]*std::exp(-0.5*arg*arg);
33 return fitval;
34}
35void myfit()
36{
37 TString dir = gROOT->GetTutorialDir();
38 dir.Append("/hsimple.C");
39 dir.ReplaceAll("/./","/");
40 if (!gInterpreter->IsLoaded(dir.Data())) gInterpreter->LoadMacro(dir.Data());
41 TFile *hsimpleFile = (TFile*)gROOT->ProcessLineFast("hsimple(1)");
42 if (!hsimpleFile) return;
43
44 TCanvas *c1 = new TCanvas("c1","the fit canvas",500,400);
45
46 TH1F *hpx = (TH1F*)hsimpleFile->Get("hpx");
47
48// Creates a Root function based on function fitf above
49 TF1 *func = new TF1("fitf",fitf,-2,2,3);
50
51// Sets initial values and parameter names
52 func->SetParameters(100,0,1);
53 func->SetParNames("Constant","Mean_value","Sigma");
54
55// Fit histogram in range defined by function
56 hpx->Fit(func,"r");
57
58// Gets integral of function between fit limits
59 printf("Integral of function = %g\n",func->Integral(-2,2));
60}
#define gInterpreter
#define gROOT
Definition TROOT.h:406
The Canvas class.
Definition TCanvas.h:23
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
1-Dim function class
Definition TF1.h:233
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition TF1.cxx:2531
virtual void SetParNames(const char *name0="", const char *name1="", const char *name2="", const char *name3="", const char *name4="", const char *name5="", const char *name6="", const char *name7="", const char *name8="", const char *name9="", const char *name10="")
Set up to 10 parameter names.
Definition TF1.cxx:3463
virtual void SetParameters(const Double_t *params)
Definition TF1.h:677
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:623
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:3906
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
TString & Append(const char *cs)
Definition TString.h:572
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17