]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliHFMassFitterTest.C
Update and addition of LS analysis (Renu, Giacomo, Francesco)
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliHFMassFitterTest.C
CommitLineData
380e2bf9 1// hinvMass -> invariant mass histogram to fit
2// typeb -> type of fit function for background. Can be 0, 1, 2, 3 (see AliHFMassFitter.cxx for details)
3// types -> type of fit function for signal. Can be 0, 1 (see AliHFMassFitter.cxx for details)
4// factor4refl -> sigmaRefl=factor4refl*sigmaSgn (have a look to AliHFMassFitter.cxx for details). Set it if types=1
5
dd5201c4 6void AliHFMassFitterTest(TH1F *hinvMass, Int_t typeb, Int_t types, Int_t factor4refl=1){
7
8 Bool_t useParFiles=kFALSE;
9 Int_t load=gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/LoadLibraries.C");
10 LoadLibraries(useParFiles);
11
380e2bf9 12 Int_t nbin=hinvMass->GetNbinsX();
13 Double_t min=hinvMass->GetBinLowEdge(1),max=hinvMass->GetBinLowEdge(nbin)+hinvMass->GetBinWidth(nbin);
14 //TH1F *hMass=new TH1F("hMass","Invariant Mass",nbin,min,max);
15
16 AliHFMassFitter *fitter=new AliHFMassFitter(hinvMass,min, max,1,typeb,types);
17
18 fitter->SetReflectionSigmaFactor(factor4refl);
19 //fitter->SetRangeFit(min,max);
a60f573d 20 fitter->MassFitter(kTRUE); //kFALSE do not draw the histogram
380e2bf9 21
22 TNtuple *ntu=fitter->NtuParamOneShot("ntuInvMass");
a60f573d 23
24 Double_t mD0pdg=TDatabasePDG::Instance()->GetParticle(421)->Mass();
25 Double_t sigmaval= 0.027;
26 Double_t min=mD0pdg-sigmaval, max=mD0pdg+sigmaval;
27 Double_t significance3sigma,significanceRange,signal3sigma,signalRange,bkg3sigma,bkgRange;
28 Double_t err3,errRange,errS3,errSRange,errB3,errBRange;
29
30 fitter->Signal(3,signal3sigma,errS3);
31 fitter->Signal(min, max,signalRange,errSRange);
32 fitter->Background(3,bkg3sigma,errB3);
33 fitter->Background(min, max,bkgRange,errBRange);
380e2bf9 34 fitter->Significance(3,significance3sigma,err3);
a60f573d 35 fitter->Significance(min, max,significanceRange,errRange);
36
37 cout<<"# 3 sigma: \n signal = "<<signal3sigma<<" +/- "<<errS3<<"\t backgdround = "<<bkg3sigma<<" +/- "<<errB3<<"\tsignificance = "<<significance3sigma<<" +/- "<<err3<<"\n";
38
39 cout<<"# range ("<<min<<", "<<max<<"): \n signal = "<<signalRange<<" +/- "<<errSRange<<"\t backgdround = "<<bkgRange<<" +/- "<<errBRange<<"\tsignificance = "<<significanceRange<<" +/- "<<errRange<<"\n";
380e2bf9 40
380e2bf9 41
42}