Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TLinearFitter.cxx
Go to the documentation of this file.
1// @(#)root/minuit:$Id$
2// Author: Anna Kreshuk 04/03/2005
3
4/*************************************************************************
5 * Copyright (C) 1995-2005, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "TLinearFitter.h"
13#include "TMath.h"
14#include "TDecompChol.h"
15#include "TGraph.h"
16#include "TGraph2D.h"
17#include "TMultiGraph.h"
18#include "TRandom2.h"
19#include "TObjString.h"
20#include "TF2.h"
21#include "TH1.h"
22#include "TList.h"
23#include "TClass.h"
24#include "TBuffer.h"
25#include "TROOT.h"
26#include "strlcpy.h"
27#include "snprintf.h"
28
30
31
32std::map<TString,TFormula*> TLinearFitter::fgFormulaMap;
33
34
35
36/**
37
38\class TLinearFitter
39
40\note An alternative to this class is to use ROOT::Fit::Fitter, calling the LinearFit() method.
41
42\ingroup MinuitOld
43
44The Linear Fitter - For fitting functions that are LINEAR IN PARAMETERS
45
46## The Linear Fitter
47
48Linear fitter is used to fit a set of data points with a linear
49combination of specified functions. Note, that "linear" in the name
50stands only for the model dependency on parameters, the specified
51functions can be nonlinear.
52The general form of this kind of model is
53~~~~
54 y(x) = a[0] + a[1]*f[1](x)+...a[n]*f[n](x)
55~~~~
56
57Functions f are fixed functions of x. For example, fitting with a
58polynomial is linear fitting in this sense.
59
60### Introduction
61
62#### The fitting method
63
64The fit is performed using the Normal Equations method with Cholesky
65decomposition.
66
67#### Why should it be used?
68
69The linear fitter is considerably faster than general non-linear
70fitters and doesn't require to set the initial values of parameters.
71
72### Using the fitter:
73
74### 1.Adding the data points:
75
76#### 1.1 To store or not to store the input data?
77 - There are 2 options in the constructor - to store or not
78 store the input data. The advantages of storing the data
79 are that you'll be able to reset the fitting model without
80 adding all the points again, and that for very large sets
81 of points the chisquare is calculated more precisely.
82 The obvious disadvantage is the amount of memory used to
83 keep all the points.
84 - Before you start adding the points, you can change the
85 store/not store option by StoreData() method.
86
87#### 1.2 The data can be added:
88 - simply point by point - AddPoint() method
89 - an array of points at once:
90 If the data is already stored in some arrays, this data
91 can be assigned to the linear fitter without physically
92 coping bytes, thanks to the Use() method of
93 TVector and TMatrix classes - AssignData() method
94
95### 2.Setting the formula
96
97#### 2.1 The linear formula syntax:
98 -Additive parts are separated by 2 plus signs "++"
99 --for example "1 ++ x" - for fitting a straight line
100 -All standard functions, undrestood by TFormula, can be used
101 as additive parts
102 --TMath functions can be used too
103 -Functions, used as additive parts, shouldn't have any parameters,
104 even if those parameters are set.
105 --for example, if normalizing a sum of a gaus(0, 1) and a
106 gaus(0, 2), don't use the built-in "gaus" of TFormula,
107 because it has parameters, take TMath::Gaus(x, 0, 1) instead.
108 -Polynomials can be used like "pol3", .."polN"
109 -If fitting a more than 3-dimensional formula, variables should
110 be numbered as follows:
111 -- x[0], x[1], x[2]... For example, to fit "1 ++ x[0] ++ x[1] ++ x[2] ++ x[3]*x[3]"
112
113#### 2.2 Setting the formula:
114
115##### 2.2.1 If fitting a 1-2-3-dimensional formula, one can create a
116 TF123 based on a linear expression and pass this function
117 to the fitter:
118 --Example:
119~~~~
120 TLinearFitter *lf = new TLinearFitter();
121 TF2 *f2 = new TF2("f2", "x ++ y ++ x*x*y*y", -2, 2, -2, 2);
122 lf->SetFormula(f2);
123~~~~
124 --The results of the fit are then stored in the function,
125 just like when the TH1::Fit or TGraph::Fit is used
126 --A linear function of this kind is by no means different
127 from any other function, it can be drawn, evaluated, etc.
128
129 --For multidimensional fitting, TFormulas of the form:
130 x[0]++...++x[n] can be used
131##### 2.2.2 There is no need to create the function if you don't want to,
132 the formula can be set by expression:
133 --Example:
134~~~~
135 // 2 is the number of dimensions
136 TLinearFitter *lf = new TLinearFitter(2);
137 lf->SetFormula("x ++ y ++ x*x*y*y");
138~~~~
139
140##### 2.2.3 The fastest functions to compute are polynomials and hyperplanes.
141 --Polynomials are set the usual way: "pol1", "pol2",...
142 --Hyperplanes are set by expression "hyp3", "hyp4", ...
143 ---The "hypN" expressions only work when the linear fitter
144 is used directly, not through TH1::Fit or TGraph::Fit.
145 To fit a graph or a histogram with a hyperplane, define
146 the function as "1++x++y".
147 ---A constant term is assumed for a hyperplane, when using
148 the "hypN" expression, so "hyp3" is in fact fitting with
149 "1++x++y++z" function.
150 --Fitting hyperplanes is much faster than fitting other
151 expressions so if performance is vital, calculate the
152 function values beforehand and give them to the fitter
153 as variables
154 --Example:
155 You want to fit "sin(x)|cos(2*x)" very fast. Calculate
156 sin(x) and cos(2*x) beforehand and store them in array *data.
157 Then:
158 TLinearFitter *lf=new TLinearFitter(2, "hyp2");
159 lf->AssignData(npoint, 2, data, y);
160
161#### 2.3 Resetting the formula
162
163##### 2.3.1 If the input data is stored (or added via AssignData() function),
164 the fitting formula can be reset without re-adding all the points.
165 --Example:
166~~~~
167 TLinearFitter *lf=new TLinearFitter("1++x++x*x");
168 lf->AssignData(n, 1, x, y, e);
169 lf->Eval()
170 //looking at the parameter significance, you see,
171 // that maybe the fit will improve, if you take out
172 // the constant term
173 lf->SetFormula("x++x*x");
174 lf->Eval();
175 ...
176~~~~
177
178##### 2.3.2 If the input data is not stored, the fitter will have to be
179 cleared and the data will have to be added again to try a
180 different formula.
181
182### 3.Accessing the fit results
183
184#### 3.1 There are methods in the fitter to access all relevant information:
185 --GetParameters, GetCovarianceMatrix, etc
186 --the t-values of parameters and their significance can be reached by
187 GetParTValue() and GetParSignificance() methods
188
189#### 3.2 If fitting with a pre-defined TF123, the fit results are also
190 written into this function.
191
192
193### 4.Robust fitting - Least Trimmed Squares regression (LTS)
194 Outliers are atypical(by definition), infrequant observations; data points
195 which do not appear to follow the characteristic distribution of the rest
196 of the data. These may reflect genuine properties of the underlying
197 phenomenon(variable), or be due to measurement errors or anomalies which
198 shouldn't be modelled. (StatSoft electronic textbook)
199
200 Even a single gross outlier can greatly influence the results of least-
201 squares fitting procedure, and in this case use of robust(resistant) methods
202 is recommended.
203
204 The method implemented here is based on the article and algorithm:
205 "Computing LTS Regression for Large Data Sets" by
206 P.J.Rousseeuw and Katrien Van Driessen
207 The idea of the method is to find the fitting coefficients for a subset
208 of h observations (out of n) with the smallest sum of squared residuals.
209 The size of the subset h should lie between (npoints + nparameters +1)/2
210 and n, and represents the minimal number of good points in the dataset.
211 The default value is set to (npoints + nparameters +1)/2, but of course
212 if you are sure that the data contains less outliers it's better to change
213 h according to your data.
214
215 To perform a robust fit, call EvalRobust() function instead of Eval() after
216 adding the points and setting the fitting function.
217 Note, that standard errors on parameters are not computed!
218
219*/
220
221
222
223////////////////////////////////////////////////////////////////////////////////
224///default c-tor, input data is stored
225///If you don't want to store the input data,
226///run the function StoreData(kFALSE) after constructor
227
229 fVal()
230{
231 fChisquare =0;
232 fNpoints =0;
233 fNdim =0;
234 fY2 =0;
235 fY2Temp =0;
236 fNfixed =0;
237 fIsSet =kFALSE;
238 fFormula =nullptr;
239 fFixedParams=nullptr;
240 fSpecial =0;
241 fInputFunction=nullptr;
244 fNfunctions = 0;
245 fFormulaSize = 0;
246 fH = 0;
247}
248
249////////////////////////////////////////////////////////////////////////////////
250///The parameter stands for number of dimensions in the fitting formula
251///The input data is stored. If you don't want to store the input data,
252///run the function StoreData(kFALSE) after constructor
253
255 fVal()
256{
257 fNdim =ndim;
258 fNpoints =0;
259 fY2 =0;
260 fY2Temp =0;
261 fNfixed =0;
262 fFixedParams=nullptr;
263 fFormula =nullptr;
264 fIsSet =kFALSE;
265 fChisquare=0;
266 fSpecial =0;
267 fInputFunction=nullptr;
269 fRobust = kFALSE;
270 fNfunctions = 0;
271 fFormulaSize = 0;
272 fH = 0;
273}
274
275////////////////////////////////////////////////////////////////////////////////
276///First parameter stands for number of dimensions in the fitting formula
277///Second parameter is the fitting formula: see class description for formula syntax
278///Options:
279///The option is to store or not to store the data
280///If you don't want to store the data, choose "" for the option, or run
281///StoreData(kFalse) member function after the constructor
282
283TLinearFitter::TLinearFitter(Int_t ndim, const char *formula, Option_t *opt)
284{
285 fNdim=ndim;
286 fNpoints=0;
287 fChisquare=0;
288 fY2=0;
289 fNfixed=0;
290 fFixedParams=nullptr;
291 fSpecial=0;
292 fInputFunction=nullptr;
293 fFormula = nullptr;
294 TString option=opt;
295 option.ToUpper();
296 if (option.Contains("D"))
298 else
301 SetFormula(formula);
302}
303
304////////////////////////////////////////////////////////////////////////////////
305///This constructor uses a linear function. How to create it?
306///TFormula now accepts formulas of the following kind:
307///TFormula("f", "x++y++z++x*x") or
308///TFormula("f", "x[0]++x[1]++x[2]*x[2]");
309///Other than the look, it's in no
310///way different from the regular formula, it can be evaluated,
311///drawn, etc.
312///The option is to store or not to store the data
313///If you don't want to store the data, choose "" for the option, or run
314///StoreData(kFalse) member function after the constructor
315
317{
318 fNdim=function->GetNdim();
319 if (!function->IsLinear()){
320 Int_t number=function->GetNumber();
321 if (number<299 || number>310){
322 Error("TLinearFitter", "Trying to fit with a nonlinear function");
323 return;
324 }
325 }
326 fNpoints=0;
327 fChisquare=0;
328 fY2=0;
329 fNfixed=0;
330 fFixedParams=nullptr;
331 fSpecial=0;
332 fFormula = nullptr;
333 TString option=opt;
334 option.ToUpper();
335 if (option.Contains("D"))
337 else
341 fInputFunction=nullptr;
342
343 SetFormula(function);
344}
345
346////////////////////////////////////////////////////////////////////////////////
347/// Copy ctor
348
350 TVirtualFitter(tlf),
351 fParams(tlf.fParams),
352 fParCovar(tlf.fParCovar),
353 fTValues(tlf.fTValues),
354 fParSign(tlf.fParSign),
355 fDesign(tlf.fDesign),
356 fDesignTemp(tlf.fDesignTemp),
357 fDesignTemp2(tlf.fDesignTemp2),
358 fDesignTemp3(tlf.fDesignTemp3),
359 fAtb(tlf.fAtb),
360 fAtbTemp(tlf.fAtbTemp),
361 fAtbTemp2(tlf.fAtbTemp2),
362 fAtbTemp3(tlf.fAtbTemp3),
363 fFunctions( * (TObjArray *)tlf.fFunctions.Clone()),
364 fY(tlf.fY),
365 fY2(tlf.fY2),
366 fY2Temp(tlf.fY2Temp),
367 fX(tlf.fX),
368 fE(tlf.fE),
369 fInputFunction(tlf.fInputFunction),
370 fVal(),
371 fNpoints(tlf.fNpoints),
372 fNfunctions(tlf.fNfunctions),
373 fFormulaSize(tlf.fFormulaSize),
374 fNdim(tlf.fNdim),
375 fNfixed(tlf.fNfixed),
376 fSpecial(tlf.fSpecial),
377 fFormula(nullptr),
378 fIsSet(tlf.fIsSet),
379 fStoreData(tlf.fStoreData),
380 fChisquare(tlf.fChisquare),
381 fH(tlf.fH),
382 fRobust(tlf.fRobust),
383 fFitsample(tlf.fFitsample),
384 fFixedParams(nullptr)
385{
386 // make a deep copy of managed objects
387 // fFormula, fFixedParams and fFunctions
388
389 if ( tlf.fFixedParams && fNfixed > 0 ) {
391 for(Int_t i=0; i<fNfixed; ++i)
392 fFixedParams[i]=tlf.fFixedParams[i];
393 }
394 if (tlf.fFormula) {
395 fFormula = new char[fFormulaSize+1];
396 strlcpy(fFormula,tlf.fFormula,fFormulaSize+1);
397 }
398
399}
400
401
402////////////////////////////////////////////////////////////////////////////////
403/// Linear fitter cleanup.
404
406{
407 if (fFormula) {
408 delete [] fFormula;
409 fFormula = nullptr;
410 }
411 if (fFixedParams) {
412 delete [] fFixedParams;
413 fFixedParams = nullptr;
414 }
415 fInputFunction = nullptr;
416
417 //fFunctions.Delete();
419
420}
421
422////////////////////////////////////////////////////////////////////////////////
423/// Assignment operator
424
426{
427 if(this!=&tlf) {
428
438
439 fAtb.ResizeTo(tlf.fAtb); fAtb=tlf.fAtb;
443
444 // use clear instead of delete
447
448 fY.ResizeTo(tlf.fY); fY = tlf.fY;
449 fX.ResizeTo(tlf.fX); fX = tlf.fX;
450 fE.ResizeTo(tlf.fE); fE = tlf.fE;
451
452 fY2 = tlf.fY2;
453 fY2Temp = tlf.fY2Temp;
454 for(Int_t i = 0; i < 1000; i++) fVal[i] = tlf.fVal[i];
455
456 if(fInputFunction) { delete fInputFunction; fInputFunction = nullptr; }
458
459 fNpoints=tlf.fNpoints;
462 fNdim=tlf.fNdim;
463 fNfixed=tlf.fNfixed;
464 fSpecial=tlf.fSpecial;
465
466 if(fFormula) { delete [] fFormula; fFormula = nullptr; }
467 if (tlf.fFormula) {
468 fFormula = new char[fFormulaSize+1];
469 strlcpy(fFormula,tlf.fFormula,fFormulaSize+1);
470 }
471
472 fIsSet=tlf.fIsSet;
475
476 fH=tlf.fH;
477 fRobust=tlf.fRobust;
479
480 if(fFixedParams) { delete [] fFixedParams; fFixedParams = nullptr; }
481 if ( tlf.fFixedParams && fNfixed > 0 ) {
483 for(Int_t i=0; i< fNfixed; ++i)
484 fFixedParams[i]=tlf.fFixedParams[i];
485 }
486
487 }
488 return *this;
489}
490
491////////////////////////////////////////////////////////////////////////////////
492///Add another linear fitter to this linear fitter. Points and Design matrices
493///are added, but the previous fitting results (if any) are deleted.
494///Fitters must have same formulas (this is not checked). Fixed parameters are not changed
495
497{
498 fParams.Zero();
499 fParCovar.Zero();
500 fTValues.Zero();
501 fParSign.Zero();
502
503 fDesign += tlf->fDesign;
504
505 fDesignTemp += tlf->fDesignTemp;
508 fAtb += tlf->fAtb;
509 fAtbTemp += tlf->fAtbTemp;
510 fAtbTemp2 += tlf->fAtbTemp2;
511 fAtbTemp3 += tlf->fAtbTemp3;
512
513 if (fStoreData){
515 Int_t newsize = fNpoints+tlf->fNpoints;
516 if (size<newsize){
517 fY.ResizeTo(newsize);
518 fE.ResizeTo(newsize);
519 fX.ResizeTo(newsize, fNdim);
520 }
521 for (Int_t i=fNpoints; i<newsize; i++){
522 fY(i)=tlf->fY(i-fNpoints);
523 fE(i)=tlf->fE(i-fNpoints);
524 for (Int_t j=0; j<fNdim; j++){
525 fX(i,j)=tlf->fX(i-fNpoints, j);
526 }
527 }
528 }
529 fY2 += tlf->fY2;
530 fY2Temp += tlf->fY2Temp;
531 fNpoints += tlf->fNpoints;
532 //fInputFunction=(TFormula*)tlf.fInputFunction->Clone();
533
534 fChisquare=0;
535 fH=0;
536 fRobust=false;
537}
538
539////////////////////////////////////////////////////////////////////////////////
540///Adds 1 point to the fitter.
541///First parameter stands for the coordinates of the point, where the function is measured
542///Second parameter - the value being fitted
543///Third parameter - weight(measurement error) of this point (=1 by default)
544
546{
547 Int_t size;
548 fNpoints++;
549 if (fStoreData){
551 if (size<fNpoints){
555 }
556
557 Int_t j=fNpoints-1;
558 fY(j)=y;
559 fE(j)=e;
560 for (Int_t i=0; i<fNdim; i++)
561 fX(j,i)=x[i];
562 }
563 //add the point to the design matrix, if the formula has been set
564 // (LM: why condition !fRobust ?? - remove it to fix Coverity 11602)
565 // when 100< fSpecial < 200 (Polynomial) fFunctions is not empty
566 // (see implementation of SetFormula method)
567 if (fFunctions.IsEmpty()&&(!fInputFunction)&&(fSpecial<=200)){
568 Error("AddPoint", "Point can't be added, because the formula hasn't been set");
569 return;
570 }
571 if (!fRobust) AddToDesign(x, y, e);
572}
573
574////////////////////////////////////////////////////////////////////////////////
575///This function is to use when you already have all the data in arrays
576///and don't want to copy them into the fitter. In this function, the Use() method
577///of TVectorD and TMatrixD is used, so no bytes are physically moved around.
578///First parameter - number of points to fit
579///Second parameter - number of variables in the model
580///Third parameter - the variables of the model, stored in the following way:
581///(x0(0), x1(0), x2(0), x3(0), x0(1), x1(1), x2(1), x3(1),...
582
584{
585 if (npoints<fNpoints){
586 Error("AddData", "Those points are already added");
587 return;
588 }
589 Bool_t same=kFALSE;
590 if (fX.GetMatrixArray()==x && fY.GetMatrixArray()==y){
591 if (e){
592 if (fE.GetMatrixArray()==e)
593 same=kTRUE;
594 }
595 }
596
597 fX.Use(npoints, xncols, x);
598 fY.Use(npoints, y);
599 if (e)
600 fE.Use(npoints, e);
601 else {
602 fE.ResizeTo(npoints);
603 fE=1;
604 }
605 Int_t xfirst;
606 if (!fFunctions.IsEmpty() || fInputFunction || fSpecial>200) {
607 if (same)
608 xfirst=fNpoints;
609
610 else
611 xfirst=0;
612 for (Int_t i=xfirst; i<npoints; i++)
613 AddToDesign(TMatrixDRow(fX, i).GetPtr(), fY(i), fE(i));
614 }
615 fNpoints=npoints;
616}
617
618////////////////////////////////////////////////////////////////////////////////
619///Add a point to the AtA matrix and to the Atb vector.
620
622{
623
624
625 Int_t i, j, ii;
626 y/=e;
627
628// Double_t fVal[1000];
629
630 if ((fSpecial>100)&&(fSpecial<200)){
631 //polynomial fitting
632 Int_t npar=fSpecial-100;
633 fVal[0]=1;
634 for (i=1; i<npar; i++)
635 fVal[i]=fVal[i-1]*x[0];
636 for (i=0; i<npar; i++)
637 fVal[i]/=e;
638 } else if (fSpecial>200){
639 //Hyperplane fitting. Constant term is added
640 Int_t npar=fSpecial-201;
641 fVal[0]=1./e;
642 for (i=0; i<npar; i++)
643 fVal[i+1]=x[i]/e;
644 } else {
645 //general case
646 for (ii=0; ii<fNfunctions; ii++){
647 if (!fFunctions.IsEmpty()){
648 // ffunctions can be TF1 or TFormula depending on how they are created
649 TObject * obj = fFunctions.UncheckedAt(ii);
650 if (obj->IsA() == TFormula::Class() ) {
651 TFormula *f1 = (TFormula*)(obj);
652 fVal[ii]=f1->EvalPar(x)/e;
653 }
654 else if (obj->IsA() == TF1::Class() ) {
655 TF1 *f1 = (TF1*)(obj);
656 fVal[ii]=f1->EvalPar(x)/e;
657 }
658 else {
659 Error("AddToDesign","Basis Function %s is of an invalid type %s",obj->GetName(),obj->IsA()->GetName());
660 return;
661 }
662 } else {
664 if (!f){
665 Error("AddToDesign","Function %s has no linear parts - maybe missing a ++ in the formula expression",fInputFunction->GetName());
666 return;
667 }
668 fVal[ii]=f->EvalPar(x)/e;
669 }
670 }
671 }
672 //additional matrices for numerical stability
673 for (i=0; i<fNfunctions; i++){
674 for (j=0; j<i; j++)
675 fDesignTemp3(j, i)+=fVal[i]*fVal[j];
676 fDesignTemp3(i, i)+=fVal[i]*fVal[i];
677 fAtbTemp3(i)+=fVal[i]*y;
678
679 }
680 fY2Temp+=y*y;
682
683 if (fNpoints % 100 == 0 && fNpoints>100){
687 fAtbTemp3.Zero();
688 if (fNpoints % 10000 == 0 && fNpoints>10000){
692 fAtbTemp2.Zero();
693 fY2+=fY2Temp;
694 fY2Temp=0;
695 if (fNpoints % 1000000 == 0 && fNpoints>1000000){
698 fAtb+=fAtbTemp;
699 fAtbTemp.Zero();
700 }
701 }
702 }
703}
704
705////////////////////////////////////////////////////////////////////////////////
706
708{
709 if (fDesignTemp3.GetNrows()>0){
718 fAtb+=fAtbTemp;
719 fAtbTemp3.Zero();
720 fAtbTemp2.Zero();
721 fAtbTemp.Zero();
722
723 fY2+=fY2Temp;
724 fY2Temp=0;
725 }
726}
727
728////////////////////////////////////////////////////////////////////////////////
729///Clears everything. Used in TH1::Fit and TGraph::Fit().
730
732{
733 fParams.Clear();
735 fTValues.Clear();
736 fParSign.Clear();
737 fDesign.Clear();
741 fAtb.Clear();
742 fAtbTemp.Clear();
746 fInputFunction=nullptr;
747 fY.Clear();
748 fX.Clear();
749 fE.Clear();
750
751 fNpoints=0;
752 fNfunctions=0;
753 fFormulaSize=0;
754 fNdim=0;
755 if (fFormula) delete [] fFormula;
756 fFormula=nullptr;
757 fIsSet=false;
758 if (fFixedParams) delete [] fFixedParams;
759 fFixedParams=nullptr;
760
761 fChisquare=0;
762 fY2=0;
763 fSpecial=0;
766}
767
768////////////////////////////////////////////////////////////////////////////////
769///To be used when different sets of points are fitted with the same formula.
770
772{
773 fDesign.Zero();
774 fAtb.Zero();
778 fAtbTemp.Zero();
779 fAtbTemp2.Zero();
780 fAtbTemp3.Zero();
781
782 fParams.Zero();
783 fParCovar.Zero();
784 fTValues.Zero();
785 fParSign.Zero();
786
787 for (Int_t i=0; i<fNfunctions; i++)
788 fFixedParams[i]=false;
789 fChisquare=0;
790 fNpoints=0;
791
792}
793
794////////////////////////////////////////////////////////////////////////////////
795///Calculates the chisquare.
796
798{
799 Int_t i, j;
800 Double_t sumtotal2;
801 Double_t temp, temp2;
802
803 if (!fStoreData){
804 sumtotal2 = 0;
805 for (i=0; i<fNfunctions; i++){
806 for (j=0; j<i; j++){
807 sumtotal2 += 2*fParams(i)*fParams(j)*fDesign(j, i);
808 }
809 sumtotal2 += fParams(i)*fParams(i)*fDesign(i, i);
810 sumtotal2 -= 2*fParams(i)*fAtb(i);
811 }
812 sumtotal2 += fY2;
813 } else {
814 sumtotal2 = 0;
815 if (fInputFunction){
816 for (i=0; i<fNpoints; i++){
817 temp = fInputFunction->EvalPar(TMatrixDRow(fX, i).GetPtr());
818 temp2 = (fY(i)-temp)*(fY(i)-temp);
819 temp2 /= fE(i)*fE(i);
820 sumtotal2 += temp2;
821 }
822 } else {
823 sumtotal2 = 0;
824 Double_t val[100];
825 for (Int_t point=0; point<fNpoints; point++){
826 temp = 0;
827 if ((fSpecial>100)&&(fSpecial<200)){
828 Int_t npar = fSpecial-100;
829 val[0] = 1;
830 for (i=1; i<npar; i++)
831 val[i] = val[i-1]*fX(point, 0);
832 for (i=0; i<npar; i++)
833 temp += fParams(i)*val[i];
834 } else {
835 if (fSpecial>200) {
836 //hyperplane case
837 Int_t npar = fSpecial-201;
838 temp+=fParams(0);
839 for (i=0; i<npar; i++)
840 temp += fParams(i+1)*fX(point, i);
841 } else {
842 for (j=0; j<fNfunctions; j++) {
844 val[j] = f1->EvalPar(TMatrixDRow(fX, point).GetPtr());
845 temp += fParams(j)*val[j];
846 }
847 }
848 }
849 temp2 = (fY(point)-temp)*(fY(point)-temp);
850 temp2 /= fE(point)*fE(point);
851 sumtotal2 += temp2;
852 }
853 }
854 }
855 fChisquare = sumtotal2;
856
857}
858
859////////////////////////////////////////////////////////////////////////////////
860/// Computes parameters' t-values and significance
861
863{
864 for (Int_t i=0; i<fNfunctions; i++){
865 fTValues(i) = fParams(i)/(TMath::Sqrt(fParCovar(i, i)));
867 }
868}
869
870////////////////////////////////////////////////////////////////////////////////
871/// Perform the fit and evaluate the parameters
872/// Returns 0 if the fit is ok, 1 if there are errors
873
875{
876 Double_t e;
877 if (fFunctions.IsEmpty()&&(!fInputFunction)&&(fSpecial<=200)){
878 Error("TLinearFitter::Eval", "The formula hasn't been set");
879 return 1;
880 }
881 //
886
887 fChisquare=0;
888
889 if (!fIsSet){
891 if (!update){
892 //no points to fit
893 fParams.Zero();
894 fParCovar.Zero();
895 fTValues.Zero();
896 fParSign.Zero();
897 fChisquare=0;
898 if (fInputFunction){
900 for (Int_t i=0; i<fNfunctions; i++)
901 ((TF1*)fInputFunction)->SetParError(i, 0);
902 ((TF1*)fInputFunction)->SetChisquare(0);
903 ((TF1*)fInputFunction)->SetNDF(0);
904 ((TF1*)fInputFunction)->SetNumberFitPoints(0);
905 }
906 return 1;
907 }
908 }
909
911
912//fixing fixed parameters, if there are any
913 Int_t i, ii, j=0;
914 if (fNfixed>0){
915 for (ii=0; ii<fNfunctions; ii++)
916 fDesignTemp(ii, fNfixed) = fAtb(ii);
917 for (i=0; i<fNfunctions; i++){
918 if (fFixedParams[i]){
919 for (ii=0; ii<i; ii++)
920 fDesignTemp(ii, j) = fDesign(ii, i);
921 for (ii=i; ii<fNfunctions; ii++)
922 fDesignTemp(ii, j) = fDesign(i, ii);
923 j++;
924 for (ii=0; ii<fNfunctions; ii++){
925 fAtb(ii)-=fParams(i)*(fDesignTemp(ii, j-1));
926 }
927 }
928 }
929 for (i=0; i<fNfunctions; i++){
930 if (fFixedParams[i]){
931 for (ii=0; ii<fNfunctions; ii++){
932 fDesign(ii, i) = 0;
933 fDesign(i, ii) = 0;
934 }
935 fDesign (i, i) = 1;
936 fAtb(i) = fParams(i);
937 }
938 }
939 }
940
941 TDecompChol chol(fDesign);
942 Bool_t ok;
943 TVectorD coef(fNfunctions);
944 coef=chol.Solve(fAtb, ok);
945 if (!ok){
946 Error("Eval","Matrix inversion failed");
947 fParams.Zero();
948 fParCovar.Zero();
949 fTValues.Zero();
950 fParSign.Zero();
951 if (fInputFunction){
954 ((TF1*)fInputFunction)->SetChisquare(0);
956 ((TF1*)fInputFunction)->SetNumberFitPoints(fNpoints);
957 }
958 }
959 return 1;
960 }
961 fParams=coef;
962 fParCovar=chol.Invert();
963
964 if (fInputFunction){
967 for (i=0; i<fNfunctions; i++){
968 e = TMath::Sqrt(fParCovar(i, i));
969 ((TF1*)fInputFunction)->SetParError(i, e);
970 }
971 if (!fObjectFit)
972 ((TF1*)fInputFunction)->SetChisquare(GetChisquare());
974 ((TF1*)fInputFunction)->SetNumberFitPoints(fNpoints);
975 }
976 }
977
978 //if parameters were fixed, change the design matrix back as it was before fixing
979 j = 0;
980 if (fNfixed>0){
981 for (i=0; i<fNfunctions; i++){
982 if (fFixedParams[i]){
983 for (ii=0; ii<i; ii++){
984 fDesign(ii, i) = fDesignTemp(ii, j);
985 fAtb(ii) = fDesignTemp(ii, fNfixed);
986 }
987 for (ii=i; ii<fNfunctions; ii++){
988 fDesign(i, ii) = fDesignTemp(ii, j);
989 fAtb(ii) = fDesignTemp(ii, fNfixed);
990 }
991 j++;
992 }
993 }
994 }
995 return 0;
996}
997
998////////////////////////////////////////////////////////////////////////////////
999///Fixes paramter `#ipar` at its current value.
1000
1002{
1003 if (fParams.NonZeros()<1){
1004 Error("FixParameter", "no value available to fix the parameter");
1005 return;
1006 }
1007 if (ipar>fNfunctions || ipar<0){
1008 Error("FixParameter", "illegal parameter value");
1009 return;
1010 }
1011 if (fNfixed==fNfunctions) {
1012 Error("FixParameter", "no free parameters left");
1013 return;
1014 }
1015 if (!fFixedParams)
1017 fFixedParams[ipar] = true;
1018 fNfixed++;
1019}
1020
1021////////////////////////////////////////////////////////////////////////////////
1022///Fixes parameter `#ipar` at value `parvalue`.
1023
1025{
1026 if (ipar>fNfunctions || ipar<0){
1027 Error("FixParameter", "illegal parameter value");
1028 return;
1029 }
1030 if (fNfixed==fNfunctions) {
1031 Error("FixParameter", "no free parameters left");
1032 return;
1033 }
1034 if(!fFixedParams)
1036 fFixedParams[ipar] = true;
1039 fParams(ipar) = parvalue;
1040 fNfixed++;
1041}
1042
1043////////////////////////////////////////////////////////////////////////////////
1044///Releases parameter `#ipar`.
1045
1047{
1048 if (ipar>fNfunctions || ipar<0){
1049 Error("ReleaseParameter", "illegal parameter value");
1050 return;
1051 }
1052 if (!fFixedParams[ipar]){
1053 Warning("ReleaseParameter","This parameter is not fixed\n");
1054 return;
1055 } else {
1056 fFixedParams[ipar] = false;
1057 fNfixed--;
1058 }
1059}
1060
1061////////////////////////////////////////////////////////////////////////////////
1062///Get the Atb vector - a vector, used for internal computations
1063
1065{
1066 if (v.GetNoElements()!=fAtb.GetNoElements())
1067 v.ResizeTo(fAtb.GetNoElements());
1068 v = fAtb;
1069}
1070
1071////////////////////////////////////////////////////////////////////////////////
1072/// Get the Chisquare.
1073
1075{
1076 if (fChisquare > 1e-16)
1077 return fChisquare;
1078 else {
1079 Chisquare();
1080 return fChisquare;
1081 }
1082}
1083
1084////////////////////////////////////////////////////////////////////////////////
1085///Computes point-by-point confidence intervals for the fitted function
1086///Parameters:
1087///n - number of points
1088///ndim - dimensions of points
1089///x - points, at which to compute the intervals, for ndim > 1
1090/// should be in order: (x0,y0, x1, y1, ... xn, yn)
1091///ci - computed intervals are returned in this array
1092///cl - confidence level, default=0.95
1093///
1094///NOTE, that this method can only be used when the fitting function inherits from a TF1,
1095///so it's not possible when the fitting function was set as a string or as a pure TFormula
1096
1098{
1099 if (fInputFunction){
1100 Double_t *grad = new Double_t[fNfunctions];
1101 Double_t *sum_vector = new Double_t[fNfunctions];
1102 Double_t c=0;
1104 Double_t t = TMath::StudentQuantile(0.5 + cl/2, df);
1105 Double_t chidf = TMath::Sqrt(fChisquare/df);
1106
1107 for (Int_t ipoint=0; ipoint<n; ipoint++){
1108 c=0;
1109 ((TF1*)(fInputFunction))->GradientPar(x+ndim*ipoint, grad);
1110 //multiply the covariance matrix by gradient
1111 for (Int_t irow=0; irow<fNfunctions; irow++){
1112 sum_vector[irow]=0;
1113 for (Int_t icol=0; icol<fNfunctions; icol++){
1114 sum_vector[irow]+=fParCovar(irow,icol)*grad[icol];
1115 }
1116 }
1117 for (Int_t i=0; i<fNfunctions; i++)
1118 c+=grad[i]*sum_vector[i];
1119 c=TMath::Sqrt(c);
1120 ci[ipoint]=c*t*chidf;
1121 }
1122
1123 delete [] grad;
1124 delete [] sum_vector;
1125 }
1126}
1127
1128////////////////////////////////////////////////////////////////////////////////
1129///Computes confidence intervals at level cl. Default is 0.95
1130///The TObject parameter can be a TGraphErrors, a TGraph2DErrors or a TH123.
1131///For Graphs, confidence intervals are computed for each point,
1132///the value of the graph at that point is set to the function value at that
1133///point, and the graph y-errors (or z-errors) are set to the value of
1134///the confidence interval at that point
1135///For Histograms, confidence intervals are computed for each bin center
1136///The bin content of this bin is then set to the function value at the bin
1137///center, and the bin error is set to the confidence interval value.
1138///Allowed combinations:
1139///Fitted object Passed object
1140///TGraph TGraphErrors, TH1
1141///TGraphErrors, AsymmErrors TGraphErrors, TH1
1142///TH1 TGraphErrors, TH1
1143///TGraph2D TGraph2DErrors, TH2
1144///TGraph2DErrors TGraph2DErrors, TH2
1145///TH2 TGraph2DErrors, TH2
1146///TH3 TH3
1147
1149{
1150 if (!fInputFunction) {
1151 Error("GetConfidenceIntervals", "The case of fitting not with a TFormula is not yet implemented");
1152 return;
1153 }
1154
1155 //TGraph//////////////////
1156
1157 if (obj->InheritsFrom(TGraph::Class())) {
1158 TGraph *gr = (TGraph*)obj;
1159 if (!gr->GetEY()){
1160 Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a graph");
1161 return;
1162 }
1164 Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a graph");
1165 return;
1166 }
1168 if (((TH1*)(fObjectFit))->GetDimension()>1){
1169 Error("GetConfidenceIntervals", "A TGraph2DErrors or a TH23 should be passed instead of a graph");
1170 return;
1171 }
1172 }
1173
1174 GetConfidenceIntervals(gr->GetN(),1,gr->GetX(), gr->GetEY(), cl);
1175 for (Int_t i=0; i<gr->GetN(); i++)
1176 gr->SetPoint(i, gr->GetX()[i], fInputFunction->Eval(gr->GetX()[i]));
1177 }
1178
1179 //TGraph2D///////////////
1180 else if (obj->InheritsFrom(TGraph2D::Class())) {
1181 TGraph2D *gr2 = (TGraph2D*)obj;
1182 if (!gr2->GetEZ()){
1183 Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a TGraph2D");
1184 return;
1185 }
1187 Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a TGraph2D");
1188 return;
1189 }
1191 if (((TH1*)(fObjectFit))->GetDimension()==1){
1192 Error("GetConfidenceIntervals", "A TGraphErrors or a TH1 should be passed instead of a graph");
1193 return;
1194 }
1195 }
1196 Double_t xy[2];
1197 Int_t np = gr2->GetN();
1198 Double_t *grad = new Double_t[fNfunctions];
1199 Double_t *sum_vector = new Double_t[fNfunctions];
1200 Double_t *x = gr2->GetX();
1201 Double_t *y = gr2->GetY();
1202 Double_t t = TMath::StudentQuantile(0.5 + cl/2, ((TF1*)(fInputFunction))->GetNDF());
1203 Double_t chidf = TMath::Sqrt(fChisquare/((TF1*)(fInputFunction))->GetNDF());
1204 Double_t c = 0;
1205 for (Int_t ipoint=0; ipoint<np; ipoint++){
1206 c=0;
1207 xy[0]=x[ipoint];
1208 xy[1]=y[ipoint];
1209 ((TF1*)(fInputFunction))->GradientPar(xy, grad);
1210 for (Int_t irow=0; irow<fNfunctions; irow++){
1211 sum_vector[irow]=0;
1212 for (Int_t icol=0; icol<fNfunctions; icol++)
1213 sum_vector[irow]+=fParCovar(irow, icol)*grad[icol];
1214 }
1215 for (Int_t i=0; i<fNfunctions; i++)
1216 c+=grad[i]*sum_vector[i];
1217 c=TMath::Sqrt(c);
1218 gr2->SetPoint(ipoint, xy[0], xy[1], fInputFunction->EvalPar(xy));
1219 gr2->GetEZ()[ipoint]=c*t*chidf;
1220 }
1221 delete [] grad;
1222 delete [] sum_vector;
1223 }
1224
1225 //TH1////////////////////////
1226 else if (obj->InheritsFrom(TH1::Class())) {
1228 if (((TH1*)obj)->GetDimension()>1){
1229 Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
1230 return;
1231 }
1232 }
1234 if (((TH1*)obj)->GetDimension()!=2){
1235 Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
1236 return;
1237 }
1238 }
1240 if (((TH1*)(fObjectFit))->GetDimension()!=((TH1*)(obj))->GetDimension()){
1241 Error("GetConfidenceIntervals", "Fitted and passed histograms have different number of dimensions");
1242 return;
1243 }
1244 }
1245
1246
1247 TH1 *hfit = (TH1*)obj;
1248 Double_t *grad = new Double_t[fNfunctions];
1249 Double_t *sum_vector = new Double_t[fNfunctions];
1250 Double_t x[3];
1251 Int_t hxfirst = hfit->GetXaxis()->GetFirst();
1252 Int_t hxlast = hfit->GetXaxis()->GetLast();
1253 Int_t hyfirst = hfit->GetYaxis()->GetFirst();
1254 Int_t hylast = hfit->GetYaxis()->GetLast();
1255 Int_t hzfirst = hfit->GetZaxis()->GetFirst();
1256 Int_t hzlast = hfit->GetZaxis()->GetLast();
1257
1258 TAxis *xaxis = hfit->GetXaxis();
1259 TAxis *yaxis = hfit->GetYaxis();
1260 TAxis *zaxis = hfit->GetZaxis();
1261 Double_t t = TMath::StudentQuantile(0.5 + cl/2, ((TF1*)(fInputFunction))->GetNDF());
1262 Double_t chidf = TMath::Sqrt(fChisquare/((TF1*)(fInputFunction))->GetNDF());
1263 Double_t c=0;
1264 for (Int_t binz=hzfirst; binz<=hzlast; binz++){
1265 x[2]=zaxis->GetBinCenter(binz);
1266 for (Int_t biny=hyfirst; biny<=hylast; biny++) {
1267 x[1]=yaxis->GetBinCenter(biny);
1268 for (Int_t binx=hxfirst; binx<=hxlast; binx++) {
1269 x[0]=xaxis->GetBinCenter(binx);
1270 ((TF1*)(fInputFunction))->GradientPar(x, grad);
1271 c=0;
1272 for (Int_t irow=0; irow<fNfunctions; irow++){
1273 sum_vector[irow]=0;
1274 for (Int_t icol=0; icol<fNfunctions; icol++)
1275 sum_vector[irow]+=fParCovar(irow, icol)*grad[icol];
1276 }
1277 for (Int_t i=0; i<fNfunctions; i++)
1278 c+=grad[i]*sum_vector[i];
1279 c=TMath::Sqrt(c);
1280 hfit->SetBinContent(binx, biny, binz, fInputFunction->EvalPar(x));
1281 hfit->SetBinError(binx, biny, binz, c*t*chidf);
1282 }
1283 }
1284 }
1285 delete [] grad;
1286 delete [] sum_vector;
1287 }
1288 else {
1289 Error("GetConfidenceIntervals", "This object type is not supported");
1290 return;
1291 }
1292}
1293
1294////////////////////////////////////////////////////////////////////////////////
1295///Returns covariance matrix
1296
1298{
1299 Double_t *p = const_cast<Double_t*>(fParCovar.GetMatrixArray());
1300 return p;
1301}
1302
1303////////////////////////////////////////////////////////////////////////////////
1304///Returns covariance matrix
1305
1307{
1308 if (matr.GetNrows()!=fNfunctions || matr.GetNcols()!=fNfunctions){
1310 }
1311 matr = fParCovar;
1312}
1313
1314////////////////////////////////////////////////////////////////////////////////
1315///Returns the internal design matrix
1316
1318{
1319 if (matr.GetNrows()!=fNfunctions || matr.GetNcols()!=fNfunctions){
1321 }
1322 matr = fDesign;
1323}
1324
1325////////////////////////////////////////////////////////////////////////////////
1326///Returns parameter errors
1327
1329{
1330 if (vpar.GetNoElements()!=fNfunctions) {
1331 vpar.ResizeTo(fNfunctions);
1332 }
1333 for (Int_t i=0; i<fNfunctions; i++)
1334 vpar(i) = TMath::Sqrt(fParCovar(i, i));
1335
1336}
1337
1338////////////////////////////////////////////////////////////////////////////////
1339///Returns parameter values
1340
1342{
1343 if (vpar.GetNoElements()!=fNfunctions) {
1344 vpar.ResizeTo(fNfunctions);
1345 }
1346 vpar=fParams;
1347}
1348
1349////////////////////////////////////////////////////////////////////////////////
1350///Returns the value and the name of the parameter `#ipar`
1351///NB: In the calling function the argument name must be set large enough
1352
1353Int_t TLinearFitter::GetParameter(Int_t ipar,char* name,Double_t& value,Double_t& /*verr*/,Double_t& /*vlow*/, Double_t& /*vhigh*/) const
1354{
1355 if (ipar<0 || ipar>fNfunctions) {
1356 Error("GetParError", "illegal value of parameter");
1357 return 0;
1358 }
1359 value = fParams(ipar);
1360 if (fInputFunction)
1361 strcpy(name, fInputFunction->GetParName(ipar));
1362 else
1363 name = nullptr;
1364 return 1;
1365}
1366
1367
1368////////////////////////////////////////////////////////////////////////////////
1369///Returns the error of parameter `#ipar`
1370
1372{
1373 if (ipar<0 || ipar>fNfunctions) {
1374 Error("GetParError", "illegal value of parameter");
1375 return 0;
1376 }
1377
1378 return TMath::Sqrt(fParCovar(ipar, ipar));
1379}
1380
1381
1382////////////////////////////////////////////////////////////////////////////////
1383///Returns name of parameter `#ipar`
1384
1385const char *TLinearFitter::GetParName(Int_t ipar) const
1386{
1387 if (ipar<0 || ipar>fNfunctions) {
1388 Error("GetParError", "illegal value of parameter");
1389 return nullptr;
1390 }
1391 if (fInputFunction)
1392 return fInputFunction->GetParName(ipar);
1393 return "";
1394}
1395
1396////////////////////////////////////////////////////////////////////////////////
1397///Returns the t-value for parameter `#ipar`
1398
1400{
1401 if (ipar<0 || ipar>fNfunctions) {
1402 Error("GetParTValue", "illegal value of parameter");
1403 return 0;
1404 }
1405 if (!fTValues.NonZeros())
1407 return fTValues(ipar);
1408}
1409
1410////////////////////////////////////////////////////////////////////////////////
1411///Returns the significance of parameter `#ipar`
1412
1414{
1415 if (ipar<0 || ipar>fNfunctions) {
1416 Error("GetParSignificance", "illegal value of parameter");
1417 return 0;
1418 }
1419 if (!fParSign.NonZeros())
1421 return fParSign(ipar);
1422}
1423
1424////////////////////////////////////////////////////////////////////////////////
1425///For robust lts fitting, returns the sample, on which the best fit was based
1426
1428{
1429 if (!fRobust){
1430 Error("GetFitSample", "there is no fit sample in ordinary least-squares fit");
1431 return;
1432 }
1433 for (Int_t i=0; i<fNpoints; i++)
1435
1436}
1437
1438////////////////////////////////////////////////////////////////////////////////
1439///Merge objects in list
1440
1442{
1443 if (!list) return -1;
1444 TIter next(list);
1445 TLinearFitter *lfit = nullptr;
1446 while ((lfit = (TLinearFitter*)next())) {
1447 if (!lfit->InheritsFrom(TLinearFitter::Class())) {
1448 Error("Add","Attempt to add object of class: %s to a %s",lfit->ClassName(),this->ClassName());
1449 return -1;
1450 }
1451 Add(lfit);
1452 }
1453 return 0;
1454}
1455////////////////////////////////////////////////////////////////////////////////
1456///set the basis functions in case the fitting function is not
1457/// set directly
1458/// The TLinearFitter will manage and delete the functions contained in the list
1459
1461{
1462 fFunctions = *(functions);
1464 int size = fFunctions.GetEntries();
1465
1467 //change the size of design matrix
1476 if (fFixedParams)
1477 delete [] fFixedParams;
1479 fDesign.Zero();
1480 fAtb.Zero();
1481 fDesignTemp.Zero();
1484 fAtbTemp.Zero();
1485 fAtbTemp2.Zero();
1486 fAtbTemp3.Zero();
1487 fY2Temp=0;
1488 fY2=0;
1489 for (int i=0; i<size; i++)
1490 fFixedParams[i]=false;
1491 fIsSet=kFALSE;
1492 fChisquare=0;
1493
1494}
1495
1496
1497////////////////////////////////////////////////////////////////////////////////
1498///set the number of dimensions
1499
1501{
1502 fNdim=ndim;
1503 fY.ResizeTo(ndim+1);
1504 fX.ResizeTo(ndim+1, ndim);
1505 fE.ResizeTo(ndim+1);
1506
1507 fNpoints=0;
1508 fIsSet=kFALSE;
1509}
1510
1511////////////////////////////////////////////////////////////////////////////////
1512///Additive parts should be separated by "++".
1513///Examples (ai are parameters to fit):
1514///1.fitting function: a0*x0 + a1*x1 + a2*x2
1515/// input formula "x[0]++x[1]++x[2]"
1516///2.TMath functions can be used:
1517/// fitting function: a0*TMath::Gaus(x, 0, 1) + a1*y
1518/// input formula: "TMath::Gaus(x, 0, 1)++y"
1519///fills the array of functions
1520
1521void TLinearFitter::SetFormula(const char *formula)
1522{
1523 Int_t size = 0, special = 0;
1524 Int_t i;
1525 //Int_t len = strlen(formula);
1526 if (fInputFunction)
1527 fInputFunction = nullptr;
1528 if (fFormula)
1529 delete [] fFormula;
1530 fFormulaSize = strlen(formula);
1531 fFormula = new char[fFormulaSize+1];
1532 strlcpy(fFormula, formula,fFormulaSize+1);
1533 fSpecial = 0;
1534 //in case of a hyperplane:
1535 char *fstring;
1536 fstring = (char *)strstr(fFormula, "hyp");
1537 if (fstring!=nullptr){
1538 // isHyper = kTRUE;
1539 fstring+=3;
1540 sscanf(fstring, "%d", &size);
1541 //+1 for the constant term
1542 size++;
1543 fSpecial=200+size;
1544 }
1545
1546 if (fSpecial==0) {
1547 //in case it's not a hyperplane
1548 TString sstring(fFormula);
1549 sstring = sstring.ReplaceAll("++", 2, "|", 1);
1550 TString replaceformula;
1551
1552 //count the number of functions
1553
1554 TObjArray *oa = sstring.Tokenize("|");
1555
1556 //change the size of functions array and clear it
1557 if (!fFunctions.IsEmpty())
1558 fFunctions.Clear();
1559
1560 // do not own the functions in this case
1562
1565
1566 //check if the old notation of xi is used somewhere instead of x[i]
1567 char pattern[12];
1568 char replacement[14];
1569 for (i=0; i<fNdim; i++){
1570 snprintf(pattern,sizeof(pattern), "x%d", i);
1571 snprintf(replacement,sizeof(replacement), "x[%d]", i);
1572 sstring = sstring.ReplaceAll(pattern, Int_t(i/10)+2, replacement, Int_t(i/10)+4);
1573 }
1574
1575 //fill the array of functions
1576 oa = sstring.Tokenize("|");
1577 for (i=0; i<fNfunctions; i++) {
1578 replaceformula = ((TObjString *)oa->UncheckedAt(i))->GetString();
1579 // look first if exists in the map
1580 TFormula * f = nullptr;
1581 auto elem = fgFormulaMap.find(replaceformula );
1582 if (elem != fgFormulaMap.end() )
1583 f = elem->second;
1584 else {
1585 // create a new formula and add in the static map
1586 f = new TFormula("f", replaceformula.Data());
1587 {
1589 fgFormulaMap[replaceformula]=f;
1590 }
1591 }
1592 if (!f) {
1593 Error("TLinearFitter", "f_linear not allocated");
1594 return;
1595 }
1596 special=f->GetNumber();
1597 fFunctions.Add(f);
1598 }
1599
1600 if ((fNfunctions==1)&&(special>299)&&(special<310)){
1601 //if fitting a polynomial
1602 size=special-299;
1603 fSpecial=100+size;
1604 } else
1606 oa->Delete();
1607 delete oa;
1608 }
1610 //change the size of design matrix
1619 if (fFixedParams)
1620 delete [] fFixedParams;
1622 fDesign.Zero();
1623 fAtb.Zero();
1624 fDesignTemp.Zero();
1627 fAtbTemp.Zero();
1628 fAtbTemp2.Zero();
1629 fAtbTemp3.Zero();
1630 fY2Temp=0;
1631 fY2=0;
1632 for (i=0; i<size; i++)
1633 fFixedParams[i]=false;
1634 fIsSet=kFALSE;
1635 fChisquare=0;
1636
1637}
1638
1639////////////////////////////////////////////////////////////////////////////////
1640///Set the fitting function.
1641
1643{
1644 Int_t special, size;
1645 fInputFunction=function;
1647 fSpecial = 0;
1648 special=fInputFunction->GetNumber();
1649 if (!fFunctions.IsEmpty())
1651
1652 if ((special>299)&&(special<310)){
1653 //if fitting a polynomial
1654 size=special-299;
1655 fSpecial=100+size;
1656 } else
1658
1660 //change the size of design matrix
1665
1668
1671 //
1672 if (fFixedParams)
1673 delete [] fFixedParams;
1675 fDesign.Zero();
1676 fAtb.Zero();
1677 fDesignTemp.Zero();
1678 fAtbTemp.Zero();
1679
1682
1683 fAtbTemp2.Zero();
1684 fAtbTemp3.Zero();
1685 fY2Temp=0;
1686 fY2=0;
1687 for (Int_t i=0; i<size; i++)
1688 fFixedParams[i]=false;
1689 //check if any parameters are fixed (not for pure TFormula)
1690
1691 if (function->InheritsFrom(TF1::Class())){
1692 Double_t al,bl;
1693 for (Int_t i=0;i<fNfunctions;i++) {
1694 ((TF1*)function)->GetParLimits(i,al,bl);
1695 if (al*bl !=0 && al >= bl) {
1696 FixParameter(i, function->GetParameter(i));
1697 }
1698 }
1699 }
1700
1701 fIsSet=kFALSE;
1702 fChisquare=0;
1703
1704}
1705
1706////////////////////////////////////////////////////////////////////////////////
1707///Update the design matrix after the formula has been changed.
1708
1710{
1711 if (fStoreData) {
1712 for (Int_t i=0; i<fNpoints; i++) {
1713 AddToDesign(TMatrixDRow(fX, i).GetPtr(), fY(i), fE(i));
1714 }
1715 return true;
1716 } else
1717 return false;
1718
1719}
1720
1721////////////////////////////////////////////////////////////////////////////////
1722///To use in TGraph::Fit and TH1::Fit().
1723
1724Int_t TLinearFitter::ExecuteCommand(const char *command, Double_t *args, Int_t /*nargs*/)
1725{
1726 if (!strcmp(command, "FitGraph")){
1727 if (args) return GraphLinearFitter(args[0]);
1728 else return GraphLinearFitter(0);
1729 }
1730 if (!strcmp(command, "FitGraph2D")){
1731 if (args) return Graph2DLinearFitter(args[0]);
1732 else return Graph2DLinearFitter(0);
1733 }
1734 if (!strcmp(command, "FitMultiGraph")){
1735 if (args) return MultiGraphLinearFitter(args[0]);
1736 else return MultiGraphLinearFitter(0);
1737 }
1738 if (!strcmp(command, "FitHist")) return HistLinearFitter();
1739// if (!strcmp(command, "FitMultiGraph")) MultiGraphLinearFitter();
1740
1741 return 0;
1742}
1743
1744////////////////////////////////////////////////////////////////////////////////
1745/// Level = 3 (to be consistent with minuit) prints parameters and parameter
1746/// errors.
1747
1748void TLinearFitter::PrintResults(Int_t level, Double_t /*amin*/) const
1749{
1750 if (level==3){
1751 if (!fRobust){
1752 printf("Fitting results:\nParameters:\nNO.\t\tVALUE\t\tERROR\n");
1753 for (Int_t i=0; i<fNfunctions; i++){
1754 printf("%d\t%e\t%e\n", i, fParams(i), TMath::Sqrt(fParCovar(i, i)));
1755 }
1756 } else {
1757 printf("Fitting results:\nParameters:\nNO.\t\tVALUE\n");
1758 for (Int_t i=0; i<fNfunctions; i++){
1759 printf("%d\t%e\n", i, fParams(i));
1760 }
1761 }
1762 }
1763}
1764
1765////////////////////////////////////////////////////////////////////////////////
1766///Used in TGraph::Fit().
1767
1769{
1771 TGraph *grr=(TGraph*)GetObjectFit();
1772 TF1 *f1=(TF1*)GetUserFunc();
1773 Foption_t fitOption=GetFitOption();
1774
1775 //Int_t np=0;
1776 Double_t *x=grr->GetX();
1777 Double_t *y=grr->GetY();
1778 Double_t e;
1779
1780 Int_t fitResult = 0;
1781 //set the fitting formula
1782 SetDim(1);
1784
1785 if (fitOption.Robust){
1786 fRobust=kTRUE;
1788 }
1789 //put the points into the fitter
1790 Int_t n=grr->GetN();
1791 for (Int_t i=0; i<n; i++){
1792 if (!f1->IsInside(&x[i])) continue;
1793 e=grr->GetErrorY(i);
1794 if (e<0 || fitOption.W1)
1795 e=1;
1796 AddPoint(&x[i], y[i], e);
1797 }
1798
1799 if (fitOption.Robust){
1800 return EvalRobust(h);
1801 }
1802
1803 fitResult = Eval();
1804
1805 //calculate the precise chisquare
1806 if (!fitResult){
1807 if (!fitOption.Nochisq){
1808 Double_t temp, temp2, sumtotal=0;
1809 for (Int_t i=0; i<n; i++){
1810 if (!f1->IsInside(&x[i])) continue;
1811 temp=f1->Eval(x[i]);
1812 temp2=(y[i]-temp)*(y[i]-temp);
1813 e=grr->GetErrorY(i);
1814 if (e<0 || fitOption.W1)
1815 e=1;
1816 temp2/=(e*e);
1817
1818 sumtotal+=temp2;
1819 }
1820 fChisquare=sumtotal;
1822 }
1823 }
1824 return fitResult;
1825}
1826
1827////////////////////////////////////////////////////////////////////////////////
1828///Minimisation function for a TGraph2D
1829
1831{
1833
1835 TF2 *f2=(TF2*)GetUserFunc();
1836
1837 Foption_t fitOption=GetFitOption();
1838 Int_t n = gr->GetN();
1839 Double_t *gx = gr->GetX();
1840 Double_t *gy = gr->GetY();
1841 Double_t *gz = gr->GetZ();
1842 Double_t x[2];
1843 Double_t z, e;
1844 Int_t fitResult=0;
1845 SetDim(2);
1846 SetFormula(f2->GetFormula());
1847
1848 if (fitOption.Robust){
1849 fRobust=kTRUE;
1851 }
1852
1853 for (Int_t bin=0;bin<n;bin++) {
1854 x[0] = gx[bin];
1855 x[1] = gy[bin];
1856 if (!f2->IsInside(x)) {
1857 continue;
1858 }
1859 z = gz[bin];
1860 e=gr->GetErrorZ(bin);
1861 if (e<0 || fitOption.W1)
1862 e=1;
1863 AddPoint(x, z, e);
1864 }
1865
1866 if (fitOption.Robust){
1867 return EvalRobust(h);
1868 }
1869
1870 fitResult = Eval();
1871
1872 if (!fitResult){
1873 if (!fitOption.Nochisq){
1874 //calculate the precise chisquare
1875 Double_t temp, temp2, sumtotal=0;
1876 for (Int_t bin=0; bin<n; bin++){
1877 x[0] = gx[bin];
1878 x[1] = gy[bin];
1879 if (!f2->IsInside(x)) {
1880 continue;
1881 }
1882 z = gz[bin];
1883
1884 temp=f2->Eval(x[0], x[1]);
1885 temp2=(z-temp)*(z-temp);
1886 e=gr->GetErrorZ(bin);
1887 if (e<0 || fitOption.W1)
1888 e=1;
1889 temp2/=(e*e);
1890
1891 sumtotal+=temp2;
1892 }
1893 fChisquare=sumtotal;
1895 }
1896 }
1897 return fitResult;
1898}
1899
1900////////////////////////////////////////////////////////////////////////////////
1901///Minimisation function for a TMultiGraph
1902
1904{
1905 Int_t n, i;
1906 Double_t *gx, *gy;
1907 Double_t e;
1909 TMultiGraph *mg = (TMultiGraph*)grFitter->GetObjectFit();
1910 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1911 Foption_t fitOption = grFitter->GetFitOption();
1912 Int_t fitResult=0;
1913 SetDim(1);
1914
1915 if (fitOption.Robust){
1916 fRobust=kTRUE;
1918 }
1920
1921 TGraph *gr;
1922 TIter next(mg->GetListOfGraphs());
1923 while ((gr = (TGraph*) next())) {
1924 n = gr->GetN();
1925 gx = gr->GetX();
1926 gy = gr->GetY();
1927 for (i=0; i<n; i++){
1928 if (!f1->IsInside(&gx[i])) continue;
1929 e=gr->GetErrorY(i);
1930 if (e<0 || fitOption.W1)
1931 e=1;
1932 AddPoint(&gx[i], gy[i], e);
1933 }
1934 }
1935
1936 if (fitOption.Robust){
1937 return EvalRobust(h);
1938 }
1939
1940 fitResult = Eval();
1941
1942 //calculate the chisquare
1943 if (!fitResult){
1944 if (!fitOption.Nochisq){
1945 Double_t temp, temp2, sumtotal=0;
1946 next.Reset();
1947 while((gr = (TGraph*)next())) {
1948 n = gr->GetN();
1949 gx = gr->GetX();
1950 gy = gr->GetY();
1951 for (i=0; i<n; i++){
1952 if (!f1->IsInside(&gx[i])) continue;
1953 temp=f1->Eval(gx[i]);
1954 temp2=(gy[i]-temp)*(gy[i]-temp);
1955 e=gr->GetErrorY(i);
1956 if (e<0 || fitOption.W1)
1957 e=1;
1958 temp2/=(e*e);
1959
1960 sumtotal+=temp2;
1961 }
1962
1963 }
1964 fChisquare=sumtotal;
1966 }
1967 }
1968 return fitResult;
1969}
1970
1971////////////////////////////////////////////////////////////////////////////////
1972/// Minimization function for H1s using a Chisquare method.
1973
1975{
1977 Double_t cu,eu;
1978 // Double_t dersum[100], grad[100];
1979 Double_t x[3];
1980 Int_t bin,binx,biny,binz;
1981 // Axis_t binlow, binup, binsize;
1982
1983 TH1 *hfit = (TH1*)GetObjectFit();
1984 TF1 *f1 = (TF1*)GetUserFunc();
1985 Int_t fitResult;
1986 Foption_t fitOption = GetFitOption();
1987 //SetDim(hfit->GetDimension());
1988 SetDim(f1->GetNdim());
1990
1991 Int_t hxfirst = GetXfirst();
1992 Int_t hxlast = GetXlast();
1993 Int_t hyfirst = GetYfirst();
1994 Int_t hylast = GetYlast();
1995 Int_t hzfirst = GetZfirst();
1996 Int_t hzlast = GetZlast();
1997 TAxis *xaxis = hfit->GetXaxis();
1998 TAxis *yaxis = hfit->GetYaxis();
1999 TAxis *zaxis = hfit->GetZaxis();
2000
2001 for (binz=hzfirst;binz<=hzlast;binz++) {
2002 x[2] = zaxis->GetBinCenter(binz);
2003 for (biny=hyfirst;biny<=hylast;biny++) {
2004 x[1] = yaxis->GetBinCenter(biny);
2005 for (binx=hxfirst;binx<=hxlast;binx++) {
2006 x[0] = xaxis->GetBinCenter(binx);
2007 if (!f1->IsInside(x)) continue;
2008 bin = hfit->GetBin(binx,biny,binz);
2009 cu = hfit->GetBinContent(bin);
2010 if (f1->GetNdim() < hfit->GetDimension()) {
2011 if (hfit->GetDimension() > 2) cu = x[2];
2012 else cu = x[1];
2013 }
2014 if (fitOption.W1) {
2015 if (fitOption.W1==1 && cu == 0) continue;
2016 eu = 1;
2017 } else {
2018 eu = hfit->GetBinError(bin);
2019 if (eu <= 0) continue;
2020 }
2021 AddPoint(x, cu, eu);
2022 }
2023 }
2024 }
2025
2026 fitResult = Eval();
2027
2028 if (!fitResult){
2029 if (!fitOption.Nochisq){
2030 Double_t temp, temp2, sumtotal=0;
2031 for (binz=hzfirst;binz<=hzlast;binz++) {
2032 x[2] = zaxis->GetBinCenter(binz);
2033 for (biny=hyfirst;biny<=hylast;biny++) {
2034 x[1] = yaxis->GetBinCenter(biny);
2035 for (binx=hxfirst;binx<=hxlast;binx++) {
2036 x[0] = xaxis->GetBinCenter(binx);
2037 if (!f1->IsInside(x)) continue;
2038 bin = hfit->GetBin(binx,biny,binz);
2039 cu = hfit->GetBinContent(bin);
2040
2041 if (fitOption.W1) {
2042 if (fitOption.W1==1 && cu == 0) continue;
2043 eu = 1;
2044 } else {
2045 eu = hfit->GetBinError(bin);
2046 if (eu <= 0) continue;
2047 }
2048 temp=f1->EvalPar(x);
2049 temp2=(cu-temp)*(cu-temp);
2050 temp2/=(eu*eu);
2051 sumtotal+=temp2;
2052 }
2053 }
2054 }
2055
2056 fChisquare=sumtotal;
2058 }
2059 }
2060 return fitResult;
2061}
2062
2063////////////////////////////////////////////////////////////////////////////////
2064
2066{
2067 if (R__b.IsReading()) {
2068 Int_t old_matr_size = fNfunctions;
2070 if (old_matr_size < fNfunctions) {
2073
2076
2079 }
2080 } else {
2081 if (fAtb.NonZeros()==0) AddTempMatrices();
2083 }
2084}
2085
2086////////////////////////////////////////////////////////////////////////////////
2087///Finds the parameters of the fitted function in case data contains
2088///outliers.
2089///Parameter h stands for the minimal fraction of good points in the
2090///dataset (h < 1, i.e. for 70% of good points take h=0.7).
2091///The default value of h*Npoints is (Npoints + Nparameters+1)/2
2092///If the user provides a value of h smaller than above, default is taken
2093///See class description for the algorithm details
2094
2096{
2097 fRobust = kTRUE;
2098 Double_t kEps = 1e-13;
2099 Int_t nmini = 300;
2100 Int_t i, j, maxind=0, k, k1 = 500;
2101 Int_t nbest = 10;
2102 Double_t chi2 = -1;
2103
2104 if (fFunctions.IsEmpty()&&(!fInputFunction)&&(fSpecial<=200)){
2105 Error("TLinearFitter::EvalRobust", "The formula hasn't been set");
2106 return 1;
2107 }
2108
2109
2110 Double_t *bestchi2 = new Double_t[nbest];
2111 for (i=0; i<nbest; i++)
2112 bestchi2[i]=1e30;
2113
2114 Int_t hdef=Int_t((fNpoints+fNfunctions+1)/2);
2115
2116 if (h>0.000001 && h<1 && fNpoints*h > hdef)
2117 fH = Int_t(fNpoints*h);
2118 else {
2119 // print message only when h is not zero
2120 if (h>0) Warning("Fitting:", "illegal value of H, default is taken, h = %3.2f",double(hdef)/fNpoints);
2121 fH=hdef;
2122 }
2123
2127
2128 Int_t *index = new Int_t[fNpoints];
2129 Double_t *residuals = new Double_t[fNpoints];
2130
2131 if (fNpoints < 2*nmini) {
2132 //when number of cases is small
2133
2134 //to store the best coefficients (columnwise)
2135 TMatrixD cstock(fNfunctions, nbest);
2136 for (k = 0; k < k1; k++) {
2138 chi2 = CStep(1, fH, residuals,index, index, -1, -1);
2139 chi2 = CStep(2, fH, residuals,index, index, -1, -1);
2140 maxind = TMath::LocMax(nbest, bestchi2);
2141 if (chi2 < bestchi2[maxind]) {
2142 bestchi2[maxind] = chi2;
2143 for (i=0; i<fNfunctions; i++)
2144 cstock(i, maxind) = fParams(i);
2145 }
2146 }
2147
2148 //for the nbest best results, perform CSteps until convergence
2149 Int_t *bestindex = new Int_t[fH];
2150 Double_t currentbest;
2151 for (i=0; i<nbest; i++) {
2152 for (j=0; j<fNfunctions; j++)
2153 fParams(j) = cstock(j, i);
2154 chi2 = 1;
2155 while (chi2 > kEps) {
2156 chi2 = CStep(2, fH, residuals,index, index, -1, -1);
2157 if (TMath::Abs(chi2 - bestchi2[i]) < kEps)
2158 break;
2159 else
2160 bestchi2[i] = chi2;
2161 }
2162 currentbest = TMath::MinElement(nbest, bestchi2);
2163 if (chi2 <= currentbest + kEps) {
2164 for (j=0; j<fH; j++){
2165 bestindex[j]=index[j];
2166 }
2167 maxind = i;
2168 }
2169 for (j=0; j<fNfunctions; j++)
2170 cstock(j, i) = fParams(j);
2171 }
2172 //report the result with the lowest chisquare
2173 for (j=0; j<fNfunctions; j++)
2174 fParams(j) = cstock(j, maxind);
2176 for (j=0; j<fH; j++){
2177 //printf("bestindex[%d]=%d\n", j, bestindex[j]);
2178 fFitsample.SetBitNumber(bestindex[j]);
2179 }
2181 ((TF1*)fInputFunction)->SetChisquare(bestchi2[maxind]);
2182 ((TF1*)fInputFunction)->SetNumberFitPoints(fH);
2183 ((TF1*)fInputFunction)->SetNDF(fH-fNfunctions);
2184 }
2185 delete [] index;
2186 delete [] bestindex;
2187 delete [] residuals;
2188 delete [] bestchi2;
2189 return 0;
2190 }
2191 //if n is large, the dataset should be partitioned
2192 Int_t indsubdat[5];
2193 for (i=0; i<5; i++)
2194 indsubdat[i] = 0;
2195
2196 Int_t nsub = Partition(nmini, indsubdat);
2197 Int_t hsub;
2198
2199 Int_t sum = TMath::Min(nmini*5, fNpoints);
2200
2201 Int_t *subdat = new Int_t[sum]; //to store the indices of selected cases
2202 RDraw(subdat, indsubdat);
2203
2204 TMatrixD cstockbig(fNfunctions, nbest*5);
2205 Int_t *beststock = new Int_t[nbest];
2206 Int_t i_start = 0;
2207 Int_t i_end = indsubdat[0];
2208 Int_t k2 = Int_t(k1/nsub);
2209 for (Int_t kgroup = 0; kgroup < nsub; kgroup++) {
2210
2211 hsub = Int_t(fH * indsubdat[kgroup]/fNpoints);
2212 for (i=0; i<nbest; i++)
2213 bestchi2[i] = 1e16;
2214 for (k=0; k<k2; k++) {
2215 CreateSubset(indsubdat[kgroup], hsub, index);
2216 chi2 = CStep(1, hsub, residuals, index, subdat, i_start, i_end);
2217 chi2 = CStep(2, hsub, residuals, index, subdat, i_start, i_end);
2218 maxind = TMath::LocMax(nbest, bestchi2);
2219 if (chi2 < bestchi2[maxind]){
2220 for (i=0; i<fNfunctions; i++)
2221 cstockbig(i, nbest*kgroup + maxind) = fParams(i);
2222 bestchi2[maxind] = chi2;
2223 }
2224 }
2225 if (kgroup != nsub - 1){
2226 i_start += indsubdat[kgroup];
2227 i_end += indsubdat[kgroup+1];
2228 }
2229 }
2230
2231 for (i=0; i<nbest; i++)
2232 bestchi2[i] = 1e30;
2233 //on the pooled subset
2234 Int_t hsub2 = Int_t(fH*sum/fNpoints);
2235 for (k=0; k<nbest*5; k++) {
2236 for (i=0; i<fNfunctions; i++)
2237 fParams(i)=cstockbig(i, k);
2238 chi2 = CStep(1, hsub2, residuals, index, subdat, 0, sum);
2239 chi2 = CStep(2, hsub2, residuals, index, subdat, 0, sum);
2240 maxind = TMath::LocMax(nbest, bestchi2);
2241 if (chi2 < bestchi2[maxind]){
2242 beststock[maxind] = k;
2243 bestchi2[maxind] = chi2;
2244 }
2245 }
2246
2247 //now the array beststock keeps indices of 10 best candidates in cstockbig matrix
2248 for (k=0; k<nbest; k++) {
2249 for (i=0; i<fNfunctions; i++)
2250 fParams(i) = cstockbig(i, beststock[k]);
2251 chi2 = CStep(1, fH, residuals, index, index, -1, -1);
2252 chi2 = CStep(2, fH, residuals, index, index, -1, -1);
2253 bestchi2[k] = chi2;
2254 }
2255
2256 maxind = TMath::LocMin(nbest, bestchi2);
2257 for (i=0; i<fNfunctions; i++)
2258 fParams(i)=cstockbig(i, beststock[maxind]);
2259
2260 chi2 = 1;
2261 while (chi2 > kEps) {
2262 chi2 = CStep(2, fH, residuals, index, index, -1, -1);
2263 if (TMath::Abs(chi2 - bestchi2[maxind]) < kEps)
2264 break;
2265 else
2266 bestchi2[maxind] = chi2;
2267 }
2268
2270 for (j=0; j<fH; j++)
2272 if (fInputFunction){
2273 ((TF1*)fInputFunction)->SetChisquare(bestchi2[maxind]);
2274 ((TF1*)fInputFunction)->SetNumberFitPoints(fH);
2275 ((TF1*)fInputFunction)->SetNDF(fH-fNfunctions);
2276 }
2277
2278 delete [] subdat;
2279 delete [] beststock;
2280 delete [] bestchi2;
2281 delete [] residuals;
2282 delete [] index;
2283
2284 return 0;
2285}
2286
2287////////////////////////////////////////////////////////////////////////////////
2288///Creates a p-subset to start
2289///ntotal - total number of points from which the subset is chosen
2290
2292{
2293 Int_t i, j;
2294 Bool_t repeat=kFALSE;
2295 Int_t nindex=0;
2296 Int_t num;
2297 for(i=0; i<ntotal; i++)
2298 index[i] = ntotal+1;
2299
2300 TRandom2 r;
2301 //create a p-subset
2302 for (i=0; i<fNfunctions; i++) {
2303 num=Int_t(r.Uniform(0, 1)*(ntotal-1));
2304 if (i>0){
2305 for(j=0; j<=i-1; j++) {
2306 if(index[j]==num)
2307 repeat = kTRUE;
2308 }
2309 }
2310 if(repeat==kTRUE) {
2311 i--;
2312 repeat = kFALSE;
2313 } else {
2314 index[i] = num;
2315 nindex++;
2316 }
2317 }
2318
2319 //compute the coefficients of a hyperplane through the p-subset
2320 fDesign.Zero();
2321 fAtb.Zero();
2322 for (i=0; i<fNfunctions; i++){
2323 AddToDesign(TMatrixDRow(fX, index[i]).GetPtr(), fY(index[i]), fE(index[i]));
2324 }
2325 Bool_t ok;
2326
2327 ok = Linf();
2328
2329 //if the chosen points don't define a hyperplane, add more
2330 while (!ok && (nindex < h)) {
2331 repeat=kFALSE;
2332 do{
2333 num=Int_t(r.Uniform(0,1)*(ntotal-1));
2334 repeat=kFALSE;
2335 for(i=0; i<nindex; i++) {
2336 if(index[i]==num) {
2337 repeat=kTRUE;
2338 break;
2339 }
2340 }
2341 } while(repeat==kTRUE);
2342
2343 index[nindex] = num;
2344 nindex++;
2345 //check if the system is of full rank now
2346 AddToDesign(TMatrixDRow(fX, index[nindex-1]).GetPtr(), fY(index[nindex-1]), fE(index[nindex-1]));
2347 ok = Linf();
2348 }
2349}
2350
2351////////////////////////////////////////////////////////////////////////////////
2352///The CStep procedure, as described in the article
2353
2354Double_t TLinearFitter::CStep(Int_t step, Int_t h, Double_t *residuals, Int_t *index, Int_t *subdat, Int_t start, Int_t end)
2355{
2357
2358 Int_t i, j, itemp, n;
2359 Double_t func;
2360 Double_t val[100];
2361 Int_t npar;
2362 if (start > -1) {
2363 n = end - start;
2364 for (i=0; i<n; i++) {
2365 func = 0;
2366 itemp = subdat[start+i];
2367 if (fInputFunction){
2369 func=fInputFunction->EvalPar(TMatrixDRow(fX, itemp).GetPtr());
2370 } else {
2371 func=0;
2372 if ((fSpecial>100)&&(fSpecial<200)){
2373 npar = fSpecial-100;
2374 val[0] = 1;
2375 for (j=1; j<npar; j++)
2376 val[j] = val[j-1]*fX(itemp, 0);
2377 for (j=0; j<npar; j++)
2378 func += fParams(j)*val[j];
2379 } else if (fSpecial>200) {
2380 //hyperplane case
2381 npar = fSpecial-201;
2382 func+=fParams(0);
2383 for (j=0; j<npar; j++)
2384 func += fParams(j+1)*fX(itemp, j);
2385 } else {
2386 // general case
2387 for (j=0; j<fNfunctions; j++) {
2388 TF1 *f1 = (TF1*)(fFunctions.UncheckedAt(j));
2389 val[j] = f1->EvalPar(TMatrixDRow(fX, itemp).GetPtr());
2390 func += fParams(j)*val[j];
2391 }
2392 }
2393 }
2394 residuals[i] = (fY(itemp) - func)*(fY(itemp) - func)/(fE(i)*fE(i));
2395 }
2396 } else {
2397 n=fNpoints;
2398 for (i=0; i<fNpoints; i++) {
2399 func = 0;
2400 if (fInputFunction){
2402 func=fInputFunction->EvalPar(TMatrixDRow(fX, i).GetPtr());
2403 } else {
2404 func=0;
2405 if ((fSpecial>100)&&(fSpecial<200)){
2406 npar = fSpecial-100;
2407 val[0] = 1;
2408 for (j=1; j<npar; j++)
2409 val[j] = val[j-1]*fX(i, 0);
2410 for (j=0; j<npar; j++)
2411 func += fParams(j)*val[j];
2412 } else if (fSpecial>200) {
2413 //hyperplane case
2414 npar = fSpecial-201;
2415 func+=fParams(0);
2416 for (j=0; j<npar; j++)
2417 func += fParams(j+1)*fX(i, j);
2418 } else {
2419 for (j=0; j<fNfunctions; j++) {
2420 TF1 *f1 = (TF1*)(fFunctions.UncheckedAt(j));
2421 val[j] = f1->EvalPar(TMatrixDRow(fX, i).GetPtr());
2422 func += fParams(j)*val[j];
2423 }
2424 }
2425 }
2426 residuals[i] = (fY(i) - func)*(fY(i) - func)/(fE(i)*fE(i));
2427 }
2428 }
2429 //take h with smallest residuals
2430 TMath::KOrdStat(n, residuals, h-1, index);
2431 //add them to the design matrix
2432 fDesign.Zero();
2433 fAtb.Zero();
2434 for (i=0; i<h; i++)
2435 AddToDesign(TMatrixDRow(fX, index[i]).GetPtr(), fY(index[i]), fE(index[i]));
2436
2437 Linf();
2438
2439 //don't calculate the chisquare at the 1st cstep
2440 if (step==1) return 0;
2441 Double_t sum=0;
2442
2443
2444 if (start > -1) {
2445 for (i=0; i<h; i++) {
2446 itemp = subdat[start+index[i]];
2447 if (fInputFunction){
2449 func=fInputFunction->EvalPar(TMatrixDRow(fX, itemp).GetPtr());
2450 } else {
2451 func=0;
2452 if ((fSpecial>100)&&(fSpecial<200)){
2453 npar = fSpecial-100;
2454 val[0] = 1;
2455 for (j=1; j<npar; j++)
2456 val[j] = val[j-1]*fX(itemp, 0);
2457 for (j=0; j<npar; j++)
2458 func += fParams(j)*val[j];
2459 } else if (fSpecial>200) {
2460 //hyperplane case
2461 npar = fSpecial-201;
2462 func+=fParams(0);
2463 for (j=0; j<npar; j++)
2464 func += fParams(j+1)*fX(itemp, j);
2465 } else {
2466 for (j=0; j<fNfunctions; j++) {
2467 TF1 *f1 = (TF1*)(fFunctions.UncheckedAt(j));
2468 val[j] = f1->EvalPar(TMatrixDRow(fX, itemp).GetPtr());
2469 func += fParams(j)*val[j];
2470 }
2471 }
2472 }
2473 sum+=(fY(itemp)-func)*(fY(itemp)-func)/(fE(itemp)*fE(itemp));
2474 }
2475 } else {
2476 for (i=0; i<h; i++) {
2477 if (fInputFunction){
2479 func=fInputFunction->EvalPar(TMatrixDRow(fX, index[i]).GetPtr());
2480 } else {
2481 func=0;
2482 if ((fSpecial>100)&&(fSpecial<200)){
2483 npar = fSpecial-100;
2484 val[0] = 1;
2485 for (j=1; j<npar; j++)
2486 val[j] = val[j-1]*fX(index[i], 0);
2487 for (j=0; j<npar; j++)
2488 func += fParams(j)*val[j];
2489 } else {
2490 if (fSpecial>200) {
2491 //hyperplane case
2492 npar = fSpecial-201;
2493 func+=fParams(0);
2494 for (j=0; j<npar; j++)
2495 func += fParams(j+1)*fX(index[i], j);
2496 } else {
2497 for (j=0; j<fNfunctions; j++) {
2498 TF1 *f1 = (TF1*)(fFunctions.UncheckedAt(j));
2499 val[j] = f1->EvalPar(TMatrixDRow(fX, index[i]).GetPtr());
2500 func += fParams(j)*val[j];
2501 }
2502 }
2503 }
2504 }
2505
2506 sum+=(fY(index[i])-func)*(fY(index[i])-func)/(fE(index[i])*fE(index[i]));
2507 }
2508 }
2509
2510 return sum;
2511}
2512
2513////////////////////////////////////////////////////////////////////////////////
2514
2516{
2517 //currently without the intercept term
2523 fDesignTemp.Zero();
2526 fAtb+=fAtbTemp;
2527 fAtbTemp3.Zero();
2528 fAtbTemp2.Zero();
2529 fAtbTemp.Zero();
2530
2531 fY2+=fY2Temp;
2532 fY2Temp=0;
2533
2534
2535 TDecompChol chol(fDesign);
2536 TVectorD temp(fNfunctions);
2537 Bool_t ok;
2538 temp = chol.Solve(fAtb, ok);
2539 if (!ok){
2540 Error("Linf","Matrix inversion failed");
2541 //fDesign.Print();
2542 fParams.Zero();
2543 return kFALSE;
2544 }
2545 fParams = temp;
2546 return ok;
2547}
2548
2549////////////////////////////////////////////////////////////////////////////////
2550///divides the elements into approximately equal subgroups
2551///number of elements in each subgroup is stored in indsubdat
2552///number of subgroups is returned
2553
2555{
2556 Int_t nsub;
2557
2558 if ((fNpoints>=2*nmini) && (fNpoints<=(3*nmini-1))) {
2559 if (fNpoints%2==1){
2560 indsubdat[0]=Int_t(fNpoints*0.5);
2561 indsubdat[1]=Int_t(fNpoints*0.5)+1;
2562 } else
2563 indsubdat[0]=indsubdat[1]=Int_t(fNpoints/2);
2564 nsub=2;
2565 }
2566 else{
2567 if((fNpoints>=3*nmini) && (fNpoints<(4*nmini -1))) {
2568 if(fNpoints%3==0){
2569 indsubdat[0]=indsubdat[1]=indsubdat[2]=Int_t(fNpoints/3);
2570 } else {
2571 indsubdat[0]=Int_t(fNpoints/3);
2572 indsubdat[1]=Int_t(fNpoints/3)+1;
2573 if (fNpoints%3==1) indsubdat[2]=Int_t(fNpoints/3);
2574 else indsubdat[2]=Int_t(fNpoints/3)+1;
2575 }
2576 nsub=3;
2577 }
2578 else{
2579 if((fNpoints>=4*nmini)&&(fNpoints<=(5*nmini-1))){
2580 if (fNpoints%4==0) indsubdat[0]=indsubdat[1]=indsubdat[2]=indsubdat[3]=Int_t(fNpoints/4);
2581 else {
2582 indsubdat[0]=Int_t(fNpoints/4);
2583 indsubdat[1]=Int_t(fNpoints/4)+1;
2584 if(fNpoints%4==1) indsubdat[2]=indsubdat[3]=Int_t(fNpoints/4);
2585 if(fNpoints%4==2) {
2586 indsubdat[2]=Int_t(fNpoints/4)+1;
2587 indsubdat[3]=Int_t(fNpoints/4);
2588 }
2589 if(fNpoints%4==3) indsubdat[2]=indsubdat[3]=Int_t(fNpoints/4)+1;
2590 }
2591 nsub=4;
2592 } else {
2593 for(Int_t i=0; i<5; i++)
2594 indsubdat[i]=nmini;
2595 nsub=5;
2596 }
2597 }
2598 }
2599 return nsub;
2600}
2601
2602////////////////////////////////////////////////////////////////////////////////
2603///Draws ngroup nonoverlapping subdatasets out of a dataset of size n
2604///such that the selected case numbers are uniformly distributed from 1 to n
2605
2606void TLinearFitter::RDraw(Int_t *subdat, Int_t *indsubdat)
2607{
2608 Int_t jndex = 0;
2609 Int_t nrand;
2610 Int_t i, k, m, j;
2611 Int_t ngroup=0;
2612 for (i=0; i<5; i++) {
2613 if (indsubdat[i]!=0)
2614 ngroup++;
2615 }
2616 TRandom r;
2617 for (k=1; k<=ngroup; k++) {
2618 for (m=1; m<=indsubdat[k-1]; m++) {
2619 nrand = Int_t(r.Uniform(0, 1) * (fNpoints-jndex)) + 1;
2620 jndex++;
2621 if (jndex==1) {
2622 subdat[0] = nrand;
2623 } else {
2624 subdat[jndex-1] = nrand + jndex - 2;
2625 for (i=1; i<=jndex-1; i++) {
2626 if(subdat[i-1] > nrand+i-2) {
2627 for(j=jndex; j>=i+1; j--) {
2628 subdat[j-1] = subdat[j-2];
2629 }
2630 subdat[i-1] = nrand+i-2;
2631 break; //breaking the loop for(i=1...
2632 }
2633 }
2634 }
2635 }
2636 }
2637}
2638
2639
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
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 np
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 r
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 value
Option_t Option_t TPoint xy
char name[80]
Definition TGX11.cxx:110
TMatrixTRow< Double_t > TMatrixDRow
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define R__LOCKGUARD(mutex)
#define snprintf
Definition civetweb.c:1540
Class to manage histogram axis.
Definition TAxis.h:31
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:469
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:458
Container of bits.
Definition TBits.h:26
void Clear(Option_t *option="") override
Clear the value.
Definition TBits.cxx:84
Bool_t TestBitNumber(UInt_t bitnumber) const
Definition TBits.h:222
void SetBitNumber(UInt_t bitnumber, Bool_t value=kTRUE)
Definition TBits.h:206
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Collection abstract base class.
Definition TCollection.h:65
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TObject * Clone(const char *newname="") const override
Make a clone of an collection using the Streamer facility.
Cholesky Decomposition class.
Definition TDecompChol.h:25
Bool_t Solve(TVectorD &b) override
Solve equations Ax=b assuming A has been factored by Cholesky.
Bool_t Invert(TMatrixDSym &inv)
For a symmetric matrix A(m,m), its inverse A_inv(m,m) is returned .
1-Dim function class
Definition TF1.h:233
static TClass * Class()
virtual void SetChisquare(Double_t chi2)
Definition TF1.h:640
virtual TFormula * GetFormula()
Definition TF1.h:481
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1468
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition TF1.cxx:1439
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition TF1.h:626
virtual Int_t GetNdim() const
Definition TF1.h:513
A 2-Dim function with parameters.
Definition TF2.h:29
Bool_t IsInside(const Double_t *x) const override
Return kTRUE is the point is inside the function range.
Definition TF2.cxx:685
The Formula class.
Definition TFormula.h:89
const TObject * GetLinearPart(Int_t i) const
Return linear part.
Double_t Eval(Args... args) const
Set first 1, 2, 3 or 4 variables (e.g.
Definition TFormula.h:324
static TClass * Class()
const char * GetParName(Int_t ipar) const
Return parameter name given by integer.
Int_t GetNumber() const
Definition TFormula.h:261
void SetParameters(const Double_t *params)
Set a vector of parameters value.
Int_t GetNpar() const
Definition TFormula.h:260
Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr) const
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition TGraph2D.h:41
static TClass * Class()
Double_t * GetY() const
Definition TGraph2D.h:124
Double_t * GetX() const
Definition TGraph2D.h:123
virtual Double_t * GetEZ() const
Definition TGraph2D.h:128
Int_t GetN() const
Definition TGraph2D.h:122
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Double_t GetErrorY(Int_t bin) const override
It returns the error along Y at point i.
Double_t * GetEY() const override
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
static TClass * Class()
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2342
Double_t * GetY() const
Definition TGraph.h:140
Int_t GetN() const
Definition TGraph.h:132
Double_t * GetX() const
Definition TGraph.h:139
virtual Double_t GetErrorY(Int_t bin) const
It always returns a negative value. Real implementation in TGraphErrors.
Definition TGraph.cxx:1369
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
TAxis * GetZaxis()
Definition TH1.h:326
static TClass * Class()
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:9054
virtual Int_t GetDimension() const
Definition TH1.h:283
TAxis * GetXaxis()
Definition TH1.h:324
virtual Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition TH1.cxx:4952
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
Definition TH1.cxx:9197
TAxis * GetYaxis()
Definition TH1.h:325
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9213
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5052
void Reset()
virtual void AddTempMatrices()
TMatrixDSym fDesignTemp2
temporary matrix, used for num.stability
Int_t GraphLinearFitter(Double_t h)
Used in TGraph::Fit().
~TLinearFitter() override
Linear fitter cleanup.
TMatrixDSym fDesignTemp
virtual Double_t GetChisquare()
Get the Chisquare.
Int_t Partition(Int_t nmini, Int_t *indsubdat)
divides the elements into approximately equal subgroups number of elements in each subgroup is stored...
virtual void GetErrors(TVectorD &vpar)
Returns parameter errors.
Double_t CStep(Int_t step, Int_t h, Double_t *residuals, Int_t *index, Int_t *subdat, Int_t start, Int_t end)
The CStep procedure, as described in the article.
virtual Int_t Merge(TCollection *list)
Merge objects in list.
Double_t fChisquare
void ComputeTValues()
Computes parameters' t-values and significance.
Bool_t * fFixedParams
TVectorD fParSign
void Clear(Option_t *option="") override
Clears everything. Used in TH1::Fit and TGraph::Fit().
TLinearFitter()
default c-tor, input data is stored If you don't want to store the input data, run the function Store...
const char * GetParName(Int_t ipar) const override
Returns name of parameter #ipar
Int_t MultiGraphLinearFitter(Double_t h)
Minimisation function for a TMultiGraph.
virtual Double_t GetParSignificance(Int_t ipar)
Returns the significance of parameter #ipar
TMatrixDSym fDesign
virtual Int_t Eval()
Perform the fit and evaluate the parameters Returns 0 if the fit is ok, 1 if there are errors.
TVectorD fAtbTemp2
temporary vector, used for num.stability
void PrintResults(Int_t level, Double_t amin=0) const override
Level = 3 (to be consistent with minuit) prints parameters and parameter errors.
Int_t HistLinearFitter()
Minimization function for H1s using a Chisquare method.
Int_t Graph2DLinearFitter(Double_t h)
Minimisation function for a TGraph2D.
static TClass * Class()
TVectorD fAtbTemp
virtual void ClearPoints()
To be used when different sets of points are fitted with the same formula.
TMatrixDSym fDesignTemp3
TObjArray fFunctions
map of basis functions and formula
Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs) override
To use in TGraph::Fit and TH1::Fit().
virtual void GetFitSample(TBits &bits)
For robust lts fitting, returns the sample, on which the best fit was based.
virtual void Add(TLinearFitter *tlf)
Add another linear fitter to this linear fitter.
virtual void GetDesignMatrix(TMatrixD &matr)
Returns the internal design matrix.
void Streamer(TBuffer &) override
Stream an object of class TObject.
virtual void AssignData(Int_t npoints, Int_t xncols, Double_t *x, Double_t *y, Double_t *e=nullptr)
This function is to use when you already have all the data in arrays and don't want to copy them into...
virtual void GetParameters(TVectorD &vpar)
Returns parameter values.
void RDraw(Int_t *subdat, Int_t *indsubdat)
Draws ngroup nonoverlapping subdatasets out of a dataset of size n such that the selected case number...
static std::map< TString, TFormula * > fgFormulaMap
virtual void SetDim(Int_t n)
set the number of dimensions
TFormula * fInputFunction
void GetConfidenceIntervals(Int_t n, Int_t ndim, const Double_t *x, Double_t *ci, Double_t cl=0.95) override
Computes point-by-point confidence intervals for the fitted function Parameters: n - number of points...
TMatrixD fX
temporary variable used for num.stability
virtual Bool_t UpdateMatrix()
Update the design matrix after the formula has been changed.
virtual void GetAtbVector(TVectorD &v)
Get the Atb vector - a vector, used for internal computations.
virtual void Chisquare()
Calculates the chisquare.
virtual void SetBasisFunctions(TObjArray *functions)
set the basis functions in case the fitting function is not set directly The TLinearFitter will manag...
Int_t fNpoints
temporary
Double_t fVal[1000]
virtual Int_t EvalRobust(Double_t h=-1)
Finds the parameters of the fitted function in case data contains outliers.
TVectorD fTValues
Double_t GetParError(Int_t ipar) const override
Returns the error of parameter #ipar
void AddToDesign(Double_t *x, Double_t y, Double_t e)
Add a point to the AtA matrix and to the Atb vector.
TVectorD fAtbTemp3
void FixParameter(Int_t ipar) override
Fixes paramter #ipar at its current value.
TLinearFitter & operator=(const TLinearFitter &tlf)
Assignment operator.
void CreateSubset(Int_t ntotal, Int_t h, Int_t *index)
Creates a p-subset to start ntotal - total number of points from which the subset is chosen.
virtual void SetFormula(const char *formula)
Additive parts should be separated by "++".
void ReleaseParameter(Int_t ipar) override
Releases parameter #ipar.
virtual void AddPoint(Double_t *x, Double_t y, Double_t e=1)
Adds 1 point to the fitter.
TMatrixDSym fParCovar
virtual Double_t GetParTValue(Int_t ipar)
Returns the t-value for parameter #ipar
Double_t * GetCovarianceMatrix() const override
Returns covariance matrix.
Double_t GetParameter(Int_t ipar) const override
virtual void StoreData(Bool_t store)
Int_t GetNrows() const
virtual TMatrixTBase< Element > & Zero()
Set matrix elements to zero.
Int_t GetNcols() const
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
void Clear(Option_t *="") override
Definition TMatrixTSym.h:96
const Element * GetMatrixArray() const override
TMatrixT< Element > & Use(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Element *data)
Use the array data to fill the matrix ([row_lwb..row_upb] x [col_lwb..col_upb])
void Clear(Option_t *="") override
Definition TMatrixT.h:123
const Element * GetMatrixArray() const override
Definition TMatrixT.h:228
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
TList * GetListOfGraphs() const
Definition TMultiGraph.h:68
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
virtual void Expand(Int_t newSize)
Expand or shrink the array to newSize elements.
void Clear(Option_t *option="") override
Remove all objects from the array.
Int_t GetEntries() const override
Return the number of objects in array (i.e.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
Bool_t IsEmpty() const override
Definition TObjArray.h:65
void Add(TObject *obj) override
Definition TObjArray.h:68
Collectable string class.
Definition TObjString.h:28
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:438
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:524
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
virtual TClass * IsA() const
Definition TObject.h:243
Random number generator class based on the maximally quidistributed combined Tausworthe generator by ...
Definition TRandom2.h:27
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
void ToUpper()
Change string to upper case.
Definition TString.cxx:1195
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition TString.cxx:2264
TVectorT< Element > & Zero()
Set vector elements to zero.
Definition TVectorT.cxx:453
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition TVectorT.cxx:294
TVectorT< Element > & Use(Int_t lwb, Int_t upb, Element *data)
Use the array data to fill the vector lwb..upb].
Definition TVectorT.cxx:349
void Clear(Option_t *="") override
Definition TVectorT.h:172
Int_t NonZeros() const
Compute the number of elements != 0.0.
Definition TVectorT.cxx:620
Int_t GetNoElements() const
Definition TVectorT.h:74
Element * GetMatrixArray()
Definition TVectorT.h:76
Abstract Base Class for Fitting.
virtual Int_t GetXlast() const
virtual Int_t GetYfirst() const
virtual TObject * GetObjectFit() const
virtual Foption_t GetFitOption() const
virtual Int_t GetZfirst() const
virtual Int_t GetZlast() const
virtual Int_t GetXfirst() const
TObject * fObjectFit
Pointer to object being fitted.
virtual Int_t GetYlast() const
TVirtualFitter & operator=(const TVirtualFitter &tvf)
assignment operator
static TVirtualFitter * GetFitter()
static: return the current Fitter
virtual TObject * GetUserFunc() const
Double_t y[n]
Definition legend1.C:17
Double_t x[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
Element KOrdStat(Size n, const Element *a, Size k, Size *work=0)
Returns k_th order statistic of the array a of size n (k_th smallest element out of n elements).
Definition TMath.h:1359
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:982
T MinElement(Long64_t n, const T *a)
Returns minimum of array a of length n.
Definition TMath.h:960
Long64_t LocMax(Long64_t n, const T *a)
Returns index of array with the maximum element.
Definition TMath.h:1012
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t StudentI(Double_t T, Double_t ndf)
Calculates the cumulative distribution function of Student's t-distribution second parameter stands f...
Definition TMath.cxx:2646
Double_t StudentQuantile(Double_t p, Double_t ndf, Bool_t lower_tail=kTRUE)
Computes quantiles of the Student's t-distribution 1st argument is the probability,...
Definition TMath.cxx:2674
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
int W1
Definition Foption.h:36
int Nochisq
Definition Foption.h:45
int Robust
Definition Foption.h:48
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345