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

Detailed Description

View in nbviewer Open in SWAN Get in memory an histogram from a root file and fit a user defined function.

Note that a user defined function must always be defined as in this example:

  • first parameter: array of variables (in this example only 1-dimension)
  • second parameter: array of parameters Note also that in case of user defined functions, one must set an initial value for each parameter.
FCN=36.7428 FROM MIGRAD STATUS=CONVERGED 91 CALLS 92 TOTAL
EDM=5.8659e-11 STRATEGY= 1 ERROR MATRIX UNCERTAINTY 2.1 per cent
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 Constant 7.97964e+02 6.80343e+00 -5.06883e-03 -2.15486e-06
2 Mean_value -8.07109e-05 7.34948e-03 -6.41911e-06 -1.19538e-03
3 Sigma 9.98753e-01 7.14444e-03 -1.19305e-06 -1.68591e-03
Integral of function = 1907.35
{
Double_t arg = 0;
if (par[2] != 0) arg = (x[0] - par[1])/par[2];
Double_t fitval = par[0]*TMath::Exp(-0.5*arg*arg);
return fitval;
}
void myfit()
{
TString dir = gROOT->GetTutorialDir();
dir.Append("/hsimple.C");
dir.ReplaceAll("/./","/");
if (!gInterpreter->IsLoaded(dir.Data())) gInterpreter->LoadMacro(dir.Data());
TFile *hsimpleFile = (TFile*)gROOT->ProcessLineFast("hsimple(1)");
if (!hsimpleFile) return;
TCanvas *c1 = new TCanvas("c1","the fit canvas",500,400);
TH1F *hpx = (TH1F*)hsimpleFile->Get("hpx");
// Creates a Root function based on function fitf above
TF1 *func = new TF1("fitf",fitf,-2,2,3);
// Sets initial values and parameter names
func->SetParameters(100,0,1);
func->SetParNames("Constant","Mean_value","Sigma");
// Fit histogram in range defined by function
hpx->Fit(func,"r");
// Gets integral of function between fit limits
printf("Integral of function = %g\n",func->Integral(-2,2));
}
double Double_t
Definition RtypesCore.h:59
#define gInterpreter
#define gROOT
Definition TROOT.h:404
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:213
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition TF1.cxx:2521
virtual void SetParameters(const Double_t *params)
Definition TF1.h:644
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:3482
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
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
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:692
TString & Append(const char *cs)
Definition TString.h:564
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
Double_t Exp(Double_t x)
Definition TMath.h:677
Author
Rene Brun

Definition in file myfit.C.