FUMILI Based on ideas, proposed by I.N. Silin [See NIM A440, 2000 (p431)] converted from FORTRAN to C by Sergey Yaschenko <s.yaschenko@fz-juelich.de>
FUMILI is used to minimize Chi-square function or to search maximum of likelihood function.
Experimentally measured values $F_i$ are fitted with theoretical functions $f_i({\vec x}_i,\vec\theta\,\,)$, where ${\vec x}_i$ are coordinates, and $\vec\theta$ -- vector of parameters.
For better convergence Chi-square function has to be the following form
$$ {\chi^2\over2}={1\over2}\sum^n_{i=1}\left(f_i(\vec x_i,\vec\theta\,\,)-F_i\over\sigma_i\right)^2 \eqno(1) $$
where $\sigma_i$ are errors of measured function.
The minimum condition is
$$ {\partial\chi^2\over\partial\theta_i}=\sum^n_{j=1}{1\over\sigma^2_j}\cdot {\partial f_j\over\partial\theta_i}\left[f_j(\vec x_j,\vec\theta\,\,)-F_j\right]=0,\qquad i=1\ldots m\eqno(2) $$
where m is the quantity of parameters.
Expanding left part of (2) over parameter increments and retaining only linear terms one gets
$$ \left(\partial\chi^2\over\theta_i\right)_{\vec\theta={\vec\theta}^0} +\sum_k\left(\partial^2\chi^2\over\partial\theta_i\partial\theta_k\right)_{ \vec\theta={\vec\theta}^0}\cdot(\theta_k-\theta_k^0) = 0\eqno(3) $$
Here ${\vec\theta}_0$ is some initial value of parameters. In general case:
$$ {\partial^2\chi^2\over\partial\theta_i\partial\theta_k}= \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i} {\partial f_j\over\theta_k} + \sum^n_{j=1}{(f_j - F_j)\over\sigma^2_j}\cdot {\partial^2f_j\over\partial\theta_i\partial\theta_k}\eqno(4) $$
In FUMILI algorithm for second derivatives of Chi-square approximate expression is used when last term in (4) is discarded. It is often done, not always wittingly, and sometimes causes troubles, for example, if user wants to limit parameters with positive values by writing down $\theta_i^2$ instead of $\theta_i$. FUMILI will fail if one tries minimize $\chi^2 = g^2(\vec\theta)$ where g is arbitrary function.
Approximate value is:
$${\partial^2\chi^2\over\partial\theta_i\partial\theta_k}\approx Z_{ik}= \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i} {\partial f_j\over\theta_k}\eqno(5) $$
Then the equations for parameter increments are
$$\left(\partial\chi^2\over\partial\theta_i\right)_{\vec\theta={\vec\theta}^0} +\sum_k Z_{ik}\cdot(\theta_k-\theta^0_k) = 0, \qquad i=1\ldots m\eqno(6) $$
Remarkable feature of algorithm is the technique for step restriction. For an initial value of parameter ${\vec\theta}^0$ a parallelepiped $P_0$ is built with the center at ${\vec\theta}^0$ and axes parallel to coordinate axes $\theta_i$. The lengths of parallelepiped sides along i-th axis is $2b_i$, where $b_i$ is such a value that the functions $f_j(\vec\theta)$ are quasi-linear all over the parallelepiped.
FUMILI takes into account simple linear inequalities in the form: $$ \theta_i^{\rm min}\le\theta_i\le\theta^{\rm max}_i\eqno(7) $$
They form parallelepiped $P$ ($P_0$ may be deformed by $P$). Very similar step formulae are used in FUMILI for negative logarithm of the likelihood function with the same idea - linearization of function argument.
virtual void | TObject::DoError(int level, const char* location, const char* fmt, va_list va) const |
void | TObject::MakeZombie() |
static TObject::(anonymous) | TObject::kBitMask | |
static TObject::EStatusBits | TObject::kCanDelete | |
static TObject::EStatusBits | TObject::kCannotPick | |
static TObject::EStatusBits | TObject::kHasUUID | |
static TObject::EStatusBits | TObject::kInvalidObject | |
static TObject::(anonymous) | TObject::kIsOnHeap | |
static TObject::EStatusBits | TObject::kIsReferenced | |
static TObject::EStatusBits | TObject::kMustCleanup | |
static TObject::EStatusBits | TObject::kNoContextMenu | |
static TObject::(anonymous) | TObject::kNotDeleted | |
static TObject::EStatusBits | TObject::kObjInCanvas | |
static TObject::(anonymous) | TObject::kOverwrite | |
static TObject::(anonymous) | TObject::kSingleKey | |
static TObject::(anonymous) | TObject::kWriteDelete | |
static TObject::(anonymous) | TObject::kZombie |
Double_t* | TVirtualFitter::fCache | [fCacheSize] array of points data (fNpoints*fPointSize < fCacheSize words) |
Int_t | TVirtualFitter::fCacheSize | Size of the fCache array |
void(*)(Int_t&,Double_t*,Double_t&,Double_t*,Int_t) | TVirtualFitter::fFCN | |
TMethodCall* | TVirtualFitter::fMethodCall | Pointer to MethodCall in case of interpreted function |
TString | TNamed::fName | object identifier |
Int_t | TVirtualFitter::fNpoints | Number of points to fit |
TObject* | TVirtualFitter::fObjectFit | pointer to object being fitted |
Foption_t | TVirtualFitter::fOption | struct with the fit options |
Int_t | TVirtualFitter::fPointSize | Number of words per point in the cache |
TString | TNamed::fTitle | object title |
TObject* | TVirtualFitter::fUserFunc | pointer to user theoretical function (a TF1*) |
Int_t | TVirtualFitter::fXfirst | first bin on X axis |
Int_t | TVirtualFitter::fXlast | last bin on X axis |
Int_t | TVirtualFitter::fYfirst | first bin on Y axis |
Int_t | TVirtualFitter::fYlast | last bin on Y axis |
Int_t | TVirtualFitter::fZfirst | first bin on Z axis |
Int_t | TVirtualFitter::fZlast | last bin on Z axis |
static TString | TVirtualFitter::fgDefault | name of the default fitter ("Minuit","Fumili",etc) |
static Double_t | TVirtualFitter::fgErrorDef | Error definition (default=1) |
static TVirtualFitter* | TVirtualFitter::fgFitter | Current fitter (default TFitter) |
static Int_t | TVirtualFitter::fgMaxiter | Maximum number of iterations |
static Int_t | TVirtualFitter::fgMaxpar | Maximum number of fit parameters for current fitter |
static Double_t | TVirtualFitter::fgPrecision | maximum precision |
Double_t* | fA | [fMaxParam] Fit parameter array |
Double_t | fAKAPPA | |
Double_t* | fAMN | [fMaxParam] Minimum param value |
Double_t* | fAMX | [fMaxParam] Maximum param value |
TString* | fANames | [fMaxParam] Parameter names |
Double_t* | fCmPar | [fMaxParam] parameters of commands |
TString | fCword | Command string |
Double_t* | fDA | [fMaxParam] Parameter step |
Bool_t | fDEBUG | debug info |
Double_t* | fDF | [fMaxParam] // First derivatives of theoretical function |
Int_t | fENDFLG | End flag of fit |
Double_t | fEPS | fEPS - required precision of parameters. If fEPS<0 then |
Double_t* | fEXDA | [fNED12] experimental data poInt_ter |
Bool_t | fGRAD | user calculated gradients |
Double_t | fGT | Expected function change in next iteration |
Double_t* | fGr | [fMaxParam] Gradients of objective function |
Int_t | fINDFLG[5] | internal flags; |
Int_t | fLastFixed | Last fixed parameter number |
Bool_t | fLogLike | LogLikelihood flag |
Int_t | fMaxParam | |
Int_t | fNED1 | Number of experimental vectors X=(x1,x2,...xK) |
Int_t | fNED12 | fNED1+fNED2 |
Int_t | fNED2 | K - Length of vector X plus 2 (for chi2) |
Int_t | fNfcn | Number of FCN calls; |
Int_t | fNlimMul | fNlimMul - after fNlimMul successful iterations permits four-fold increasing of fPL |
Int_t | fNlog | |
Int_t | fNmaxIter | fNmaxIter - maximum number of iterations |
Int_t | fNpar | fNpar - number of parameters |
Int_t | fNstepDec | fNstepDec - maximum number of step decreasing counter |
Bool_t | fNumericDerivatives | |
Double_t* | fPL | [fMaxParam] Limits for parameters step. If <0, then parameter is fixed |
Double_t* | fPL0 | [fMaxParam] Step initial bounds |
Double_t* | fParamError | [fMaxParam] Parameter errors |
Double_t* | fR | [fMaxParam] Correlation factors |
Double_t | fRP | Precision of fit ( machine zero on CDC 6000) quite old yeh? |
Double_t | fS | fS - objective function value (return) |
Double_t* | fSumLog | [fNlog] |
Bool_t | fWARN | warnings |
Double_t* | fZ | [fMaxParam2] Invers fZ0 matrix - covariance matrix |
Double_t* | fZ0 | [fMaxParam2] Matrix of approximate second derivatives of objective function |
Resets all parameter names, values and errors to zero Argument opt is ignored NB: this procedure doesn't reset parameter limits
Calculates partial derivatives of theoretical function Input: fX - vector of data point Output: DF - array of derivatives ARITHM.F Converted from CERNLIB
Evaluate the minimisation function Input parameters: npar: number of currently variable parameters par: array of (constant and variable) parameters flag: Indicates what is to be calculated grad: array of gradients Output parameters: fval: The calculated function value. grad: The vector of first derivatives. The meaning of the parameters par is of course defined by the user, who uses the values of those parameters to calculate their function value. The starting values must be specified by the user. Inside FCN user has to define Z-matrix by means TFumili::GetZ and TFumili::Derivatives, set theoretical function by means of TFumili::SetUserFunc, but first - pass number of parameters by TFumili::SetParNumber Later values are determined by Fumili as it searches for the minimum or performs whatever analysis is requested by the user. The default function calls the function specified in SetFCN
Evaluate theoretical function df: array of partial derivatives X: vector of theoretical function argument
Execute MINUIT commands. MINImize, SIMplex, MIGrad and FUMili all will call TFumili::Minimize method. For full command list see MINUIT. Reference Manual. CERN Program Library Long Writeup D506. Improvement and errors calculation are not yet implemented as well as Monte-Carlo seeking and minimization. Contour commands are also unsupported. command : command string args : array of arguments nargs : number of arguments
Called from TFumili::ExecuteCommand in case of "SET xxx" and "SHOW xxx".
return element i,j from the covariance matrix
Get various ipar parameter attributs: cname: parameter name value: parameter value verr: parameter error vlow: lower limit vhigh: upper limit WARNING! parname must be suitably dimensionned in the calling function.
Return errors after MINOs not implemented
return global fit parameters amin : chisquare edm : estimated distance to minimum errdef nvpar : number of variable parameters nparx : total number of parameters
Inverts packed diagonal matrix Z by square-root method. Matrix elements corresponding to fix parameters are removed. n: number of variable parameters
* FUMILI Based on ideas, proposed by I.N. Silin [See NIM A440, 2000 (p431)] converted from FORTRAN to C by Sergey Yaschenko <s.yaschenko@fz-juelich.de> * This function is called after setting theoretical function by means of TFumili::SetUserFunc and initializing parameters. Optionally one can set FCN function (see TFumili::SetFCN and TFumili::Eval) If FCN is undefined then user has to provide data arrays by calling TFumili::SetData procedure. TFumili::Minimize return following values: 0 - fit is converged -2 - function is not decreasing (or bad derivatives) -3 - error estimations are infinite -4 - maximum number of iterations is exceeded
Prints fit results. ikode is the type of printing parameters p is function value ikode = 1 - print values, errors and limits ikode = 2 - print values, errors and steps ikode = 3 - print values, errors, steps and derivatives ikode = 4 - print only values and errors
Sets pointer to data array provided by user. Necessary if SetFCN is not called. numpoints: number of experimental points vecsize: size of data point vector + 2 (for N-dimensional fit vecsize=N+2) exdata: data array with following format exdata[0] = ExpValue_0 - experimental data value number 0 exdata[1] = ExpSigma_0 - error of value number 0 exdata[2] = X_0[0] exdata[3] = X_0[1] exdata[vecsize-1] = X_0[vecsize-3] exdata[vecsize] = ExpValue_1 exdata[vecsize+1] = ExpSigma_1 exdata[vecsize+2] = X_1[0] exdata[vecsize*(numpoints-1)] = ExpValue_(numpoints-1) exdata[vecsize*numpoints-1] = X_(numpoints-1)[vecsize-3]
Sets for prameter number ipar initial parameter value, name parname, initial error verr and limits vlow and vhigh If vlow = vhigh but not equil to zero, parameter will be fixed. If vlow = vhigh = 0, parameter is released and its limits are discarded
Minimization function for H1s using a Chisquare method Default method (function evaluated at center of bin) for each point the cache contains the following info -1D : bc,e,xc (bin content, error, x of center of bin) -2D : bc,e,xc,yc -3D : bc,e,xc,yc,zc
Minimization function for H1s using a Chisquare method The "I"ntegral method is used for each point the cache contains the following info -1D : bc,e,xc,xw (bin content, error, x of center of bin, x bin width of bin) -2D : bc,e,xc,xw,yc,yw -3D : bc,e,xc,xw,yc,yw,zc,zw
Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood by determining the Poisson probability that given a number of entries in a particular bin, the fit would predict it's value. This is then done for each bin, and the sum of the logs is taken as the likelihood. Default method (function evaluated at center of bin) for each point the cache contains the following info -1D : bc,e,xc (bin content, error, x of center of bin) -2D : bc,e,xc,yc -3D : bc,e,xc,yc,zc
Minimization function for H1s using a Likelihood method*-*-*-*-*-*
Basically, it forms the likelihood by determining the Poisson
probability that given a number of entries in a particular bin,
the fit would predict it's value. This is then done for each bin,
and the sum of the logs is taken as the likelihood.
The "I"ntegral method is used
for each point the cache contains the following info
-1D : bc,e,xc,xw (bin content, error, x of center of bin, x bin width of bin)
-2D : bc,e,xc,xw,yc,yw
-3D : bc,e,xc,xw,yc,yw,zc,zw