/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /////////////////////////////////////////////////////////////////////////// // Class TStatToolkit // // Subset of matheamtical functions not included in the TMath // /////////////////////////////////////////////////////////////////////////// #include "TMath.h" #include "Riostream.h" #include "TH1F.h" #include "TH3.h" #include "TF1.h" #include "TTree.h" #include "TChain.h" #include "TObjString.h" #include "TLinearFitter.h" // // includes neccessary for test functions // #include "TSystem.h" #include "TRandom.h" #include "TStopwatch.h" #include "TTreeStream.h" #include "TStatToolkit.h" ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O TStatToolkit::TStatToolkit() : TObject() { // // Default constructor // } /////////////////////////////////////////////////////////////////////////// TStatToolkit::~TStatToolkit() { // // Destructor // } //_____________________________________________________________________________ void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean , Double_t &sigma, Int_t hh) { // // Robust estimator in 1D case MI version - (faster than ROOT version) // // For the univariate case // estimates of location and scatter are returned in mean and sigma parameters // the algorithm works on the same principle as in multivariate case - // it finds a subset of size hh with smallest sigma, and then returns mean and // sigma of this subset // if (hh==0) hh=(nvectors+2)/2; 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}; Int_t *index=new Int_t[nvectors]; TMath::Sort(nvectors, data, index, kFALSE); Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11); Double_t factor = faclts[TMath::Max(0,nquant-1)]; Double_t sumx =0; Double_t sumx2 =0; Int_t bestindex = -1; Double_t bestmean = 0; Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.); // maximal possible sigma bestsigma *=bestsigma; for (Int_t i=0; i0){ // fix proper normalization - Anja factor = faclts[nquant-1]; } // // Double_t sumx =0; Double_t sumx2 =0; Int_t bestindex = -1; Double_t bestmean = 0; Double_t bestsigma = -1; for (Int_t i=0; iGetNbinsX(); Float_t nentries = his->GetEntries(); Float_t sum =0; Float_t mean = 0; Float_t sigma2 = 0; Float_t ncumul=0; for (Int_t ibin=1;ibinGetBinContent(ibin); Float_t fraction = Float_t(ncumul)/Float_t(nentries); if (fraction>down && fractionGetBinContent(ibin); mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin); sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin); } } mean/=sum; sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean)); if (param){ (*param)[0] = his->GetMaximum(); (*param)[1] = mean; (*param)[2] = sigma2; } if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2); } void TStatToolkit::LTM(TH1F * his, TVectorD *param , Float_t fraction, Bool_t verbose){ // // LTM // Int_t nbins = his->GetNbinsX(); Int_t nentries = (Int_t)his->GetEntries(); Double_t *data = new Double_t[nentries]; Int_t npoints=0; for (Int_t ibin=1;ibinGetBinContent(ibin); Float_t xcenter= his->GetBinCenter(ibin); for (Int_t ic=0; icGetMaximum(); (*param)[1] = mean; (*param)[2] = sigma; } } Double_t TStatToolkit::FitGaus(TH1F* his, TVectorD *param, TMatrixD *matrix, Float_t xmin, Float_t xmax, Bool_t verbose){ // // Fit histogram with gaussian function // // Prameters: // return value- chi2 - if negative ( not enough points) // his - input histogram // param - vector with parameters // xmin, xmax - range to fit - if xmin=xmax=0 - the full histogram range used // Fitting: // 1. Step - make logarithm // 2. Linear fit (parabola) - more robust - always converge // 3. In case of small statistic bins are averaged // static TLinearFitter fitter(3,"pol2"); TVectorD par(3); TVectorD sigma(3); TMatrixD mat(3,3); if (his->GetMaximum()<4) return -1; if (his->GetEntries()<12) return -1; if (his->GetRMS()GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS())); Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate)); if (maxEstimate<1) return -1; Int_t nbins = his->GetNbinsX(); Int_t npoints=0; // if (xmin>=xmax){ xmin = his->GetXaxis()->GetXmin(); xmax = his->GetXaxis()->GetXmax(); } for (Int_t iter=0; iter<2; iter++){ fitter.ClearPoints(); npoints=0; for (Int_t ibin=1;ibinGetBinContent(ibin); for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){ if (ibin+delta>1 &&ibin+deltaGetBinContent(ibin+delta); countB++; } } entriesI/=countB; Double_t xcenter= his->GetBinCenter(ibin); if (xcenterxmax) continue; Double_t error=1./TMath::Sqrt(countB); Float_t cont=2; if (iter>0){ if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0; cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter); if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB)); } if (entriesI>1&&cont>1){ fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error); npoints++; } } if (npoints>3){ fitter.Eval(); fitter.GetParameters(par); }else{ break; } } if (npoints<=3){ return -1; } fitter.GetParameters(par); fitter.GetCovarianceMatrix(mat); if (TMath::Abs(par[1])Print(); printf("Chi2=%f\n",chi2); TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax()); f1->SetParameter(0, (*param)[0]); f1->SetParameter(1, (*param)[1]); f1->SetParameter(2, (*param)[2]); f1->Draw("same"); } return chi2; } Double_t TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD *matrix, Bool_t verbose){ // // Fit histogram with gaussian function // // Prameters: // nbins: size of the array and number of histogram bins // xMin, xMax: histogram range // param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma) // matrix: covariance matrix -- not implemented yet, pass dummy matrix!!! // // Return values: // >0: the chi2 returned by TLinearFitter // -3: only three points have been used for the calculation - no fitter was used // -2: only two points have been used for the calculation - center of gravity was uesed for calculation // -1: only one point has been used for the calculation - center of gravity was uesed for calculation // -4: invalid result!! // // Fitting: // 1. Step - make logarithm // 2. Linear fit (parabola) - more robust - always converge // static TLinearFitter fitter(3,"pol2"); static TMatrixD mat(3,3); static Double_t kTol = mat.GetTol(); fitter.StoreData(kFALSE); fitter.ClearPoints(); TVectorD par(3); TVectorD sigma(3); TMatrixD A(3,3); TMatrixD b(3,1); Float_t rms = TMath::RMS(nBins,arr); Float_t max = TMath::MaxElement(nBins,arr); Float_t binWidth = (xMax-xMin)/(Float_t)nBins; Float_t meanCOG = 0; Float_t rms2COG = 0; Float_t sumCOG = 0; Float_t entries = 0; Int_t nfilled=0; for (Int_t i=0; i0) nfilled++; } if (max<4) return -4; if (entries<12) return -4; if (rms1){ Double_t xcenter = xMin+(ibin+0.5)*binWidth; Float_t error = 1./TMath::Sqrt(entriesI); Float_t val = TMath::Log(Float_t(entriesI)); fitter.AddPoint(&xcenter,val,error); if (npoints<3){ A(npoints,0)=1; A(npoints,1)=xcenter; A(npoints,2)=xcenter*xcenter; b(npoints,0)=val; meanCOG+=xcenter*entriesI; rms2COG +=xcenter*entriesI*xcenter; sumCOG +=entriesI; } npoints++; } } Double_t chi2 = 0; if (npoints>=3){ if ( npoints == 3 ){ //analytic calculation of the parameters for three points A.Invert(); TMatrixD res(1,3); res.Mult(A,b); par[0]=res(0,0); par[1]=res(0,1); par[2]=res(0,2); chi2 = -3.; } else { // use fitter for more than three points fitter.Eval(); fitter.GetParameters(par); fitter.GetCovarianceMatrix(mat); chi2 = fitter.GetChisquare()/Float_t(npoints); } if (TMath::Abs(par[1])307 ) return -4; (*param)[0] = TMath::Exp(lnparam0); if (verbose){ par.Print(); mat.Print(); param->Print(); printf("Chi2=%f\n",chi2); TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax); f1->SetParameter(0, (*param)[0]); f1->SetParameter(1, (*param)[1]); f1->SetParameter(2, (*param)[2]); f1->Draw("same"); } return chi2; } if (npoints == 2){ //use center of gravity for 2 points meanCOG/=sumCOG; rms2COG /=sumCOG; (*param)[0] = max; (*param)[1] = meanCOG; (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG)); chi2=-2.; } if ( npoints == 1 ){ meanCOG/=sumCOG; (*param)[0] = max; (*param)[1] = meanCOG; (*param)[2] = binWidth/TMath::Sqrt(12); chi2=-1.; } return chi2; } Float_t TStatToolkit::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum) { // // calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax // return COG; in case of failure return xMin // Float_t meanCOG = 0; Float_t rms2COG = 0; Float_t sumCOG = 0; Int_t npoints = 0; Float_t binWidth = (xMax-xMin)/(Float_t)nBins; for (Int_t ibin=0; ibin0 ){ meanCOG += xcenter*entriesI; rms2COG += xcenter*entriesI*xcenter; sumCOG += entriesI; npoints++; } } if ( sumCOG == 0 ) return xMin; meanCOG/=sumCOG; if ( rms ){ rms2COG /=sumCOG; (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG)); if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12); } if ( sum ) (*sum) = sumCOG; return meanCOG; } /////////////////////////////////////////////////////////////// ////////////// TEST functions ///////////////////////// /////////////////////////////////////////////////////////////// void TStatToolkit::TestGausFit(Int_t nhistos){ // // Test performance of the parabolic - gaussian fit - compare it with // ROOT gauss fit // nhistos - number of histograms to be used for test // TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root"); Float_t *xTrue = new Float_t[nhistos]; Float_t *sTrue = new Float_t[nhistos]; TVectorD **par1 = new TVectorD*[nhistos]; TVectorD **par2 = new TVectorD*[nhistos]; TMatrixD dummy(3,3); TH1F **h1f = new TH1F*[nhistos]; TF1 *myg = new TF1("myg","gaus"); TF1 *fit = new TF1("fit","gaus"); gRandom->SetSeed(0); //init for (Int_t i=0;iRndm(); gSystem->Sleep(2); sTrue[i]= .75+gRandom->Rndm()*.5; myg->SetParameters(1,xTrue[i],sTrue[i]); h1f[i]->FillRandom("myg"); } TStopwatch s; s.Start(); //standard gaus fit for (Int_t i=0; iFit(fit,"0q"); (*par1[i])(0) = fit->GetParameter(0); (*par1[i])(1) = fit->GetParameter(1); (*par1[i])(2) = fit->GetParameter(2); } s.Stop(); printf("Gaussian fit\t"); s.Print(); s.Start(); //TStatToolkit gaus fit for (Int_t i=0; iGetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy); } s.Stop(); printf("Parabolic fit\t"); s.Print(); //write stream for (Int_t i=0;iGetXaxis(); TAxis * yaxis = his->GetYaxis(); // TAxis * zaxis = his->GetZaxis(); Int_t nbinx = xaxis->GetNbins(); Int_t nbiny = yaxis->GetNbins(); char name[1000]; Int_t icount=0; TGraph2D *graph = new TGraph2D(nbinx*nbiny); TF1 f1("f1","gaus"); for (Int_t ix=0; ixGetBinCenter(ix); Float_t ycenter = yaxis->GetBinCenter(iy); sprintf(name,"%s_%d_%d",his->GetName(), ix,iy); TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1); Float_t stat= 0; if (type==0) stat = projection->GetMean(); if (type==1) stat = projection->GetRMS(); if (type==2 || type==3){ TVectorD vec(3); TStatToolkit::LTM((TH1F*)projection,&vec,0.7); if (type==2) stat= vec[1]; if (type==3) stat= vec[0]; } if (type==4|| type==5){ projection->Fit(&f1); if (type==4) stat= f1.GetParameter(1); if (type==5) stat= f1.GetParameter(2); } //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat); graph->SetPoint(icount,xcenter, ycenter, stat); icount++; } return graph; } TGraph * TStatToolkit::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){ // // // // delta - number of bins to integrate // type - 0 - mean value TAxis * xaxis = his->GetXaxis(); TAxis * yaxis = his->GetYaxis(); // TAxis * zaxis = his->GetZaxis(); Int_t nbinx = xaxis->GetNbins(); Int_t nbiny = yaxis->GetNbins(); char name[1000]; Int_t icount=0; TGraph *graph = new TGraph(nbinx); TF1 f1("f1","gaus"); for (Int_t ix=0; ixGetBinCenter(ix); // Float_t ycenter = yaxis->GetBinCenter(iy); sprintf(name,"%s_%d",his->GetName(), ix); TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny); Float_t stat= 0; if (type==0) stat = projection->GetMean(); if (type==1) stat = projection->GetRMS(); if (type==2 || type==3){ TVectorD vec(3); TStatToolkit::LTM((TH1F*)projection,&vec,0.7); if (type==2) stat= vec[1]; if (type==3) stat= vec[0]; } if (type==4|| type==5){ projection->Fit(&f1); if (type==4) stat= f1.GetParameter(1); if (type==5) stat= f1.GetParameter(2); } //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat); graph->SetPoint(icount,xcenter, stat); icount++; } return graph; } 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){ // // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts // returns chi2, fitParam and covMatrix // returns TString with fitted formula // TString formulaStr(formula); TString drawStr(drawCommand); TString cutStr(cuts); TString ferr("1"); TString strVal(drawCommand); if (strVal.Contains(":")){ TObjArray* valTokens = strVal.Tokenize(":"); drawStr = valTokens->At(0)->GetName(); ferr = valTokens->At(1)->GetName(); } formulaStr.ReplaceAll("++", "~"); TObjArray* formulaTokens = formulaStr.Tokenize("~"); Int_t dim = formulaTokens->GetEntriesFast(); fitParam.ResizeTo(dim); covMatrix.ResizeTo(dim,dim); TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim)); fitter->StoreData(kTRUE); fitter->ClearPoints(); Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start); if (entries == -1) return new TString("An ERROR has occured during fitting!"); Double_t **values = new Double_t*[dim+1] ; // entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start); if (entries == -1) return new TString("An ERROR has occured during fitting!"); Double_t *errors = new Double_t[entries]; memcpy(errors, tree->GetV1(), entries*sizeof(Double_t)); for (Int_t i = 0; i < dim + 1; i++){ Int_t centries = 0; if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start); else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start); if (entries != centries) return new TString("An ERROR has occured during fitting!"); values[i] = new Double_t[entries]; memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t)); } // add points to the fitter for (Int_t i = 0; i < entries; i++){ Double_t x[1000]; for (Int_t j=0; jAddPoint(x, values[dim][i], errors[i]); } fitter->Eval(); if (frac>0.5 && frac<1){ fitter->EvalRobust(frac); } fitter->GetParameters(fitParam); fitter->GetCovarianceMatrix(covMatrix); chi2 = fitter->GetChisquare(); chi2 = chi2; npoints = entries; // TString *preturnFormula = new TString(Form("%f*(",fitParam[0])), &returnFormula = *preturnFormula; // for (Int_t iparam = 0; iparam < dim; iparam++) { // returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]/fitParam[0])); // if (iparam < dim-1) returnFormula.Append("+"); // } // returnFormula.Append(" )"); TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; for (Int_t iparam = 0; iparam < dim; iparam++) { returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1])); if (iparam < dim-1) returnFormula.Append("+"); } returnFormula.Append(" )"); delete formulaTokens; delete fitter; delete[] values; return preturnFormula; }