ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
myfit.C
Go to the documentation of this file.
1 // Get in memory an histogram from a root file and fit a user defined function.
2 // Note that a user defined function must always be defined
3 // as in this example:
4 // - first parameter: array of variables (in this example only 1-dimension)
5 // - second parameter: array of parameters
6 // Note also that in case of user defined functions, one must set
7 // an initial value for each parameter.
8 //Author: Rene Brun
9 
11 {
12  Double_t arg = 0;
13  if (par[2] != 0) arg = (x[0] - par[1])/par[2];
14 
15  Double_t fitval = par[0]*TMath::Exp(-0.5*arg*arg);
16  return fitval;
17 }
18 void myfit()
19 {
20  TString dir = gSystem->UnixPathName(__FILE__);
21  dir.ReplaceAll("myfit.C","../hsimple.C");
22  dir.ReplaceAll("/./","/");
23  if (!gInterpreter->IsLoaded(dir.Data())) gInterpreter->LoadMacro(dir.Data());
24  TFile *hsimple = (TFile*)gROOT->ProcessLineFast("hsimple(1)");
25  if (!hsimple) return;
26 
27  TCanvas *c1 = new TCanvas("c1","the fit canvas",500,400);
28 
29  TH1F *hpx = (TH1F*)hsimple->Get("hpx");
30 
31 // Creates a Root function based on function fitf above
32  TF1 *func = new TF1("fitf",fitf,-2,2,3);
33 
34 // Sets initial values and parameter names
35  func->SetParameters(100,0,1);
36  func->SetParNames("Constant","Mean_value","Sigma");
37 
38 // Fit histogram in range defined by function
39  hpx->Fit(func,"r");
40 
41 // Gets integral of function between fit limits
42  printf("Integral of function = %g\n",func->Integral(-2,2));
43 }
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
virtual void SetParNames(const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
Set up to 10 parameter names.
Definition: TF1.cxx:3147
TCanvas * c1
Definition: legend1.C:2
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
#define gROOT
Definition: TROOT.h:344
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
#define gInterpreter
Definition: TInterpreter.h:502
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
virtual const char * UnixPathName(const char *unixpathname)
Convert from a Unix pathname to a local pathname.
Definition: TSystem.cxx:1020
const char * Data() const
Definition: TString.h:349
Double_t x[n]
Definition: legend1.C:17
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
The Canvas class.
Definition: TCanvas.h:48
Double_t Exp(Double_t x)
Definition: TMath.h:495
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
void dir(char *path=0)
Definition: rootalias.C:30
double func(double *x, double *p)
Definition: stressTF1.cxx:213
Double_t fitf(Double_t *x, Double_t *par)
Definition: myfit.C:10
TFile * hsimple(Int_t get=0)
Definition: hsimple.C:13
1-Dim function class
Definition: TF1.h:149
void myfit()
Definition: myfit.C:18
TH1F * hpx
Definition: hcons.C:32
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:3607