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