// Author: Daniel.Lohner@cern.ch #include "AliTRDNDFast.h" #include "TCanvas.h" extern Double_t langaufunc(Double_t *x, Double_t *par) { // needs protection, parameter [3]>0!!!!! // needs protection, parameter [4]>0!!!!! //Fit parameters: //par[0]=Width (scale) parameter of Landau density //par[1]=Most Probable (MP, location) parameter of Landau density //par[2]=Total area (integral -inf to inf, normalization constant) //par[3]=Width (sigma) of convoluted Gaussian function //par[4]=Exponential Slope Parameter // //In the Landau distribution (represented by the CERNLIB approximation), //the maximum is located at x=-0.22278298 with the location parameter=0. //This shift is corrected within this function, so that the actual //maximum is identical to the MP parameter. // Numeric constants Double_t invsq2pi = 0.3989422804014; // (2 pi)^(-1/2) Double_t mpshift = -0.22278298; // Landau maximum location // Control constants Double_t np = 100.0; // number of convolution stpdf Double_t sc = 5.0; // convolution extends to +-sc Gaussian sigmas // Variables Double_t xx; Double_t mpc; Double_t fland; Double_t sum = 0.0; Double_t xlow,xupp; Double_t step; Double_t i; // MP shift correction mpc = par[1] - mpshift * par[0]; // Range of convolution integral xlow = x[0] - sc * par[3]; xupp = x[0] + sc * par[3]; if(xlow<0)xlow=0; if(xuppSetDirectory(0); } } AliTRDNDFast::AliTRDNDFast(const AliTRDNDFast &ref) : TObject(ref), fNDim(ref.fNDim), fTitle(ref.fTitle), fFunc(NULL), fHistos(NULL), fPars() { Cleanup(); Init(); for(Int_t idim=0;idimClone(Form("%s_axis_%d",GetName(),idim)); fHistos[idim]->SetDirectory(0); for(Int_t ipar=0;iparGetParameter(ii); func->GetParLimits(ii,low[ii],high[ii]); } Double_t scalefactor=mpv/start[1]; for(Int_t ii=0;iiSetParLimits(ii, low[ii], high[ii]); } func->SetParameters(start); } //--------------------------------------------------------------- //--------------------------------------------------------------- TF1 *AliTRDNDFast::GetLangauFun(TString funcname,Double_t range[2],Double_t scalefactor){ Double_t start[kNpar] = {120, 1000, 1, 100, 1e-5}; Double_t low[kNpar] = {10, 300, 0.01, 1, 0.}; Double_t high[kNpar] = {1000, 3000, 100, 1000, 1.}; TF1 *hlandauE=new TF1(funcname.Data(),langaufunc,0,range[1],kNpar); hlandauE->SetParameters(start); hlandauE->SetParNames("Width","MP","Area","GSigma","delta"); for (int i=0; iSetParLimits(i, low[i], high[i]); } if(scalefactor!=1){ScaleLangauFun(hlandauE,scalefactor*start[1]);} return hlandauE; } TF1 *AliTRDNDFast::FitLandau(TString name,TH1F *htemp,Double_t range[2],TString option){ TF1 *fitlandau1D=GetLangauFun(name,range); TF1 fit("land","landau"); Double_t max=htemp->GetXaxis()->GetBinCenter(htemp->GetMaximumBin()); fit.SetParameter(1,max); fit.SetParLimits(1,0,htemp->GetXaxis()->GetXmax()); fit.SetParameter(2,0.3*max); // MPV may never become negative!!!!! htemp->Fit("land",option.Data(),"",range[0],range[1]); ScaleLangauFun(fitlandau1D,fit.GetParameter(1)); htemp->Fit(fitlandau1D,option.Data(),"",range[0],range[1]); // range for references return fitlandau1D; } void AliTRDNDFast::BuildHistos(){ for(Int_t idim=0;idimReset(); // Fill Histo Double_t pars[kNpar]; for(Int_t ipar=0;iparGetNbinsX()+1;ii++){ Double_t xval=fHistos[idim]->GetXaxis()->GetBinCenter(ii); Double_t val=langaufunc(&xval,&pars[0]); //Double_t val=fFunc[idim]->Eval(xval); fHistos[idim]->SetBinContent(ii,val); } // Normalization fHistos[idim]->Scale(1./fHistos[idim]->Integral(1,fHistos[idim]->GetNbinsX(),"width")); } } void AliTRDNDFast::Build(Double_t **pars){ // pars[ndim][npar] for(Int_t idim=0;idimDivide(2,1); if(fNDim==3)CANV->Divide(2,2); if(fNDim>6)CANV->Divide(3,3); // Fit and Extract Parameters for(Int_t idim=0;idimcd(idim+1); gPad->SetLogy(); range[0]=hdEdx[idim]->GetXaxis()->GetXmin(); range[1]=hdEdx[idim]->GetXaxis()->GetXmax(); // Norm Histogram hdEdx[idim]->Scale(1./hdEdx[idim]->Integral(1,hdEdx[idim]->GetNbinsX(),"width")); // Fit Histogram fFunc[idim]=FitLandau(Form("fit%d",idim),hdEdx[idim],range,"RMB0"); // Norm Landau fFunc[idim]->SetParameter(2,fFunc[idim]->GetParameter(2)/fFunc[idim]->Integral(range[0],range[1])); hdEdx[idim]->DrawCopy(); fFunc[idim]->DrawCopy("same"); // Set Pars for(Int_t ipar=0;iparGetParameter(ipar),idim); } if(path!="")CANV->Print(Form("%s/%s_Build.pdf",path.Data(),fTitle.Data())); delete CANV; BuildHistos(); } Double_t AliTRDNDFast::Eval(Double_t *point) const{ Double_t val=1; for(Int_t idim=0;idimGetXaxis()->FindBin(point[idim]); val*=fHistos[idim]->GetBinContent(bin); } return val; } void AliTRDNDFast::Random(Double_t *point) const{ for(Int_t idim=0;idimGetRandom(); } } void AliTRDNDFast::Random(Double_t *point,AliTRDNDFast *nd0,AliTRDNDFast *nd1,Double_t w0,Double_t w1){ for(Int_t idim=0;idimfNDim;idim++){ point[idim]=GetRandomInterpolation(nd0->fHistos[idim],nd1->fHistos[idim],w0,w1); } } Int_t AliTRDNDFast::BinarySearchInterpolation(Int_t start,Int_t end,Double_t *a0,Double_t *a1,Double_t w0,Double_t w1,Double_t val){ if((end-start)==1)return start; Int_t mid=Int_t((end+start)/2); Double_t valmid=(w0*a0[mid]+w1*a1[mid])/(w0+w1); if(val>=valmid)return BinarySearchInterpolation(mid,end,a0,a1,w0,w1,val); return BinarySearchInterpolation(start,mid,a0,a1,w0,w1,val); } Double_t AliTRDNDFast::GetRandomInterpolation(TH1F *hist0,TH1F *hist1,Double_t w0,Double_t w1){ // Draw Random Variable Double_t rand=gRandom->Rndm(); // Get Integral Arrays Double_t *integral0=hist0->GetIntegral(); Double_t *integral1=hist1->GetIntegral(); Int_t nbinsX=hist0->GetNbinsX(); Int_t ibin=BinarySearchInterpolation(1,nbinsX+1,integral0,integral1,w0,w1,rand); return hist0->GetXaxis()->GetBinCenter(ibin); }