Logo ROOT  
Reference Guide
Ifit.C File Reference

Detailed Description

View in nbviewer Open in SWAN Example of a program to fit non-equidistant data points

The fitting function fcn is a simple chisquare function The data consists of 5 data points (arrays x,y,z) + the errors in errorsz More details on the various functions or parameters for these functions can be obtained in an interactive ROOT session with:

Root > TMinuit *minuit = new TMinuit(10);
Implementation in C++ of the Minuit package written by Fred James.
Definition: TMinuit.h:27
const char * Root
Definition: TXMLSetup.cxx:46
Root > minuit->mnhelp("*") to see the list of possible keywords
Root > minuit->mnhelp("SET") explains most parameters
**********
** 1 **SET ERR 1
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 a1 3.00000e+00 1.00000e-01 no limits
2 a2 1.00000e+00 1.00000e-01 no limits
3 a3 1.00000e-01 1.00000e-02 no limits
4 a4 1.00000e-02 1.00000e-03 no limits
**********
** 2 **MIGRAD 500 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=90047.1 FROM MIGRAD STATUS=INITIATE 14 CALLS 15 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a1 3.00000e+00 1.00000e-01 1.00000e-01 2.81614e+05
2 a2 1.00000e+00 1.00000e-01 1.00000e-01 -2.73395e+05
3 a3 1.00000e-01 1.00000e-02 1.00000e-02 -3.08505e+05
4 a4 1.00000e-02 1.00000e-03 1.00000e-03 3.53925e+05
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
EIGENVALUES OF SECOND-DERIVATIVE MATRIX:
-6.7573e-05 7.2964e-05 4.1376e-02 3.9586e+00
MINUIT WARNING IN HESSE
============== MATRIX FORCED POS-DEF BY ADDING 0.004026 TO DIAGONAL.
FCN=10.3986 FROM HESSE STATUS=NOT POSDEF 23 CALLS 112 TOTAL
EDM=0.11462 STRATEGY= 1 ERR MATRIX NOT POS-DEF
EXT PARAMETER APPROXIMATE STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a1 2.51776e+00 4.70562e-02 6.10330e-06 1.74962e+00
2 a2 1.39551e+00 8.92750e-02 1.31410e-05 6.70438e+00
3 a3 2.53639e-01 9.47120e-02 1.16109e-05 4.39325e+00
4 a4 6.07417e-02 4.27988e-02 9.90333e-06 -3.17219e+00
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=5.83333 FROM MIGRAD STATUS=CONVERGED 232 CALLS 233 TOTAL
EDM=0.000805786 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a1 2.15261e+00 1.03157e-01 3.03942e-06 4.51168e+00
2 a2 8.11244e-01 2.53808e-01 5.62214e-06 -2.40912e+00
3 a3 1.71311e-01 4.00479e-01 4.93858e-06 -2.98641e+00
4 a4 1.01582e-01 1.60309e-01 4.18094e-06 3.54102e+00
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 4 ERR DEF=1
1.064e-02 1.057e-02 1.068e-02 2.336e-03
1.057e-02 6.442e-02 -7.900e-02 -3.414e-02
1.068e-02 -7.900e-02 1.604e-01 6.356e-02
2.336e-03 -3.414e-02 6.356e-02 2.570e-02
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3 4
1 0.99974 1.000 0.404 0.259 0.141
2 0.99985 0.404 1.000 -0.777 -0.839
3 0.99995 0.259 -0.777 1.000 0.990
4 0.99980 0.141 -0.839 0.990 1.000
#include "TMinuit.h"
Float_t z[5],x[5],y[5],errorz[5];
//______________________________________________________________________________
Double_t func(float x,float y,Double_t *par)
{
Double_t value=( (par[0]*par[0])/(x*x)-1)/ ( par[1]+par[2]*y-par[3]*y*y);
return value;
}
//______________________________________________________________________________
void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
{
const Int_t nbins = 5;
Int_t i;
//calculate chisquare
Double_t chisq = 0;
Double_t delta;
for (i=0;i<nbins; i++) {
delta = (z[i]-func(x[i],y[i],par))/errorz[i];
chisq += delta*delta;
}
f = chisq;
}
//______________________________________________________________________________
void Ifit()
{
// The z values
z[0]=1;
z[1]=0.96;
z[2]=0.89;
z[3]=0.85;
z[4]=0.78;
// The errors on z values
Float_t error = 0.01;
errorz[0]=error;
errorz[1]=error;
errorz[2]=error;
errorz[3]=error;
errorz[4]=error;
// the x values
x[0]=1.5751;
x[1]=1.5825;
x[2]=1.6069;
x[3]=1.6339;
x[4]=1.6706;
// the y values
y[0]=1.0642;
y[1]=0.97685;
y[2]=1.13168;
y[3]=1.128654;
y[4]=1.44016;
TMinuit *gMinuit = new TMinuit(5); //initialize TMinuit with a maximum of 5 params
gMinuit->SetFCN(fcn);
Double_t arglist[10];
Int_t ierflg = 0;
arglist[0] = 1;
gMinuit->mnexcm("SET ERR", arglist ,1,ierflg);
// Set starting values and step sizes for parameters
static Double_t vstart[4] = {3, 1 , 0.1 , 0.01};
static Double_t step[4] = {0.1 , 0.1 , 0.01 , 0.001};
gMinuit->mnparm(0, "a1", vstart[0], step[0], 0,0,ierflg);
gMinuit->mnparm(1, "a2", vstart[1], step[1], 0,0,ierflg);
gMinuit->mnparm(2, "a3", vstart[2], step[2], 0,0,ierflg);
gMinuit->mnparm(3, "a4", vstart[3], step[3], 0,0,ierflg);
// Now ready for minimization step
arglist[0] = 500;
arglist[1] = 1.;
gMinuit->mnexcm("MIGRAD", arglist ,2,ierflg);
// Print results
Double_t amin,edm,errdef;
Int_t nvpar,nparx,icstat;
gMinuit->mnstat(amin,edm,errdef,nvpar,nparx,icstat);
//gMinuit->mnprin(3,amin);
}
#define f(i)
Definition: RSha256.hxx:104
int Int_t
Definition: RtypesCore.h:43
double Double_t
Definition: RtypesCore.h:57
float Float_t
Definition: RtypesCore.h:55
R__EXTERN TMinuit * gMinuit
Definition: TMinuit.h:271
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
To set the address of the minimization function.
Definition: TMinuit.cxx:929
virtual void mnexcm(const char *comand, Double_t *plist, Int_t llist, Int_t &ierflg)
Interprets a command and takes appropriate action.
Definition: TMinuit.cxx:2673
virtual void mnstat(Double_t &fmin, Double_t &fedm, Double_t &errdef, Int_t &npari, Int_t &nparx, Int_t &istat)
Returns concerning the current status of the minimization.
Definition: TMinuit.cxx:7647
virtual void mnparm(Int_t k, TString cnamj, Double_t uk, Double_t wk, Double_t a, Double_t b, Int_t &ierflg)
Implements one parameter definition.
Definition: TMinuit.cxx:5676
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
Author
Rene Brun

Definition in file Ifit.C.