TFumili.cxx

Go to the documentation of this file.
00001 // @(#)root/fumili:$Id: TFumili.cxx 34991 2010-08-25 10:32:26Z brun $
00002 // Author: Stanislav Nesterov  07/05/2003
00003 
00004 //______________________________________________________________________________
00005 //         FUMILI
00006 //  Based on ideas, proposed by I.N. Silin
00007 //    [See NIM A440, 2000 (p431)]
00008 // converted from FORTRAN to C  by
00009 //     Sergey Yaschenko <s.yaschenko@fz-juelich.de>
00010 //
00011 //______________________________________________________________________________
00012 //BEGIN_HTML <!--
00013 /* -->
00014 <H2>FUMILI minimization package</H2>
00015 <p>FUMILI is used to minimize Chi-square function or to search maximum of
00016 likelihood function.
00017 
00018 <p>Experimentally measured values $F_i$ are fitted with theoretical
00019 functions $f_i({\vec x}_i,\vec\theta\,\,)$, where ${\vec x}_i$ are
00020 coordinates, and $\vec\theta$ -- vector of parameters.
00021 
00022 <p>For better convergence Chi-square function has to be the following form
00023 
00024 <p>$$
00025 {\chi^2\over2}={1\over2}\sum^n_{i=1}\left(f_i(\vec
00026 x_i,\vec\theta\,\,)-F_i\over\sigma_i\right)^2 \eqno(1)
00027 $$
00028 <p>where $\sigma_i$ are errors of measured function.
00029 
00030 <p>The minimum condition is
00031 <p>$$
00032 {\partial\chi^2\over\partial\theta_i}=\sum^n_{j=1}{1\over\sigma^2_j}\cdot
00033 {\partial f_j\over\partial\theta_i}\left[f_j(\vec
00034 x_j,\vec\theta\,\,)-F_j\right]=0,\qquad i=1\ldots m\eqno(2)
00035 $$
00036 <p>where m is the quantity of parameters.
00037 
00038 <p>Expanding left part of (2) over parameter increments and
00039 retaining only linear terms one gets
00040 <p>$$
00041 \left(\partial\chi^2\over\theta_i\right)_{\vec\theta={\vec\theta}^0}
00042 +\sum_k\left(\partial^2\chi^2\over\partial\theta_i\partial\theta_k\right)_{
00043 \vec\theta={\vec\theta}^0}\cdot(\theta_k-\theta_k^0)
00044 = 0\eqno(3)
00045 $$
00046 
00047  <p>Here ${\vec\theta}_0$ is some initial value of parameters. In general
00048 case:
00049 <p>$$
00050 {\partial^2\chi^2\over\partial\theta_i\partial\theta_k}=
00051 \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i}
00052 {\partial f_j\over\theta_k} +
00053 \sum^n_{j=1}{(f_j - F_j)\over\sigma^2_j}\cdot
00054 {\partial^2f_j\over\partial\theta_i\partial\theta_k}\eqno(4)
00055 $$
00056 
00057 <p>In FUMILI algorithm for second derivatives of Chi-square approximate
00058 expression is used when last term in (4) is discarded. It is often
00059 done, not always wittingly, and sometimes causes troubles, for example,
00060 if user wants to limit parameters with positive values by writing down
00061 $\theta_i^2$ instead of $\theta_i$. FUMILI will fail if one tries
00062 minimize $\chi^2 = g^2(\vec\theta)$ where g is arbitrary function.
00063 
00064 <p>Approximate value is:
00065 <p>$${\partial^2\chi^2\over\partial\theta_i\partial\theta_k}\approx
00066 Z_{ik}=
00067 \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i}
00068 {\partial f_j\over\theta_k}\eqno(5)
00069 $$
00070 
00071 <p>Then the equations for parameter increments are
00072 <p>$$\left(\partial\chi^2\over\partial\theta_i\right)_{\vec\theta={\vec\theta}^0}
00073 +\sum_k Z_{ik}\cdot(\theta_k-\theta^0_k) = 0,
00074 \qquad i=1\ldots m\eqno(6)
00075 $$
00076 
00077 <p>Remarkable feature of algorithm is the technique for step
00078 restriction. For an initial value of parameter ${\vec\theta}^0$ a
00079 parallelepiped $P_0$ is built with the center at ${\vec\theta}^0$ and
00080 axes parallel to coordinate axes $\theta_i$. The lengths of
00081 parallelepiped sides along i-th axis is $2b_i$, where $b_i$ is such a
00082 value that the functions $f_j(\vec\theta)$ are quasi-linear all over
00083 the parallelepiped.
00084 
00085 <p>FUMILI takes into account simple linear inequalities in the form:
00086 $$
00087 \theta_i^{\rm min}\le\theta_i\le\theta^{\rm max}_i\eqno(7)
00088 $$
00089 
00090 <p>They form parallelepiped $P$ ($P_0$ may be deformed by $P$).
00091 Very similar step formulae are used in FUMILI for negative logarithm
00092 of the likelihood function with the same idea - linearization of
00093 function argument.
00094 
00095 <!--*/
00096 // -->END_HTML
00097 //______________________________________________________________________________
00098 
00099 
00100 
00101 #include "TROOT.h"
00102 #include "TFumili.h"
00103 #include "TMath.h"
00104 #include "TF1.h"
00105 #include "TH1.h"
00106 #include "TGraphAsymmErrors.h"
00107 
00108 #include "Riostream.h"
00109 
00110 
00111 extern void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00112 extern void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00113 extern void GraphFitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00114 
00115 
00116 ClassImp(TFumili);
00117 
00118 TFumili *gFumili=0;
00119 // Machine dependent values  fiXME!!
00120 // But don't set min=max=0 if param is unlimited
00121 static const Double_t gMAXDOUBLE=1e300;
00122 static const Double_t gMINDOUBLE=-1e300;
00123 
00124 
00125 //______________________________________________________________________________
00126 TFumili::TFumili(Int_t maxpar)
00127 {//----------- FUMILI constructor ---------
00128    // maxpar is the maximum number of parameters used with TFumili object
00129    //
00130    fMaxParam = TMath::Max(maxpar,25);
00131    if (fMaxParam>200) fMaxParam=200;
00132    fMaxParam2 *= fMaxParam;
00133    BuildArrays();
00134 
00135    fNumericDerivatives = true;
00136    fLogLike = false;
00137    fNpar    = fMaxParam;
00138    fGRAD    = false;
00139    fWARN    = true;
00140    fDEBUG   = false;
00141    fNlog    = 0;
00142    fSumLog  = 0;
00143    fNED1    = 0;
00144    fNED2    = 0;
00145    fNED12   = fNED1+fNED2;
00146    fEXDA    = 0;
00147    fFCN     = 0;
00148    fNfcn    = 0;
00149    fRP      = 1.e-15; //precision
00150    fS       = 1e10;
00151    fEPS     =0.01;
00152    fENDFLG  = 0;
00153    fNlimMul = 2;
00154    fNmaxIter= 150;
00155    fNstepDec= 3;
00156    fLastFixed = -1;
00157 
00158    SetName("Fumili");
00159    gFumili = this;
00160    gROOT->GetListOfSpecials()->Add(gFumili);
00161 }
00162 
00163 //______________________________________________________________________________
00164 void TFumili::BuildArrays(){
00165    //
00166    //   Allocates memory for internal arrays. Called by TFumili::TFumili
00167    //
00168    fCmPar      = new Double_t[fMaxParam];
00169    fA          = new Double_t[fMaxParam];
00170    fPL0        = new Double_t[fMaxParam];
00171    fPL         = new Double_t[fMaxParam];
00172    fParamError = new Double_t[fMaxParam];
00173    fDA         = new Double_t[fMaxParam];
00174    fAMX        = new Double_t[fMaxParam];
00175    fAMN        = new Double_t[fMaxParam];
00176    fR          = new Double_t[fMaxParam];
00177    fDF         = new Double_t[fMaxParam];
00178    fGr         = new Double_t[fMaxParam];
00179    fANames     = new TString[fMaxParam];
00180 
00181    //   fX = new Double_t[10];
00182 
00183    Int_t zSize = fMaxParam*(fMaxParam+1)/2;
00184    fZ0 = new Double_t[zSize];
00185    fZ  = new Double_t[zSize];
00186 
00187    for (Int_t i=0;i<fMaxParam;i++) {
00188       fA[i] =0.;
00189       fDF[i]=0.;
00190       fAMN[i]=gMINDOUBLE;
00191       fAMX[i]=gMAXDOUBLE;
00192       fPL0[i]=.1;
00193       fPL[i] =.1;
00194       fParamError[i]=0.;
00195       fANames[i]=Form("%d",i);
00196    }
00197 }
00198 
00199 
00200 //______________________________________________________________________________
00201 TFumili::~TFumili() {
00202    //
00203    // TFumili destructor
00204    //
00205    DeleteArrays();
00206    gROOT->GetListOfSpecials()->Remove(this);
00207    if (gFumili == this) gFumili = 0;
00208 }
00209 
00210 //______________________________________________________________________________
00211 Double_t TFumili::Chisquare(Int_t npar, Double_t *params) const
00212 {
00213    // return a chisquare equivalent
00214 
00215    Double_t amin = 0;
00216    H1FitChisquareFumili(npar,params,amin,params,1);
00217    return 2*amin;
00218 }
00219 
00220 
00221 //______________________________________________________________________________
00222 void TFumili::Clear(Option_t *)
00223 {
00224    //
00225    // Resets all parameter names, values and errors to zero
00226    //
00227    // Argument opt is ignored
00228    //
00229    // NB: this procedure doesn't reset parameter limits
00230    //
00231    fNpar = fMaxParam;
00232    fNfcn = 0;
00233    for (Int_t i=0;i<fNpar;i++) {
00234       fA[i]   =0.;
00235       fDF[i]  =0.;
00236       fPL0[i] =.1;
00237       fPL[i]  =.1;
00238       fAMN[i] = gMINDOUBLE;
00239       fAMX[i] = gMAXDOUBLE;
00240       fParamError[i]=0.;
00241       fANames[i]=Form("%d",i);
00242    }
00243 }
00244 
00245 
00246 //______________________________________________________________________________
00247 void TFumili::DeleteArrays(){
00248    //
00249    // Deallocates memory. Called from destructor TFumili::~TFumili
00250    //
00251    delete[] fCmPar;
00252    delete[] fANames;
00253    delete[] fDF;
00254    // delete[] fX;
00255    delete[] fZ0;
00256    delete[] fZ;
00257    delete[] fGr;
00258    delete[] fA;
00259    delete[] fPL0;
00260    delete[] fPL;
00261    delete[] fDA;
00262    delete[] fAMN;
00263    delete[] fAMX;
00264    delete[] fParamError;
00265    delete[] fR;
00266 }
00267 
00268 
00269 //______________________________________________________________________________
00270 void TFumili::Derivatives(Double_t *df,Double_t *fX){
00271    //
00272    // Calculates partial derivatives of theoretical function
00273    //
00274    // Input:
00275    //    fX  - vector of data point
00276    // Output:
00277    //    DF - array of derivatives
00278    //
00279    // ARITHM.F
00280    // Converted from CERNLIB
00281    //
00282    Double_t ff,ai,hi,y,pi;
00283    y = EvalTFN(df,fX);
00284    for (Int_t i=0;i<fNpar;i++) {
00285       df[i]=0;
00286       if(fPL0[i]>0.) {
00287          ai = fA[i]; // save current parameter value
00288          hi = 0.01*fPL0[i]; // diff step
00289          pi = fRP*TMath::Abs(ai);
00290          if (hi<pi) hi = pi; // if diff step is less than precision
00291          fA[i] = ai+hi;
00292 
00293          if (fA[i]>fAMX[i]) { // if param is out of limits
00294             fA[i] = ai-hi;
00295             hi = -hi;
00296             if (fA[i]<fAMN[i]) { // again out of bounds
00297                fA[i] = fAMX[i];   // set param to high limit
00298                hi = fAMX[i]-ai;
00299                if (fAMN[i]-ai+hi<0) { // if hi < (ai-fAMN)
00300                   fA[i]=fAMN[i];
00301                   hi=fAMN[i]-ai;
00302                }
00303             }
00304          }
00305          ff = EvalTFN(df,fX);
00306          df[i] = (ff-y)/hi;
00307          fA[i] = ai;
00308       }
00309    }
00310 }
00311 
00312 
00313 
00314 //______________________________________________________________________________
00315 Int_t TFumili::Eval(Int_t& npar, Double_t *grad, Double_t &fval, Double_t *par, Int_t flag)
00316 {
00317    // Evaluate the minimisation function
00318    //  Input parameters:
00319    //    npar:    number of currently variable parameters
00320    //    par:     array of (constant and variable) parameters
00321    //    flag:    Indicates what is to be calculated
00322    //    grad:    array of gradients
00323    //  Output parameters:
00324    //    fval:    The calculated function value.
00325    //    grad:    The vector of first derivatives.
00326    //
00327    // The meaning of the parameters par is of course defined by the user,
00328    // who uses the values of those parameters to calculate his function value.
00329    // The starting values must be specified by the user.
00330    // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00331    // Inside FCN user has to define Z-matrix by means TFumili::GetZ
00332    //  and TFumili::Derivatives,
00333    // set theoretical function by means of TFumili::SetUserFunc,
00334    // but first - pass number of parameters by TFumili::SetParNumber
00335    // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00336    // Later values are determined by Fumili as it searches for the minimum
00337    // or performs whatever analysis is requested by the user.
00338    //
00339    // The default function calls the function specified in SetFCN
00340    //
00341 
00342    if (fFCN) (*fFCN)(npar,grad,fval,par,flag);
00343    return npar;
00344 }
00345 
00346 
00347 //______________________________________________________________________________
00348 Double_t TFumili::EvalTFN(Double_t * /*df*/, Double_t *X)
00349 {
00350    // Evaluate theoretical function
00351    // df: array of partial derivatives
00352    // X:  vector of theoretical function argument
00353 
00354    // for the time being disable possibility to compute derivatives
00355    //if(fTFN)
00356    //  return (*fTFN)(df,X,fA);
00357    //else if(fTFNF1) {
00358 
00359    TF1 *f1 = (TF1*)fUserFunc;
00360    return f1->EvalPar(X,fA);
00361    //}
00362    //return 0.;
00363 }
00364 
00365 //______________________________________________________________________________
00366 Int_t TFumili::ExecuteCommand(const char *command, Double_t *args, Int_t nargs){
00367    //
00368    //  Execute MINUIT commands. MINImize, SIMplex, MIGrad and FUMili all
00369    //  will call TFumili::Minimize method.
00370    //
00371    //  For full command list see
00372    //  MINUIT. Reference Manual. CERN Program Library Long Writeup D506.
00373    //
00374    //  Improvement and errors calculation are not yet implemented as well
00375    //  as Monte-Carlo seeking and minimization.
00376    //  Contour commands are also unsupported.
00377    //
00378    //  command   : command string
00379    //  args      : array of arguments
00380    //  nargs     : number of arguments
00381    //
00382    TString comand = command;
00383    static TString clower = "abcdefghijklmnopqrstuvwxyz";
00384    static TString cupper = "ABCDEFGHIJKLMNOPQRSTUVWXYZ";
00385    const Int_t nntot = 40;
00386    const char *cname[nntot] = {
00387       "MINImize  ",    //  0    checked
00388       "SEEk      ",    //  1    none
00389       "SIMplex   ",    //  2    checked same as 0
00390       "MIGrad    ",    //  3    checked  same as 0
00391       "MINOs     ",    //  4    none
00392       "SET xxx   ",    //  5 lot of stuff
00393       "SHOw xxx  ",    //  6 -----------
00394       "TOP of pag",    //  7 .
00395       "fiX       ",   //   8 .
00396       "REStore   ",   //   9 .
00397       "RELease   ",   //   10 .
00398       "SCAn      ",   //   11  not yet implemented
00399       "CONtour   ",   //   12  not yet implemented
00400       "HESse     ",   //   13  not yet implemented
00401       "SAVe      ",   //   14  obsolete
00402       "IMProve   ",   //   15  not yet implemented
00403       "CALl fcn  ",   //   16 .
00404       "STAndard  ",   //   17 .
00405       "END       ",   //   18 .
00406       "EXIt      ",   //   19 .
00407       "RETurn    ",   //   20 .
00408       "CLEar     ",   //   21 .
00409       "HELP      ",   //   22 not yet implemented
00410       "MNContour ",   //   23 not yet implemented
00411       "STOp      ",   //   24 .
00412       "JUMp      ",   //   25 not yet implemented
00413       "          ",   //
00414       "          ",   //
00415       "FUMili    ",   //    28 checked same as 0
00416       "          ",   //
00417       "          ",  //
00418       "          ",  //
00419       "          ",  //
00420       "COVARIANCE",  // 33
00421       "PRINTOUT  ",  // 34
00422       "GRADIENT  ",  // 35
00423       "MATOUT    ",  // 36
00424       "ERROR DEF ",  // 37
00425       "LIMITS    ",  // 38
00426       "PUNCH     "}; // 39
00427 
00428 
00429    fCword = comand;
00430    fCword.ToUpper();
00431    if (nargs<=0) fCmPar[0] = 0;
00432    Int_t i;
00433    for(i=0;i<fMaxParam;i++) {
00434       if(i<=nargs) fCmPar[i] = args[i];
00435    }
00436    /*
00437    fNmaxIter = int(fCmPar[0]);
00438    if (fNmaxIter <= 0) {
00439        fNmaxIter = fNpar*10 + 20 + fNpar*M*5;
00440    }
00441    fEPS = fCmPar[1];
00442    */
00443    //*-*-               look for command in list CNAME . . . . . . . . . .
00444    TString ctemp = fCword(0,3);
00445    Int_t ind;
00446    for (ind = 0; ind < nntot; ++ind) {
00447       if (strncmp(ctemp.Data(),cname[ind],3) == 0) break;
00448    }
00449    if (ind==nntot) return -3; // Unknow command - input ignored
00450    if (fCword(0,4) == "MINO") ind=3;
00451    switch (ind) {
00452       case 0:  case 3: case 2: case 28:
00453          // MINImize [maxcalls] [tolerance]
00454          // also SIMplex, MIGrad  and  FUMili
00455          if(nargs>=1)
00456          fNmaxIter=TMath::Max(Int_t(fCmPar[0]),fNmaxIter); // fiXME!!
00457          if(nargs==2)
00458             fEPS=fCmPar[1];
00459          return Minimize();
00460       case 1:
00461          // SEEk not implemented in this package
00462          return -10;
00463 
00464       case 4: // MINos errors analysis not implemented
00465          return -10;
00466 
00467       case 5: case 6: // SET xxx & SHOW xxx
00468          return ExecuteSetCommand(nargs);
00469 
00470       case 7: // Obsolete command
00471          Printf("1");
00472          return 0;
00473       case 8: // fiX <parno> ....
00474          if (nargs<1) return -1; // No parameters specified
00475          for (i=0;i<nargs;i++) {
00476             Int_t parnum = Int_t(fCmPar[i])-1;
00477             FixParameter(parnum);
00478          }
00479          return 0;
00480       case 9: // REStore <code>
00481          if (nargs<1) return 0;
00482          if(fCmPar[0]==0.)
00483             for (i=0;i<fNpar;i++)
00484                ReleaseParameter(i);
00485          else
00486             if(fCmPar[0]==1.) {
00487                ReleaseParameter(fLastFixed);
00488                cout <<fLastFixed<<endl;
00489             }
00490          return 0;
00491       case 10: // RELease <parno> ...
00492          if (nargs<1) return -1; // No parameters specified
00493          for (i=0;i<nargs;i++) {
00494             Int_t parnum = Int_t(fCmPar[i])-1;
00495             ReleaseParameter(parnum);
00496          }
00497          return 0;
00498       case 11: // SCAn not implemented
00499          return -10;
00500       case 12: // CONt not implemented
00501          return -10;
00502 
00503       case 13: // HESSe not implemented
00504          return -10;
00505       case 14: // SAVe
00506          Printf("SAVe command is obsolete");
00507          return -10;
00508       case 15: // IMProve not implemented
00509          return -10;
00510       case 16: // CALl fcn <iflag>
00511          {if(nargs<1) return -1;
00512          Int_t flag = Int_t(fCmPar[0]);
00513          Double_t fval;
00514          Eval(fNpar,fGr,fval,fA,flag);
00515          return 0;}
00516       case 17: // STAndard must call function STAND
00517          return 0;
00518       case 18:   case 19:
00519       case 20:  case 24: {
00520          Double_t fval;
00521          Int_t flag = 3;
00522          Eval(fNpar,fGr,fval,fA,flag);
00523          return 0;
00524       }
00525       case 21:
00526          Clear();
00527          return 0;
00528       case 22: //HELp not implemented
00529       case 23: //MNContour not implemented
00530       case 25: // JUMp not implemented
00531          return -10;
00532       case 26:   case 27:  case 29:  case 30:  case 31:  case 32:
00533          return 0; // blank commands
00534       case 33:   case 34:   case 35:  case 36:   case 37:  case 38:
00535       case 39:
00536          Printf("Obsolete command. Use corresponding SET command instead");
00537          return -10;
00538       default:
00539          break;
00540    }
00541    return 0;
00542 }
00543 
00544 
00545 
00546 //______________________________________________________________________________
00547 Int_t TFumili::ExecuteSetCommand(Int_t nargs){
00548    //
00549    // Called from TFumili::ExecuteCommand in case
00550    // of "SET xxx" and "SHOW xxx".
00551    //
00552    static Int_t nntot = 30;
00553    static const char *cname[30] = {
00554       "FCN value ", // 0 .
00555       "PARameters", // 1 .
00556       "LIMits    ", // 2 .
00557       "COVariance", // 3 .
00558       "CORrelatio", // 4 .
00559       "PRInt levl", // 5 not implemented yet
00560       "NOGradient", // 6 .
00561       "GRAdient  ", // 7 .
00562       "ERRor def ", // 8 not sure how to implement - by time being ignored
00563       "INPut file", // 9 not implemented
00564       "WIDth page", // 10 not implemented yet
00565       "LINes page", // 11 not implemented yet
00566       "NOWarnings", // 12 .
00567       "WARnings  ", // 13 .
00568       "RANdom gen", // 14 not implemented
00569       "TITle     ", // 15 ignored
00570       "STRategy  ", // 16 ignored
00571       "EIGenvalue", // 17 not implemented yet
00572       "PAGe throw", // 18 ignored
00573       "MINos errs", // 19 not implemented yet
00574       "EPSmachine", // 20 .
00575       "OUTputfile", // 21 not implemented
00576       "BATch     ", // 22 ignored
00577       "INTeractiv", // 23 ignored
00578       "VERsion   ", // 24 .
00579       "reserve   ", // 25 .
00580       "NODebug   ", // 26 .
00581       "DEBug     ", // 27 .
00582       "SHOw      ", // 28 err
00583       "SET       "};// 29 err
00584 
00585    TString  cfname, cmode, ckind,  cwarn, copt, ctemp, ctemp2;
00586    Int_t i, ind;
00587    Bool_t setCommand=kFALSE;
00588    for (ind = 0; ind < nntot; ++ind) {
00589       ctemp  = cname[ind];
00590       ckind  = ctemp(0,3);
00591       ctemp2 = fCword(4,6);
00592       if (strstr(ctemp2.Data(),ckind.Data())) break;
00593    }
00594    ctemp2 = fCword(0,3);
00595    if(ctemp2.Contains("SET")) setCommand=true;
00596    if(ctemp2.Contains("HEL") || ctemp2.Contains("SHO")) setCommand=false;
00597 
00598    if (ind>=nntot) return -3;
00599 
00600    switch(ind) {
00601       case 0: // SET FCN value illegial // SHOw only
00602          if(!setCommand) Printf("FCN=%f",fS);
00603          return 0;
00604       case 1: // PARameter <parno> <value>
00605          {
00606          if (nargs<2 && setCommand) return -1;
00607          Int_t parnum;
00608          Double_t val;
00609          if(setCommand) {
00610             parnum = Int_t(fCmPar[0])-1;
00611             val= fCmPar[1];
00612             if(parnum<0 || parnum>=fNpar) return -2; //no such parameter
00613             fA[parnum] = val;
00614          } else {
00615             if (nargs>0) {
00616                parnum = Int_t(fCmPar[0])-1;
00617                if(parnum<0 || parnum>=fNpar) return -2; //no such parameter
00618                Printf("Parameter %s = %E",fANames[parnum].Data(),fA[parnum]);
00619             } else
00620                for (i=0;i<fNpar;i++)
00621                Printf("Parameter %s = %E",fANames[i].Data(),fA[i]);
00622 
00623          }
00624          return 0;
00625          }
00626       case 2: // LIMits [parno] [ <lolim> <uplim> ]
00627          {
00628          Int_t parnum;
00629          Double_t lolim,uplim;
00630          if (nargs<1) {
00631             for(i=0;i<fNpar;i++)
00632                if(setCommand) {
00633                   fAMN[i] = gMINDOUBLE;
00634                   fAMX[i] = gMAXDOUBLE;
00635                } else
00636                   Printf("Limits for param %s: Low=%E, High=%E",
00637                          fANames[i].Data(),fAMN[i],fAMX[i]);
00638          } else {
00639             parnum = Int_t(fCmPar[0])-1;
00640             if(parnum<0 || parnum>=fNpar)return -1;
00641             if(setCommand) {
00642                if(nargs>2) {
00643                   lolim = fCmPar[1];
00644                   uplim = fCmPar[2];
00645                   if(uplim==lolim) return -1;
00646                   if(lolim>uplim) {
00647                      Double_t tmp = lolim;
00648                      lolim = uplim;
00649                      uplim = tmp;
00650                   }
00651                } else {
00652                   lolim = gMINDOUBLE;
00653                   uplim = gMAXDOUBLE;
00654                }
00655                fAMN[parnum] = lolim;
00656                fAMX[parnum] = uplim;
00657             } else
00658                Printf("Limits for param %s Low=%E, High=%E",
00659                     fANames[parnum].Data(),fAMN[parnum],fAMX[parnum]);
00660             }
00661             return 0;
00662          }
00663       case 3:
00664          {
00665          if(setCommand) return 0;
00666          Printf("\nCovariant matrix ");
00667          Int_t l = 0,nn=0,nnn=0;
00668          for (i=0;i<fNpar;i++) if(fPL0[i]>0.) nn++;
00669          for (i=0;i<nn;i++) {
00670             for(;fPL0[nnn]<=0.;nnn++) { }
00671             printf("%5s: ",fANames[nnn++].Data());
00672             for (Int_t j=0;j<=i;j++)
00673                printf("%11.2E",fZ[l++]);
00674             cout<<endl;
00675          }
00676          cout<<endl;
00677          return 0;
00678          }
00679       case 4:
00680          if(setCommand) return 0;
00681          Printf("\nGlobal correlation factors (maximum correlation of the parameter\n  with arbitrary linear combination of other parameters)");
00682          for(i=0;i<fNpar;i++) {
00683             printf("%5s: ",fANames[i].Data());
00684             printf("%11.3E\n",TMath::Sqrt(1-1/((fR[i]!=0.)?fR[i]:1.)) );
00685          }
00686          cout<<endl;
00687          return 0;
00688       case 5:   // PRIntout not implemented
00689          return -10;
00690       case 6: // NOGradient
00691          if(!setCommand) return 0;
00692          fGRAD = false;
00693          return 0;
00694       case 7: // GRAdient
00695          if(!setCommand) return 0;
00696          fGRAD = true;
00697          return 0;
00698       case 8: // ERRordef - now ignored
00699          return 0;
00700       case 9: // INPut - not implemented
00701          return -10;
00702       case 10: // WIDthpage - not implemented
00703          return -10;
00704       case 11: // LINesperpage - not implemented
00705          return -10;
00706       case 12: //NOWarnings
00707          if(!setCommand) return 0;
00708          fWARN = false;
00709          return 0;
00710       case 13: // WARnings
00711          if(!setCommand) return 0;
00712          fWARN = true;
00713          return 0;
00714       case 14: // RANdomgenerator - not implemented
00715          return -10;
00716       case 15: // TITle - ignored
00717          return 0;
00718       case 16: // STRategy - ignored
00719          return 0;
00720       case 17: // EIGenvalues - not implemented
00721          return -10;
00722       case 18: // PAGethrow - ignored
00723          return 0;
00724       case 19: // MINos errors - not implemented
00725          return -10;
00726       case 20: //EPSmachine
00727          if(!setCommand) {
00728             Printf("Relative floating point presicion RP=%E",fRP);
00729          } else
00730             if (nargs>0) {
00731                Double_t pres=fCmPar[0];
00732                if (pres<1e-5 && pres>1e-34) fRP=pres;
00733             }
00734          return 0;
00735       case 21: // OUTputfile - not implemented
00736          return -10;
00737       case 22: // BATch - ignored
00738          return 0;
00739       case 23: // INTerative - ignored
00740          return 0;
00741       case 24: // VERsion
00742          if(setCommand) return 0;
00743          Printf("FUMILI-ROOT version 0.1");
00744          return 0;
00745       case 25: // reserved
00746          return 0;
00747       case 26: // NODebug
00748          if(!setCommand) return 0;
00749          fDEBUG = false;
00750          return 0;
00751       case 27: // DEBug
00752          if(!setCommand) return 0;
00753          fDEBUG = true;
00754          return 0;
00755       case 28:
00756       case 29:
00757          return -3;
00758       default:
00759          break;
00760    }
00761    return -3;
00762 }
00763 
00764 //______________________________________________________________________________
00765 void TFumili::FixParameter(Int_t ipar) {
00766    // Fixes parameter number ipar
00767 
00768    if(ipar>=0 && ipar<fNpar && fPL0[ipar]>0.) {
00769       fPL0[ipar] = -fPL0[ipar];
00770       fLastFixed = ipar;
00771    }
00772 }
00773 
00774 //______________________________________________________________________________
00775 Double_t *TFumili::GetCovarianceMatrix() const
00776 {
00777    // return a pointer to the covariance matrix
00778 
00779    return fZ;
00780 
00781 }
00782 
00783 //______________________________________________________________________________
00784 Double_t TFumili::GetCovarianceMatrixElement(Int_t i, Int_t j) const
00785 {
00786    // return element i,j from the covariance matrix
00787 
00788    if (!fZ) return 0;
00789    if (i < 0 || i >= fNpar || j < 0 || j >= fNpar) {
00790       Error("GetCovarianceMatrixElement","Illegal arguments i=%d, j=%d",i,j);
00791       return 0;
00792    }
00793    return fZ[j+fNpar*i];
00794 }
00795 
00796 
00797 //______________________________________________________________________________
00798 Int_t TFumili::GetNumberTotalParameters() const
00799 {
00800    // return the total number of parameters (free + fixed)
00801 
00802    return fNpar;
00803 }
00804 
00805 //______________________________________________________________________________
00806 Int_t TFumili::GetNumberFreeParameters() const
00807 {
00808    // return the number of free parameters
00809 
00810    Int_t nfree = fNpar;
00811    for (Int_t i=0;i<fNpar;i++) {
00812       if (IsFixed(i)) nfree--;
00813    }
00814    return nfree;
00815 }
00816 
00817 //______________________________________________________________________________
00818 Double_t TFumili::GetParError(Int_t ipar) const
00819 {
00820    // return error of parameter ipar
00821 
00822    if (ipar<0 || ipar>=fNpar) return 0;
00823    else return fParamError[ipar];
00824 }
00825 
00826 //______________________________________________________________________________
00827 Double_t TFumili::GetParameter(Int_t ipar) const
00828 {
00829    // return current value of parameter ipar
00830 
00831    if (ipar<0 || ipar>=fNpar) return 0;
00832    else return fA[ipar];
00833 }
00834 
00835 
00836 //______________________________________________________________________________
00837 Int_t TFumili::GetParameter(Int_t ipar,char *cname,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const
00838 {
00839    // Get various ipar parameter attributs:
00840    //
00841    // cname:    parameter name
00842    // value:    parameter value
00843    // verr:     parameter error
00844    // vlow:     lower limit
00845    // vhigh:    upper limit
00846    // WARNING! parname must be suitably dimensionned in the calling function.
00847 
00848    if (ipar<0 || ipar>=fNpar) {
00849       value = 0;
00850       verr  = 0;
00851       vlow  = 0;
00852       vhigh = 0;
00853       return -1;
00854    }
00855    strcpy(cname,fANames[ipar].Data());
00856    value = fA[ipar];
00857    verr  = fParamError[ipar];
00858    vlow  = fAMN[ipar];
00859    vhigh = fAMX[ipar];
00860    return 0;
00861 }
00862 
00863 //______________________________________________________________________________
00864 const char *TFumili::GetParName(Int_t ipar) const
00865 {
00866    // return name of parameter ipar
00867 
00868    if (ipar < 0 || ipar > fNpar) return "";
00869    return fANames[ipar].Data();
00870 }
00871 
00872 //______________________________________________________________________________
00873 Int_t TFumili::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
00874 {
00875    // Return errors after MINOs
00876    // not implemented
00877    eparab = 0;
00878    globcc = 0;
00879    if (ipar<0 || ipar>=fNpar) {
00880       eplus  = 0;
00881       eminus = 0;
00882       return -1;
00883    }
00884    eplus=fParamError[ipar];
00885    eminus=-eplus;
00886    return 0;
00887 }
00888 
00889 //______________________________________________________________________________
00890 Int_t TFumili::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
00891 {
00892    // return global fit parameters
00893    //   amin     : chisquare
00894    //   edm      : estimated distance to minimum
00895    //   errdef
00896    //   nvpar    : number of variable parameters
00897    //   nparx    : total number of parameters
00898    amin   = 2*fS;
00899    edm    = fGT; //
00900    errdef = 0; // ??
00901    nparx  = fNpar;
00902    nvpar  = 0;
00903    for(Int_t ii=0; ii<fNpar; ii++) {
00904       if(fPL0[ii]>0.) nvpar++;
00905    }
00906    return 0;
00907 }
00908 
00909 
00910 
00911 //______________________________________________________________________________
00912 Double_t TFumili::GetSumLog(Int_t n)
00913 {
00914    // return Sum(log(i) i=0,n
00915    // used by log likelihood fits
00916 
00917    if (n < 0) return 0;
00918    if (n > fNlog) {
00919       if (fSumLog) delete [] fSumLog;
00920       fNlog = 2*n+1000;
00921       fSumLog = new Double_t[fNlog+1];
00922       Double_t fobs = 0;
00923       for (Int_t j=0;j<=fNlog;j++) {
00924          if (j > 1) fobs += TMath::Log(j);
00925          fSumLog[j] = fobs;
00926       }
00927    }
00928    if (fSumLog) return fSumLog[n];
00929    return 0;
00930 }
00931 
00932 
00933 
00934 //______________________________________________________________________________
00935 void TFumili::InvertZ(Int_t n)
00936 {
00937    // Inverts packed diagonal matrix Z by square-root method.
00938    //  Matrix elements corresponding to
00939    // fix parameters are removed.
00940    //
00941    // n: number of variable parameters
00942    //
00943    static Double_t am = 3.4e138;
00944    static Double_t rp = 5.0e-14;
00945    Double_t  ap, aps, c, d;
00946    Double_t *r_1=fR;
00947    Double_t *pl_1=fPL;
00948    Double_t *z_1=fZ;
00949    Int_t i, k, l, ii, ki, li, kk, ni, ll, nk, nl, ir, lk;
00950    if (n < 1) {
00951       return;
00952    }
00953    --pl_1;
00954    --r_1;
00955    --z_1;
00956    aps = am / n;
00957    aps = sqrt(aps);
00958    ap = 1.0e0 / (aps * aps);
00959    ir = 0;
00960    for (i = 1; i <= n; ++i) {
00961       L1:
00962          ++ir;
00963          if (pl_1[ir] <= 0.0e0) goto L1;
00964          else                   goto L2;
00965       L2:
00966          ni = i * (i - 1) / 2;
00967          ii = ni + i;
00968          k = n + 1;
00969          if (z_1[ii] <= rp * TMath::Abs(r_1[ir]) || z_1[ii] <= ap) {
00970             goto L19;
00971          }
00972          z_1[ii] = 1.0e0 / sqrt(z_1[ii]);
00973          nl = ii - 1;
00974       L3:
00975          if (nl - ni <= 0) goto L5;
00976          else              goto L4;
00977          L4:
00978          z_1[nl] *= z_1[ii];
00979          if (TMath::Abs(z_1[nl]) >= aps) {
00980             goto L16;
00981          }
00982          --nl;
00983          goto L3;
00984       L5:
00985          if (i - n >= 0) goto L12;
00986          else            goto L6;
00987       L6:
00988          --k;
00989          nk = k * (k - 1) / 2;
00990          nl = nk;
00991          kk = nk + i;
00992          d = z_1[kk] * z_1[ii];
00993          c = d * z_1[ii];
00994          l = k;
00995       L7:
00996          ll = nk + l;
00997          li = nl + i;
00998          z_1[ll] -= z_1[li] * c;
00999          --l;
01000          nl -= l;
01001          if (l - i <= 0) goto L9;
01002          else            goto L7;
01003       L8:
01004          ll = nk + l;
01005          li = ni + l;
01006          z_1[ll] -= z_1[li] * d;
01007       L9:
01008          --l;
01009          if (l <= 0) goto L10;
01010          else        goto L8;
01011       L10:
01012          z_1[kk] = -c;
01013          if (k - i - 1 <= 0) goto L11;
01014          else                goto L6;
01015       L11:
01016          ;
01017    }
01018    L12:
01019       for (i = 1; i <= n; ++i) {
01020          for (k = i; k <= n; ++k) {
01021             nl = k * (k - 1) / 2;
01022             ki = nl + i;
01023             d = 0.0e0;
01024             for (l = k; l <= n; ++l) {
01025                li = nl + i;
01026                lk = nl + k;
01027                d += z_1[li] * z_1[lk];
01028                nl += l;
01029             }
01030             ki = k * (k - 1) / 2 + i;
01031             z_1[ki] = d;
01032          }
01033       }
01034    L15:
01035       return;
01036    L16:
01037       k = i + nl - ii;
01038       ir = 0;
01039       for (i = 1; i <= k; ++i) {
01040          L17:
01041             ++ir;
01042             if (pl_1[ir] <= 0.0e0) {
01043                goto L17;
01044             }
01045       }
01046    L19:
01047       pl_1[ir] = -2.0e0;
01048       r_1[ir] = 0.0e0;
01049       fINDFLG[0] = ir - 1;
01050       goto L15;
01051 }
01052 
01053 
01054 
01055 
01056 //______________________________________________________________________________
01057 Bool_t TFumili::IsFixed(Int_t ipar) const
01058 {
01059    //return kTRUE if parameter ipar is fixed, kFALSE othersise)
01060 
01061    if(ipar < 0 || ipar >= fNpar) {
01062       Warning("IsFixed","Illegal parameter number :%d",ipar);
01063       return kFALSE;
01064    }
01065    if (fPL0[ipar] < 0) return kTRUE;
01066    else                return kFALSE;
01067 }
01068 
01069 
01070 //______________________________________________________________________________
01071 Int_t TFumili::Minimize()
01072 {// Main minimization procedure
01073    //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*//
01074    //         FUMILI
01075    //  Based on ideas, proposed by I.N. Silin
01076    //    [See NIM A440, 2000 (p431)]
01077    // converted from FORTRAN to C  by
01078    //     Sergey Yaschenko <s.yaschenko@fz-juelich.de>
01079    //
01080    //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*//
01081    //
01082    // This function is called after setting theoretical function
01083    // by means of TFumili::SetUserFunc and initializing parameters.
01084    // Optionally one can set FCN function (see TFumili::SetFCN and TFumili::Eval)
01085    // If FCN is undefined then user has to provide data arrays by calling
01086    //  TFumili::SetData procedure.
01087    //
01088    // TFumili::Minimize return following values:
01089    //    0  - fit is converged
01090    //   -2  - function is not decreasing (or bad derivatives)
01091    //   -3  - error estimations are infinite
01092    //   -4  - maximum number of iterations is exceeded
01093    //
01094    Int_t i;
01095    // Flag3 - is fit is chi2 or likelihood? 0 - chi2, 1 - likelihood
01096    fINDFLG[2]=0;
01097    //
01098    // Are the parameters outside of the boundaries ?
01099    //
01100    Int_t parn;
01101 
01102    if(fFCN) {
01103       Eval(parn,fGr,fS,fA,9); fNfcn++;
01104    }
01105    for( i = 0; i < fNpar; i++) {
01106       if(fA[i] > fAMX[i]) fA[i] = fAMX[i];
01107       if(fA[i] < fAMN[i]) fA[i] = fAMN[i];
01108    }
01109 
01110    Int_t nn2, n, fixFLG,  ifix1, fi, nn3, nn1, n0;
01111    Double_t t1;
01112    Double_t sp, t, olds=0;
01113    Double_t bi, aiMAX=0, amb;
01114    Double_t afix, sigi, akap;
01115    Double_t alambd, al, bm, abi, abm;
01116    Int_t l1, k, ifix;
01117 
01118    nn2=0;
01119 
01120    // Number of parameters;
01121    n=fNpar;
01122    fixFLG=0;
01123 
01124    // Exit flag
01125    fENDFLG=0;
01126 
01127    // Flag2
01128    fINDFLG[1] = 0;
01129    ifix1=-1;
01130    fi=0;
01131    nn3=0;
01132 
01133    // Initialize param.step limits
01134    for( i=0; i < n; i++) {
01135       fR[i]=0.;
01136       if ( fEPS > 0.) fParamError[i] = 0.;
01137       fPL[i] = fPL0[i];
01138    }
01139 
01140 L3: // Start Iteration
01141 
01142    nn1 = 1;
01143    t1 = 1.;
01144 
01145 L4: // New iteration
01146 
01147    // fS - objective function value - zero first
01148    fS = 0.;
01149    // n0 - number of variable parameters in fit
01150    n0 = 0;
01151    for( i = 0; i < n; i++) {
01152       fGr[i]=0.; // zero gradients
01153       if (fPL0[i] > .0) {
01154          n0=n0+1;
01155          // new iteration - new parallelepiped
01156          if (fPL[i] > .0) fPL0[i]=fPL[i];
01157       }
01158    }
01159    Int_t nn0;
01160    // Calculate number of fZ-matrix elements as nn0=1+2+..+n0
01161    nn0 = n0*(n0+1)/2;
01162    // if (nn0 >= 1) ????
01163    // fZ-matrix is initialized
01164    for( i=0; i < nn0; i++) fZ[i]=0.;
01165 
01166    // Flag1
01167    fINDFLG[0] = 0;
01168    Int_t ijkl=1;
01169 
01170    // Calculate fS - objective function, fGr - gradients, fZ - fZ-matrix
01171    if(fFCN) {
01172       Eval(parn,fGr,fS,fA,2);
01173       fNfcn++;
01174    } else
01175       ijkl = SGZ();
01176    if(!ijkl) return 10;
01177    if (ijkl == -1) fINDFLG[0]=1;
01178 
01179    // sp - scaled on fS machine precision
01180    sp=fRP*TMath::Abs(fS);
01181 
01182    // save fZ-matrix
01183    for( i=0; i < nn0; i++) fZ0[i] = fZ[i];
01184    if (nn3 > 0) {
01185       if (nn1 <= fNstepDec) {
01186          t=2.*(fS-olds-fGT);
01187          if (fINDFLG[0] == 0) {
01188             if (TMath::Abs(fS-olds) <= sp && -fGT <= sp) goto L19;
01189             if(        0.59*t < -fGT) goto L19;
01190             t = -fGT/t;
01191             if (t < 0.25 ) t = 0.25;
01192          }
01193          else   t = 0.25;
01194          fGT = fGT*t;
01195          t1 = t1*t;
01196          nn2=0;
01197          for( i = 0; i < n; i++) {
01198             if (fPL[i] > 0.) {
01199                fA[i]=fA[i]-fDA[i];
01200                fPL[i]=fPL[i]*t;
01201                fDA[i]=fDA[i]*t;
01202                fA[i]=fA[i]+fDA[i];
01203             }
01204          }
01205          nn1=nn1+1;
01206          goto L4;
01207       }
01208    }
01209 
01210 L19:
01211 
01212    if(fINDFLG[0] != 0) {
01213       fENDFLG=-4;
01214       printf("trying to execute an illegal junp at L85\n");
01215       //goto L85;
01216    }
01217 
01218 
01219    Int_t k1, k2, i1, j, l;
01220    k1 = 1;
01221    k2 = 1;
01222    i1 = 1;
01223    // In this cycle we removed from fZ contributions from fixed parameters
01224    // We'll get fixed parameters after boudary check
01225    for( i = 0; i < n; i++) {
01226       if (fPL0[i] > .0) {
01227          // if parameter was fixed - release it
01228          if (fPL[i] == 0.) fPL[i]=fPL0[i];
01229          if (fPL[i] > .0) { // ??? it is already non-zero
01230             // if derivative is negative and we above maximum
01231             // or vice versa then fix parameter again and increment k1 by i1
01232             if ((fA[i] >= fAMX[i] && fGr[i] < 0.) ||
01233                    (fA[i] <= fAMN[i] && fGr[i] > 0.)) {
01234                fPL[i] = 0.;
01235                k1 = k1 + i1; // i1 stands for fZ-matrix row-number multiplier
01236                ///  - skip this row
01237                //  in case we are fixing parameter number i
01238             } else {
01239                for( j=0; j <= i; j++) {// cycle on columns of fZ-matrix
01240                   if (fPL0[j] > .0) {
01241                      // if parameter is not fixed then fZ = fZ0
01242                      // Now matrix fZ of other dimension
01243                      if (fPL[j] > .0) {
01244                         fZ[k2 -1] = fZ0[k1 -1];
01245                         k2=k2+1;
01246                      }
01247                      k1=k1+1;
01248                   }
01249                }
01250             }
01251          }
01252          else k1 = k1 + i1; // In case of negative fPL[i] - after mconvd
01253          i1=i1+1;  // Next row of fZ0
01254       }
01255    }
01256 
01257    // INVERT fZ-matrix (mconvd() procedure)
01258    i1 = 1;
01259    l  = 1;
01260    for( i = 0; i < n; i++) {// extract diagonal elements to fR-vector
01261       if (fPL[i] > .0) {
01262          fR[i] = fZ[l - 1];
01263          i1 = i1+1;
01264          l = l + i1;
01265       }
01266    }
01267 
01268    n0 = i1 - 1;
01269    InvertZ(n0);
01270 
01271    // fZ matrix now is inversed
01272    if (fINDFLG[0] != 0) { // problems
01273       // some PLs now have negative values, try to reduce fZ-matrix again
01274       fINDFLG[0] = 0;
01275       fINDFLG[1] = 1; // errors can be infinite
01276       fixFLG = fixFLG + 1;
01277       fi = 0;
01278       goto L19;
01279    }
01280 
01281    // ... CALCULATE THEORETICAL STEP TO MINIMUM
01282    i1 = 1;
01283    for( i = 0; i < n; i++) {
01284       fDA[i]=0.; // initial step is zero
01285       if (fPL[i] > .0) {   // for non-fixed parameters
01286          l1=1;
01287          for( l = 0; l < n; l++) {
01288             if (fPL[l] > .0) {
01289                // Caluclate offset of Z^-1(i1,l1) element in packed matrix
01290                // because we skip fixed param numbers we need also i,l
01291                if (i1 <= l1 ) k=l1*(l1-1)/2+i1;
01292                else k=i1*(i1-1)/2+l1;
01293                // dA_i = \sum (-Z^{-1}_{il}*grad(fS)_l)
01294                fDA[i]=fDA[i]-fGr[l]*fZ[k - 1];
01295                l1=l1+1;
01296             }
01297          }
01298          i1=i1+1;
01299       }
01300    }
01301    //          ... CHECK FOR PARAMETERS ON BOUNDARY
01302 
01303    afix=0.;
01304    ifix = -1;
01305    i1 = 1;
01306    l = i1;
01307    for( i = 0; i < n; i++)
01308       if (fPL[i] > .0) {
01309          sigi = TMath::Sqrt(TMath::Abs(fZ[l - 1])); // calculate \sqrt{Z^{-1}_{ii}}
01310          fR[i] = fR[i]*fZ[l - 1];      // Z_ii * Z^-1_ii
01311          if (fEPS > .0) fParamError[i]=sigi;
01312          if ((fA[i] >= fAMX[i] && fDA[i] > 0.) || (fA[i] <= fAMN[i]
01313                                                && fDA[i] < .0)) {
01314             // if parameter out of bounds and if step is making things worse
01315 
01316             akap = TMath::Abs(fDA[i]/sigi);
01317             // let's found maximum of dA/sigi - the worst of parameter steps
01318             if (akap > afix) {
01319                afix=akap;
01320                ifix=i;
01321                ifix1=i;
01322             }
01323          }
01324          i1=i1+1;
01325          l=l+i1;
01326       }
01327    if (ifix != -1) {
01328       // so the worst parameter is found - fix it and exclude,
01329       //  reduce fZ-matrix again
01330       fPL[ifix] = -1.;
01331       fixFLG = fixFLG + 1;
01332       fi = 0;
01333       //.. REPEAT CALCULATION OF THEORETICAL STEP AFTER fiXING EACH PARAMETER
01334       goto L19;
01335    }
01336 
01337    //... CALCULATE STEP CORRECTION FACTOR
01338 
01339    alambd = 1.;
01340    fAKAPPA = 0.;
01341    Int_t imax;
01342    imax = -1;
01343 
01344 
01345    for( i = 0; i < n; i++) {
01346       if (fPL[i] > .0) {
01347          bm = fAMX[i] - fA[i];
01348          abi = fA[i] + fPL[i]; // upper  parameter limit
01349          abm = fAMX[i];
01350          if (fDA[i] <= .0) {
01351             bm = fA[i] - fAMN[i];
01352             abi = fA[i] - fPL[i]; // lower parameter limit
01353             abm = fAMN[i];
01354          }
01355          bi = fPL[i];
01356          // if parallelepiped boundary is crossing limits
01357          // then reduce it (deforming)
01358          if ( bi > bm) {
01359             bi = bm;
01360             abi = abm;
01361          }
01362          // if calculated step is out of bounds
01363          if ( TMath::Abs(fDA[i]) > bi) {
01364             // derease step splitter alambdA if needed
01365             al = TMath::Abs(bi/fDA[i]);
01366             if (alambd > al) {
01367                imax=i;
01368                aiMAX=abi;
01369                alambd=al;
01370             }
01371          }
01372          // fAKAPPA - parameter will be <fEPS if fit is converged
01373          akap = TMath::Abs(fDA[i]/fParamError[i]);
01374          if (akap > fAKAPPA) fAKAPPA=akap;
01375       }
01376    }
01377    //... CALCULATE NEW CORRECTED STEP
01378    fGT = 0.;
01379    amb = 1.e18;
01380    // alambd - multiplier to split teoretical step dA
01381    if (alambd > .0) amb = 0.25/alambd;
01382    for( i = 0; i < n; i++) {
01383       if (fPL[i] > .0) {
01384          if (nn2 > fNlimMul ) {
01385             if (TMath::Abs(fDA[i]/fPL[i]) > amb ) {
01386                fPL[i] = 4.*fPL[i]; // increase parallelepiped
01387                t1=4.; // flag - that fPL was increased
01388             }
01389          }
01390          // cut step
01391          fDA[i] = fDA[i]*alambd;
01392          // expected function value change in next iteration
01393          fGT = fGT + fDA[i]*fGr[i];
01394       }
01395    }
01396 
01397    //.. CHECK IF MINIMUM ATTAINED AND SET EXIT MODE
01398    // if expected fGT smaller than precision
01399    // and other stuff
01400    if (-fGT <= sp && t1 < 1. && alambd < 1.)fENDFLG = -1; // function is not decreasing
01401 
01402    if (fENDFLG >= 0) {
01403       if (fAKAPPA < TMath::Abs(fEPS)) { // fit is converging
01404          if (fixFLG == 0)
01405             fENDFLG=1; // successful fit
01406          else {// we have fixed parameters
01407             if (fENDFLG == 0) {
01408                //... CHECK IF fiXING ON BOUND IS CORRECT
01409                fENDFLG = 1;
01410                fixFLG = 0;
01411                ifix1=-1;
01412                // release fixed parameters
01413                for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
01414                fINDFLG[1] = 0;
01415                // and repeat iteration
01416                goto L19;
01417             } else {
01418                if( ifix1 >= 0) {
01419                   fi = fi + 1;
01420                   fENDFLG = 0;
01421                }
01422             }
01423          }
01424       } else { // fit does not converge
01425          if( fixFLG != 0) {
01426             if( fi > fixFLG ) {
01427                //... CHECK IF fiXING ON BOUND IS CORRECT
01428                fENDFLG = 1;
01429                fixFLG = 0;
01430                ifix1=-1;
01431                for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
01432                fINDFLG[1] = 0;
01433                goto L19;
01434             } else {
01435                fi = fi + 1;
01436                fENDFLG = 0;
01437             }
01438          } else {
01439             fi = fi + 1;
01440             fENDFLG = 0;
01441          }
01442       }
01443    }
01444 
01445 // L85:
01446    // iteration number limit is exceeded
01447    if(fENDFLG == 0 && nn3 >= fNmaxIter) fENDFLG=-3;
01448 
01449    // fit errors are infinite;
01450    if(fENDFLG > 0 && fINDFLG[1] > 0) fENDFLG=-2;
01451 
01452    //MONITO (fS,fNpar,nn3,IT,fEPS,fGT,fAKAPPA,alambd);
01453    if (fENDFLG == 0) { // make step
01454       for ( i = 0; i < n; i++) fA[i] = fA[i] + fDA[i];
01455       if (imax >= 0) fA[imax] = aiMAX;
01456       olds=fS;
01457       nn2=nn2+1;
01458       nn3=nn3+1;
01459    } else {
01460       // fill covariant matrix VL
01461       // fill parameter error matrix up
01462       Int_t il;
01463       il = 0;
01464       for( Int_t ip = 0; ip < fNpar; ip++) {
01465          if( fPL0[ip] > .0) {
01466             for( Int_t jp = 0; jp <= ip; jp++) {
01467                if(fPL0[jp] > .0) {
01468                   //         VL[ind(ip,jp)] = fZ[il];
01469                   il = il + 1;
01470                }
01471             }
01472          }
01473       }
01474       return fENDFLG - 1;
01475    }
01476    goto L3;
01477 }
01478 
01479 //______________________________________________________________________________
01480 void TFumili::PrintResults(Int_t ikode,Double_t p) const
01481 {
01482    // Prints fit results.
01483    //
01484    // ikode is the type of printing parameters
01485    // p is function value
01486    //
01487    //  ikode = 1   - print values, errors and limits
01488    //  ikode = 2   - print values, errors and steps
01489    //  ikode = 3   - print values, errors, steps and derivatives
01490    //  ikode = 4   - print only values and errors
01491    //
01492    TString exitStatus="";
01493    TString xsexpl="";
01494    TString colhdu[3],colhdl[3],cx2,cx3;
01495    switch (fENDFLG) {
01496    case 1:
01497       exitStatus="CONVERGED";
01498       break;
01499    case -1:
01500       exitStatus="CONST FCN";
01501       xsexpl="****\n* FUNCTION IS NOT DECREASING OR BAD DERIVATIVES\n****";
01502       break;
01503    case -2:
01504       exitStatus="ERRORS INF";
01505       xsexpl="****\n* ESTIMATED ERRORS ARE INfiNITE\n****";
01506       break;
01507    case -3:
01508       exitStatus="MAX ITER.";
01509       xsexpl="****\n* MAXIMUM NUMBER OF ITERATIONS IS EXCEEDED\n****";
01510       break;
01511    case -4:
01512       exitStatus="ZERO PROBAB";
01513       xsexpl="****\n* PROBABILITY OF LIKLIHOOD FUNCTION IS NEGATIVE OR ZERO\n****";
01514       break;
01515    default:
01516       exitStatus="UNDEfiNED";
01517       xsexpl="****\n* fiT IS IN PROGRESS\n****";
01518       break;
01519    }
01520    if (ikode == 1) {
01521       colhdu[0] = "              ";
01522       colhdl[0] = "      ERROR   ";
01523       colhdu[1] = "      PHYSICAL";
01524       colhdu[2] = " LIMITS       ";
01525       colhdl[1] = "    NEGATIVE  ";
01526       colhdl[2] = "    POSITIVE  ";
01527    }
01528    if (ikode == 2) {
01529       colhdu[0] = "              ";
01530       colhdl[0] = "      ERROR   ";
01531       colhdu[1] = "    INTERNAL  ";
01532       colhdl[1] = "    STEP SIZE ";
01533       colhdu[2] = "    INTERNAL  ";
01534       colhdl[2] = "      VALUE   ";
01535    }
01536    if (ikode == 3) {
01537       colhdu[0] = "              ";
01538       colhdl[0] = "      ERROR   ";
01539       colhdu[1] = "       STEP   ";
01540       colhdl[1] = "       SIZE   ";
01541       colhdu[2] = "       fiRST  ";
01542       colhdl[2] = "    DERIVATIVE";
01543    }
01544    if (ikode == 4) {
01545       colhdu[0] = "    PARABOLIC ";
01546       colhdl[0] = "      ERROR   ";
01547       colhdu[1] = "        MINOS ";
01548       colhdu[2] = "ERRORS        ";
01549       colhdl[1] = "   NEGATIVE   ";
01550       colhdl[2] = "   POSITIVE   ";
01551    }
01552    if(fENDFLG<1)Printf("%s", xsexpl.Data());
01553    Printf(" FCN=%g FROM FUMILI  STATUS=%-10s %9d CALLS OF FCN",
01554          p,exitStatus.Data(),fNfcn);
01555    Printf(" EDM=%g ",-fGT);
01556    Printf("  EXT PARAMETER              %-14s%-14s%-14s",
01557          (const char*)colhdu[0].Data()
01558          ,(const char*)colhdu[1].Data()
01559          ,(const char*)colhdu[2].Data());
01560    Printf("  NO.   NAME          VALUE  %-14s%-14s%-14s",
01561          (const char*)colhdl[0].Data()
01562          ,(const char*)colhdl[1].Data()
01563          ,(const char*)colhdl[2].Data());
01564 
01565    for (Int_t i=0;i<fNpar;i++) {
01566       if (ikode==3) {
01567          cx2 = Form("%14.5e",fDA[i]);
01568          cx3 = Form("%14.5e",fGr[i]);
01569 
01570       }
01571       if (ikode==1) {
01572          cx2 = Form("%14.5e",fAMN[i]);
01573          cx3 = Form("%14.5e",fAMX[i]);
01574       }
01575       if (ikode==2) {
01576          cx2 = Form("%14.5e",fDA[i]);
01577          cx3 = Form("%14.5e",fA[i]);
01578       }
01579       if(ikode==4){
01580          cx2 = " *undefined*  ";
01581          cx3 = " *undefined*  ";
01582       }
01583       if(fPL0[i]<=0.) { cx2="    *fixed*   ";cx3=""; }
01584       Printf("%4d %-11s%14.5e%14.5e%-14s%-14s",i+1
01585            ,fANames[i].Data(),fA[i],fParamError[i]
01586            ,cx2.Data(),cx3.Data());
01587    }
01588 }
01589 
01590 
01591 //______________________________________________________________________________
01592 void TFumili::ReleaseParameter(Int_t ipar) {
01593    // Releases parameter number ipar
01594 
01595    if(ipar>=0 && ipar<fNpar && fPL0[ipar]<=0.) {
01596       fPL0[ipar] = -fPL0[ipar];
01597       if (fPL0[ipar] == 0. || fPL0[ipar]>=1.) fPL0[ipar]=.1;
01598    }
01599 }
01600 
01601 
01602 //______________________________________________________________________________
01603 void TFumili::SetData(Double_t *exdata,Int_t numpoints,Int_t vecsize){
01604    // Sets pointer to data array provided by user.
01605    // Necessary if SetFCN is not called.
01606    //
01607    // numpoints:    number of experimental points
01608    // vecsize:      size of data point vector + 2
01609    //               (for N-dimensional fit vecsize=N+2)
01610    // exdata:       data array with following format
01611    //
01612    //   exdata[0] = ExpValue_0     - experimental data value number 0
01613    //   exdata[1] = ExpSigma_0     - error of value number 0
01614    //   exdata[2] = X_0[0]
01615    //   exdata[3] = X_0[1]
01616    //       .........
01617    //   exdata[vecsize-1] = X_0[vecsize-3]
01618    //   exdata[vecsize]   = ExpValue_1
01619    //   exdata[vecsize+1] = ExpSigma_1
01620    //   exdata[vecsize+2] = X_1[0]
01621    //       .........
01622    //   exdata[vecsize*(numpoints-1)] = ExpValue_(numpoints-1)
01623    //       .........
01624    //   exdata[vecsize*numpoints-1] = X_(numpoints-1)[vecsize-3]
01625    //
01626 
01627    if(exdata){
01628       fNED1 = numpoints;
01629       fNED2 = vecsize;
01630       fEXDA = exdata;
01631    }
01632 }
01633 
01634 
01635 //______________________________________________________________________________
01636 void TFumili::SetFitMethod(const char *name)
01637 {
01638    // ret fit method (chisquare or loglikelihood)
01639 
01640    if (!strcmp(name,"H1FitChisquare"))    SetFCN(H1FitChisquareFumili);
01641    if (!strcmp(name,"H1FitLikelihood"))   SetFCN(H1FitLikelihoodFumili);
01642    if (!strcmp(name,"GraphFitChisquare")) SetFCN(GraphFitChisquareFumili);
01643 }
01644 
01645 
01646 //______________________________________________________________________________
01647 Int_t TFumili::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh) {
01648    // Sets for prameter number ipar initial parameter value,
01649    // name parname, initial error verr and limits vlow and vhigh
01650    // If vlow = vhigh but not equil to zero, parameter will be fixed.
01651    // If vlow = vhigh = 0, parameter is released and its limits are discarded
01652    //
01653    if (ipar<0 || ipar>=fNpar) return -1;
01654    fANames[ipar] = parname;
01655    fA[ipar] = value;
01656    fParamError[ipar] = verr;
01657    if(vlow<vhigh) {
01658       fAMN[ipar] = vlow;
01659       fAMX[ipar] = vhigh;
01660    } else {
01661       if(vhigh<vlow) {
01662          fAMN[ipar] = vhigh;
01663          fAMX[ipar] = vlow;
01664       }
01665       if(vhigh==vlow) {
01666          if(vhigh==0.) {
01667             ReleaseParameter(ipar);
01668             fAMN[ipar] = gMINDOUBLE;
01669             fAMX[ipar] = gMAXDOUBLE;
01670          }
01671          if(vlow!=0) FixParameter(ipar);
01672       }
01673    }
01674    return 0;
01675 }
01676 
01677 //______________________________________________________________________________
01678 Int_t TFumili::SGZ()
01679 {
01680    //  Evaluates objective function ( chi-square ), gradients and
01681    //  Z-matrix using data provided by user via TFumili::SetData
01682    //
01683    fS = 0.;
01684    Int_t i,j,l,k2=1,k1,ki=0;
01685    Double_t *x  = new Double_t[fNED2];
01686    Double_t *df = new Double_t[fNpar];
01687    Int_t nx = fNED2-2;
01688    for (l=0;l<fNED1;l++) { // cycle on all exp. points
01689       k1 = k2;
01690       if (fLogLike) {
01691          fNumericDerivatives = kTRUE;
01692          nx  = fNED2;
01693          k1 -= 2;
01694       }
01695 
01696       for (i=0;i<nx;i++){
01697          ki  += 1+i;
01698          x[i] = fEXDA[ki];
01699       }
01700       //  Double_t y = ARITHM(df,x);
01701       Double_t y = EvalTFN(df,x);
01702       if(fNumericDerivatives) Derivatives(df,x);
01703       Double_t sig=1.;
01704       if(fLogLike) { // Likelihood method
01705          if(y>0.) {
01706             fS = fS - log(y);
01707             y  = -y;
01708             sig= y;
01709          } else { //
01710             delete [] x;
01711             delete [] df;
01712             fS = 1e10;
01713             return -1; // indflg[0] = 1;
01714          }
01715       } else { // Chi2 method
01716          sig = fEXDA[k2]; // sigma of experimental point
01717          y = y - fEXDA[k1-1]; // f(x_i) - F_i
01718          fS = fS + (y*y/(sig*sig))*.5; // simple chi2/2
01719       }
01720       Int_t n = 0;
01721       for (i=0;i<fNpar;i++) {
01722          if (fPL0[i]>0){
01723             df[n]   = df[i]/sig; // left only non-fixed param derivatives div by Sig
01724             fGr[i] += df[n]*(y/sig);
01725             n++;
01726          }
01727       }
01728       l = 0;
01729       for (i=0;i<n;i++)
01730          for (j=0;j<=i;j++)
01731             fZ[l++] += df[i]*df[j];
01732       k2 += fNED2;
01733    }
01734 
01735    delete[] df;
01736    delete[] x;
01737    return 1;
01738 }
01739 
01740 
01741 //______________________________________________________________________________
01742 void TFumili::FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01743 {
01744    //  Minimization function for H1s using a Chisquare method
01745    //  Default method (function evaluated at center of bin)
01746    //  for each point the cache contains the following info
01747    //    -1D : bc,e,xc  (bin content, error, x of center of bin)
01748    //    -2D : bc,e,xc,yc
01749    //    -3D : bc,e,xc,yc,zc
01750 
01751    Foption_t fitOption = GetFitOption();
01752    if (fitOption.Integral) {
01753       FitChisquareI(npar,gin,f,u,flag);
01754       return;
01755    }
01756    Double_t cu,eu,fu,fsum;
01757    Double_t x[3];
01758    Double_t *zik=0;
01759    Double_t *pl0=0;
01760 
01761    TH1 *hfit = (TH1*)GetObjectFit();
01762    TF1 *f1   = (TF1*)GetUserFunc();
01763    Int_t nd  = hfit->GetDimension();
01764    Int_t j;
01765 
01766    npar = f1->GetNpar();
01767    SetParNumber(npar);
01768    if(flag == 9) return;
01769    zik = GetZ();
01770    pl0 = GetPL0();
01771 
01772    Double_t *df = new Double_t[npar];
01773    f1->InitArgs(x,u);
01774    f = 0;
01775 
01776    Int_t npfit = 0;
01777    Double_t *cache = fCache;
01778    for (Int_t i=0;i<fNpoints;i++) {
01779       if (nd > 2) x[2]     = cache[4];
01780       if (nd > 1) x[1]     = cache[3];
01781       x[0]     = cache[2];
01782       cu  = cache[0];
01783       TF1::RejectPoint(kFALSE);
01784       fu = f1->EvalPar(x,u);
01785       if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01786       eu = cache[1];
01787       Derivatives(df,x);
01788       Int_t n = 0;
01789       fsum = (fu-cu)/eu;
01790       if (flag!=1) {
01791          for (j=0;j<npar;j++) {
01792             if (pl0[j]>0) {
01793                df[n] = df[j]/eu;
01794                // left only non-fixed param derivatives / by Sigma
01795                gin[j] += df[n]*fsum;
01796                n++;
01797             }
01798          }
01799          Int_t l = 0;
01800          for (j=0;j<n;j++)
01801             for (Int_t k=0;k<=j;k++)
01802                zik[l++] += df[j]*df[k];
01803       }
01804       f += .5*fsum*fsum;
01805       npfit++;
01806       cache += fPointSize;
01807    }
01808    f1->SetNumberFitPoints(npfit);
01809    delete [] df;
01810 }
01811 
01812 
01813 //______________________________________________________________________________
01814 void TFumili::FitChisquareI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01815 {
01816    //  Minimization function for H1s using a Chisquare method
01817    //  The "I"ntegral method is used
01818    //  for each point the cache contains the following info
01819    //    -1D : bc,e,xc,xw  (bin content, error, x of center of bin, x bin width of bin)
01820    //    -2D : bc,e,xc,xw,yc,yw
01821    //    -3D : bc,e,xc,xw,yc,yw,zc,zw
01822 
01823    Double_t cu,eu,fu,fsum;
01824    Double_t x[3];
01825    Double_t *zik=0;
01826    Double_t *pl0=0;
01827 
01828    TH1 *hfit = (TH1*)GetObjectFit();
01829    TF1 *f1   = (TF1*)GetUserFunc();
01830    Int_t nd  = hfit->GetDimension();
01831    Int_t j;
01832 
01833    f1->InitArgs(x,u);
01834    npar = f1->GetNpar();
01835    SetParNumber(npar);
01836    if(flag == 9) return;
01837    zik = GetZ();
01838    pl0 = GetPL0();
01839 
01840    Double_t *df=new Double_t[npar];
01841    f = 0;
01842 
01843    Int_t npfit = 0;
01844    Double_t *cache = fCache;
01845    for (Int_t i=0;i<fNpoints;i++) {
01846       cu  = cache[0];
01847       TF1::RejectPoint(kFALSE);
01848       f1->SetParameters(u);
01849       if (nd < 2) {
01850          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
01851       } else if (nd < 3) {
01852          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
01853       } else {
01854          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
01855       }
01856       if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01857       eu = cache[1];
01858       Derivatives(df,x);
01859       Int_t n = 0;
01860       fsum = (fu-cu)/eu;
01861       if (flag!=1) {
01862          for (j=0;j<npar;j++) {
01863             if (pl0[j]>0){
01864                df[n] = df[j]/eu;
01865                // left only non-fixed param derivatives / by Sigma
01866                gin[j] += df[n]*fsum;
01867                n++;
01868             }
01869          }
01870          Int_t l = 0;
01871          for (j=0;j<n;j++)
01872             for (Int_t k=0;k<=j;k++)
01873                zik[l++] += df[j]*df[k];
01874       }
01875       f += .5*fsum*fsum;
01876       npfit++;
01877       cache += fPointSize;
01878    }
01879    f1->SetNumberFitPoints(npfit);
01880    delete[] df;
01881 }
01882 
01883 
01884 //______________________________________________________________________________
01885 void TFumili::FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01886 {
01887    //  Minimization function for H1s using a Likelihood method*-*-*-*-*-*
01888    //     Basically, it forms the likelihood by determining the Poisson
01889    //     probability that given a number of entries in a particular bin,
01890    //     the fit would predict it's value.  This is then done for each bin,
01891    //     and the sum of the logs is taken as the likelihood.
01892    //  Default method (function evaluated at center of bin)
01893    //  for each point the cache contains the following info
01894    //    -1D : bc,e,xc  (bin content, error, x of center of bin)
01895    //    -2D : bc,e,xc,yc
01896    //    -3D : bc,e,xc,yc,zc
01897 
01898    Foption_t fitOption = GetFitOption();
01899    if (fitOption.Integral) {
01900       FitLikelihoodI(npar,gin,f,u,flag);
01901       return;
01902    }
01903    Double_t cu,fu,fobs,fsub;
01904    Double_t dersum[100];
01905    Double_t x[3];
01906    Int_t icu;
01907 
01908    TH1 *hfit = (TH1*)GetObjectFit();
01909    TF1 *f1   = (TF1*)GetUserFunc();
01910    Int_t nd = hfit->GetDimension();
01911    Int_t j;
01912    Double_t *zik = GetZ();
01913    Double_t *pl0 = GetPL0();
01914 
01915    npar = f1->GetNpar();
01916    SetParNumber(npar);
01917    if(flag == 9) return;
01918    Double_t *df=new Double_t[npar];
01919    if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
01920    f1->InitArgs(x,u);
01921    f = 0;
01922 
01923    Int_t npfit = 0;
01924    Double_t *cache = fCache;
01925    for (Int_t i=0;i<fNpoints;i++) {
01926       if (nd > 2) x[2] = cache[4];
01927       if (nd > 1) x[1] = cache[3];
01928       x[0]     = cache[2];
01929       cu  = cache[0];
01930       TF1::RejectPoint(kFALSE);
01931       fu = f1->EvalPar(x,u);
01932       if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01933       if (flag == 2) {
01934          for (j=0;j<npar;j++) {
01935             dersum[j] += 1; //should be the derivative
01936             //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
01937          }
01938       }
01939       if (fu < 1.e-9) fu = 1.e-9;
01940       icu  = Int_t(cu);
01941       fsub = -fu +icu*TMath::Log(fu);
01942       fobs = GetSumLog(icu);
01943       fsub -= fobs;
01944       Derivatives(df,x);
01945       int n=0;
01946       // Here we need gradients of Log likelihood function
01947       //
01948       for (j=0;j<npar;j++) {
01949          if (pl0[j]>0){
01950             df[n]   = df[j]*(icu/fu-1);
01951             gin[j] -= df[n];
01952             n++;
01953          }
01954       }
01955       Int_t l = 0;
01956       // Z-matrix here - production of first derivatives
01957       //  of log-likelihood function
01958       for (j=0;j<n;j++)
01959          for (Int_t k=0;k<=j;k++)
01960             zik[l++] += df[j]*df[k];
01961 
01962       f -= fsub;
01963       npfit++;
01964       cache += fPointSize;
01965    }
01966    f *= 2;
01967    f1->SetNumberFitPoints(npfit);
01968    delete[] df;
01969 }
01970 
01971 
01972 //______________________________________________________________________________
01973 void TFumili::FitLikelihoodI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01974 {
01975    //  Minimization function for H1s using a Likelihood method*-*-*-*-*-*
01976    //     Basically, it forms the likelihood by determining the Poisson
01977    //     probability that given a number of entries in a particular bin,
01978    //     the fit would predict it's value.  This is then done for each bin,
01979    //     and the sum of the logs is taken as the likelihood.
01980    //  The "I"ntegral method is used
01981    //  for each point the cache contains the following info
01982    //    -1D : bc,e,xc,xw  (bin content, error, x of center of bin, x bin width of bin)
01983    //    -2D : bc,e,xc,xw,yc,yw
01984    //    -3D : bc,e,xc,xw,yc,yw,zc,zw
01985 
01986    Double_t cu,fu,fobs,fsub;
01987    Double_t dersum[100];
01988    Double_t x[3];
01989    Int_t icu;
01990 
01991    TH1 *hfit = (TH1*)GetObjectFit();
01992    TF1 *f1   = (TF1*)GetUserFunc();
01993    Int_t nd = hfit->GetDimension();
01994    Int_t j;
01995    Double_t *zik = GetZ();
01996    Double_t *pl0 = GetPL0();
01997 
01998    Double_t *df=new Double_t[npar];
01999 
02000    npar = f1->GetNpar();
02001    SetParNumber(npar);
02002    if(flag == 9) {delete [] df; return;}
02003    if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
02004    f1->InitArgs(x,u);
02005    f = 0;
02006 
02007    Int_t npfit = 0;
02008    Double_t *cache = fCache;
02009    for (Int_t i=0;i<fNpoints;i++) {
02010       if (nd > 2) x[2] = cache[4];
02011       if (nd > 1) x[1] = cache[3];
02012       x[0]     = cache[2];
02013       cu  = cache[0];
02014       TF1::RejectPoint(kFALSE);
02015       if (nd < 2) {
02016          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
02017       } else if (nd < 3) {
02018          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
02019       } else {
02020          fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
02021       }
02022       if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
02023       if (flag == 2) {
02024          for (j=0;j<npar;j++) {
02025             dersum[j] += 1; //should be the derivative
02026             //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
02027          }
02028       }
02029       if (fu < 1.e-9) fu = 1.e-9;
02030       icu  = Int_t(cu);
02031       fsub = -fu +icu*TMath::Log(fu);
02032       fobs = GetSumLog(icu);
02033       fsub -= fobs;
02034       Derivatives(df,x);
02035       int n=0;
02036       // Here we need gradients of Log likelihood function
02037       //
02038       for (j=0;j<npar;j++) {
02039          if (pl0[j]>0){
02040             df[n]   = df[j]*(icu/fu-1);
02041             gin[j] -= df[n];
02042             n++;
02043          }
02044       }
02045       Int_t l = 0;
02046       // Z-matrix here - production of first derivatives
02047       //  of log-likelihood function
02048       for (j=0;j<n;j++)
02049          for (Int_t k=0;k<=j;k++)
02050             zik[l++] += df[j]*df[k];
02051 
02052       f -= fsub;
02053       npfit++;
02054       cache += fPointSize;
02055    }
02056    f *= 2;
02057    f1->SetNumberFitPoints(npfit);
02058    delete[] df;
02059 }
02060 
02061 
02062 //______________________________________________________________________________
02063 //
02064 //  STATIC functions
02065 //______________________________________________________________________________
02066 
02067 //______________________________________________________________________________
02068 void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
02069 {
02070 //           Minimization function for H1s using a Chisquare method
02071 //           ======================================================
02072 
02073    TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
02074    hFitter->FitChisquare(npar, gin, f, u, flag);
02075 }
02076 
02077 //______________________________________________________________________________
02078 void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
02079 {
02080 //   -*-*-*-*Minimization function for H1s using a Likelihood method*-*-*-*-*-*
02081 //           =======================================================
02082 //     Basically, it forms the likelihood by determining the Poisson
02083 //     probability that given a number of entries in a particular bin,
02084 //     the fit would predict it's value.  This is then done for each bin,
02085 //     and the sum of the logs is taken as the likelihood.
02086 //     PDF:  P=exp(-f(x_i))/[F_i]!*(f(x_i))^[F_i]
02087 //    where F_i - experimental value, f(x_i) - expected theoretical value
02088 //    [F_i] - integer part of F_i.
02089 //    drawback is that if F_i>Int_t - GetSumLog will fail
02090 //    for big F_i is faster to use Euler's Gamma-function
02091 
02092 
02093    TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
02094    hFitter->FitLikelihood(npar, gin, f, u, flag);
02095 }
02096 
02097 
02098 
02099 //______________________________________________________________________________
02100 void GraphFitChisquareFumili(Int_t &npar, Double_t * gin, Double_t &f,
02101                        Double_t *u, Int_t flag)
02102 {
02103 //*-*-*-*-*-*Minimization function for Graphs using a Chisquare method*-*-*-*-*
02104 //*-*        =========================================================
02105 //
02106 // In case of a TGraphErrors object, ex, the error along x,  is projected
02107 // along the y-direction by calculating the function at the points x-exlow and
02108 // x+exhigh.
02109 //
02110 // The chisquare is computed as the sum of the quantity below at each point:
02111 //
02112 //                     (y - f(x))**2
02113 //         -----------------------------------
02114 //         ey**2 + (0.5*(exl + exh)*f'(x))**2
02115 //
02116 // where x and y are the point coordinates and f'(x) is the derivative of function f(x).
02117 // This method to approximate the uncertainty in y because of the errors in x, is called
02118 // "effective variance" method.
02119 // The improvement, compared to the previously used  method (f(x+ exhigh) - f(x-exlow))/2
02120 // is of (error of x)**2 order.
02121 //  NOTE:
02122 //  1) By using the "effective variance" method a simple linear regression
02123 //      becomes a non-linear case , which takes several iterations
02124 //      instead of 0 as in the linear case .
02125 //
02126 //  2) The effective variance technique assumes that there is no correlation
02127 //      between the x and y coordinate .
02128 //
02129 // In case the function lies below (above) the data point, ey is ey_low (ey_high).
02130 
02131    Double_t cu,eu,exl,exh,ey,eux,fu,fsum;
02132    Double_t x[1];
02133    Int_t i, bin, npfits=0;
02134 
02135    TFumili *grFitter = (TFumili*)TVirtualFitter::GetFitter();
02136    TGraph *gr     = (TGraph*)grFitter->GetObjectFit();
02137    TF1 *f1   = (TF1*)grFitter->GetUserFunc();
02138    Foption_t fitOption = grFitter->GetFitOption();
02139 
02140    Int_t n        = gr->GetN();
02141    Double_t *gx   = gr->GetX();
02142    Double_t *gy   = gr->GetY();
02143    npar           = f1->GetNpar();
02144 
02145    grFitter->SetParNumber(npar);
02146 
02147    if(flag == 9) return;
02148    Double_t *zik = grFitter->GetZ();
02149    Double_t *pl0 = grFitter->GetPL0();
02150    Double_t *df  = new Double_t[npar];
02151 
02152 
02153    f1->InitArgs(x,u);
02154    f      = 0;
02155    for (bin=0;bin<n;bin++) {
02156       x[0] = gx[bin];
02157       if (!f1->IsInside(x)) continue;
02158       cu   = gy[bin];
02159       TF1::RejectPoint(kFALSE);
02160       fu   = f1->EvalPar(x,u);
02161       if (TF1::RejectedPoint()) continue;
02162       npfits++;
02163       Double_t eusq=1.;
02164       if (fitOption.W1) {
02165          eu = 1.;
02166       } else {
02167          exh  = gr->GetErrorXhigh(bin);
02168          exl  = gr->GetErrorXlow(bin);
02169          ey  = gr->GetErrorY(bin);
02170          if (exl < 0) exl = 0;
02171          if (exh < 0) exh = 0;
02172          if (ey < 0)  ey  = 0;
02173          if (exh > 0 && exl > 0) {
02174 //          "Effective variance" method for projecting errors
02175             eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
02176          } else
02177             eux = 0.;
02178          eu = ey*ey+eux*eux;
02179          if (eu <= 0) eu = 1;
02180          eusq = TMath::Sqrt(eu);
02181       }
02182       grFitter->Derivatives(df,x);
02183       n = 0;
02184       fsum = (fu-cu)/eusq;
02185       for (i=0;i<npar;i++) {
02186          if (pl0[i]>0){
02187             df[n] = df[i]/eusq;
02188             // left only non-fixed param derivatives / by Sigma
02189             gin[i] += df[n]*fsum;
02190             n++;
02191          }
02192       }
02193       Int_t l = 0;
02194       for (i=0;i<n;i++)
02195          for (Int_t j=0;j<=i;j++)
02196             zik[l++] += df[i]*df[j];
02197       f += .5*fsum*fsum;
02198 
02199    }
02200    delete [] df;
02201    f1->SetNumberFitPoints(npfits);
02202 }
02203 

Generated on Tue Jul 5 14:30:50 2011 for ROOT_528-00b_version by  doxygen 1.5.1