ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TFumili.cxx
Go to the documentation of this file.
1 // @(#)root/fumili:$Id$
2 // Author: Stanislav Nesterov 07/05/2003
3 
4 //______________________________________________________________________________
5 // FUMILI
6 // Based on ideas, proposed by I.N. Silin
7 // [See NIM A440, 2000 (p431)]
8 // converted from FORTRAN to C by
9 // Sergey Yaschenko <s.yaschenko@fz-juelich.de>
10 //
11 //______________________________________________________________________________
12 //BEGIN_HTML <!--
13 /* -->
14 <H2>FUMILI minimization package</H2>
15 <p>FUMILI is used to minimize Chi-square function or to search maximum of
16 likelihood function.
17 
18 <p>Experimentally measured values $F_i$ are fitted with theoretical
19 functions $f_i({\vec x}_i,\vec\theta\,\,)$, where ${\vec x}_i$ are
20 coordinates, and $\vec\theta$ -- vector of parameters.
21 
22 <p>For better convergence Chi-square function has to be the following form
23 
24 <p>$$
25 {\chi^2\over2}={1\over2}\sum^n_{i=1}\left(f_i(\vec
26 x_i,\vec\theta\,\,)-F_i\over\sigma_i\right)^2 \eqno(1)
27 $$
28 <p>where $\sigma_i$ are errors of measured function.
29 
30 <p>The minimum condition is
31 <p>$$
32 {\partial\chi^2\over\partial\theta_i}=\sum^n_{j=1}{1\over\sigma^2_j}\cdot
33 {\partial f_j\over\partial\theta_i}\left[f_j(\vec
34 x_j,\vec\theta\,\,)-F_j\right]=0,\qquad i=1\ldots m\eqno(2)
35 $$
36 <p>where m is the quantity of parameters.
37 
38 <p>Expanding left part of (2) over parameter increments and
39 retaining only linear terms one gets
40 <p>$$
41 \left(\partial\chi^2\over\theta_i\right)_{\vec\theta={\vec\theta}^0}
42 +\sum_k\left(\partial^2\chi^2\over\partial\theta_i\partial\theta_k\right)_{
43 \vec\theta={\vec\theta}^0}\cdot(\theta_k-\theta_k^0)
44 = 0\eqno(3)
45 $$
46 
47  <p>Here ${\vec\theta}_0$ is some initial value of parameters. In general
48 case:
49 <p>$$
50 {\partial^2\chi^2\over\partial\theta_i\partial\theta_k}=
51 \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i}
52 {\partial f_j\over\theta_k} +
53 \sum^n_{j=1}{(f_j - F_j)\over\sigma^2_j}\cdot
54 {\partial^2f_j\over\partial\theta_i\partial\theta_k}\eqno(4)
55 $$
56 
57 <p>In FUMILI algorithm for second derivatives of Chi-square approximate
58 expression is used when last term in (4) is discarded. It is often
59 done, not always wittingly, and sometimes causes troubles, for example,
60 if user wants to limit parameters with positive values by writing down
61 $\theta_i^2$ instead of $\theta_i$. FUMILI will fail if one tries
62 minimize $\chi^2 = g^2(\vec\theta)$ where g is arbitrary function.
63 
64 <p>Approximate value is:
65 <p>$${\partial^2\chi^2\over\partial\theta_i\partial\theta_k}\approx
66 Z_{ik}=
67 \sum^n_{j=1}{1\over\sigma^2_j}{\partial f_j\over\theta_i}
68 {\partial f_j\over\theta_k}\eqno(5)
69 $$
70 
71 <p>Then the equations for parameter increments are
72 <p>$$\left(\partial\chi^2\over\partial\theta_i\right)_{\vec\theta={\vec\theta}^0}
73 +\sum_k Z_{ik}\cdot(\theta_k-\theta^0_k) = 0,
74 \qquad i=1\ldots m\eqno(6)
75 $$
76 
77 <p>Remarkable feature of algorithm is the technique for step
78 restriction. For an initial value of parameter ${\vec\theta}^0$ a
79 parallelepiped $P_0$ is built with the center at ${\vec\theta}^0$ and
80 axes parallel to coordinate axes $\theta_i$. The lengths of
81 parallelepiped sides along i-th axis is $2b_i$, where $b_i$ is such a
82 value that the functions $f_j(\vec\theta)$ are quasi-linear all over
83 the parallelepiped.
84 
85 <p>FUMILI takes into account simple linear inequalities in the form:
86 $$
87 \theta_i^{\rm min}\le\theta_i\le\theta^{\rm max}_i\eqno(7)
88 $$
89 
90 <p>They form parallelepiped $P$ ($P_0$ may be deformed by $P$).
91 Very similar step formulae are used in FUMILI for negative logarithm
92 of the likelihood function with the same idea - linearization of
93 function argument.
94 
95 <!--*/
96 // -->END_HTML
97 //______________________________________________________________________________
98 
99 
100 
101 #include "TROOT.h"
102 #include "TFumili.h"
103 #include "TMath.h"
104 #include "TF1.h"
105 #include "TF2.h"
106 #include "TF3.h"
107 #include "TH1.h"
108 #include "TGraphAsymmErrors.h"
109 
110 #include "Riostream.h"
111 
112 
113 extern void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
114 extern void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
115 extern void GraphFitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
116 
117 
119 
120 TFumili *gFumili=0;
121 // Machine dependent values fiXME!!
122 // But don't set min=max=0 if param is unlimited
123 static const Double_t gMAXDOUBLE=1e300;
124 static const Double_t gMINDOUBLE=-1e300;
125 
126 
127 ////////////////////////////////////////////////////////////////////////////////
128 
130 {//----------- FUMILI constructor ---------
131  // maxpar is the maximum number of parameters used with TFumili object
132  //
133  fMaxParam = TMath::Max(maxpar,25);
134  if (fMaxParam>200) fMaxParam=200;
135  BuildArrays();
136 
137  fNumericDerivatives = true;
138  fLogLike = false;
139  fNpar = fMaxParam;
140  fGRAD = false;
141  fWARN = true;
142  fDEBUG = false;
143  fNlog = 0;
144  fSumLog = 0;
145  fNED1 = 0;
146  fNED2 = 0;
147  fNED12 = fNED1+fNED2;
148  fEXDA = 0;
149  fFCN = 0;
150  fNfcn = 0;
151  fRP = 1.e-15; //precision
152  fS = 1e10;
153  fEPS =0.01;
154  fENDFLG = 0;
155  fNlimMul = 2;
156  fNmaxIter= 150;
157  fNstepDec= 3;
158  fLastFixed = -1;
159 
160  fAKAPPA = 0.;
161  fGT = 0.;
162  for (int i = 0; i<5; ++i) fINDFLG[i] = 0;
163 
164  SetName("Fumili");
165  gFumili = this;
166  gROOT->GetListOfSpecials()->Add(gFumili);
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 ///
171 /// Allocates memory for internal arrays. Called by TFumili::TFumili
172 ///
173 
175  fCmPar = new Double_t[fMaxParam];
176  fA = new Double_t[fMaxParam];
177  fPL0 = new Double_t[fMaxParam];
178  fPL = new Double_t[fMaxParam];
180  fDA = new Double_t[fMaxParam];
181  fAMX = new Double_t[fMaxParam];
182  fAMN = new Double_t[fMaxParam];
183  fR = new Double_t[fMaxParam];
184  fDF = new Double_t[fMaxParam];
185  fGr = new Double_t[fMaxParam];
186  fANames = new TString[fMaxParam];
187 
188  // fX = new Double_t[10];
189 
190  Int_t zSize = fMaxParam*(fMaxParam+1)/2;
191  fZ0 = new Double_t[zSize];
192  fZ = new Double_t[zSize];
193 
194  for (Int_t i=0;i<fMaxParam;i++) {
195  fA[i] =0.;
196  fDF[i]=0.;
197  fAMN[i]=gMINDOUBLE;
198  fAMX[i]=gMAXDOUBLE;
199  fPL0[i]=.1;
200  fPL[i] =.1;
201  fParamError[i]=0.;
202  fANames[i]=Form("%d",i);
203  }
204 }
205 
206 
207 ////////////////////////////////////////////////////////////////////////////////
208 ///
209 /// TFumili destructor
210 ///
211 
213  DeleteArrays();
214  gROOT->GetListOfSpecials()->Remove(this);
215  if (gFumili == this) gFumili = 0;
216 }
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 /// return a chisquare equivalent
220 
222 {
223  Double_t amin = 0;
224  H1FitChisquareFumili(npar,params,amin,params,1);
225  return 2*amin;
226 }
227 
228 
229 ////////////////////////////////////////////////////////////////////////////////
230 ///
231 /// Resets all parameter names, values and errors to zero
232 ///
233 /// Argument opt is ignored
234 ///
235 /// NB: this procedure doesn't reset parameter limits
236 ///
237 
239 {
240  fNpar = fMaxParam;
241  fNfcn = 0;
242  for (Int_t i=0;i<fNpar;i++) {
243  fA[i] =0.;
244  fDF[i] =0.;
245  fPL0[i] =.1;
246  fPL[i] =.1;
247  fAMN[i] = gMINDOUBLE;
248  fAMX[i] = gMAXDOUBLE;
249  fParamError[i]=0.;
250  fANames[i]=Form("%d",i);
251  }
252 }
253 
254 
255 ////////////////////////////////////////////////////////////////////////////////
256 ///
257 /// Deallocates memory. Called from destructor TFumili::~TFumili
258 ///
259 
261  delete[] fCmPar;
262  delete[] fANames;
263  delete[] fDF;
264  // delete[] fX;
265  delete[] fZ0;
266  delete[] fZ;
267  delete[] fGr;
268  delete[] fA;
269  delete[] fPL0;
270  delete[] fPL;
271  delete[] fDA;
272  delete[] fAMN;
273  delete[] fAMX;
274  delete[] fParamError;
275  delete[] fR;
276 }
277 
278 
279 ////////////////////////////////////////////////////////////////////////////////
280 ///
281 /// Calculates partial derivatives of theoretical function
282 ///
283 /// Input:
284 /// fX - vector of data point
285 /// Output:
286 /// DF - array of derivatives
287 ///
288 /// ARITHM.F
289 /// Converted from CERNLIB
290 ///
291 
293  Double_t ff,ai,hi,y,pi;
294  y = EvalTFN(df,fX);
295  for (Int_t i=0;i<fNpar;i++) {
296  df[i]=0;
297  if(fPL0[i]>0.) {
298  ai = fA[i]; // save current parameter value
299  hi = 0.01*fPL0[i]; // diff step
300  pi = fRP*TMath::Abs(ai);
301  if (hi<pi) hi = pi; // if diff step is less than precision
302  fA[i] = ai+hi;
303 
304  if (fA[i]>fAMX[i]) { // if param is out of limits
305  fA[i] = ai-hi;
306  hi = -hi;
307  if (fA[i]<fAMN[i]) { // again out of bounds
308  fA[i] = fAMX[i]; // set param to high limit
309  hi = fAMX[i]-ai;
310  if (fAMN[i]-ai+hi<0) { // if hi < (ai-fAMN)
311  fA[i]=fAMN[i];
312  hi=fAMN[i]-ai;
313  }
314  }
315  }
316  ff = EvalTFN(df,fX);
317  df[i] = (ff-y)/hi;
318  fA[i] = ai;
319  }
320  }
321 }
322 
323 
324 
325 ////////////////////////////////////////////////////////////////////////////////
326 /// Evaluate the minimisation function
327 /// Input parameters:
328 /// npar: number of currently variable parameters
329 /// par: array of (constant and variable) parameters
330 /// flag: Indicates what is to be calculated
331 /// grad: array of gradients
332 /// Output parameters:
333 /// fval: The calculated function value.
334 /// grad: The vector of first derivatives.
335 ///
336 /// The meaning of the parameters par is of course defined by the user,
337 /// who uses the values of those parameters to calculate their function value.
338 /// The starting values must be specified by the user.
339 /// !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
340 /// Inside FCN user has to define Z-matrix by means TFumili::GetZ
341 /// and TFumili::Derivatives,
342 /// set theoretical function by means of TFumili::SetUserFunc,
343 /// but first - pass number of parameters by TFumili::SetParNumber
344 /// !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345 /// Later values are determined by Fumili as it searches for the minimum
346 /// or performs whatever analysis is requested by the user.
347 ///
348 /// The default function calls the function specified in SetFCN
349 ///
350 
352 {
353  if (fFCN) (*fFCN)(npar,grad,fval,par,flag);
354  return npar;
355 }
356 
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 /// Evaluate theoretical function
360 /// df: array of partial derivatives
361 /// X: vector of theoretical function argument
362 
364 {
365  // for the time being disable possibility to compute derivatives
366  //if(fTFN)
367  // return (*fTFN)(df,X,fA);
368  //else if(fTFNF1) {
369 
370  TF1 *f1 = (TF1*)fUserFunc;
371  return f1->EvalPar(X,fA);
372  //}
373  //return 0.;
374 }
375 
376 ////////////////////////////////////////////////////////////////////////////////
377 ///
378 /// Execute MINUIT commands. MINImize, SIMplex, MIGrad and FUMili all
379 /// will call TFumili::Minimize method.
380 ///
381 /// For full command list see
382 /// MINUIT. Reference Manual. CERN Program Library Long Writeup D506.
383 ///
384 /// Improvement and errors calculation are not yet implemented as well
385 /// as Monte-Carlo seeking and minimization.
386 /// Contour commands are also unsupported.
387 ///
388 /// command : command string
389 /// args : array of arguments
390 /// nargs : number of arguments
391 ///
392 
393 Int_t TFumili::ExecuteCommand(const char *command, Double_t *args, Int_t nargs){
394  TString comand = command;
395  static TString clower = "abcdefghijklmnopqrstuvwxyz";
396  static TString cupper = "ABCDEFGHIJKLMNOPQRSTUVWXYZ";
397  const Int_t nntot = 40;
398  const char *cname[nntot] = {
399  "MINImize ", // 0 checked
400  "SEEk ", // 1 none
401  "SIMplex ", // 2 checked same as 0
402  "MIGrad ", // 3 checked same as 0
403  "MINOs ", // 4 none
404  "SET xxx ", // 5 lot of stuff
405  "SHOw xxx ", // 6 -----------
406  "TOP of pag", // 7 .
407  "fiX ", // 8 .
408  "REStore ", // 9 .
409  "RELease ", // 10 .
410  "SCAn ", // 11 not yet implemented
411  "CONtour ", // 12 not yet implemented
412  "HESse ", // 13 not yet implemented
413  "SAVe ", // 14 obsolete
414  "IMProve ", // 15 not yet implemented
415  "CALl fcn ", // 16 .
416  "STAndard ", // 17 .
417  "END ", // 18 .
418  "EXIt ", // 19 .
419  "RETurn ", // 20 .
420  "CLEar ", // 21 .
421  "HELP ", // 22 not yet implemented
422  "MNContour ", // 23 not yet implemented
423  "STOp ", // 24 .
424  "JUMp ", // 25 not yet implemented
425  " ", //
426  " ", //
427  "FUMili ", // 28 checked same as 0
428  " ", //
429  " ", //
430  " ", //
431  " ", //
432  "COVARIANCE", // 33
433  "PRINTOUT ", // 34
434  "GRADIENT ", // 35
435  "MATOUT ", // 36
436  "ERROR DEF ", // 37
437  "LIMITS ", // 38
438  "PUNCH "}; // 39
439 
440 
441  fCword = comand;
442  fCword.ToUpper();
443  if (nargs<=0) fCmPar[0] = 0;
444  Int_t i;
445  for(i=0;i<fMaxParam;i++) {
446  if(i<nargs) fCmPar[i] = args[i];
447  }
448  /*
449  fNmaxIter = int(fCmPar[0]);
450  if (fNmaxIter <= 0) {
451  fNmaxIter = fNpar*10 + 20 + fNpar*M*5;
452  }
453  fEPS = fCmPar[1];
454  */
455  //*-*- look for command in list CNAME . . . . . . . . . .
456  TString ctemp = fCword(0,3);
457  Int_t ind;
458  for (ind = 0; ind < nntot; ++ind) {
459  if (strncmp(ctemp.Data(),cname[ind],3) == 0) break;
460  }
461  if (ind==nntot) return -3; // Unknow command - input ignored
462  if (fCword(0,4) == "MINO") ind=3;
463  switch (ind) {
464  case 0: case 3: case 2: case 28:
465  // MINImize [maxcalls] [tolerance]
466  // also SIMplex, MIGrad and FUMili
467  if(nargs>=1)
468  fNmaxIter=TMath::Max(Int_t(fCmPar[0]),fNmaxIter); // fiXME!!
469  if(nargs==2)
470  fEPS=fCmPar[1];
471  return Minimize();
472  case 1:
473  // SEEk not implemented in this package
474  return -10;
475 
476  case 4: // MINos errors analysis not implemented
477  return -10;
478 
479  case 5: case 6: // SET xxx & SHOW xxx
480  return ExecuteSetCommand(nargs);
481 
482  case 7: // Obsolete command
483  Printf("1");
484  return 0;
485  case 8: // fiX <parno> ....
486  if (nargs<1) return -1; // No parameters specified
487  for (i=0;i<nargs;i++) {
488  Int_t parnum = Int_t(fCmPar[i])-1;
489  FixParameter(parnum);
490  }
491  return 0;
492  case 9: // REStore <code>
493  if (nargs<1) return 0;
494  if(fCmPar[0]==0.)
495  for (i=0;i<fNpar;i++)
496  ReleaseParameter(i);
497  else
498  if(fCmPar[0]==1.) {
500  std::cout <<fLastFixed<<std::endl;
501  }
502  return 0;
503  case 10: // RELease <parno> ...
504  if (nargs<1) return -1; // No parameters specified
505  for (i=0;i<nargs;i++) {
506  Int_t parnum = Int_t(fCmPar[i])-1;
507  ReleaseParameter(parnum);
508  }
509  return 0;
510  case 11: // SCAn not implemented
511  return -10;
512  case 12: // CONt not implemented
513  return -10;
514 
515  case 13: // HESSe not implemented
516  return -10;
517  case 14: // SAVe
518  Printf("SAVe command is obsolete");
519  return -10;
520  case 15: // IMProve not implemented
521  return -10;
522  case 16: // CALl fcn <iflag>
523  {if(nargs<1) return -1;
524  Int_t flag = Int_t(fCmPar[0]);
525  Double_t fval;
526  Eval(fNpar,fGr,fval,fA,flag);
527  return 0;}
528  case 17: // STAndard must call function STAND
529  return 0;
530  case 18: case 19:
531  case 20: case 24: {
532  Double_t fval;
533  Int_t flag = 3;
534  Eval(fNpar,fGr,fval,fA,flag);
535  return 0;
536  }
537  case 21:
538  Clear();
539  return 0;
540  case 22: //HELp not implemented
541  case 23: //MNContour not implemented
542  case 25: // JUMp not implemented
543  return -10;
544  case 26: case 27: case 29: case 30: case 31: case 32:
545  return 0; // blank commands
546  case 33: case 34: case 35: case 36: case 37: case 38:
547  case 39:
548  Printf("Obsolete command. Use corresponding SET command instead");
549  return -10;
550  default:
551  break;
552  }
553  return 0;
554 }
555 
556 
557 
558 ////////////////////////////////////////////////////////////////////////////////
559 ///
560 /// Called from TFumili::ExecuteCommand in case
561 /// of "SET xxx" and "SHOW xxx".
562 ///
563 
565  static Int_t nntot = 30;
566  static const char *cname[30] = {
567  "FCN value ", // 0 .
568  "PARameters", // 1 .
569  "LIMits ", // 2 .
570  "COVariance", // 3 .
571  "CORrelatio", // 4 .
572  "PRInt levl", // 5 not implemented yet
573  "NOGradient", // 6 .
574  "GRAdient ", // 7 .
575  "ERRor def ", // 8 not sure how to implement - by time being ignored
576  "INPut file", // 9 not implemented
577  "WIDth page", // 10 not implemented yet
578  "LINes page", // 11 not implemented yet
579  "NOWarnings", // 12 .
580  "WARnings ", // 13 .
581  "RANdom gen", // 14 not implemented
582  "TITle ", // 15 ignored
583  "STRategy ", // 16 ignored
584  "EIGenvalue", // 17 not implemented yet
585  "PAGe throw", // 18 ignored
586  "MINos errs", // 19 not implemented yet
587  "EPSmachine", // 20 .
588  "OUTputfile", // 21 not implemented
589  "BATch ", // 22 ignored
590  "INTeractiv", // 23 ignored
591  "VERsion ", // 24 .
592  "reserve ", // 25 .
593  "NODebug ", // 26 .
594  "DEBug ", // 27 .
595  "SHOw ", // 28 err
596  "SET "};// 29 err
597 
598  TString cfname, cmode, ckind, cwarn, copt, ctemp, ctemp2;
599  Int_t i, ind;
600  Bool_t setCommand=kFALSE;
601  for (ind = 0; ind < nntot; ++ind) {
602  ctemp = cname[ind];
603  ckind = ctemp(0,3);
604  ctemp2 = fCword(4,6);
605  if (strstr(ctemp2.Data(),ckind.Data())) break;
606  }
607  ctemp2 = fCword(0,3);
608  if(ctemp2.Contains("SET")) setCommand=true;
609  if(ctemp2.Contains("HEL") || ctemp2.Contains("SHO")) setCommand=false;
610 
611  if (ind>=nntot) return -3;
612 
613  switch(ind) {
614  case 0: // SET FCN value illegial // SHOw only
615  if(!setCommand) Printf("FCN=%f",fS);
616  return 0;
617  case 1: // PARameter <parno> <value>
618  {
619  if (nargs<2 && setCommand) return -1;
620  Int_t parnum;
621  Double_t val;
622  if(setCommand) {
623  parnum = Int_t(fCmPar[0])-1;
624  val= fCmPar[1];
625  if(parnum<0 || parnum>=fNpar) return -2; //no such parameter
626  fA[parnum] = val;
627  } else {
628  if (nargs>0) {
629  parnum = Int_t(fCmPar[0])-1;
630  if(parnum<0 || parnum>=fNpar) return -2; //no such parameter
631  Printf("Parameter %s = %E",fANames[parnum].Data(),fA[parnum]);
632  } else
633  for (i=0;i<fNpar;i++)
634  Printf("Parameter %s = %E",fANames[i].Data(),fA[i]);
635 
636  }
637  return 0;
638  }
639  case 2: // LIMits [parno] [ <lolim> <uplim> ]
640  {
641  Int_t parnum;
642  Double_t lolim,uplim;
643  if (nargs<1) {
644  for(i=0;i<fNpar;i++)
645  if(setCommand) {
646  fAMN[i] = gMINDOUBLE;
647  fAMX[i] = gMAXDOUBLE;
648  } else
649  Printf("Limits for param %s: Low=%E, High=%E",
650  fANames[i].Data(),fAMN[i],fAMX[i]);
651  } else {
652  parnum = Int_t(fCmPar[0])-1;
653  if(parnum<0 || parnum>=fNpar)return -1;
654  if(setCommand) {
655  if(nargs>2) {
656  lolim = fCmPar[1];
657  uplim = fCmPar[2];
658  if(uplim==lolim) return -1;
659  if(lolim>uplim) {
660  Double_t tmp = lolim;
661  lolim = uplim;
662  uplim = tmp;
663  }
664  } else {
665  lolim = gMINDOUBLE;
666  uplim = gMAXDOUBLE;
667  }
668  fAMN[parnum] = lolim;
669  fAMX[parnum] = uplim;
670  } else
671  Printf("Limits for param %s Low=%E, High=%E",
672  fANames[parnum].Data(),fAMN[parnum],fAMX[parnum]);
673  }
674  return 0;
675  }
676  case 3:
677  {
678  if(setCommand) return 0;
679  Printf("\nCovariant matrix ");
680  Int_t l = 0,nn=0,nnn=0;
681  for (i=0;i<fNpar;i++) if(fPL0[i]>0.) nn++;
682  for (i=0;i<nn;i++) {
683  for(;fPL0[nnn]<=0.;nnn++) { }
684  printf("%5s: ",fANames[nnn++].Data());
685  for (Int_t j=0;j<=i;j++)
686  printf("%11.2E",fZ[l++]);
687  std::cout<<std::endl;
688  }
689  std::cout<<std::endl;
690  return 0;
691  }
692  case 4:
693  if(setCommand) return 0;
694  Printf("\nGlobal correlation factors (maximum correlation of the parameter\n with arbitrary linear combination of other parameters)");
695  for(i=0;i<fNpar;i++) {
696  printf("%5s: ",fANames[i].Data());
697  printf("%11.3E\n",TMath::Sqrt(1-1/((fR[i]!=0.)?fR[i]:1.)) );
698  }
699  std::cout<<std::endl;
700  return 0;
701  case 5: // PRIntout not implemented
702  return -10;
703  case 6: // NOGradient
704  if(!setCommand) return 0;
705  fGRAD = false;
706  return 0;
707  case 7: // GRAdient
708  if(!setCommand) return 0;
709  fGRAD = true;
710  return 0;
711  case 8: // ERRordef - now ignored
712  return 0;
713  case 9: // INPut - not implemented
714  return -10;
715  case 10: // WIDthpage - not implemented
716  return -10;
717  case 11: // LINesperpage - not implemented
718  return -10;
719  case 12: //NOWarnings
720  if(!setCommand) return 0;
721  fWARN = false;
722  return 0;
723  case 13: // WARnings
724  if(!setCommand) return 0;
725  fWARN = true;
726  return 0;
727  case 14: // RANdomgenerator - not implemented
728  return -10;
729  case 15: // TITle - ignored
730  return 0;
731  case 16: // STRategy - ignored
732  return 0;
733  case 17: // EIGenvalues - not implemented
734  return -10;
735  case 18: // PAGethrow - ignored
736  return 0;
737  case 19: // MINos errors - not implemented
738  return -10;
739  case 20: //EPSmachine
740  if(!setCommand) {
741  Printf("Relative floating point presicion RP=%E",fRP);
742  } else
743  if (nargs>0) {
744  Double_t pres=fCmPar[0];
745  if (pres<1e-5 && pres>1e-34) fRP=pres;
746  }
747  return 0;
748  case 21: // OUTputfile - not implemented
749  return -10;
750  case 22: // BATch - ignored
751  return 0;
752  case 23: // INTerative - ignored
753  return 0;
754  case 24: // VERsion
755  if(setCommand) return 0;
756  Printf("FUMILI-ROOT version 0.1");
757  return 0;
758  case 25: // reserved
759  return 0;
760  case 26: // NODebug
761  if(!setCommand) return 0;
762  fDEBUG = false;
763  return 0;
764  case 27: // DEBug
765  if(!setCommand) return 0;
766  fDEBUG = true;
767  return 0;
768  case 28:
769  case 29:
770  return -3;
771  default:
772  break;
773  }
774  return -3;
775 }
776 
777 ////////////////////////////////////////////////////////////////////////////////
778 /// Fixes parameter number ipar
779 
781  if(ipar>=0 && ipar<fNpar && fPL0[ipar]>0.) {
782  fPL0[ipar] = -fPL0[ipar];
783  fLastFixed = ipar;
784  }
785 }
786 
787 ////////////////////////////////////////////////////////////////////////////////
788 /// return a pointer to the covariance matrix
789 
791 {
792  return fZ;
793 
794 }
795 
796 ////////////////////////////////////////////////////////////////////////////////
797 /// return element i,j from the covariance matrix
798 
800 {
801  if (!fZ) return 0;
802  if (i < 0 || i >= fNpar || j < 0 || j >= fNpar) {
803  Error("GetCovarianceMatrixElement","Illegal arguments i=%d, j=%d",i,j);
804  return 0;
805  }
806  return fZ[j+fNpar*i];
807 }
808 
809 
810 ////////////////////////////////////////////////////////////////////////////////
811 /// return the total number of parameters (free + fixed)
812 
814 {
815  return fNpar;
816 }
817 
818 ////////////////////////////////////////////////////////////////////////////////
819 /// return the number of free parameters
820 
822 {
823  Int_t nfree = fNpar;
824  for (Int_t i=0;i<fNpar;i++) {
825  if (IsFixed(i)) nfree--;
826  }
827  return nfree;
828 }
829 
830 ////////////////////////////////////////////////////////////////////////////////
831 /// return error of parameter ipar
832 
834 {
835  if (ipar<0 || ipar>=fNpar) return 0;
836  else return fParamError[ipar];
837 }
838 
839 ////////////////////////////////////////////////////////////////////////////////
840 /// return current value of parameter ipar
841 
843 {
844  if (ipar<0 || ipar>=fNpar) return 0;
845  else return fA[ipar];
846 }
847 
848 
849 ////////////////////////////////////////////////////////////////////////////////
850 /// Get various ipar parameter attributs:
851 ///
852 /// cname: parameter name
853 /// value: parameter value
854 /// verr: parameter error
855 /// vlow: lower limit
856 /// vhigh: upper limit
857 /// WARNING! parname must be suitably dimensionned in the calling function.
858 
859 Int_t TFumili::GetParameter(Int_t ipar,char *cname,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const
860 {
861  if (ipar<0 || ipar>=fNpar) {
862  value = 0;
863  verr = 0;
864  vlow = 0;
865  vhigh = 0;
866  return -1;
867  }
868  strcpy(cname,fANames[ipar].Data());
869  value = fA[ipar];
870  verr = fParamError[ipar];
871  vlow = fAMN[ipar];
872  vhigh = fAMX[ipar];
873  return 0;
874 }
875 
876 ////////////////////////////////////////////////////////////////////////////////
877 /// return name of parameter ipar
878 
879 const char *TFumili::GetParName(Int_t ipar) const
880 {
881  if (ipar < 0 || ipar > fNpar) return "";
882  return fANames[ipar].Data();
883 }
884 
885 ////////////////////////////////////////////////////////////////////////////////
886 /// Return errors after MINOs
887 /// not implemented
888 
889 Int_t TFumili::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
890 {
891  eparab = 0;
892  globcc = 0;
893  if (ipar<0 || ipar>=fNpar) {
894  eplus = 0;
895  eminus = 0;
896  return -1;
897  }
898  eplus=fParamError[ipar];
899  eminus=-eplus;
900  return 0;
901 }
902 
903 ////////////////////////////////////////////////////////////////////////////////
904 /// return global fit parameters
905 /// amin : chisquare
906 /// edm : estimated distance to minimum
907 /// errdef
908 /// nvpar : number of variable parameters
909 /// nparx : total number of parameters
910 
911 Int_t TFumili::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
912 {
913  amin = 2*fS;
914  edm = fGT; //
915  errdef = 0; // ??
916  nparx = fNpar;
917  nvpar = 0;
918  for(Int_t ii=0; ii<fNpar; ii++) {
919  if(fPL0[ii]>0.) nvpar++;
920  }
921  return 0;
922 }
923 
924 
925 
926 ////////////////////////////////////////////////////////////////////////////////
927 /// return Sum(log(i) i=0,n
928 /// used by log likelihood fits
929 
931 {
932  if (n < 0) return 0;
933  if (n > fNlog) {
934  if (fSumLog) delete [] fSumLog;
935  fNlog = 2*n+1000;
936  fSumLog = new Double_t[fNlog+1];
937  Double_t fobs = 0;
938  for (Int_t j=0;j<=fNlog;j++) {
939  if (j > 1) fobs += TMath::Log(j);
940  fSumLog[j] = fobs;
941  }
942  }
943  if (fSumLog) return fSumLog[n];
944  return 0;
945 }
946 
947 
948 
949 ////////////////////////////////////////////////////////////////////////////////
950 /// Inverts packed diagonal matrix Z by square-root method.
951 /// Matrix elements corresponding to
952 /// fix parameters are removed.
953 ///
954 /// n: number of variable parameters
955 ///
956 
958 {
959  static Double_t am = 3.4e138;
960  static Double_t rp = 5.0e-14;
961  Double_t ap, aps, c, d;
962  Double_t *r_1=fR;
963  Double_t *pl_1=fPL;
964  Double_t *z_1=fZ;
965  Int_t i, k, l, ii, ki, li, kk, ni, ll, nk, nl, ir, lk;
966  if (n < 1) {
967  return;
968  }
969  --pl_1;
970  --r_1;
971  --z_1;
972  aps = am / n;
973  aps = sqrt(aps);
974  ap = 1.0e0 / (aps * aps);
975  ir = 0;
976  for (i = 1; i <= n; ++i) {
977  L1:
978  ++ir;
979  if (pl_1[ir] <= 0.0e0) goto L1;
980  else goto L2;
981  L2:
982  ni = i * (i - 1) / 2;
983  ii = ni + i;
984  k = n + 1;
985  if (z_1[ii] <= rp * TMath::Abs(r_1[ir]) || z_1[ii] <= ap) {
986  goto L19;
987  }
988  z_1[ii] = 1.0e0 / sqrt(z_1[ii]);
989  nl = ii - 1;
990  L3:
991  if (nl - ni <= 0) goto L5;
992  else goto L4;
993  L4:
994  z_1[nl] *= z_1[ii];
995  if (TMath::Abs(z_1[nl]) >= aps) {
996  goto L16;
997  }
998  --nl;
999  goto L3;
1000  L5:
1001  if (i - n >= 0) goto L12;
1002  else goto L6;
1003  L6:
1004  --k;
1005  nk = k * (k - 1) / 2;
1006  nl = nk;
1007  kk = nk + i;
1008  d = z_1[kk] * z_1[ii];
1009  c = d * z_1[ii];
1010  l = k;
1011  L7:
1012  ll = nk + l;
1013  li = nl + i;
1014  z_1[ll] -= z_1[li] * c;
1015  --l;
1016  nl -= l;
1017  if (l - i <= 0) goto L9;
1018  else goto L7;
1019  L8:
1020  ll = nk + l;
1021  li = ni + l;
1022  z_1[ll] -= z_1[li] * d;
1023  L9:
1024  --l;
1025  if (l <= 0) goto L10;
1026  else goto L8;
1027  L10:
1028  z_1[kk] = -c;
1029  if (k - i - 1 <= 0) goto L11;
1030  else goto L6;
1031  L11:
1032  ;
1033  }
1034  L12:
1035  for (i = 1; i <= n; ++i) {
1036  for (k = i; k <= n; ++k) {
1037  nl = k * (k - 1) / 2;
1038  ki = nl + i;
1039  d = 0.0e0;
1040  for (l = k; l <= n; ++l) {
1041  li = nl + i;
1042  lk = nl + k;
1043  d += z_1[li] * z_1[lk];
1044  nl += l;
1045  }
1046  ki = k * (k - 1) / 2 + i;
1047  z_1[ki] = d;
1048  }
1049  }
1050  L15:
1051  return;
1052  L16:
1053  k = i + nl - ii;
1054  ir = 0;
1055  for (i = 1; i <= k; ++i) {
1056  L17:
1057  ++ir;
1058  if (pl_1[ir] <= 0.0e0) {
1059  goto L17;
1060  }
1061  }
1062  L19:
1063  pl_1[ir] = -2.0e0;
1064  r_1[ir] = 0.0e0;
1065  fINDFLG[0] = ir - 1;
1066  goto L15;
1067 }
1068 
1069 
1070 
1071 
1072 ////////////////////////////////////////////////////////////////////////////////
1073 ///return kTRUE if parameter ipar is fixed, kFALSE othersise)
1074 
1076 {
1077  if(ipar < 0 || ipar >= fNpar) {
1078  Warning("IsFixed","Illegal parameter number :%d",ipar);
1079  return kFALSE;
1080  }
1081  if (fPL0[ipar] < 0) return kTRUE;
1082  else return kFALSE;
1083 }
1084 
1085 
1086 ////////////////////////////////////////////////////////////////////////////////
1087 
1089 {// Main minimization procedure
1090  //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*//
1091  // FUMILI
1092  // Based on ideas, proposed by I.N. Silin
1093  // [See NIM A440, 2000 (p431)]
1094  // converted from FORTRAN to C by
1095  // Sergey Yaschenko <s.yaschenko@fz-juelich.de>
1096  //
1097  //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*//
1098  //
1099  // This function is called after setting theoretical function
1100  // by means of TFumili::SetUserFunc and initializing parameters.
1101  // Optionally one can set FCN function (see TFumili::SetFCN and TFumili::Eval)
1102  // If FCN is undefined then user has to provide data arrays by calling
1103  // TFumili::SetData procedure.
1104  //
1105  // TFumili::Minimize return following values:
1106  // 0 - fit is converged
1107  // -2 - function is not decreasing (or bad derivatives)
1108  // -3 - error estimations are infinite
1109  // -4 - maximum number of iterations is exceeded
1110  //
1111  Int_t i;
1112  // Flag3 - is fit is chi2 or likelihood? 0 - chi2, 1 - likelihood
1113  fINDFLG[2]=0;
1114  //
1115  // Are the parameters outside of the boundaries ?
1116  //
1117  Int_t parn;
1118 
1119  if(fFCN) {
1120  Eval(parn,fGr,fS,fA,9); fNfcn++;
1121  }
1122  for( i = 0; i < fNpar; i++) {
1123  if(fA[i] > fAMX[i]) fA[i] = fAMX[i];
1124  if(fA[i] < fAMN[i]) fA[i] = fAMN[i];
1125  }
1126 
1127  Int_t nn2, n, fixFLG, ifix1, fi, nn3, nn1, n0;
1128  Double_t t1;
1129  Double_t sp, t, olds=0;
1130  Double_t bi, aiMAX=0, amb;
1131  Double_t afix, sigi, akap;
1132  Double_t alambd, al, bm, abi, abm;
1133  Int_t l1, k, ifix;
1134 
1135  nn2=0;
1136 
1137  // Number of parameters;
1138  n=fNpar;
1139  fixFLG=0;
1140 
1141  // Exit flag
1142  fENDFLG=0;
1143 
1144  // Flag2
1145  fINDFLG[1] = 0;
1146  ifix1=-1;
1147  fi=0;
1148  nn3=0;
1149 
1150  // Initialize param.step limits
1151  for( i=0; i < n; i++) {
1152  fR[i]=0.;
1153  if ( fEPS > 0.) fParamError[i] = 0.;
1154  fPL[i] = fPL0[i];
1155  }
1156 
1157 L3: // Start Iteration
1158 
1159  nn1 = 1;
1160  t1 = 1.;
1161 
1162 L4: // New iteration
1163 
1164  // fS - objective function value - zero first
1165  fS = 0.;
1166  // n0 - number of variable parameters in fit
1167  n0 = 0;
1168  for( i = 0; i < n; i++) {
1169  fGr[i]=0.; // zero gradients
1170  if (fPL0[i] > .0) {
1171  n0=n0+1;
1172  // new iteration - new parallelepiped
1173  if (fPL[i] > .0) fPL0[i]=fPL[i];
1174  }
1175  }
1176  Int_t nn0;
1177  // Calculate number of fZ-matrix elements as nn0=1+2+..+n0
1178  nn0 = n0*(n0+1)/2;
1179  // if (nn0 >= 1) ????
1180  // fZ-matrix is initialized
1181  for( i=0; i < nn0; i++) fZ[i]=0.;
1182 
1183  // Flag1
1184  fINDFLG[0] = 0;
1185  Int_t ijkl=1;
1186 
1187  // Calculate fS - objective function, fGr - gradients, fZ - fZ-matrix
1188  if(fFCN) {
1189  Eval(parn,fGr,fS,fA,2);
1190  fNfcn++;
1191  } else
1192  ijkl = SGZ();
1193  if(!ijkl) return 10;
1194  if (ijkl == -1) fINDFLG[0]=1;
1195 
1196  // sp - scaled on fS machine precision
1197  sp=fRP*TMath::Abs(fS);
1198 
1199  // save fZ-matrix
1200  for( i=0; i < nn0; i++) fZ0[i] = fZ[i];
1201  if (nn3 > 0) {
1202  if (nn1 <= fNstepDec) {
1203  t=2.*(fS-olds-fGT);
1204  if (fINDFLG[0] == 0) {
1205  if (TMath::Abs(fS-olds) <= sp && -fGT <= sp) goto L19;
1206  if( 0.59*t < -fGT) goto L19;
1207  t = -fGT/t;
1208  if (t < 0.25 ) t = 0.25;
1209  }
1210  else t = 0.25;
1211  fGT = fGT*t;
1212  t1 = t1*t;
1213  nn2=0;
1214  for( i = 0; i < n; i++) {
1215  if (fPL[i] > 0.) {
1216  fA[i]=fA[i]-fDA[i];
1217  fPL[i]=fPL[i]*t;
1218  fDA[i]=fDA[i]*t;
1219  fA[i]=fA[i]+fDA[i];
1220  }
1221  }
1222  nn1=nn1+1;
1223  goto L4;
1224  }
1225  }
1226 
1227 L19:
1228 
1229  if(fINDFLG[0] != 0) {
1230  fENDFLG=-4;
1231  printf("trying to execute an illegal junp at L85\n");
1232  //goto L85;
1233  }
1234 
1235 
1236  Int_t k1, k2, i1, j, l;
1237  k1 = 1;
1238  k2 = 1;
1239  i1 = 1;
1240  // In this cycle we removed from fZ contributions from fixed parameters
1241  // We'll get fixed parameters after boudary check
1242  for( i = 0; i < n; i++) {
1243  if (fPL0[i] > .0) {
1244  // if parameter was fixed - release it
1245  if (fPL[i] == 0.) fPL[i]=fPL0[i];
1246  if (fPL[i] > .0) { // ??? it is already non-zero
1247  // if derivative is negative and we above maximum
1248  // or vice versa then fix parameter again and increment k1 by i1
1249  if ((fA[i] >= fAMX[i] && fGr[i] < 0.) ||
1250  (fA[i] <= fAMN[i] && fGr[i] > 0.)) {
1251  fPL[i] = 0.;
1252  k1 = k1 + i1; // i1 stands for fZ-matrix row-number multiplier
1253  /// - skip this row
1254  // in case we are fixing parameter number i
1255  } else {
1256  for( j=0; j <= i; j++) {// cycle on columns of fZ-matrix
1257  if (fPL0[j] > .0) {
1258  // if parameter is not fixed then fZ = fZ0
1259  // Now matrix fZ of other dimension
1260  if (fPL[j] > .0) {
1261  fZ[k2 -1] = fZ0[k1 -1];
1262  k2=k2+1;
1263  }
1264  k1=k1+1;
1265  }
1266  }
1267  }
1268  }
1269  else k1 = k1 + i1; // In case of negative fPL[i] - after mconvd
1270  i1=i1+1; // Next row of fZ0
1271  }
1272  }
1273 
1274  // INVERT fZ-matrix (mconvd() procedure)
1275  i1 = 1;
1276  l = 1;
1277  for( i = 0; i < n; i++) {// extract diagonal elements to fR-vector
1278  if (fPL[i] > .0) {
1279  fR[i] = fZ[l - 1];
1280  i1 = i1+1;
1281  l = l + i1;
1282  }
1283  }
1284 
1285  n0 = i1 - 1;
1286  InvertZ(n0);
1287 
1288  // fZ matrix now is inversed
1289  if (fINDFLG[0] != 0) { // problems
1290  // some PLs now have negative values, try to reduce fZ-matrix again
1291  fINDFLG[0] = 0;
1292  fINDFLG[1] = 1; // errors can be infinite
1293  fixFLG = fixFLG + 1;
1294  fi = 0;
1295  goto L19;
1296  }
1297 
1298  // ... CALCULATE THEORETICAL STEP TO MINIMUM
1299  i1 = 1;
1300  for( i = 0; i < n; i++) {
1301  fDA[i]=0.; // initial step is zero
1302  if (fPL[i] > .0) { // for non-fixed parameters
1303  l1=1;
1304  for( l = 0; l < n; l++) {
1305  if (fPL[l] > .0) {
1306  // Caluclate offset of Z^-1(i1,l1) element in packed matrix
1307  // because we skip fixed param numbers we need also i,l
1308  if (i1 <= l1 ) k=l1*(l1-1)/2+i1;
1309  else k=i1*(i1-1)/2+l1;
1310  // dA_i = \sum (-Z^{-1}_{il}*grad(fS)_l)
1311  fDA[i]=fDA[i]-fGr[l]*fZ[k - 1];
1312  l1=l1+1;
1313  }
1314  }
1315  i1=i1+1;
1316  }
1317  }
1318  // ... CHECK FOR PARAMETERS ON BOUNDARY
1319 
1320  afix=0.;
1321  ifix = -1;
1322  i1 = 1;
1323  l = i1;
1324  for( i = 0; i < n; i++)
1325  if (fPL[i] > .0) {
1326  sigi = TMath::Sqrt(TMath::Abs(fZ[l - 1])); // calculate \sqrt{Z^{-1}_{ii}}
1327  fR[i] = fR[i]*fZ[l - 1]; // Z_ii * Z^-1_ii
1328  if (fEPS > .0) fParamError[i]=sigi;
1329  if ((fA[i] >= fAMX[i] && fDA[i] > 0.) || (fA[i] <= fAMN[i]
1330  && fDA[i] < .0)) {
1331  // if parameter out of bounds and if step is making things worse
1332 
1333  akap = TMath::Abs(fDA[i]/sigi);
1334  // let's found maximum of dA/sigi - the worst of parameter steps
1335  if (akap > afix) {
1336  afix=akap;
1337  ifix=i;
1338  ifix1=i;
1339  }
1340  }
1341  i1=i1+1;
1342  l=l+i1;
1343  }
1344  if (ifix != -1) {
1345  // so the worst parameter is found - fix it and exclude,
1346  // reduce fZ-matrix again
1347  fPL[ifix] = -1.;
1348  fixFLG = fixFLG + 1;
1349  fi = 0;
1350  //.. REPEAT CALCULATION OF THEORETICAL STEP AFTER fiXING EACH PARAMETER
1351  goto L19;
1352  }
1353 
1354  //... CALCULATE STEP CORRECTION FACTOR
1355 
1356  alambd = 1.;
1357  fAKAPPA = 0.;
1358  Int_t imax;
1359  imax = -1;
1360 
1361 
1362  for( i = 0; i < n; i++) {
1363  if (fPL[i] > .0) {
1364  bm = fAMX[i] - fA[i];
1365  abi = fA[i] + fPL[i]; // upper parameter limit
1366  abm = fAMX[i];
1367  if (fDA[i] <= .0) {
1368  bm = fA[i] - fAMN[i];
1369  abi = fA[i] - fPL[i]; // lower parameter limit
1370  abm = fAMN[i];
1371  }
1372  bi = fPL[i];
1373  // if parallelepiped boundary is crossing limits
1374  // then reduce it (deforming)
1375  if ( bi > bm) {
1376  bi = bm;
1377  abi = abm;
1378  }
1379  // if calculated step is out of bounds
1380  if ( TMath::Abs(fDA[i]) > bi) {
1381  // derease step splitter alambdA if needed
1382  al = TMath::Abs(bi/fDA[i]);
1383  if (alambd > al) {
1384  imax=i;
1385  aiMAX=abi;
1386  alambd=al;
1387  }
1388  }
1389  // fAKAPPA - parameter will be <fEPS if fit is converged
1390  akap = TMath::Abs(fDA[i]/fParamError[i]);
1391  if (akap > fAKAPPA) fAKAPPA=akap;
1392  }
1393  }
1394  //... CALCULATE NEW CORRECTED STEP
1395  fGT = 0.;
1396  amb = 1.e18;
1397  // alambd - multiplier to split teoretical step dA
1398  if (alambd > .0) amb = 0.25/alambd;
1399  for( i = 0; i < n; i++) {
1400  if (fPL[i] > .0) {
1401  if (nn2 > fNlimMul ) {
1402  if (TMath::Abs(fDA[i]/fPL[i]) > amb ) {
1403  fPL[i] = 4.*fPL[i]; // increase parallelepiped
1404  t1=4.; // flag - that fPL was increased
1405  }
1406  }
1407  // cut step
1408  fDA[i] = fDA[i]*alambd;
1409  // expected function value change in next iteration
1410  fGT = fGT + fDA[i]*fGr[i];
1411  }
1412  }
1413 
1414  //.. CHECK IF MINIMUM ATTAINED AND SET EXIT MODE
1415  // if expected fGT smaller than precision
1416  // and other stuff
1417  if (-fGT <= sp && t1 < 1. && alambd < 1.)fENDFLG = -1; // function is not decreasing
1418 
1419  if (fENDFLG >= 0) {
1420  if (fAKAPPA < TMath::Abs(fEPS)) { // fit is converging
1421  if (fixFLG == 0)
1422  fENDFLG=1; // successful fit
1423  else {// we have fixed parameters
1424  if (fENDFLG == 0) {
1425  //... CHECK IF fiXING ON BOUND IS CORRECT
1426  fENDFLG = 1;
1427  fixFLG = 0;
1428  ifix1=-1;
1429  // release fixed parameters
1430  for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
1431  fINDFLG[1] = 0;
1432  // and repeat iteration
1433  goto L19;
1434  } else {
1435  if( ifix1 >= 0) {
1436  fi = fi + 1;
1437  fENDFLG = 0;
1438  }
1439  }
1440  }
1441  } else { // fit does not converge
1442  if( fixFLG != 0) {
1443  if( fi > fixFLG ) {
1444  //... CHECK IF fiXING ON BOUND IS CORRECT
1445  fENDFLG = 1;
1446  fixFLG = 0;
1447  ifix1=-1;
1448  for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
1449  fINDFLG[1] = 0;
1450  goto L19;
1451  } else {
1452  fi = fi + 1;
1453  fENDFLG = 0;
1454  }
1455  } else {
1456  fi = fi + 1;
1457  fENDFLG = 0;
1458  }
1459  }
1460  }
1461 
1462 // L85:
1463  // iteration number limit is exceeded
1464  if(fENDFLG == 0 && nn3 >= fNmaxIter) fENDFLG=-3;
1465 
1466  // fit errors are infinite;
1467  if(fENDFLG > 0 && fINDFLG[1] > 0) fENDFLG=-2;
1468 
1469  //MONITO (fS,fNpar,nn3,IT,fEPS,fGT,fAKAPPA,alambd);
1470  if (fENDFLG == 0) { // make step
1471  for ( i = 0; i < n; i++) fA[i] = fA[i] + fDA[i];
1472  if (imax >= 0) fA[imax] = aiMAX;
1473  olds=fS;
1474  nn2=nn2+1;
1475  nn3=nn3+1;
1476  } else {
1477  // fill covariant matrix VL
1478  // fill parameter error matrix up
1479  Int_t il;
1480  il = 0;
1481  for( Int_t ip = 0; ip < fNpar; ip++) {
1482  if( fPL0[ip] > .0) {
1483  for( Int_t jp = 0; jp <= ip; jp++) {
1484  if(fPL0[jp] > .0) {
1485  // VL[ind(ip,jp)] = fZ[il];
1486  il = il + 1;
1487  }
1488  }
1489  }
1490  }
1491  return fENDFLG - 1;
1492  }
1493  goto L3;
1494 }
1495 
1496 ////////////////////////////////////////////////////////////////////////////////
1497 /// Prints fit results.
1498 ///
1499 /// ikode is the type of printing parameters
1500 /// p is function value
1501 ///
1502 /// ikode = 1 - print values, errors and limits
1503 /// ikode = 2 - print values, errors and steps
1504 /// ikode = 3 - print values, errors, steps and derivatives
1505 /// ikode = 4 - print only values and errors
1506 ///
1507 
1509 {
1510  TString exitStatus="";
1511  TString xsexpl="";
1512  TString colhdu[3],colhdl[3],cx2,cx3;
1513  switch (fENDFLG) {
1514  case 1:
1515  exitStatus="CONVERGED";
1516  break;
1517  case -1:
1518  exitStatus="CONST FCN";
1519  xsexpl="****\n* FUNCTION IS NOT DECREASING OR BAD DERIVATIVES\n****";
1520  break;
1521  case -2:
1522  exitStatus="ERRORS INF";
1523  xsexpl="****\n* ESTIMATED ERRORS ARE INfiNITE\n****";
1524  break;
1525  case -3:
1526  exitStatus="MAX ITER.";
1527  xsexpl="****\n* MAXIMUM NUMBER OF ITERATIONS IS EXCEEDED\n****";
1528  break;
1529  case -4:
1530  exitStatus="ZERO PROBAB";
1531  xsexpl="****\n* PROBABILITY OF LIKLIHOOD FUNCTION IS NEGATIVE OR ZERO\n****";
1532  break;
1533  default:
1534  exitStatus="UNDEfiNED";
1535  xsexpl="****\n* fiT IS IN PROGRESS\n****";
1536  break;
1537  }
1538  if (ikode == 1) {
1539  colhdu[0] = " ";
1540  colhdl[0] = " ERROR ";
1541  colhdu[1] = " PHYSICAL";
1542  colhdu[2] = " LIMITS ";
1543  colhdl[1] = " NEGATIVE ";
1544  colhdl[2] = " POSITIVE ";
1545  }
1546  if (ikode == 2) {
1547  colhdu[0] = " ";
1548  colhdl[0] = " ERROR ";
1549  colhdu[1] = " INTERNAL ";
1550  colhdl[1] = " STEP SIZE ";
1551  colhdu[2] = " INTERNAL ";
1552  colhdl[2] = " VALUE ";
1553  }
1554  if (ikode == 3) {
1555  colhdu[0] = " ";
1556  colhdl[0] = " ERROR ";
1557  colhdu[1] = " STEP ";
1558  colhdl[1] = " SIZE ";
1559  colhdu[2] = " fiRST ";
1560  colhdl[2] = " DERIVATIVE";
1561  }
1562  if (ikode == 4) {
1563  colhdu[0] = " PARABOLIC ";
1564  colhdl[0] = " ERROR ";
1565  colhdu[1] = " MINOS ";
1566  colhdu[2] = "ERRORS ";
1567  colhdl[1] = " NEGATIVE ";
1568  colhdl[2] = " POSITIVE ";
1569  }
1570  if(fENDFLG<1)Printf("%s", xsexpl.Data());
1571  Printf(" FCN=%g FROM FUMILI STATUS=%-10s %9d CALLS OF FCN",
1572  p,exitStatus.Data(),fNfcn);
1573  Printf(" EDM=%g ",-fGT);
1574  Printf(" EXT PARAMETER %-14s%-14s%-14s",
1575  (const char*)colhdu[0].Data()
1576  ,(const char*)colhdu[1].Data()
1577  ,(const char*)colhdu[2].Data());
1578  Printf(" NO. NAME VALUE %-14s%-14s%-14s",
1579  (const char*)colhdl[0].Data()
1580  ,(const char*)colhdl[1].Data()
1581  ,(const char*)colhdl[2].Data());
1582 
1583  for (Int_t i=0;i<fNpar;i++) {
1584  if (ikode==3) {
1585  cx2 = Form("%14.5e",fDA[i]);
1586  cx3 = Form("%14.5e",fGr[i]);
1587 
1588  }
1589  if (ikode==1) {
1590  cx2 = Form("%14.5e",fAMN[i]);
1591  cx3 = Form("%14.5e",fAMX[i]);
1592  }
1593  if (ikode==2) {
1594  cx2 = Form("%14.5e",fDA[i]);
1595  cx3 = Form("%14.5e",fA[i]);
1596  }
1597  if(ikode==4){
1598  cx2 = " *undefined* ";
1599  cx3 = " *undefined* ";
1600  }
1601  if(fPL0[i]<=0.) { cx2=" *fixed* ";cx3=""; }
1602  Printf("%4d %-11s%14.5e%14.5e%-14s%-14s",i+1
1603  ,fANames[i].Data(),fA[i],fParamError[i]
1604  ,cx2.Data(),cx3.Data());
1605  }
1606 }
1607 
1608 
1609 ////////////////////////////////////////////////////////////////////////////////
1610 /// Releases parameter number ipar
1611 
1613  if(ipar>=0 && ipar<fNpar && fPL0[ipar]<=0.) {
1614  fPL0[ipar] = -fPL0[ipar];
1615  if (fPL0[ipar] == 0. || fPL0[ipar]>=1.) fPL0[ipar]=.1;
1616  }
1617 }
1618 
1619 
1620 ////////////////////////////////////////////////////////////////////////////////
1621 /// Sets pointer to data array provided by user.
1622 /// Necessary if SetFCN is not called.
1623 ///
1624 /// numpoints: number of experimental points
1625 /// vecsize: size of data point vector + 2
1626 /// (for N-dimensional fit vecsize=N+2)
1627 /// exdata: data array with following format
1628 ///
1629 /// exdata[0] = ExpValue_0 - experimental data value number 0
1630 /// exdata[1] = ExpSigma_0 - error of value number 0
1631 /// exdata[2] = X_0[0]
1632 /// exdata[3] = X_0[1]
1633 /// .........
1634 /// exdata[vecsize-1] = X_0[vecsize-3]
1635 /// exdata[vecsize] = ExpValue_1
1636 /// exdata[vecsize+1] = ExpSigma_1
1637 /// exdata[vecsize+2] = X_1[0]
1638 /// .........
1639 /// exdata[vecsize*(numpoints-1)] = ExpValue_(numpoints-1)
1640 /// .........
1641 /// exdata[vecsize*numpoints-1] = X_(numpoints-1)[vecsize-3]
1642 ///
1643 
1644 void TFumili::SetData(Double_t *exdata,Int_t numpoints,Int_t vecsize){
1645  if(exdata){
1646  fNED1 = numpoints;
1647  fNED2 = vecsize;
1648  fEXDA = exdata;
1649  }
1650 }
1651 
1652 
1653 ////////////////////////////////////////////////////////////////////////////////
1654 /// ret fit method (chisquare or loglikelihood)
1655 
1656 void TFumili::SetFitMethod(const char *name)
1657 {
1658  if (!strcmp(name,"H1FitChisquare")) SetFCN(H1FitChisquareFumili);
1659  if (!strcmp(name,"H1FitLikelihood")) SetFCN(H1FitLikelihoodFumili);
1660  if (!strcmp(name,"GraphFitChisquare")) SetFCN(GraphFitChisquareFumili);
1661 }
1662 
1663 
1664 ////////////////////////////////////////////////////////////////////////////////
1665 /// Sets for prameter number ipar initial parameter value,
1666 /// name parname, initial error verr and limits vlow and vhigh
1667 /// If vlow = vhigh but not equil to zero, parameter will be fixed.
1668 /// If vlow = vhigh = 0, parameter is released and its limits are discarded
1669 ///
1670 
1671 Int_t TFumili::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh) {
1672  if (ipar<0 || ipar>=fNpar) return -1;
1673  fANames[ipar] = parname;
1674  fA[ipar] = value;
1675  fParamError[ipar] = verr;
1676  if(vlow<vhigh) {
1677  fAMN[ipar] = vlow;
1678  fAMX[ipar] = vhigh;
1679  } else {
1680  if(vhigh<vlow) {
1681  fAMN[ipar] = vhigh;
1682  fAMX[ipar] = vlow;
1683  }
1684  if(vhigh==vlow) {
1685  if(vhigh==0.) {
1686  ReleaseParameter(ipar);
1687  fAMN[ipar] = gMINDOUBLE;
1688  fAMX[ipar] = gMAXDOUBLE;
1689  }
1690  if(vlow!=0) FixParameter(ipar);
1691  }
1692  }
1693  return 0;
1694 }
1695 
1696 ////////////////////////////////////////////////////////////////////////////////
1697 /// Evaluates objective function ( chi-square ), gradients and
1698 /// Z-matrix using data provided by user via TFumili::SetData
1699 ///
1700 
1702 {
1703  fS = 0.;
1704  Int_t i,j,l,k2=1,k1,ki=0;
1705  Double_t *x = new Double_t[fNED2];
1706  Double_t *df = new Double_t[fNpar];
1707  Int_t nx = fNED2-2;
1708  for (l=0;l<fNED1;l++) { // cycle on all exp. points
1709  k1 = k2;
1710  if (fLogLike) {
1712  nx = fNED2;
1713  k1 -= 2;
1714  }
1715 
1716  for (i=0;i<nx;i++){
1717  ki += 1+i;
1718  x[i] = fEXDA[ki];
1719  }
1720  // Double_t y = ARITHM(df,x);
1721  Double_t y = EvalTFN(df,x);
1722  if(fNumericDerivatives) Derivatives(df,x);
1723  Double_t sig=1.;
1724  if(fLogLike) { // Likelihood method
1725  if(y>0.) {
1726  fS = fS - log(y);
1727  y = -y;
1728  sig= y;
1729  } else { //
1730  delete [] x;
1731  delete [] df;
1732  fS = 1e10;
1733  return -1; // indflg[0] = 1;
1734  }
1735  } else { // Chi2 method
1736  sig = fEXDA[k2]; // sigma of experimental point
1737  y = y - fEXDA[k1-1]; // f(x_i) - F_i
1738  fS = fS + (y*y/(sig*sig))*.5; // simple chi2/2
1739  }
1740  Int_t n = 0;
1741  for (i=0;i<fNpar;i++) {
1742  if (fPL0[i]>0){
1743  df[n] = df[i]/sig; // left only non-fixed param derivatives div by Sig
1744  fGr[i] += df[n]*(y/sig);
1745  n++;
1746  }
1747  }
1748  l = 0;
1749  for (i=0;i<n;i++)
1750  for (j=0;j<=i;j++)
1751  fZ[l++] += df[i]*df[j];
1752  k2 += fNED2;
1753  }
1754 
1755  delete[] df;
1756  delete[] x;
1757  return 1;
1758 }
1759 
1760 
1761 ////////////////////////////////////////////////////////////////////////////////
1762 /// Minimization function for H1s using a Chisquare method
1763 /// Default method (function evaluated at center of bin)
1764 /// for each point the cache contains the following info
1765 /// -1D : bc,e,xc (bin content, error, x of center of bin)
1766 /// -2D : bc,e,xc,yc
1767 /// -3D : bc,e,xc,yc,zc
1768 
1770 {
1771  Foption_t fitOption = GetFitOption();
1772  if (fitOption.Integral) {
1773  FitChisquareI(npar,gin,f,u,flag);
1774  return;
1775  }
1776  Double_t cu,eu,fu,fsum;
1777  Double_t x[3];
1778  Double_t *zik=0;
1779  Double_t *pl0=0;
1780 
1781  TH1 *hfit = (TH1*)GetObjectFit();
1782  TF1 *f1 = (TF1*)GetUserFunc();
1783  Int_t nd = hfit->GetDimension();
1784  Int_t j;
1785 
1786  npar = f1->GetNpar();
1787  SetParNumber(npar);
1788  if(flag == 9) return;
1789  zik = GetZ();
1790  pl0 = GetPL0();
1791 
1792  Double_t *df = new Double_t[npar];
1793  f1->InitArgs(x,u);
1794  f = 0;
1795 
1796  Int_t npfit = 0;
1797  Double_t *cache = fCache;
1798  for (Int_t i=0;i<fNpoints;i++) {
1799  if (nd > 2) x[2] = cache[4];
1800  if (nd > 1) x[1] = cache[3];
1801  x[0] = cache[2];
1802  cu = cache[0];
1804  fu = f1->EvalPar(x,u);
1805  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
1806  eu = cache[1];
1807  Derivatives(df,x);
1808  Int_t n = 0;
1809  fsum = (fu-cu)/eu;
1810  if (flag!=1) {
1811  for (j=0;j<npar;j++) {
1812  if (pl0[j]>0) {
1813  df[n] = df[j]/eu;
1814  // left only non-fixed param derivatives / by Sigma
1815  gin[j] += df[n]*fsum;
1816  n++;
1817  }
1818  }
1819  Int_t l = 0;
1820  for (j=0;j<n;j++)
1821  for (Int_t k=0;k<=j;k++)
1822  zik[l++] += df[j]*df[k];
1823  }
1824  f += .5*fsum*fsum;
1825  npfit++;
1826  cache += fPointSize;
1827  }
1828  f1->SetNumberFitPoints(npfit);
1829  delete [] df;
1830 }
1831 
1832 
1833 ////////////////////////////////////////////////////////////////////////////////
1834 /// Minimization function for H1s using a Chisquare method
1835 /// The "I"ntegral method is used
1836 /// for each point the cache contains the following info
1837 /// -1D : bc,e,xc,xw (bin content, error, x of center of bin, x bin width of bin)
1838 /// -2D : bc,e,xc,xw,yc,yw
1839 /// -3D : bc,e,xc,xw,yc,yw,zc,zw
1840 
1842 {
1843  Double_t cu,eu,fu,fsum;
1844  Double_t x[3];
1845  Double_t *zik=0;
1846  Double_t *pl0=0;
1847 
1848  TH1 *hfit = (TH1*)GetObjectFit();
1849  TF1 *f1 = (TF1*)GetUserFunc();
1850  Int_t nd = hfit->GetDimension();
1851  Int_t j;
1852 
1853  f1->InitArgs(x,u);
1854  npar = f1->GetNpar();
1855  SetParNumber(npar);
1856  if(flag == 9) return;
1857  zik = GetZ();
1858  pl0 = GetPL0();
1859 
1860  Double_t *df=new Double_t[npar];
1861  f = 0;
1862 
1863  Int_t npfit = 0;
1864  Double_t *cache = fCache;
1865  for (Int_t i=0;i<fNpoints;i++) {
1866  cu = cache[0];
1868  f1->SetParameters(u);
1869  if (nd < 2) {
1870  fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3])/cache[3];
1871  } else if (nd < 3) {
1872  fu = ((TF2*)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]);
1873  } else {
1874  fu = ((TF3*)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]);
1875  }
1876  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
1877  eu = cache[1];
1878  Derivatives(df,x);
1879  Int_t n = 0;
1880  fsum = (fu-cu)/eu;
1881  if (flag!=1) {
1882  for (j=0;j<npar;j++) {
1883  if (pl0[j]>0){
1884  df[n] = df[j]/eu;
1885  // left only non-fixed param derivatives / by Sigma
1886  gin[j] += df[n]*fsum;
1887  n++;
1888  }
1889  }
1890  Int_t l = 0;
1891  for (j=0;j<n;j++)
1892  for (Int_t k=0;k<=j;k++)
1893  zik[l++] += df[j]*df[k];
1894  }
1895  f += .5*fsum*fsum;
1896  npfit++;
1897  cache += fPointSize;
1898  }
1899  f1->SetNumberFitPoints(npfit);
1900  delete[] df;
1901 }
1902 
1903 
1904 ////////////////////////////////////////////////////////////////////////////////
1905 /// Minimization function for H1s using a Likelihood method*-*-*-*-*-*
1906 /// Basically, it forms the likelihood by determining the Poisson
1907 /// probability that given a number of entries in a particular bin,
1908 /// the fit would predict it's value. This is then done for each bin,
1909 /// and the sum of the logs is taken as the likelihood.
1910 /// Default method (function evaluated at center of bin)
1911 /// for each point the cache contains the following info
1912 /// -1D : bc,e,xc (bin content, error, x of center of bin)
1913 /// -2D : bc,e,xc,yc
1914 /// -3D : bc,e,xc,yc,zc
1915 
1917 {
1918  Foption_t fitOption = GetFitOption();
1919  if (fitOption.Integral) {
1920  FitLikelihoodI(npar,gin,f,u,flag);
1921  return;
1922  }
1923  Double_t cu,fu,fobs,fsub;
1924  Double_t dersum[100];
1925  Double_t x[3];
1926  Int_t icu;
1927 
1928  TH1 *hfit = (TH1*)GetObjectFit();
1929  TF1 *f1 = (TF1*)GetUserFunc();
1930  Int_t nd = hfit->GetDimension();
1931  Int_t j;
1932  Double_t *zik = GetZ();
1933  Double_t *pl0 = GetPL0();
1934 
1935  npar = f1->GetNpar();
1936  SetParNumber(npar);
1937  if(flag == 9) return;
1938  Double_t *df=new Double_t[npar];
1939  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
1940  f1->InitArgs(x,u);
1941  f = 0;
1942 
1943  Int_t npfit = 0;
1944  Double_t *cache = fCache;
1945  for (Int_t i=0;i<fNpoints;i++) {
1946  if (nd > 2) x[2] = cache[4];
1947  if (nd > 1) x[1] = cache[3];
1948  x[0] = cache[2];
1949  cu = cache[0];
1951  fu = f1->EvalPar(x,u);
1952  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
1953  if (flag == 2) {
1954  for (j=0;j<npar;j++) {
1955  dersum[j] += 1; //should be the derivative
1956  //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
1957  }
1958  }
1959  if (fu < 1.e-9) fu = 1.e-9;
1960  icu = Int_t(cu);
1961  fsub = -fu +icu*TMath::Log(fu);
1962  fobs = GetSumLog(icu);
1963  fsub -= fobs;
1964  Derivatives(df,x);
1965  int n=0;
1966  // Here we need gradients of Log likelihood function
1967  //
1968  for (j=0;j<npar;j++) {
1969  if (pl0[j]>0){
1970  df[n] = df[j]*(icu/fu-1);
1971  gin[j] -= df[n];
1972  n++;
1973  }
1974  }
1975  Int_t l = 0;
1976  // Z-matrix here - production of first derivatives
1977  // of log-likelihood function
1978  for (j=0;j<n;j++)
1979  for (Int_t k=0;k<=j;k++)
1980  zik[l++] += df[j]*df[k];
1981 
1982  f -= fsub;
1983  npfit++;
1984  cache += fPointSize;
1985  }
1986  f *= 2;
1987  f1->SetNumberFitPoints(npfit);
1988  delete[] df;
1989 }
1990 
1991 
1992 ////////////////////////////////////////////////////////////////////////////////
1993 /// Minimization function for H1s using a Likelihood method*-*-*-*-*-*
1994 /// Basically, it forms the likelihood by determining the Poisson
1995 /// probability that given a number of entries in a particular bin,
1996 /// the fit would predict it's value. This is then done for each bin,
1997 /// and the sum of the logs is taken as the likelihood.
1998 /// The "I"ntegral method is used
1999 /// for each point the cache contains the following info
2000 /// -1D : bc,e,xc,xw (bin content, error, x of center of bin, x bin width of bin)
2001 /// -2D : bc,e,xc,xw,yc,yw
2002 /// -3D : bc,e,xc,xw,yc,yw,zc,zw
2003 
2005 {
2006  Double_t cu,fu,fobs,fsub;
2007  Double_t dersum[100];
2008  Double_t x[3];
2009  Int_t icu;
2010 
2011  TH1 *hfit = (TH1*)GetObjectFit();
2012  TF1 *f1 = (TF1*)GetUserFunc();
2013  Int_t nd = hfit->GetDimension();
2014  Int_t j;
2015  Double_t *zik = GetZ();
2016  Double_t *pl0 = GetPL0();
2017 
2018  Double_t *df=new Double_t[npar];
2019 
2020  npar = f1->GetNpar();
2021  SetParNumber(npar);
2022  if(flag == 9) {delete [] df; return;}
2023  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
2024  f1->InitArgs(x,u);
2025  f = 0;
2026 
2027  Int_t npfit = 0;
2028  Double_t *cache = fCache;
2029  for (Int_t i=0;i<fNpoints;i++) {
2030  if (nd > 2) x[2] = cache[4];
2031  if (nd > 1) x[1] = cache[3];
2032  x[0] = cache[2];
2033  cu = cache[0];
2035  if (nd < 2) {
2036  fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3])/cache[3];
2037  } else if (nd < 3) {
2038  fu = ((TF2*)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]);
2039  } else {
2040  fu = ((TF3*)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]);
2041  }
2042  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
2043  if (flag == 2) {
2044  for (j=0;j<npar;j++) {
2045  dersum[j] += 1; //should be the derivative
2046  //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
2047  }
2048  }
2049  if (fu < 1.e-9) fu = 1.e-9;
2050  icu = Int_t(cu);
2051  fsub = -fu +icu*TMath::Log(fu);
2052  fobs = GetSumLog(icu);
2053  fsub -= fobs;
2054  Derivatives(df,x);
2055  int n=0;
2056  // Here we need gradients of Log likelihood function
2057  //
2058  for (j=0;j<npar;j++) {
2059  if (pl0[j]>0){
2060  df[n] = df[j]*(icu/fu-1);
2061  gin[j] -= df[n];
2062  n++;
2063  }
2064  }
2065  Int_t l = 0;
2066  // Z-matrix here - production of first derivatives
2067  // of log-likelihood function
2068  for (j=0;j<n;j++)
2069  for (Int_t k=0;k<=j;k++)
2070  zik[l++] += df[j]*df[k];
2071 
2072  f -= fsub;
2073  npfit++;
2074  cache += fPointSize;
2075  }
2076  f *= 2;
2077  f1->SetNumberFitPoints(npfit);
2078  delete[] df;
2079 }
2080 
2081 
2082 //______________________________________________________________________________
2083 //
2084 // STATIC functions
2085 //______________________________________________________________________________
2086 
2087 ////////////////////////////////////////////////////////////////////////////////
2088 /// Minimization function for H1s using a Chisquare method
2089 /// ======================================================
2090 
2092 {
2093  TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
2094  hFitter->FitChisquare(npar, gin, f, u, flag);
2095 }
2096 
2097 ////////////////////////////////////////////////////////////////////////////////
2098 /// -*-*-*-*Minimization function for H1s using a Likelihood method*-*-*-*-*-*
2099 /// =======================================================
2100 /// Basically, it forms the likelihood by determining the Poisson
2101 /// probability that given a number of entries in a particular bin,
2102 /// the fit would predict it's value. This is then done for each bin,
2103 /// and the sum of the logs is taken as the likelihood.
2104 /// PDF: P=exp(-f(x_i))/[F_i]!*(f(x_i))^[F_i]
2105 /// where F_i - experimental value, f(x_i) - expected theoretical value
2106 /// [F_i] - integer part of F_i.
2107 /// drawback is that if F_i>Int_t - GetSumLog will fail
2108 /// for big F_i is faster to use Euler's Gamma-function
2109 
2111 {
2112 
2113  TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
2114  hFitter->FitLikelihood(npar, gin, f, u, flag);
2115 }
2116 
2117 
2118 
2119 ////////////////////////////////////////////////////////////////////////////////
2120 ///*-*-*-*-*-*Minimization function for Graphs using a Chisquare method*-*-*-*-*
2121 ///*-* =========================================================
2122 ///
2123 /// In case of a TGraphErrors object, ex, the error along x, is projected
2124 /// along the y-direction by calculating the function at the points x-exlow and
2125 /// x+exhigh.
2126 ///
2127 /// The chisquare is computed as the sum of the quantity below at each point:
2128 ///
2129 /// (y - f(x))**2
2130 /// -----------------------------------
2131 /// ey**2 + (0.5*(exl + exh)*f'(x))**2
2132 ///
2133 /// where x and y are the point coordinates and f'(x) is the derivative of function f(x).
2134 /// This method to approximate the uncertainty in y because of the errors in x, is called
2135 /// "effective variance" method.
2136 /// The improvement, compared to the previously used method (f(x+ exhigh) - f(x-exlow))/2
2137 /// is of (error of x)**2 order.
2138 /// NOTE:
2139 /// 1) By using the "effective variance" method a simple linear regression
2140 /// becomes a non-linear case , which takes several iterations
2141 /// instead of 0 as in the linear case .
2142 ///
2143 /// 2) The effective variance technique assumes that there is no correlation
2144 /// between the x and y coordinate .
2145 ///
2146 /// In case the function lies below (above) the data point, ey is ey_low (ey_high).
2147 
2149  Double_t *u, Int_t flag)
2150 {
2151  Double_t cu,eu,exl,exh,ey,eux,fu,fsum;
2152  Double_t x[1];
2153  Int_t i, bin, npfits=0;
2154 
2155  TFumili *grFitter = (TFumili*)TVirtualFitter::GetFitter();
2156  TGraph *gr = (TGraph*)grFitter->GetObjectFit();
2157  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
2158  Foption_t fitOption = grFitter->GetFitOption();
2159 
2160  Int_t n = gr->GetN();
2161  Double_t *gx = gr->GetX();
2162  Double_t *gy = gr->GetY();
2163  npar = f1->GetNpar();
2164 
2165  grFitter->SetParNumber(npar);
2166 
2167  if(flag == 9) return;
2168  Double_t *zik = grFitter->GetZ();
2169  Double_t *pl0 = grFitter->GetPL0();
2170  Double_t *df = new Double_t[npar];
2171 
2172 
2173  f1->InitArgs(x,u);
2174  f = 0;
2175  for (bin=0;bin<n;bin++) {
2176  x[0] = gx[bin];
2177  if (!f1->IsInside(x)) continue;
2178  cu = gy[bin];
2180  fu = f1->EvalPar(x,u);
2181  if (TF1::RejectedPoint()) continue;
2182  npfits++;
2183  Double_t eusq=1.;
2184  if (fitOption.W1) {
2185  eu = 1.;
2186  } else {
2187  exh = gr->GetErrorXhigh(bin);
2188  exl = gr->GetErrorXlow(bin);
2189  ey = gr->GetErrorY(bin);
2190  if (exl < 0) exl = 0;
2191  if (exh < 0) exh = 0;
2192  if (ey < 0) ey = 0;
2193  if (exh > 0 && exl > 0) {
2194 // "Effective variance" method for projecting errors
2195  eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
2196  } else
2197  eux = 0.;
2198  eu = ey*ey+eux*eux;
2199  if (eu <= 0) eu = 1;
2200  eusq = TMath::Sqrt(eu);
2201  }
2202  grFitter->Derivatives(df,x);
2203  n = 0;
2204  fsum = (fu-cu)/eusq;
2205  for (i=0;i<npar;i++) {
2206  if (pl0[i]>0){
2207  df[n] = df[i]/eusq;
2208  // left only non-fixed param derivatives / by Sigma
2209  gin[i] += df[n]*fsum;
2210  n++;
2211  }
2212  }
2213  Int_t l = 0;
2214  for (i=0;i<n;i++)
2215  for (Int_t j=0;j<=i;j++)
2216  zik[l++] += df[i]*df[j];
2217  f += .5*fsum*fsum;
2218 
2219  }
2220  delete [] df;
2221  f1->SetNumberFitPoints(npfits);
2222 }
2223 
const int nx
Definition: kalman.C:16
virtual Int_t GetNumberTotalParameters() const
return the total number of parameters (free + fixed)
Definition: TFumili.cxx:813
double par[1]
Definition: unuranDistr.cxx:38
Double_t * fCache
Double_t * fDF
Definition: TFumili.h:64
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
void grad()
Definition: grad.C:17
static float fu(float x)
Definition: main.cpp:53
Double_t Log(Double_t x)
Definition: TMath.h:526
Double_t * fR
Definition: TFumili.h:62
Int_t Eval(Int_t &npar, Double_t *grad, Double_t &fval, Double_t *par, Int_t flag)
Evaluate the minimisation function Input parameters: npar: number of currently variable parameters pa...
Definition: TFumili.cxx:351
const double pi
return c
const char Option_t
Definition: RtypesCore.h:62
virtual Foption_t GetFitOption() const
Double_t fRP
Definition: TFumili.h:69
virtual void FixParameter(Int_t ipar)
Fixes parameter number ipar.
Definition: TFumili.cxx:780
virtual Int_t GetDimension() const
Definition: TH1.h:283
virtual Double_t Chisquare(Int_t npar, Double_t *params) const
return a chisquare equivalent
Definition: TFumili.cxx:221
Int_t fNlog
Definition: TFumili.h:24
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
TFumili(Int_t maxpar=25)
Definition: TFumili.cxx:129
virtual Double_t GetCovarianceMatrixElement(Int_t i, Int_t j) const
return element i,j from the covariance matrix
Definition: TFumili.cxx:799
virtual Int_t GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
return global fit parameters amin : chisquare edm : estimated distance to minimum errdef nvpar : numb...
Definition: TFumili.cxx:911
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1088
#define gROOT
Definition: TROOT.h:344
Basic string class.
Definition: TString.h:137
virtual TObject * GetUserFunc() const
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:421
Double_t * fDA
Definition: TFumili.h:59
Double_t * fParamError
Definition: TFumili.h:49
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t * fPL
Definition: TFumili.h:56
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SetFCN(void *fcn)
To set the address of the minimization objective function.
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
Int_t GetN() const
Definition: TGraph.h:132
TLatex * t1
Definition: textangle.C:20
Double_t * fEXDA
Definition: TFumili.h:51
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Int_t Minimize()
Definition: TFumili.cxx:1088
Bool_t fNumericDerivatives
Definition: TFumili.h:42
TFile * f
Int_t fNED1
Definition: TFumili.h:26
virtual Double_t GetErrorXhigh(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1390
virtual Double_t GetParameter(Int_t ipar) const
return current value of parameter ipar
Definition: TFumili.cxx:842
const char * Data() const
Definition: TString.h:349
Double_t * GetY() const
Definition: TGraph.h:140
virtual Int_t SetParameter(Int_t ipar, const char *parname, Double_t value, Double_t verr, Double_t vlow, Double_t vhigh)
Sets for prameter number ipar initial parameter value, name parname, initial error verr and limits vl...
Definition: TFumili.cxx:1671
double sqrt(double)
void InvertZ(Int_t)
Inverts packed diagonal matrix Z by square-root method.
Definition: TFumili.cxx:957
Double_t x[n]
Definition: legend1.C:17
virtual void PrintResults(Int_t k, Double_t p) const
Prints fit results.
Definition: TFumili.cxx:1508
void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method
Definition: TFumili.cxx:2091
virtual Double_t GetParError(Int_t ipar) const
return error of parameter ipar
Definition: TFumili.cxx:833
void BuildArrays()
Allocates memory for internal arrays.
Definition: TFumili.cxx:174
TString * fANames
Definition: TFumili.h:72
int d
Definition: tornado.py:11
TLine l1(2.5, 4.5, 15.5, 4.5)
Int_t fNpar
Definition: TFumili.h:29
static const Double_t gMINDOUBLE
Definition: TFumili.cxx:124
std::vector< std::vector< double > > Data
Double_t EvalTFN(Double_t *, Double_t *)
Evaluate theoretical function df: array of partial derivatives X: vector of theoretical function argu...
Definition: TFumili.cxx:363
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
virtual Int_t GetNumberFreeParameters() const
return the number of free parameters
Definition: TFumili.cxx:821
TString fCword
Definition: TFumili.h:73
Bool_t fDEBUG
Definition: TFumili.h:40
Bool_t fLogLike
Definition: TFumili.h:41
Double_t * fZ0
Definition: TFumili.h:44
virtual Double_t GetErrorXlow(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1400
Double_t * fCmPar
Definition: TFumili.h:65
TObject * fUserFunc
TThread * t[5]
Definition: threadsh1.C:13
Int_t fNED2
Definition: TFumili.h:27
Int_t ExecuteSetCommand(Int_t)
Called from TFumili::ExecuteCommand in case of "SET xxx" and "SHOW xxx".
Definition: TFumili.cxx:564
Double_t * GetX() const
Definition: TGraph.h:139
Double_t * GetPL0() const
Definition: TFumili.h:105
A 3-Dim function with parameters.
Definition: TF3.h:30
virtual Bool_t IsFixed(Int_t ipar) const
return kTRUE if parameter ipar is fixed, kFALSE othersise)
Definition: TFumili.cxx:1075
static const Double_t gMAXDOUBLE
Definition: TFumili.cxx:123
void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
-*-*-*-*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.
Definition: TFumili.cxx:2110
void DeleteArrays()
Deallocates memory.
Definition: TFumili.cxx:260
Double_t fAKAPPA
Definition: TFumili.h:70
char * Form(const char *fmt,...)
static TVirtualFitter * GetFitter()
static: return the current Fitter
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition: TF1.cxx:3367
TLine * l
Definition: textangle.C:4
Double_t fS
Definition: TFumili.h:67
int Integral
Definition: Foption.h:44
Int_t fMaxParam
Definition: TFumili.h:23
A 2-Dim function with parameters.
Definition: TF2.h:33
virtual Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs)
Execute MINUIT commands.
Definition: TFumili.cxx:393
Double_t * fGr
Definition: TFumili.h:48
virtual void SetFitMethod(const char *name)
ret fit method (chisquare or loglikelihood)
Definition: TFumili.cxx:1656
TGraphErrors * gr
Definition: legend1.C:25
int W1
Definition: Foption.h:36
#define Printf
Definition: TGeoToOCC.h:18
virtual void FitChisquareI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method The "I"ntegral method is used for each point t...
Definition: TFumili.cxx:1841
Int_t fNED12
Definition: TFumili.h:28
virtual Double_t Derivative(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the first derivative of the function at point x, computed by Richardson's extrapolation metho...
Definition: TF1.cxx:829
Bool_t fWARN
Definition: TFumili.h:39
Int_t fLastFixed
Definition: TFumili.h:33
Double_t * fAMX
Definition: TFumili.h:60
virtual void FitLikelihoodI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood...
Definition: TFumili.cxx:2004
Int_t fENDFLG
Definition: TFumili.h:34
double Double_t
Definition: RtypesCore.h:55
Double_t * fSumLog
Definition: TFumili.h:50
Int_t fNstepDec
Definition: TFumili.h:30
virtual Double_t * GetCovarianceMatrix() const
return a pointer to the covariance matrix
Definition: TFumili.cxx:790
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
virtual ~TFumili()
TFumili destructor.
Definition: TFumili.cxx:212
void Derivatives(Double_t *, Double_t *)
Calculates partial derivatives of theoretical function.
Definition: TFumili.cxx:292
Double_t y[n]
Definition: legend1.C:17
Double_t ey[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:80
virtual void FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method Default method (function evaluated at center o...
Definition: TFumili.cxx:1769
static const double eu
Definition: Vavilov.cxx:52
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition: TF1.cxx:2198
#define name(a, b)
Definition: linkTestLib0.cpp:5
Double_t * fAMN
Definition: TFumili.h:61
virtual Double_t GetErrorY(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1380
virtual void ReleaseParameter(Int_t ipar)
Releases parameter number ipar.
Definition: TFumili.cxx:1612
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3376
virtual Double_t GetSumLog(Int_t)
return Sum(log(i) i=0,n used by log likelihood fits
Definition: TFumili.cxx:930
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
virtual TObject * GetObjectFit() const
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
TF1 * f1
Definition: legend1.C:11
void SetData(Double_t *, Int_t, Int_t)
Sets pointer to data array provided by user.
Definition: TFumili.cxx:1644
Double_t * fZ
Definition: TFumili.h:47
void(* fFCN)(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Double_t * fPL0
Definition: TFumili.h:55
virtual void Clear(Option_t *opt="")
Resets all parameter names, values and errors to zero.
Definition: TFumili.cxx:238
Int_t fNfcn
Definition: TFumili.h:25
Double_t fEPS
Definition: TFumili.h:68
Double_t * GetZ() const
Definition: TFumili.h:112
float type_of_call hi(const int &, const int &)
void SetParNumber(Int_t ParNum)
Definition: TFumili.h:122
Double_t * fA
Definition: TFumili.h:54
virtual void FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood...
Definition: TFumili.cxx:1916
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition: TF1.h:404
Bool_t fGRAD
Definition: TFumili.h:38
const Bool_t kTRUE
Definition: Rtypes.h:91
ClassImp(TFumili)
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1192
void GraphFitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
*-*-*-*-*-*Minimization function for Graphs using a Chisquare method*-*-*-*-* *-* ===================...
Definition: TFumili.cxx:2148
virtual Int_t GetErrors(Int_t ipar, Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
Return errors after MINOs not implemented.
Definition: TFumili.cxx:889
virtual const char * GetParName(Int_t ipar) const
return name of parameter ipar
Definition: TFumili.cxx:879
float value
Definition: math.cpp:443
Int_t fINDFLG[5]
Definition: TFumili.h:35
const Int_t n
Definition: legend1.C:16
Double_t fGT
Definition: TFumili.h:71
virtual Int_t GetNpar() const
Definition: TF1.h:342
Int_t SGZ()
Evaluates objective function ( chi-square ), gradients and Z-matrix using data provided by user via T...
Definition: TFumili.cxx:1701
Int_t npfits
Definition: fit2dHist.C:46
double log(double)
Int_t fNmaxIter
Definition: TFumili.h:32
Int_t fNlimMul
Definition: TFumili.h:31
int ii
Definition: hprod.C:34
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904