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 "AliMathBase.h"
26 #include "Riostream.h"
30 #include "TLinearFitter.h"
33 // includes neccessary for test functions
38 #include "TStopwatch.h"
39 #include "TTreeStream.h"
41 ClassImp(AliMathBase) // Class implementation to enable ROOT I/O
43 AliMathBase::AliMathBase() : TObject()
46 // Default constructor
49 ///////////////////////////////////////////////////////////////////////////
50 AliMathBase::~AliMathBase()
58 //_____________________________________________________________________________
59 void AliMathBase::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
60 , Double_t &sigma, Int_t hh)
63 // Robust estimator in 1D case MI version - (faster than ROOT version)
65 // For the univariate case
66 // estimates of location and scatter are returned in mean and sigma parameters
67 // the algorithm works on the same principle as in multivariate case -
68 // it finds a subset of size hh with smallest sigma, and then returns mean and
69 // sigma of this subset
74 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};
75 Int_t *index=new Int_t[nvectors];
76 TMath::Sort(nvectors, data, index, kFALSE);
78 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
79 Double_t factor = faclts[TMath::Max(0,nquant-1)];
84 Double_t bestmean = 0;
85 Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.); // maximal possible sigma
86 bestsigma *=bestsigma;
88 for (Int_t i=0; i<hh; i++){
89 sumx += data[index[i]];
90 sumx2 += data[index[i]]*data[index[i]];
93 Double_t norm = 1./Double_t(hh);
94 Double_t norm2 = 1./Double_t(hh-1);
95 for (Int_t i=hh; i<nvectors; i++){
96 Double_t cmean = sumx*norm;
97 Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
98 if (csigma<bestsigma){
104 sumx += data[index[i]]-data[index[i-hh]];
105 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
108 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
117 void AliMathBase::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh, Float_t externalfactor)
119 // Modified version of ROOT robust EvaluateUni
120 // robust estimator in 1D case MI version
121 // added external factor to include precision of external measurement
126 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};
127 Int_t *index=new Int_t[nvectors];
128 TMath::Sort(nvectors, data, index, kFALSE);
130 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
131 Double_t factor = faclts[0];
133 // fix proper normalization - Anja
134 factor = faclts[nquant-1];
141 Int_t bestindex = -1;
142 Double_t bestmean = 0;
143 Double_t bestsigma = -1;
144 for (Int_t i=0; i<hh; i++){
145 sumx += data[index[i]];
146 sumx2 += data[index[i]]*data[index[i]];
149 Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
150 Double_t norm = 1./Double_t(hh);
151 for (Int_t i=hh; i<nvectors; i++){
152 Double_t cmean = sumx*norm;
153 Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
154 if (csigma<bestsigma || bestsigma<0){
161 sumx += data[index[i]]-data[index[i-hh]];
162 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
165 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
172 //_____________________________________________________________________________
173 Int_t AliMathBase::Freq(Int_t n, const Int_t *inlist
174 , Int_t *outlist, Bool_t down)
177 // Sort eleements according occurancy
178 // The size of output array has is 2*n
181 Int_t * sindexS = new Int_t[n]; // temp array for sorting
182 Int_t * sindexF = new Int_t[2*n];
183 for (Int_t i=0;i<n;i++) sindexF[i]=0;
185 TMath::Sort(n,inlist, sindexS, down);
186 Int_t last = inlist[sindexS[0]];
193 for(Int_t i=1;i<n; i++){
194 val = inlist[sindexS[i]];
195 if (last == val) sindexF[countPos]++;
198 sindexF[countPos+n] = val;
203 if (last==val) countPos++;
204 // sort according frequency
205 TMath::Sort(countPos, sindexF, sindexS, kTRUE);
206 for (Int_t i=0;i<countPos;i++){
207 outlist[2*i ] = sindexF[sindexS[i]+n];
208 outlist[2*i+1] = sindexF[sindexS[i]];
217 //___AliMathBase__________________________________________________________________________
218 void AliMathBase::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
222 Int_t nbins = his->GetNbinsX();
223 Float_t nentries = his->GetEntries();
228 for (Int_t ibin=1;ibin<nbins; ibin++){
229 ncumul+= his->GetBinContent(ibin);
230 Float_t fraction = Float_t(ncumul)/Float_t(nentries);
231 if (fraction>down && fraction<up){
232 sum+=his->GetBinContent(ibin);
233 mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
234 sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);
238 sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
240 (*param)[0] = his->GetMaximum();
242 (*param)[2] = sigma2;
245 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
248 void AliMathBase::LTM(TH1F * his, TVectorD *param , Float_t fraction, Bool_t verbose){
252 Int_t nbins = his->GetNbinsX();
253 Int_t nentries = (Int_t)his->GetEntries();
254 Double_t *data = new Double_t[nentries];
256 for (Int_t ibin=1;ibin<nbins; ibin++){
257 Float_t entriesI = his->GetBinContent(ibin);
258 Float_t xcenter= his->GetBinCenter(ibin);
259 for (Int_t ic=0; ic<entriesI; ic++){
260 if (npoints<nentries){
261 data[npoints]= xcenter;
266 Double_t mean, sigma;
267 Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
268 npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
269 AliMathBase::EvaluateUni(npoints, data, mean,sigma,npoints2);
271 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
272 (*param)[0] = his->GetMaximum();
278 Double_t AliMathBase::FitGaus(TH1F* his, TVectorD *param, TMatrixD *matrix, Float_t xmin, Float_t xmax, Bool_t verbose){
280 // Fit histogram with gaussian function
283 // return value- chi2 - if negative ( not enough points)
284 // his - input histogram
285 // param - vector with parameters
286 // xmin, xmax - range to fit - if xmin=xmax=0 - the full histogram range used
288 // 1. Step - make logarithm
289 // 2. Linear fit (parabola) - more robust - always converge
290 // 3. In case of small statistic bins are averaged
292 static TLinearFitter fitter(3,"pol2");
296 if (his->GetMaximum()<4) return -1;
297 if (his->GetEntries()<12) return -1;
298 if (his->GetRMS()<mat.GetTol()) return -1;
299 Float_t maxEstimate = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
300 Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
302 if (maxEstimate<1) return -1;
303 Int_t nbins = his->GetNbinsX();
309 xmin = his->GetXaxis()->GetXmin();
310 xmax = his->GetXaxis()->GetXmax();
312 for (Int_t iter=0; iter<2; iter++){
313 fitter.ClearPoints();
315 for (Int_t ibin=1;ibin<nbins+1; ibin++){
317 Float_t entriesI = his->GetBinContent(ibin);
318 for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
319 if (ibin+delta>1 &&ibin+delta<nbins-1){
320 entriesI += his->GetBinContent(ibin+delta);
325 Double_t xcenter= his->GetBinCenter(ibin);
326 if (xcenter<xmin || xcenter>xmax) continue;
327 Double_t error=1./TMath::Sqrt(countB);
330 if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
331 cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
332 if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
334 if (entriesI>1&&cont>1){
335 fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
341 fitter.GetParameters(par);
349 fitter.GetParameters(par);
350 fitter.GetCovarianceMatrix(mat);
351 if (TMath::Abs(par[1])<mat.GetTol()) return -1;
352 if (TMath::Abs(par[2])<mat.GetTol()) return -1;
353 Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
354 //fitter.GetParameters();
355 if (!param) param = new TVectorD(3);
356 if (!matrix) matrix = new TMatrixD(3,3);
357 (*param)[1] = par[1]/(-2.*par[2]);
358 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
359 (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1]);
364 printf("Chi2=%f\n",chi2);
365 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
366 f1->SetParameter(0, (*param)[0]);
367 f1->SetParameter(1, (*param)[1]);
368 f1->SetParameter(2, (*param)[2]);
374 Double_t AliMathBase::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD *matrix, Bool_t verbose){
376 // Fit histogram with gaussian function
379 // nbins: size of the array and number of histogram bins
380 // xMin, xMax: histogram range
381 // param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
382 // matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
385 // >0: the chi2 returned by TLinearFitter
386 // -3: only three points have been used for the calculation - no fitter was used
387 // -2: only two points have been used for the calculation - center of gravity was uesed for calculation
388 // -1: only one point has been used for the calculation - center of gravity was uesed for calculation
389 // -4: invalid result!!
392 // 1. Step - make logarithm
393 // 2. Linear fit (parabola) - more robust - always converge
395 static TLinearFitter fitter(3,"pol2");
396 static TMatrixD mat(3,3);
397 static Double_t kTol = mat.GetTol();
398 fitter.StoreData(kFALSE);
399 fitter.ClearPoints();
404 Float_t rms = TMath::RMS(nBins,arr);
405 Float_t max = TMath::MaxElement(nBins,arr);
406 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
415 for (Int_t i=0; i<nBins; i++){
417 if (arr[i]>0) nfilled++;
420 if (max<4) return -4;
421 if (entries<12) return -4;
422 if (rms<kTol) return -4;
428 for (Int_t ibin=0;ibin<nBins; ibin++){
429 Float_t entriesI = arr[ibin];
431 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
433 Float_t error = 1./TMath::Sqrt(entriesI);
434 Float_t val = TMath::Log(Float_t(entriesI));
435 fitter.AddPoint(&xcenter,val,error);
438 A(npoints,1)=xcenter;
439 A(npoints,2)=xcenter*xcenter;
441 meanCOG+=xcenter*entriesI;
442 rms2COG +=xcenter*entriesI*xcenter;
453 //analytic calculation of the parameters for three points
462 // use fitter for more than three points
464 fitter.GetParameters(par);
465 fitter.GetCovarianceMatrix(mat);
466 chi2 = fitter.GetChisquare()/Float_t(npoints);
468 if (TMath::Abs(par[1])<kTol) return -4;
469 if (TMath::Abs(par[2])<kTol) return -4;
471 if (!param) param = new TVectorD(3);
472 if (!matrix) matrix = new TMatrixD(3,3); // !!!!might be a memory leek. use dummy matrix pointer to call this function!
474 (*param)[1] = par[1]/(-2.*par[2]);
475 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
476 Double_t lnparam0 = par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1];
477 if ( lnparam0>307 ) return -4;
478 (*param)[0] = TMath::Exp(lnparam0);
483 printf("Chi2=%f\n",chi2);
484 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
485 f1->SetParameter(0, (*param)[0]);
486 f1->SetParameter(1, (*param)[1]);
487 f1->SetParameter(2, (*param)[2]);
494 //use center of gravity for 2 points
498 (*param)[1] = meanCOG;
499 (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
505 (*param)[1] = meanCOG;
506 (*param)[2] = binWidth/TMath::Sqrt(12);
514 Float_t AliMathBase::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
517 // calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
518 // return COG; in case of failure return xMin
525 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
527 for (Int_t ibin=0; ibin<nBins; ibin++){
528 Float_t entriesI = (Float_t)arr[ibin];
529 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
531 meanCOG += xcenter*entriesI;
532 rms2COG += xcenter*entriesI*xcenter;
537 if ( sumCOG == 0 ) return xMin;
542 (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
543 if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
554 ///////////////////////////////////////////////////////////////
555 ////////////// TEST functions /////////////////////////
556 ///////////////////////////////////////////////////////////////
562 void AliMathBase::TestGausFit(Int_t nhistos){
564 // Test performance of the parabolic - gaussian fit - compare it with
566 // nhistos - number of histograms to be used for test
568 TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
570 Float_t *xTrue = new Float_t[nhistos];
571 Float_t *sTrue = new Float_t[nhistos];
572 TVectorD **par1 = new TVectorD*[nhistos];
573 TVectorD **par2 = new TVectorD*[nhistos];
577 TH1F **h1f = new TH1F*[nhistos];
578 TF1 *myg = new TF1("myg","gaus");
579 TF1 *fit = new TF1("fit","gaus");
583 for (Int_t i=0;i<nhistos; i++){
584 par1[i] = new TVectorD(3);
585 par2[i] = new TVectorD(3);
586 h1f[i] = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
587 xTrue[i]= gRandom->Rndm();
589 sTrue[i]= .75+gRandom->Rndm()*.5;
590 myg->SetParameters(1,xTrue[i],sTrue[i]);
591 h1f[i]->FillRandom("myg");
597 for (Int_t i=0; i<nhistos; i++){
598 h1f[i]->Fit(fit,"0q");
599 (*par1[i])(0) = fit->GetParameter(0);
600 (*par1[i])(1) = fit->GetParameter(1);
601 (*par1[i])(2) = fit->GetParameter(2);
604 printf("Gaussian fit\t");
608 //AliMathBase gaus fit
609 for (Int_t i=0; i<nhistos; i++){
610 AliMathBase::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
614 printf("Parabolic fit\t");
617 for (Int_t i=0;i<nhistos; i++){
618 Float_t xt = xTrue[i];
619 Float_t st = sTrue[i];
628 for (Int_t i=0;i<nhistos; i++){
645 TGraph2D * AliMathBase::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
649 // delta - number of bins to integrate
650 // type - 0 - mean value
652 TAxis * xaxis = his->GetXaxis();
653 TAxis * yaxis = his->GetYaxis();
654 // TAxis * zaxis = his->GetZaxis();
655 Int_t nbinx = xaxis->GetNbins();
656 Int_t nbiny = yaxis->GetNbins();
659 TGraph2D *graph = new TGraph2D(nbinx*nbiny);
661 for (Int_t ix=0; ix<nbinx;ix++)
662 for (Int_t iy=0; iy<nbiny;iy++){
663 Float_t xcenter = xaxis->GetBinCenter(ix);
664 Float_t ycenter = yaxis->GetBinCenter(iy);
665 sprintf(name,"%s_%d_%d",his->GetName(), ix,iy);
666 TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
668 if (type==0) stat = projection->GetMean();
669 if (type==1) stat = projection->GetRMS();
670 if (type==2 || type==3){
672 AliMathBase::LTM((TH1F*)projection,&vec,0.7);
673 if (type==2) stat= vec[1];
674 if (type==3) stat= vec[0];
676 if (type==4|| type==5){
677 projection->Fit(&f1);
678 if (type==4) stat= f1.GetParameter(1);
679 if (type==5) stat= f1.GetParameter(2);
681 //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
682 graph->SetPoint(icount,xcenter, ycenter, stat);
688 TGraph * AliMathBase::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
692 // delta - number of bins to integrate
693 // type - 0 - mean value
695 TAxis * xaxis = his->GetXaxis();
696 TAxis * yaxis = his->GetYaxis();
697 // TAxis * zaxis = his->GetZaxis();
698 Int_t nbinx = xaxis->GetNbins();
699 Int_t nbiny = yaxis->GetNbins();
702 TGraph *graph = new TGraph(nbinx);
704 for (Int_t ix=0; ix<nbinx;ix++){
705 Float_t xcenter = xaxis->GetBinCenter(ix);
706 // Float_t ycenter = yaxis->GetBinCenter(iy);
707 sprintf(name,"%s_%d",his->GetName(), ix);
708 TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
710 if (type==0) stat = projection->GetMean();
711 if (type==1) stat = projection->GetRMS();
712 if (type==2 || type==3){
714 AliMathBase::LTM((TH1F*)projection,&vec,0.7);
715 if (type==2) stat= vec[1];
716 if (type==3) stat= vec[0];
718 if (type==4|| type==5){
719 projection->Fit(&f1);
720 if (type==4) stat= f1.GetParameter(1);
721 if (type==5) stat= f1.GetParameter(2);
723 //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
724 graph->SetPoint(icount,xcenter, stat);