1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////
20 // Subset of matheamtical functions not included in the TMath
23 ///////////////////////////////////////////////////////////////////////////
25 #include "Riostream.h"
31 #include "TObjString.h"
32 #include "TLinearFitter.h"
35 // includes neccessary for test functions
39 #include "TStopwatch.h"
40 #include "TTreeStream.h"
42 #include "TStatToolkit.h"
45 ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O
47 TStatToolkit::TStatToolkit() : TObject()
50 // Default constructor
53 ///////////////////////////////////////////////////////////////////////////
54 TStatToolkit::~TStatToolkit()
62 //_____________________________________________________________________________
63 void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
64 , Double_t &sigma, Int_t hh)
67 // Robust estimator in 1D case MI version - (faster than ROOT version)
69 // For the univariate case
70 // estimates of location and scatter are returned in mean and sigma parameters
71 // the algorithm works on the same principle as in multivariate case -
72 // it finds a subset of size hh with smallest sigma, and then returns mean and
73 // sigma of this subset
78 Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
79 Int_t *index=new Int_t[nvectors];
80 TMath::Sort(nvectors, data, index, kFALSE);
82 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
83 Double_t factor = faclts[TMath::Max(0,nquant-1)];
88 Double_t bestmean = 0;
89 Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.); // maximal possible sigma
90 bestsigma *=bestsigma;
92 for (Int_t i=0; i<hh; i++){
93 sumx += data[index[i]];
94 sumx2 += data[index[i]]*data[index[i]];
97 Double_t norm = 1./Double_t(hh);
98 Double_t norm2 = 1./Double_t(hh-1);
99 for (Int_t i=hh; i<nvectors; i++){
100 Double_t cmean = sumx*norm;
101 Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
102 if (csigma<bestsigma){
108 sumx += data[index[i]]-data[index[i-hh]];
109 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
112 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
121 void TStatToolkit::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh, Float_t externalfactor)
123 // Modified version of ROOT robust EvaluateUni
124 // robust estimator in 1D case MI version
125 // added external factor to include precision of external measurement
130 Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
131 Int_t *index=new Int_t[nvectors];
132 TMath::Sort(nvectors, data, index, kFALSE);
134 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
135 Double_t factor = faclts[0];
137 // fix proper normalization - Anja
138 factor = faclts[nquant-1];
145 Int_t bestindex = -1;
146 Double_t bestmean = 0;
147 Double_t bestsigma = -1;
148 for (Int_t i=0; i<hh; i++){
149 sumx += data[index[i]];
150 sumx2 += data[index[i]]*data[index[i]];
153 Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
154 Double_t norm = 1./Double_t(hh);
155 for (Int_t i=hh; i<nvectors; i++){
156 Double_t cmean = sumx*norm;
157 Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
158 if (csigma<bestsigma || bestsigma<0){
165 sumx += data[index[i]]-data[index[i-hh]];
166 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
169 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
176 //_____________________________________________________________________________
177 Int_t TStatToolkit::Freq(Int_t n, const Int_t *inlist
178 , Int_t *outlist, Bool_t down)
181 // Sort eleements according occurancy
182 // The size of output array has is 2*n
185 Int_t * sindexS = new Int_t[n]; // temp array for sorting
186 Int_t * sindexF = new Int_t[2*n];
187 for (Int_t i=0;i<n;i++) sindexF[i]=0;
189 TMath::Sort(n,inlist, sindexS, down);
190 Int_t last = inlist[sindexS[0]];
197 for(Int_t i=1;i<n; i++){
198 val = inlist[sindexS[i]];
199 if (last == val) sindexF[countPos]++;
202 sindexF[countPos+n] = val;
207 if (last==val) countPos++;
208 // sort according frequency
209 TMath::Sort(countPos, sindexF, sindexS, kTRUE);
210 for (Int_t i=0;i<countPos;i++){
211 outlist[2*i ] = sindexF[sindexS[i]+n];
212 outlist[2*i+1] = sindexF[sindexS[i]];
221 //___TStatToolkit__________________________________________________________________________
222 void TStatToolkit::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
226 Int_t nbins = his->GetNbinsX();
227 Float_t nentries = his->GetEntries();
232 for (Int_t ibin=1;ibin<nbins; ibin++){
233 ncumul+= his->GetBinContent(ibin);
234 Float_t fraction = Float_t(ncumul)/Float_t(nentries);
235 if (fraction>down && fraction<up){
236 sum+=his->GetBinContent(ibin);
237 mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
238 sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);
242 sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
244 (*param)[0] = his->GetMaximum();
246 (*param)[2] = sigma2;
249 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
252 void TStatToolkit::LTM(TH1F * his, TVectorD *param , Float_t fraction, Bool_t verbose){
256 Int_t nbins = his->GetNbinsX();
257 Int_t nentries = (Int_t)his->GetEntries();
258 Double_t *data = new Double_t[nentries];
260 for (Int_t ibin=1;ibin<nbins; ibin++){
261 Float_t entriesI = his->GetBinContent(ibin);
262 Float_t xcenter= his->GetBinCenter(ibin);
263 for (Int_t ic=0; ic<entriesI; ic++){
264 if (npoints<nentries){
265 data[npoints]= xcenter;
270 Double_t mean, sigma;
271 Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
272 npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
273 TStatToolkit::EvaluateUni(npoints, data, mean,sigma,npoints2);
275 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
276 (*param)[0] = his->GetMaximum();
282 Double_t TStatToolkit::FitGaus(TH1F* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
284 // Fit histogram with gaussian function
287 // return value- chi2 - if negative ( not enough points)
288 // his - input histogram
289 // param - vector with parameters
290 // xmin, xmax - range to fit - if xmin=xmax=0 - the full histogram range used
292 // 1. Step - make logarithm
293 // 2. Linear fit (parabola) - more robust - always converge
294 // 3. In case of small statistic bins are averaged
296 static TLinearFitter fitter(3,"pol2");
300 if (his->GetMaximum()<4) return -1;
301 if (his->GetEntries()<12) return -1;
302 if (his->GetRMS()<mat.GetTol()) return -1;
303 Float_t maxEstimate = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
304 Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
306 if (maxEstimate<1) return -1;
307 Int_t nbins = his->GetNbinsX();
313 xmin = his->GetXaxis()->GetXmin();
314 xmax = his->GetXaxis()->GetXmax();
316 for (Int_t iter=0; iter<2; iter++){
317 fitter.ClearPoints();
319 for (Int_t ibin=1;ibin<nbins+1; ibin++){
321 Float_t entriesI = his->GetBinContent(ibin);
322 for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
323 if (ibin+delta>1 &&ibin+delta<nbins-1){
324 entriesI += his->GetBinContent(ibin+delta);
329 Double_t xcenter= his->GetBinCenter(ibin);
330 if (xcenter<xmin || xcenter>xmax) continue;
331 Double_t error=1./TMath::Sqrt(countB);
334 if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
335 cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
336 if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
338 if (entriesI>1&&cont>1){
339 fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
345 fitter.GetParameters(par);
353 fitter.GetParameters(par);
354 fitter.GetCovarianceMatrix(mat);
355 if (TMath::Abs(par[1])<mat.GetTol()) return -1;
356 if (TMath::Abs(par[2])<mat.GetTol()) return -1;
357 Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
358 //fitter.GetParameters();
359 if (!param) param = new TVectorD(3);
360 // if (!matrix) matrix = new TMatrixD(3,3); // Covariance matrix to be implemented
361 (*param)[1] = par[1]/(-2.*par[2]);
362 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
363 (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1]);
368 printf("Chi2=%f\n",chi2);
369 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
370 f1->SetParameter(0, (*param)[0]);
371 f1->SetParameter(1, (*param)[1]);
372 f1->SetParameter(2, (*param)[2]);
378 Double_t TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
380 // Fit histogram with gaussian function
383 // nbins: size of the array and number of histogram bins
384 // xMin, xMax: histogram range
385 // param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
386 // matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
389 // >0: the chi2 returned by TLinearFitter
390 // -3: only three points have been used for the calculation - no fitter was used
391 // -2: only two points have been used for the calculation - center of gravity was uesed for calculation
392 // -1: only one point has been used for the calculation - center of gravity was uesed for calculation
393 // -4: invalid result!!
396 // 1. Step - make logarithm
397 // 2. Linear fit (parabola) - more robust - always converge
399 static TLinearFitter fitter(3,"pol2");
400 static TMatrixD mat(3,3);
401 static Double_t kTol = mat.GetTol();
402 fitter.StoreData(kFALSE);
403 fitter.ClearPoints();
408 Float_t rms = TMath::RMS(nBins,arr);
409 Float_t max = TMath::MaxElement(nBins,arr);
410 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
419 for (Int_t i=0; i<nBins; i++){
421 if (arr[i]>0) nfilled++;
424 if (max<4) return -4;
425 if (entries<12) return -4;
426 if (rms<kTol) return -4;
432 for (Int_t ibin=0;ibin<nBins; ibin++){
433 Float_t entriesI = arr[ibin];
435 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
437 Float_t error = 1./TMath::Sqrt(entriesI);
438 Float_t val = TMath::Log(Float_t(entriesI));
439 fitter.AddPoint(&xcenter,val,error);
442 A(npoints,1)=xcenter;
443 A(npoints,2)=xcenter*xcenter;
445 meanCOG+=xcenter*entriesI;
446 rms2COG +=xcenter*entriesI*xcenter;
457 //analytic calculation of the parameters for three points
466 // use fitter for more than three points
468 fitter.GetParameters(par);
469 fitter.GetCovarianceMatrix(mat);
470 chi2 = fitter.GetChisquare()/Float_t(npoints);
472 if (TMath::Abs(par[1])<kTol) return -4;
473 if (TMath::Abs(par[2])<kTol) return -4;
475 if (!param) param = new TVectorD(3);
476 //if (!matrix) matrix = new TMatrixD(3,3); // !!!!might be a memory leek. use dummy matrix pointer to call this function! // Covariance matrix to be implemented
478 (*param)[1] = par[1]/(-2.*par[2]);
479 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
480 Double_t lnparam0 = par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1];
481 if ( lnparam0>307 ) return -4;
482 (*param)[0] = TMath::Exp(lnparam0);
487 printf("Chi2=%f\n",chi2);
488 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
489 f1->SetParameter(0, (*param)[0]);
490 f1->SetParameter(1, (*param)[1]);
491 f1->SetParameter(2, (*param)[2]);
498 //use center of gravity for 2 points
502 (*param)[1] = meanCOG;
503 (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
509 (*param)[1] = meanCOG;
510 (*param)[2] = binWidth/TMath::Sqrt(12);
518 Float_t TStatToolkit::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
521 // calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
522 // return COG; in case of failure return xMin
529 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
531 for (Int_t ibin=0; ibin<nBins; ibin++){
532 Float_t entriesI = (Float_t)arr[ibin];
533 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
535 meanCOG += xcenter*entriesI;
536 rms2COG += xcenter*entriesI*xcenter;
541 if ( sumCOG == 0 ) return xMin;
546 (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
547 if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
558 ///////////////////////////////////////////////////////////////
559 ////////////// TEST functions /////////////////////////
560 ///////////////////////////////////////////////////////////////
566 void TStatToolkit::TestGausFit(Int_t nhistos){
568 // Test performance of the parabolic - gaussian fit - compare it with
570 // nhistos - number of histograms to be used for test
572 TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
574 Float_t *xTrue = new Float_t[nhistos];
575 Float_t *sTrue = new Float_t[nhistos];
576 TVectorD **par1 = new TVectorD*[nhistos];
577 TVectorD **par2 = new TVectorD*[nhistos];
581 TH1F **h1f = new TH1F*[nhistos];
582 TF1 *myg = new TF1("myg","gaus");
583 TF1 *fit = new TF1("fit","gaus");
587 for (Int_t i=0;i<nhistos; i++){
588 par1[i] = new TVectorD(3);
589 par2[i] = new TVectorD(3);
590 h1f[i] = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
591 xTrue[i]= gRandom->Rndm();
593 sTrue[i]= .75+gRandom->Rndm()*.5;
594 myg->SetParameters(1,xTrue[i],sTrue[i]);
595 h1f[i]->FillRandom("myg");
601 for (Int_t i=0; i<nhistos; i++){
602 h1f[i]->Fit(fit,"0q");
603 (*par1[i])(0) = fit->GetParameter(0);
604 (*par1[i])(1) = fit->GetParameter(1);
605 (*par1[i])(2) = fit->GetParameter(2);
608 printf("Gaussian fit\t");
612 //TStatToolkit gaus fit
613 for (Int_t i=0; i<nhistos; i++){
614 TStatToolkit::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
618 printf("Parabolic fit\t");
621 for (Int_t i=0;i<nhistos; i++){
622 Float_t xt = xTrue[i];
623 Float_t st = sTrue[i];
632 for (Int_t i=0;i<nhistos; i++){
649 TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
653 // delta - number of bins to integrate
654 // type - 0 - mean value
656 TAxis * xaxis = his->GetXaxis();
657 TAxis * yaxis = his->GetYaxis();
658 // TAxis * zaxis = his->GetZaxis();
659 Int_t nbinx = xaxis->GetNbins();
660 Int_t nbiny = yaxis->GetNbins();
663 TGraph2D *graph = new TGraph2D(nbinx*nbiny);
665 for (Int_t ix=0; ix<nbinx;ix++)
666 for (Int_t iy=0; iy<nbiny;iy++){
667 Float_t xcenter = xaxis->GetBinCenter(ix);
668 Float_t ycenter = yaxis->GetBinCenter(iy);
669 snprintf(name,1000,"%s_%d_%d",his->GetName(), ix,iy);
670 TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
672 if (type==0) stat = projection->GetMean();
673 if (type==1) stat = projection->GetRMS();
674 if (type==2 || type==3){
676 TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
677 if (type==2) stat= vec[1];
678 if (type==3) stat= vec[0];
680 if (type==4|| type==5){
681 projection->Fit(&f1);
682 if (type==4) stat= f1.GetParameter(1);
683 if (type==5) stat= f1.GetParameter(2);
685 //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
686 graph->SetPoint(icount,xcenter, ycenter, stat);
692 TGraph * TStatToolkit::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
696 // delta - number of bins to integrate
697 // type - 0 - mean value
699 TAxis * xaxis = his->GetXaxis();
700 TAxis * yaxis = his->GetYaxis();
701 // TAxis * zaxis = his->GetZaxis();
702 Int_t nbinx = xaxis->GetNbins();
703 Int_t nbiny = yaxis->GetNbins();
706 TGraph *graph = new TGraph(nbinx);
708 for (Int_t ix=0; ix<nbinx;ix++){
709 Float_t xcenter = xaxis->GetBinCenter(ix);
710 // Float_t ycenter = yaxis->GetBinCenter(iy);
711 snprintf(name,1000,"%s_%d",his->GetName(), ix);
712 TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
714 if (type==0) stat = projection->GetMean();
715 if (type==1) stat = projection->GetRMS();
716 if (type==2 || type==3){
718 TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
719 if (type==2) stat= vec[1];
720 if (type==3) stat= vec[0];
722 if (type==4|| type==5){
723 projection->Fit(&f1);
724 if (type==4) stat= f1.GetParameter(1);
725 if (type==5) stat= f1.GetParameter(2);
727 //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
728 graph->SetPoint(icount,xcenter, stat);
738 TString* TStatToolkit::FitPlane(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Bool_t fix0){
740 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
741 // returns chi2, fitParam and covMatrix
742 // returns TString with fitted formula
745 TString formulaStr(formula);
746 TString drawStr(drawCommand);
747 TString cutStr(cuts);
750 TString strVal(drawCommand);
751 if (strVal.Contains(":")){
752 TObjArray* valTokens = strVal.Tokenize(":");
753 drawStr = valTokens->At(0)->GetName();
754 ferr = valTokens->At(1)->GetName();
758 formulaStr.ReplaceAll("++", "~");
759 TObjArray* formulaTokens = formulaStr.Tokenize("~");
760 Int_t dim = formulaTokens->GetEntriesFast();
762 fitParam.ResizeTo(dim);
763 covMatrix.ResizeTo(dim,dim);
765 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
766 fitter->StoreData(kTRUE);
767 fitter->ClearPoints();
769 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
770 if (entries == -1) return new TString("An ERROR has occured during fitting!");
771 Double_t **values = new Double_t*[dim+1] ;
773 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
774 if (entries == -1) return new TString("An ERROR has occured during fitting!");
775 Double_t *errors = new Double_t[entries];
776 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
778 for (Int_t i = 0; i < dim + 1; i++){
780 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
781 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
783 if (entries != centries) return new TString("An ERROR has occured during fitting!");
784 values[i] = new Double_t[entries];
785 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
788 // add points to the fitter
789 for (Int_t i = 0; i < entries; i++){
791 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
792 fitter->AddPoint(x, values[dim][i], errors[i]);
796 if (frac>0.5 && frac<1){
797 fitter->EvalRobust(frac);
800 fitter->FixParameter(0,0);
804 fitter->GetParameters(fitParam);
805 fitter->GetCovarianceMatrix(covMatrix);
806 chi2 = fitter->GetChisquare();
808 // TString *preturnFormula = new TString(Form("%f*(",fitParam[0])), &returnFormula = *preturnFormula;
810 // for (Int_t iparam = 0; iparam < dim; iparam++) {
811 // returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]/fitParam[0]));
812 // if (iparam < dim-1) returnFormula.Append("+");
814 // returnFormula.Append(" )");
816 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
818 for (Int_t iparam = 0; iparam < dim; iparam++) {
819 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
820 if (iparam < dim-1) returnFormula.Append("+");
822 returnFormula.Append(" )");
827 delete formulaTokens;
830 return preturnFormula;
833 TString* TStatToolkit::FitPlaneConstrain(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Double_t constrain){
835 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
836 // returns chi2, fitParam and covMatrix
837 // returns TString with fitted formula
840 TString formulaStr(formula);
841 TString drawStr(drawCommand);
842 TString cutStr(cuts);
845 TString strVal(drawCommand);
846 if (strVal.Contains(":")){
847 TObjArray* valTokens = strVal.Tokenize(":");
848 drawStr = valTokens->At(0)->GetName();
849 ferr = valTokens->At(1)->GetName();
853 formulaStr.ReplaceAll("++", "~");
854 TObjArray* formulaTokens = formulaStr.Tokenize("~");
855 Int_t dim = formulaTokens->GetEntriesFast();
857 fitParam.ResizeTo(dim);
858 covMatrix.ResizeTo(dim,dim);
860 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
861 fitter->StoreData(kTRUE);
862 fitter->ClearPoints();
864 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
865 if (entries == -1) return new TString("An ERROR has occured during fitting!");
866 Double_t **values = new Double_t*[dim+1] ;
868 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
869 if (entries == -1) return new TString("An ERROR has occured during fitting!");
870 Double_t *errors = new Double_t[entries];
871 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
873 for (Int_t i = 0; i < dim + 1; i++){
875 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
876 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
878 if (entries != centries) return new TString("An ERROR has occured during fitting!");
879 values[i] = new Double_t[entries];
880 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
883 // add points to the fitter
884 for (Int_t i = 0; i < entries; i++){
886 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
887 fitter->AddPoint(x, values[dim][i], errors[i]);
890 for (Int_t i = 0; i < dim; i++){
892 for (Int_t j=0; j<dim;j++) if (i!=j) x[j]=0;
894 fitter->AddPoint(x, 0, constrain);
900 if (frac>0.5 && frac<1){
901 fitter->EvalRobust(frac);
903 fitter->GetParameters(fitParam);
904 fitter->GetCovarianceMatrix(covMatrix);
905 chi2 = fitter->GetChisquare();
907 // TString *preturnFormula = new TString(Form("%f*(",fitParam[0])), &returnFormula = *preturnFormula;
909 // for (Int_t iparam = 0; iparam < dim; iparam++) {
910 // returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]/fitParam[0]));
911 // if (iparam < dim-1) returnFormula.Append("+");
913 // returnFormula.Append(" )");
915 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
917 for (Int_t iparam = 0; iparam < dim; iparam++) {
918 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
919 if (iparam < dim-1) returnFormula.Append("+");
921 returnFormula.Append(" )");
926 delete formulaTokens;
929 return preturnFormula;
934 TString* TStatToolkit::FitPlaneFixed(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop){
936 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
937 // returns chi2, fitParam and covMatrix
938 // returns TString with fitted formula
941 TString formulaStr(formula);
942 TString drawStr(drawCommand);
943 TString cutStr(cuts);
946 TString strVal(drawCommand);
947 if (strVal.Contains(":")){
948 TObjArray* valTokens = strVal.Tokenize(":");
949 drawStr = valTokens->At(0)->GetName();
950 ferr = valTokens->At(1)->GetName();
954 formulaStr.ReplaceAll("++", "~");
955 TObjArray* formulaTokens = formulaStr.Tokenize("~");
956 Int_t dim = formulaTokens->GetEntriesFast();
958 fitParam.ResizeTo(dim);
959 covMatrix.ResizeTo(dim,dim);
960 TString fitString="x0";
961 for (Int_t i=1; i<dim; i++) fitString+=Form("++x%d",i);
962 TLinearFitter* fitter = new TLinearFitter(dim, fitString.Data());
963 fitter->StoreData(kTRUE);
964 fitter->ClearPoints();
966 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
967 if (entries == -1) return new TString("An ERROR has occured during fitting!");
968 Double_t **values = new Double_t*[dim+1] ;
970 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
971 if (entries == -1) return new TString("An ERROR has occured during fitting!");
972 Double_t *errors = new Double_t[entries];
973 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
975 for (Int_t i = 0; i < dim + 1; i++){
977 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
978 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
980 if (entries != centries) return new TString("An ERROR has occured during fitting!");
981 values[i] = new Double_t[entries];
982 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
985 // add points to the fitter
986 for (Int_t i = 0; i < entries; i++){
988 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
989 fitter->AddPoint(x, values[dim][i], errors[i]);
993 if (frac>0.5 && frac<1){
994 fitter->EvalRobust(frac);
996 fitter->GetParameters(fitParam);
997 fitter->GetCovarianceMatrix(covMatrix);
998 chi2 = fitter->GetChisquare();
1000 // TString *preturnFormula = new TString(Form("%f*(",fitParam[0])), &returnFormula = *preturnFormula;
1002 // for (Int_t iparam = 0; iparam < dim; iparam++) {
1003 // returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]/fitParam[0]));
1004 // if (iparam < dim-1) returnFormula.Append("+");
1006 // returnFormula.Append(" )");
1008 TString *preturnFormula = new TString("("), &returnFormula = *preturnFormula;
1010 for (Int_t iparam = 0; iparam < dim; iparam++) {
1011 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam]));
1012 if (iparam < dim-1) returnFormula.Append("+");
1014 returnFormula.Append(" )");
1019 delete formulaTokens;
1022 return preturnFormula;
1029 Int_t TStatToolkit::GetFitIndex(TString fString, TString subString){
1031 // fitString - ++ separated list of fits
1032 // substring - ++ separated list of the requiered substrings
1034 // return the last occurance of substring in fit string
1036 TObjArray *arrFit = fString.Tokenize("++");
1037 TObjArray *arrSub = subString.Tokenize("++");
1039 for (Int_t i=0; i<arrFit->GetEntries(); i++){
1041 TString str =arrFit->At(i)->GetName();
1042 for (Int_t isub=0; isub<arrSub->GetEntries(); isub++){
1043 if (str.Contains(arrSub->At(isub)->GetName())==0) isOK=kFALSE;
1051 TString TStatToolkit::FilterFit(TString &input, TString filter, TVectorD ¶m, TMatrixD & covar){
1053 // Filter fit expression make sub-fit
1055 TObjArray *array0= input.Tokenize("++");
1056 TObjArray *array1= filter.Tokenize("++");
1057 //TString *presult=new TString("(0");
1058 TString result="(0.0";
1059 for (Int_t i=0; i<array0->GetEntries(); i++){
1061 TString str(array0->At(i)->GetName());
1062 for (Int_t j=0; j<array1->GetEntries(); j++){
1063 if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;
1067 result+=Form("*(%f)",param[i+1]);
1068 printf("%f\t%f\t%s\n",param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());
1075 void TStatToolkit::Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &vecXk, TMatrixD &covXk){
1077 // Update parameters and covariance - with one measurement
1079 // vecXk - input vector - Updated in function
1080 // covXk - covariance matrix - Updated in function
1081 // delta, sigma, s1 - new measurement, rms of new measurement and the index of measurement
1082 const Int_t knMeas=1;
1083 Int_t knElem=vecXk.GetNrows();
1085 TMatrixD mat1(knElem,knElem); // update covariance matrix
1086 TMatrixD matHk(1,knElem); // vector to mesurement
1087 TMatrixD vecYk(knMeas,1); // Innovation or measurement residual
1088 TMatrixD matHkT(knElem,knMeas); // helper matrix Hk transpose
1089 TMatrixD matSk(knMeas,knMeas); // Innovation (or residual) covariance
1090 TMatrixD matKk(knElem,knMeas); // Optimal Kalman gain
1091 TMatrixD covXk2(knElem,knElem); // helper matrix
1092 TMatrixD covXk3(knElem,knElem); // helper matrix
1093 TMatrixD vecZk(1,1);
1094 TMatrixD measR(1,1);
1096 measR(0,0)=sigma*sigma;
1099 for (Int_t iel=0;iel<knElem;iel++)
1100 for (Int_t ip=0;ip<knMeas;ip++) matHk(ip,iel)=0;
1102 for (Int_t iel=0;iel<knElem;iel++) {
1103 for (Int_t jel=0;jel<knElem;jel++) mat1(iel,jel)=0;
1108 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
1109 matHkT=matHk.T(); matHk.T();
1110 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
1112 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
1113 vecXk += matKk*vecYk; // updated vector
1114 covXk2= (mat1-(matKk*matHk));
1115 covXk3 = covXk2*covXk;
1117 Int_t nrows=covXk3.GetNrows();
1119 for (Int_t irow=0; irow<nrows; irow++)
1120 for (Int_t icol=0; icol<nrows; icol++){
1121 // rounding problems - make matrix again symteric
1122 covXk(irow,icol)=(covXk3(irow,icol)+covXk3(icol,irow))*0.5;
1128 void TStatToolkit::Constrain1D(TString &input, TString filter, TVectorD ¶m, TMatrixD & covar, Double_t mean, Double_t sigma){
1130 // constrain linear fit
1131 // input - string description of fit function
1132 // filter - string filter to select sub fits
1133 // param,covar - parameters and covariance matrix of the fit
1134 // mean,sigma - new measurement uning which the fit is updated
1136 TObjArray *array0= input.Tokenize("++");
1137 TObjArray *array1= filter.Tokenize("++");
1138 TMatrixD paramM(param.GetNrows(),1);
1139 for (Int_t i=0; i<=array0->GetEntries(); i++){paramM(i,0)=param(i);}
1141 for (Int_t i=0; i<array0->GetEntries(); i++){
1143 TString str(array0->At(i)->GetName());
1144 for (Int_t j=0; j<array1->GetEntries(); j++){
1145 if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;
1148 TStatToolkit::Update1D(mean, sigma, i+1, paramM, covar);//
1151 for (Int_t i=0; i<=array0->GetEntries(); i++){
1152 param(i)=paramM(i,0);
1156 TString TStatToolkit::MakeFitString(TString &input, TVectorD ¶m, TMatrixD & covar){
1160 TObjArray *array0= input.Tokenize("++");
1161 TString result="(0.0";
1162 for (Int_t i=0; i<array0->GetEntries(); i++){
1163 TString str(array0->At(i)->GetName());
1165 result+=Form("*(%f)",param[i+1]);
1166 printf("%f\t%f\t%s\n", param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());