/************************************************************************** * Copyright(c) 1998-2009, 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. * **************************************************************************/ ///////////////////////////////////////////////////////////// // // AliHFMassFitter for the fit of invariant mass distribution // of charmed mesons // // Author: C.Bianchin, chiara.bianchin@pd.infn.it ///////////////////////////////////////////////////////////// #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include "AliHFMassFitter.h" ClassImp(AliHFMassFitter) //************** constructors AliHFMassFitter::AliHFMassFitter() : TNamed(), fhistoInvMass(0), fminMass(0), fmaxMass(0), fNbin(0), fParsSize(1), fNFinalPars(1), fFitPars(0), fWithBkg(0), ftypeOfFit4Bkg(0), ftypeOfFit4Sgn(0), ffactor(1), fntuParam(0), fMass(1.85), fSigmaSgn(0.012), fSideBands(0), fFixPar(0), fSideBandl(0), fSideBandr(0), fcounter(0), fContourGraph(0) { // default constructor cout<<"Default constructor:"<Clone("fhistoInvMass"); fhistoInvMass->SetDirectory(0); fminMass=minvalue; fmaxMass=maxvalue; if(rebin!=1) RebinMass(rebin); else fNbin=(Int_t)fhistoInvMass->GetNbinsX(); CheckRangeFit(); ftypeOfFit4Bkg=fittypeb; ftypeOfFit4Sgn=fittypes; if(ftypeOfFit4Bkg!=0 && ftypeOfFit4Bkg!=1 && ftypeOfFit4Bkg!=2) fWithBkg=kFALSE; else fWithBkg=kTRUE; if (!fWithBkg) cout<<"Fit Histogram of Signal only"< 0){ if(fFitPars) { delete[] fFitPars; fFitPars=NULL; } fFitPars=new Float_t[fParsSize]; memcpy(fFitPars,mfit.fFitPars,mfit.fParsSize*sizeof(Float_t)); if(fFixPar) { delete[] fFixPar; fFixPar=NULL; } fFixPar=new Bool_t[fNFinalPars]; memcpy(fFixPar,mfit.fFixPar,mfit.fNFinalPars*sizeof(Float_t)); } // fFitPars=new Float_t[fParsSize]; // for(Int_t i=0;i=fNFinalPars) { cout<<"Error! Parameter out of bounds! Max is "<=fNFinalPars) { cout<<"Error! Parameter out of bounds! Max is "<Clone(); fhistoInvMass = new TH1F(*histoToFit); fhistoInvMass->SetDirectory(0); cout<<"SetHisto pointer "<SetBranchAddress("intbkg1",&fFitPars[0]); fntuParam->SetBranchAddress("slope1",&fFitPars[1]); fntuParam->SetBranchAddress("conc1",&fFitPars[2]); fntuParam->SetBranchAddress("intGB",&fFitPars[3]); fntuParam->SetBranchAddress("meanGB",&fFitPars[4]); fntuParam->SetBranchAddress("sigmaGB",&fFitPars[5]); fntuParam->SetBranchAddress("intbkg2",&fFitPars[6]); fntuParam->SetBranchAddress("slope2",&fFitPars[7]); fntuParam->SetBranchAddress("conc2",&fFitPars[8]); fntuParam->SetBranchAddress("inttot",&fFitPars[9]); fntuParam->SetBranchAddress("slope3",&fFitPars[10]); fntuParam->SetBranchAddress("conc3",&fFitPars[11]); fntuParam->SetBranchAddress("intsgn",&fFitPars[12]); fntuParam->SetBranchAddress("meansgn",&fFitPars[13]); fntuParam->SetBranchAddress("sigmasgn",&fFitPars[14]); fntuParam->SetBranchAddress("intbkg1Err",&fFitPars[15]); fntuParam->SetBranchAddress("slope1Err",&fFitPars[16]); fntuParam->SetBranchAddress("conc1Err",&fFitPars[17]); fntuParam->SetBranchAddress("intGBErr",&fFitPars[18]); fntuParam->SetBranchAddress("meanGBErr",&fFitPars[19]); fntuParam->SetBranchAddress("sigmaGBErr",&fFitPars[20]); fntuParam->SetBranchAddress("intbkg2Err",&fFitPars[21]); fntuParam->SetBranchAddress("slope2Err",&fFitPars[22]); fntuParam->SetBranchAddress("conc2Err",&fFitPars[23]); fntuParam->SetBranchAddress("inttotErr",&fFitPars[24]); fntuParam->SetBranchAddress("slope3Err",&fFitPars[25]); fntuParam->SetBranchAddress("conc3Err",&fFitPars[26]); fntuParam->SetBranchAddress("intsgnErr",&fFitPars[27]); fntuParam->SetBranchAddress("meansgnErr",&fFitPars[28]); fntuParam->SetBranchAddress("sigmasgnErr",&fFitPars[29]); } else { if(ftypeOfFit4Bkg==3){ fntuParam->SetBranchAddress("intbkg1",&fFitPars[0]); fntuParam->SetBranchAddress("slope1",¬hing); fntuParam->SetBranchAddress("conc1",¬hing); fntuParam->SetBranchAddress("intGB",&fFitPars[1]); fntuParam->SetBranchAddress("meanGB",&fFitPars[2]); fntuParam->SetBranchAddress("sigmaGB",&fFitPars[3]); fntuParam->SetBranchAddress("intbkg2",&fFitPars[4]); fntuParam->SetBranchAddress("slope2",¬hing); fntuParam->SetBranchAddress("conc2",¬hing); fntuParam->SetBranchAddress("inttot",&fFitPars[6]); fntuParam->SetBranchAddress("slope3",¬hing); fntuParam->SetBranchAddress("conc3",¬hing); fntuParam->SetBranchAddress("intsgn",&fFitPars[6]); fntuParam->SetBranchAddress("meansgn",&fFitPars[7]); fntuParam->SetBranchAddress("sigmasgn",&fFitPars[8]); fntuParam->SetBranchAddress("intbkg1Err",&fFitPars[9]); fntuParam->SetBranchAddress("slope1Err",¬hing); fntuParam->SetBranchAddress("conc1Err",¬hing); fntuParam->SetBranchAddress("intGBErr",&fFitPars[10]); fntuParam->SetBranchAddress("meanGBErr",&fFitPars[11]); fntuParam->SetBranchAddress("sigmaGBErr",&fFitPars[12]); fntuParam->SetBranchAddress("intbkg2Err",&fFitPars[13]); fntuParam->SetBranchAddress("slope2Err",¬hing); fntuParam->SetBranchAddress("conc2Err",¬hing); fntuParam->SetBranchAddress("inttotErr",&fFitPars[15]); fntuParam->SetBranchAddress("slope3Err",¬hing); fntuParam->SetBranchAddress("conc3Err",¬hing); fntuParam->SetBranchAddress("intsgnErr",&fFitPars[15]); fntuParam->SetBranchAddress("meansgnErr",&fFitPars[16]); fntuParam->SetBranchAddress("sigmasgnErr",&fFitPars[17]); } else{ fntuParam->SetBranchAddress("intbkg1",&fFitPars[0]); fntuParam->SetBranchAddress("slope1",&fFitPars[1]); fntuParam->SetBranchAddress("conc1",¬hing); fntuParam->SetBranchAddress("intGB",&fFitPars[2]); fntuParam->SetBranchAddress("meanGB",&fFitPars[3]); fntuParam->SetBranchAddress("sigmaGB",&fFitPars[4]); fntuParam->SetBranchAddress("intbkg2",&fFitPars[5]); fntuParam->SetBranchAddress("slope2",&fFitPars[6]); fntuParam->SetBranchAddress("conc2",¬hing); fntuParam->SetBranchAddress("inttot",&fFitPars[7]); fntuParam->SetBranchAddress("slope3",&fFitPars[8]); fntuParam->SetBranchAddress("conc3",¬hing); fntuParam->SetBranchAddress("intsgn",&fFitPars[9]); fntuParam->SetBranchAddress("meansgn",&fFitPars[10]); fntuParam->SetBranchAddress("sigmasgn",&fFitPars[11]); fntuParam->SetBranchAddress("intbkg1Err",&fFitPars[12]); fntuParam->SetBranchAddress("slope1Err",&fFitPars[13]); fntuParam->SetBranchAddress("conc1Err",¬hing); fntuParam->SetBranchAddress("intGBErr",&fFitPars[14]); fntuParam->SetBranchAddress("meanGBErr",&fFitPars[15]); fntuParam->SetBranchAddress("sigmaGBErr",&fFitPars[16]); fntuParam->SetBranchAddress("intbkg2Err",&fFitPars[17]); fntuParam->SetBranchAddress("slope2Err",&fFitPars[18]); fntuParam->SetBranchAddress("conc2Err",¬hing); fntuParam->SetBranchAddress("inttotErr",&fFitPars[19]); fntuParam->SetBranchAddress("slope3Err",&fFitPars[20]); fntuParam->SetBranchAddress("conc3Err",¬hing); fntuParam->SetBranchAddress("intsgnErr",&fFitPars[21]); fntuParam->SetBranchAddress("meansgnErr",&fFitPars[22]); fntuParam->SetBranchAddress("sigmasgnErr",&fFitPars[23]); } } fntuParam->TTree::Fill(); } //_________________________________________________________________________ TNtuple* AliHFMassFitter::NtuParamOneShot(char *ntuname){ // Create, fill and return ntuple with fit parameters InitNtuParam(ntuname); FillNtuParam(); return fntuParam; } //_________________________________________________________________________ void AliHFMassFitter::RebinMass(Int_t bingroup){ // Rebin invariant mass histogram if(!fhistoInvMass){ cout<<"Histogram not set"<GetNbinsX(); if(bingroup<1){ cout<<"Error! Cannot group "<Rebin(bingroup); if(nbinshisto%bingroup != 0) CheckRangeFit(); fNbin = fhistoInvMass->GetNbinsX(); cout<<"New number of bins: "< integral(exponential)=A/B*exp(B*x)](min,max) //-> A = B*integral/(exp(B*max)-exp(B*min)) where integral can be written //as integralTot- integralGaus (=par [2]) //Par: // * [0] = integralBkg; // * [1] = B; //exponential = [1]*[0]/(exp([1]*max)-exp([1]*min))*exp([1]*x) total = par[0+firstPar]*par[1+firstPar]/(TMath::Exp(par[1+firstPar]*fmaxMass)-TMath::Exp(par[1+firstPar]*fminMass))*TMath::Exp(par[1+firstPar]*x[0]); break; case 1: //linear //y=a+b*x -> integral = a(max-min)+1/2*b*(max^2-min^2) -> a = (integral-1/2*b*(max^2-min^2))/(max-min)=integral/(max-min)-1/2*b*(max+min) // * [0] = integralBkg; // * [1] = b; total= par[0+firstPar]/(fmaxMass-fminMass)+par[1+firstPar]*(x[0]-0.5*(fmaxMass+fminMass)); break; case 2: //polynomial //y=a+b*x+c*x**2 -> integral = a(max-min) + 1/2*b*(max^2-min^2) + //+ 1/3*c*(max^3-min^3) -> //a = (integral-1/2*b*(max^2-min^2)-1/3*c*(max^3-min^3))/(max-min) // * [0] = integralBkg; // * [1] = b; // * [2] = c; total = par[0+firstPar]/(fmaxMass-fminMass)+par[1]*(x[0]-0.5*(fmaxMass+fminMass))+par[2+firstPar]*(x[0]*x[0]-1/3.*(fmaxMass*fmaxMass*fmaxMass-fminMass*fminMass*fminMass)/(fmaxMass-fminMass)); break; case 3: total=par[0+firstPar]; break; // default: // Types of Fit Functions for Background: // * 0 = exponential; // * 1 = linear; // * 2 = polynomial 2nd order // * 3 = no background"<GetBinWidth(8); if (fNbin==0) fNbin=fhistoInvMass->GetNbinsX(); Double_t minHisto=fhistoInvMass->GetBinLowEdge(1); Double_t maxHisto=fhistoInvMass->GetBinLowEdge(fNbin)+width; if(fMass-fminMass < 0 || fmaxMass-fMass <0) { cout<<"Left limit of range > mean or right limit of range < mean: change left/right limit or initial mean value"<GetNbinsX(); Double_t width=fhistoInvMass->GetBinWidth(1); Double_t minhisto=fhistoInvMass->GetBinLowEdge(1), maxhisto=fhistoInvMass->GetBinLowEdge(nbins)+width; //check if limits are inside histogram range if( fminMass-minhisto < 0. ) { cout<<"Out of histogram left bound!"< 0. ) { cout<<"Out of histogram right bound!"<FindBin(fminMass); if (fminMass > fhistoInvMass->GetBinCenter(binl)) binl++; fminMass=fhistoInvMass->GetBinLowEdge(binl); if(TMath::Abs(tmp-fminMass) > 1e-6){ cout<<"Left bound "<FindBin(fmaxMass); if (fmaxMass < fhistoInvMass->GetBinCenter(binr)) binr--; fmaxMass=fhistoInvMass->GetBinLowEdge(binr)+width; if(TMath::Abs(tmp-fmaxMass) > 1e-6){ cout<<"Right bound "<SetOwner(); } fContourGraph->SetName(listname); //function names TString bkgname="funcbkg"; TString bkg1name="funcbkg1"; TString massname="funcmass"; //Total integral Double_t totInt = fhistoInvMass->Integral(fhistoInvMass->FindBin(fminMass), fhistoInvMass->FindBin(fmaxMass), "width"); fSideBands = kTRUE; Double_t width=fhistoInvMass->GetBinWidth(8); cout<<"fNbin"<GetNbinsX(); //Double_t minHisto=fhistoInvMass->GetBinLowEdge(1); //Double_t maxHisto=fhistoInvMass->GetBinLowEdge(fNbin)+width; Bool_t ok=SideBandsBounds(); if(!ok) return kFALSE; //sidebands integral - first approx (from histo) Double_t sideBandsInt=(Double_t)fhistoInvMass->Integral(1,fSideBandl,"width") + (Double_t)fhistoInvMass->Integral(fSideBandr,fNbin,"width"); cout<<"------nbin = "<SetParameters(sideBandsInt,-2.); break; case 1: funcbkg->SetParNames("BkgInt","Slope"); funcbkg->SetParameters(sideBandsInt,-100.); break; case 2: funcbkg->SetParNames("BkgInt","Coef1","Coef2"); funcbkg->SetParameters(sideBandsInt,-10.,5); break; case 3: if(ftypeOfFit4Sgn==0){ funcbkg->SetParNames("Const"); funcbkg->SetParameter(0,0.); funcbkg->FixParameter(0,0.); } break; default: cout<<"Wrong choise of ftypeOfFit4Bkg ("<Fit(bkgname.Data(),"R,L,E,0"); for(Int_t i=0;iGetParameter(i); //cout<GetParameter(i)<<"\t"; fFitPars[nFitPars+2*bkgPar+3+i]= funcbkg->GetParError(i); //cout<GetParError(i)<GetParameter(0); funcbkg->SetRange(fminMass,fmaxMass); intbkg1 = funcbkg->Integral(fminMass,fmaxMass); if(ftypeOfFit4Bkg!=3) slope1 = funcbkg->GetParameter(1); if(ftypeOfFit4Bkg==2) conc1 = funcbkg->GetParameter(2); cout<<"First fit: \nintbkg1 = "<SetParNames("IntGB","MeanGB","SigmaGB","BkgInt","Coef1","Coef2"); funcbkg1->SetParameters(0.5*(totInt-intbkg1),fMass,ffactor*fSigmaSgn,intbkg1,slope1,conc1); //cout<<"Parameters set to: "<<0.5*(totInt-intbkg1)<<"\t"<SetParNames("IntGB","MeanGB","SigmaGB","Const"); funcbkg1->SetParameters(0.5*totInt,fMass,ffactor*fSigmaSgn,0.); funcbkg1->FixParameter(3,0.); } else{ //expo or linear if(ftypeOfFit4Bkg==0) cout<<"*** Exponential Fit ***"<SetParNames("IntGB","MeanGB","SigmaGB","BkgInt","Slope"); funcbkg1->SetParameters(0.5*(totInt-intbkg1),fMass,ffactor*fSigmaSgn,intbkg1,slope1); } } Int_t status=fhistoInvMass->Fit(bkg1name.Data(),"R,L,E,+,0"); if (status != 0){ cout<<"Minuit returned "<GetParameter(i); //cout<GetParameter(i); fFitPars[nFitPars+3*bkgPar-6+i]= funcbkg1->GetParError(i); //cout<<"\t"<GetParError(i)<GetParameter(3); if(ftypeOfFit4Bkg!=3) slope1 = funcbkg1->GetParameter(4); if(ftypeOfFit4Bkg==2) conc1 = funcbkg1->GetParameter(5); } else { bkgPar+=3; for(Int_t i=0;i<3;i++){ fFitPars[bkgPar-3+i]=0.; cout<GetParameter(i); cout<GetParameter(i)<<"\t"; fFitPars[nFitPars+3*bkgPar-3+i]= funcbkg->GetParError(i); cout<GetParError(i)<SetParNames("TotInt","Slope","SgnInt","Mean","Sigma"); funcmass->SetParameters(totInt,slope1,sgnInt,fMass,fSigmaSgn); //cout<<"Parameters set to: "<FixParameter(0,totInt); } if(fFixPar[1]){ cout<<"fix2"<FixParameter(1,slope1); } if(fFixPar[2]){ cout<<"fix3"<FixParameter(2,sgnInt); } if(fFixPar[3]){ cout<<"fix4"<FixParameter(3,fMass); } if(fFixPar[4]){ cout<<"fix5"<FixParameter(4,fSigmaSgn); } } if (nFitPars==6){ funcmass->SetParNames("TotInt","Coef1","Coef2","SgnInt","Mean","Sigma"); funcmass->SetParameters(totInt,slope1,conc1,sgnInt,fMass,fSigmaSgn); //cout<<"Parameters set to: "<SetParameters(0.,0.5*totInt,fMass,fSigmaSgn); else funcmass->SetParameters(0.,totInt,fMass,fSigmaSgn); if(fFixPar[0]) funcmass->FixParameter(0,0.); //cout<<"Parameters set to: "<<0.5*totInt<<"\t"<GetParameter(nFitPars-2); fSigmaSgn=funcmass->GetParameter(nFitPars-1); for(Int_t i=0;iGetParameter(i); fFitPars[nFitPars+4*bkgPar-6+i]= funcmass->GetParError(i); //cout<GetParameter(i)<<"\t\t"<GetParError(i)<GetName(),fhistoInvMass->GetName()); // TH1F *fhistocopy=new TH1F(*fhistoInvMass); // canvas->cd(); // fhistocopy->DrawClone(); // if(ftypeOfFit4Sgn == 1) funcbkg1->DrawClone("sames"); // else funcbkg->DrawClone("sames"); // funcmass->DrawClone("sames"); // cout<<"Drawn"<GetParameter(nFitPars-1) <0 || funcmass->GetParameter(nFitPars-2) <0 || funcmass->GetParameter(nFitPars-3) <0 ) { cout<<"IntS or mean or sigma negative. You may tray to SetInitialGaussianSigma(..) and SetInitialGaussianMean(..)"<SetErrorDef(errDef[i]); cont[i] = (TGraph*)gMinuit->Contour(80,kpar,jpar); cout<<"Minuit Status = "<GetStatus()<GetParName(kpar); TString titleY = funcmass->GetParName(jpar); for (Int_t i=0; i<2; i++) { cont[i]->SetName( Form("cperr%d_%d%d", i, kpar, jpar) ); cont[i]->SetTitle(title); cont[i]->GetXaxis()->SetTitle(titleX); cont[i]->GetYaxis()->SetTitle(titleY); cont[i]->GetYaxis()->SetLabelSize(0.033); cont[i]->GetYaxis()->SetTitleSize(0.033); cont[i]->GetYaxis()->SetTitleOffset(1.67); fContourGraph->Add(cont[i]); } // plot them TString cvname = Form("c%d%d", kpar, jpar); TCanvas *c4=new TCanvas(cvname,cvname,600,600); c4->cd(); cont[1]->SetFillColor(38); cont[1]->Draw("alf"); cont[0]->SetFillColor(9); cont[0]->Draw("lf"); } } } if (ftypeOfFit4Sgn == 1 && funcbkg1) { delete funcbkg1; funcbkg1=NULL; } if (funcbkg) { delete funcbkg; funcbkg=NULL; } if (funcmass) { delete funcmass; funcmass=NULL; } AddFunctionsToHisto(); if (draw) DrawFit(); return kTRUE; } //_________________________________________________________________________ Double_t AliHFMassFitter::GetChiSquare() const{ TF1 *funcmass=(TF1*)fhistoInvMass->GetFunction("funcmass"); return funcmass->GetChisquare(); } //_________________________________________________________________________ Double_t AliHFMassFitter::GetReducedChiSquare() const{ TF1 *funcmass=(TF1*)fhistoInvMass->GetFunction("funcmass"); return funcmass->GetChisquare()/funcmass->GetNDF(); } //*********output //_________________________________________________________________________ void AliHFMassFitter::GetFitPars(Float_t *vector) const { // Return fit parameters for(Int_t i=0;iFindObject(testname.Data()); if(testfunc){ cout<<"AddFunctionsToHisto already used: exiting...."<GetListOfFunctions(); hlist->ls(); TF1 *b=(TF1*)hlist->FindObject(bkgname.Data()); if(!b){ cout<GetName()<SetParName(i,b->GetParName(i)); bfullrange->SetParameter(i,b->GetParameter(i)); bfullrange->SetParError(i,b->GetParError(i)); } bfullrange->SetLineStyle(4); bfullrange->SetLineColor(14); bkgnamesave += "Recalc"; TF1 *blastpar=new TF1(bkgnamesave.Data(),this,&AliHFMassFitter::FitFunction4Bkg,fminMass,fmaxMass,np,"AliHFMassFitter","FitFunction4Bkg"); TF1 *mass=fhistoInvMass->GetFunction("funcmass"); if (!mass){ cout<<"funcmass doesn't exist "<SetParameter(0,mass->GetParameter(0)-mass->GetParameter(np)); blastpar->SetParError(0,mass->GetParError(np)); if (np>=2) { blastpar->SetParameter(1,mass->GetParameter(1)); blastpar->SetParError(1,mass->GetParError(1)); } if (np==3) { blastpar->SetParameter(2,mass->GetParameter(2)); blastpar->SetParError(2,mass->GetParError(2)); } blastpar->SetLineStyle(1); blastpar->SetLineColor(2); hlist->Add((TF1*)bfullrange->Clone()); hlist->Add((TF1*)blastpar->Clone()); hlist->ls(); if(bfullrange) { delete bfullrange; bfullrange=NULL; } if(blastpar) { delete blastpar; blastpar=NULL; } } //_________________________________________________________________________ TH1F* AliHFMassFitter::GetHistoClone() const{ TH1F* hout=(TH1F*)fhistoInvMass->Clone(fhistoInvMass->GetName()); return hout; } //_________________________________________________________________________ void AliHFMassFitter::WriteHisto(TString path) const { //Write the histogram in the default file HFMassFitterOutput.root if (fcounter == 0) { cout<<"Use MassFitter method before WriteHisto"<Clone(); path += "HFMassFitterOutput.root"; TFile *output; if (fcounter == 1) output = new TFile(path.Data(),"recreate"); else output = new TFile(path.Data(),"update"); output->cd(); hget->Write(); fContourGraph->Write(); output->Close(); cout<GetName()<<" written in "<Clone(); path += "HFMassFitterOutput.root"; TFile *output = new TFile(path.Data(),"update"); output->cd(); fntuParam->Write(); //nget->Write(); output->Close(); //cout<GetName()<<" written in "<GetName()<<" written in "<SetOptStat(0); gStyle->SetCanvasColor(0); gStyle->SetFrameFillColor(0); TString type=""; switch (ftypeOfFit4Bkg){ case 0: type="Exp"; //3+2 break; case 1: type="Lin"; //3+2 break; case 2: type="Pl2"; //3+3 break; case 3: type="noB"; //3+1 break; } TString filename=Form("%sMassFit.root",type.Data()); filename.Prepend(userIDstring); path.Append(filename); TFile* outputcv=new TFile(path.Data(),"update"); TCanvas* c=(TCanvas*)GetPad(nsigma,writeFitInfo); c->SetName(Form("%s%s%s",c->GetName(),userIDstring.Data(),type.Data())); if(draw)c->DrawClone(); outputcv->cd(); c->Write(); outputcv->Close(); } //_________________________________________________________________________ TVirtualPad* AliHFMassFitter::GetPad(Double_t nsigma,Int_t writeFitInfo)const{ //return a TVirtualPad with the fitted histograms and info TString cvtitle="fit of "; cvtitle+=fhistoInvMass->GetName(); TString cvname="c"; cvname+=fcounter; TCanvas *c=new TCanvas(cvname,cvtitle); PlotFit(c->cd(),nsigma,writeFitInfo); return c->cd(); } //_________________________________________________________________________ void AliHFMassFitter::PlotFit(TVirtualPad* pd,Double_t nsigma,Int_t writeFitInfo)const{ //plot histogram, fit functions and write parameters according to verbosity level (0,1,>1) gStyle->SetOptStat(0); gStyle->SetCanvasColor(0); gStyle->SetFrameFillColor(0); cout<<"nsigma = "<GetFunction("funcbkgFullRange")->DrawClone("same"); hdraw->GetFunction("funcbkgRecalc")->DrawClone("same"); hdraw->GetFunction("funcmass")->DrawClone("same"); if(writeFitInfo > 0){ TPaveText *pinfob=new TPaveText(0.6,0.86,1.,1.,"NDC"); TPaveText *pinfom=new TPaveText(0.6,0.7,1.,.87,"NDC"); pinfob->SetBorderSize(0); pinfob->SetFillStyle(0); pinfom->SetBorderSize(0); pinfom->SetFillStyle(0); TF1* ff=fhistoInvMass->GetFunction("funcmass"); for (Int_t i=fNFinalPars-3;iSetTextColor(kBlue); TString str=Form("%s = %f #pm %f",ff->GetParName(i),ff->GetParameter(i),ff->GetParError(i)); if(!(writeFitInfo==1 && i==fNFinalPars-3)) pinfom->AddText(str); } pd->cd(); pinfom->DrawClone(); TPaveText *pinfo2=new TPaveText(0.1,0.1,0.6,0.4,"NDC"); pinfo2->SetBorderSize(0); pinfo2->SetFillStyle(0); Double_t signif, signal, bkg, errsignif, errsignal, errbkg; Significance(nsigma,signif,errsignif); Signal(nsigma,signal,errsignal); Background(nsigma,bkg, errbkg); /* Significance(1.828,1.892,signif,errsignif); Signal(1.828,1.892,signal,errsignal); Background(1.828,1.892,bkg, errbkg); */ TString str=Form("Significance (%.0f#sigma) %.1f #pm %.1f ",nsigma,signif,errsignif); pinfo2->AddText(str); str=Form("S (%.0f#sigma) %.0f #pm %.0f ",nsigma,signal,errsignal); pinfo2->AddText(str); str=Form("B (%.0f#sigma) %.0f #pm %.0f",nsigma,bkg,errbkg); pinfo2->AddText(str); pd->cd(); pinfo2->Draw(); if(writeFitInfo > 1){ for (Int_t i=0;iSetTextColor(kRed); TString str=Form("%s = %f #pm %f",ff->GetParName(i),ff->GetParameter(i),ff->GetParError(i)); pinfob->AddText(str); } pd->cd(); pinfob->DrawClone(); } } return; } //_________________________________________________________________________ void AliHFMassFitter::DrawHere(TVirtualPad* pd,Double_t nsigma,Int_t writeFitInfo) const { //draws histogram together with fit functions with default nice colors in user canvas PlotFit(pd,nsigma,writeFitInfo); pd->Draw(); } //_________________________________________________________________________ void AliHFMassFitter::DrawFit(Double_t nsigma) const{ //draws histogram together with fit functions with default nice colors gStyle->SetOptStat(0); gStyle->SetCanvasColor(0); gStyle->SetFrameFillColor(0); TCanvas* c=(TCanvas*)GetPad(nsigma,1); c->Draw(); } //_________________________________________________________________________ void AliHFMassFitter::PrintParTitles() const{ //prints to screen the parameters names TF1 *f=fhistoInvMass->GetFunction("funcmass"); if(!f) { cout<<"Fit function not found"<GetParName(i)<GetFunction(massname.Data()); // if(!funcmass){ // cout<<"AliHFMassFitter::Signal() ERROR -> Mass distr function not found!"<GetFunction(bkgname.Data()); // else funcbkg=fhistoInvMass->GetFunction(bkg1name.Data()); // if(!funcbkg){ // cout<<"AliHFMassFitter::Signal() ERROR -> Bkg function not found!"<GetParameter(np); // intSerr=funcmass->GetParError(np); // cout<<"Sgn relative error evaluation from fit: "<Integral(min, max)/fhistoInvMass->GetBinWidth(4); // signal=mass - background; // errsignal=TMath::Sqrt((intSerr/intS*mass)*(intSerr/intS*mass)/*assume relative error is the same as for total integral*/ + errbackground*errbackground); return; } //_________________________________________________________________________ void AliHFMassFitter::Signal(Double_t min, Double_t max, Double_t &signal,Double_t &errsignal) const { // Return signal integral in a range if(fcounter==0) { cout<<"Use MassFitter method before Signal"<GetFunction(massname.Data()); if(!funcmass){ cout<<"AliHFMassFitter::Signal() ERROR -> Mass distr function not found!"<GetFunction(bkgname.Data()); else funcbkg=fhistoInvMass->GetFunction(bkg1name.Data()); if(!funcbkg){ cout<<"AliHFMassFitter::Signal() ERROR -> Bkg function not found!"<GetParameter(np); intSerr=funcmass->GetParError(np); cout<<"Sgn relative error evaluation from fit: "<Integral(min, max)/fhistoInvMass->GetBinWidth(4); signal=mass - background; errsignal=(intSerr/intS)*signal;/*assume relative error is the same as for total integral*/ } //_________________________________________________________________________ void AliHFMassFitter::Background(Double_t nOfSigma,Double_t &background,Double_t &errbackground) const { // Return background integral in mean+- n sigma if(fcounter==0) { cout<<"Use MassFitter method before Background"<GetFunction(bkgname.Data()); // else funcbkg=fhistoInvMass->GetFunction(bkg1name.Data()); // if(!funcbkg){ // cout<<"AliHFMassFitter::Background() ERROR -> Bkg function not found!"<GetParameter(0); // intBerr=funcbkg->GetParError(0); // cout<<"Bkg relative error evaluation: from final parameters of the fit: "<Integral(1,fSideBandl)+fhistoInvMass->Integral(fSideBandr,fNbin); // Double_t sum2=0; // for(Int_t i=1;i<=fSideBandl;i++){ // sum2+=fhistoInvMass->GetBinError(i)*fhistoInvMass->GetBinError(i); // } // for(Int_t i=fSideBandr;i<=fNbin;i++){ // sum2+=fhistoInvMass->GetBinError(i)*fhistoInvMass->GetBinError(i); // } // intBerr=TMath::Sqrt(sum2); // cout<<"Bkg relative error evaluation: from histo: "<Integral(min,max)/(Double_t)fhistoInvMass->GetBinWidth(2); // errbackground=intBerr/intB*background; // assume relative error is the same as for total integral // //cout<<"integral = "<Integral(min, max)<<"\tbinW = "<GetBinWidth(2)<GetFunction(bkgname.Data()); else funcbkg=fhistoInvMass->GetFunction(bkg1name.Data()); if(!funcbkg){ cout<<"AliHFMassFitter::Background() ERROR -> Bkg function not found!"<GetParameter(0); intBerr=funcbkg->GetParError(0); cout<<"Bkg relative error evaluation: from final parameters of the fit: "<Integral(1,fSideBandl)+fhistoInvMass->Integral(fSideBandr,fNbin); Double_t sum2=0; for(Int_t i=1;i<=fSideBandl;i++){ sum2+=fhistoInvMass->GetBinError(i)*fhistoInvMass->GetBinError(i); } for(Int_t i=fSideBandr;i<=fNbin;i++){ sum2+=fhistoInvMass->GetBinError(i)*fhistoInvMass->GetBinError(i); } intBerr=TMath::Sqrt(sum2); cout<<"Bkg relative error evaluation: from histo: "<Integral(min,max)/(Double_t)fhistoInvMass->GetBinWidth(2); errbackground=intBerr/intB*background; // assume relative error is the same as for total integral //cout<<"integral = "<Integral(min, max)<<"\tbinW = "<GetBinWidth(2)<