]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PWG3/vertexingHF/AliHFMassFitter.h
Bug fix
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliHFMassFitter.h
... / ...
CommitLineData
1#ifndef ALIHFMASSFITTER_H
2#define ALIHFMASSFITTER_H
3/* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
5
6/////////////////////////////////////////////////////////////
7//
8// AliHFMassFitter for the fit of invariant mass distribution
9// of charmed mesons
10//
11// Author: C.Bianchin, chiara.bianchin@pd.infn.it
12/////////////////////////////////////////////////////////////
13
14#include <TNamed.h>
15#include <TString.h>
16
17class TF1;
18class TNtuple;
19class TFile;
20class TList;
21class TH1F;
22
23class AliHFMassFitter : public TNamed {
24
25 public:
26 AliHFMassFitter();
27 AliHFMassFitter(const TH1F* histoToFit, Double_t minvalue, Double_t maxvalue, Int_t rebin=1,Int_t fittypeb=0,Int_t fittypes=0);
28 virtual ~AliHFMassFitter();
29
30 AliHFMassFitter(const AliHFMassFitter &mfit);
31 AliHFMassFitter& operator=(const AliHFMassFitter &mfit);
32
33 //setters
34 void SetHisto(const TH1F *histoToFit);
35 void SetRangeFit(Double_t minvalue, Double_t maxvalue){fminMass=minvalue; fmaxMass=maxvalue;}
36 void SetMinRangeFit(Double_t minvalue){fminMass=minvalue;}
37 void SetMaxRangeFit(Double_t maxvalue){fmaxMass=maxvalue;}
38 void SetBinN(Int_t newbinN){fNbin=newbinN;}
39 void SetType(Int_t fittypeb, Int_t fittypes);
40 void SetReflectionSigmaFactor(Int_t constant) {ffactor=constant;}
41 void SetInitialGaussianMean(Double_t mean) {fMass=mean;} // change the default value of the mean
42 void SetInitialGaussianSigma(Double_t sigma) {fSigmaSgn=sigma;} // change the default value of the sigma
43 void SetSideBands(Bool_t onlysidebands=kTRUE) {fSideBands=onlysidebands;} // consider only side bands
44
45 //getters
46 TH1F* GetHistoClone() const; //return the histogram
47 void GetRangeFit(Double_t &minvalue, Double_t &maxvalue) const {minvalue=fminMass; maxvalue=fmaxMass;}
48 Double_t GetMinRangeFit()const {return fminMass;}
49 Double_t GetMaxRangeFit()const {return fmaxMass;}
50 Int_t GetBinN() const {return fNbin;}
51 void GetFitPars(Float_t* pars) const;
52 void GetTypeOfFit(Bool_t &background, Int_t &typeb) const {background = fWithBkg; typeb = ftypeOfFit4Bkg;}
53 Int_t GetReflectionSigmaFactor() const {return ffactor;}
54 Double_t GetMean() const {return fMass;}
55 Double_t GetSigma()const {return fSigmaSgn;}
56 Double_t GetChiSquare() const;
57 Double_t GetReducedChiSquare() const;
58 void GetSideBandsBounds(Int_t& lb, Int_t& hb) const;
59
60 void PrintParTitles() const;
61
62 void InitNtuParam(char *ntuname="ntupar"); // initialize TNtuple to store the parameters
63 void FillNtuParam(); //Fill the TNtuple with the current parameters
64 TNtuple* GetNtuParam() const {return fntuParam;} // return the TNtuple
65 TNtuple* NtuParamOneShot(char *ntuname="ntupar"); // the three functions above all together
66 void WriteHisto(TString path="./") const; // write the histogram
67 void WriteNtuple(TString path="./") const; // write the TNtuple
68 void DrawFit() const;
69 void Reset();
70
71 void IntS(Float_t *valuewitherror) const; // integral of signal given my the fit with error
72 Double_t IntTot() const {return fhistoInvMass->Integral("width");} // return total integral of the histogram
73 void Signal(Double_t nOfSigma,Double_t &signal,Double_t &errsignal) const; // signal in nsigma with error
74 void Signal(Double_t min,Double_t max,Double_t &signal,Double_t &errsignal) const; // signal in (min, max) with error
75 void Background(Double_t nOfSigma,Double_t &background,Double_t &errbackground) const; // backgournd in nsigma with error
76 void Background(Double_t min,Double_t max,Double_t &background,Double_t &errbackground) const; // backgournd in (min, max) with error
77 void Significance(Double_t nOfSigma,Double_t &significance,Double_t &errsignificance) const; // significance in nsigma with error
78 void Significance(Double_t min,Double_t max,Double_t &significance,Double_t &errsignificance) const; // significance in (min, max) with error
79
80 Double_t FitFunction4MassDistr (Double_t* x, Double_t* par);
81 Double_t FitFunction4Sgn (Double_t* x, Double_t* par);
82 Double_t FitFunction4Bkg (Double_t* x, Double_t* par);
83 Bool_t MassFitter(Bool_t draw=kTRUE);
84 void RebinMass(Int_t binground=1);
85
86
87 private:
88
89 void ComputeParSize();
90 Bool_t SideBandsBounds();
91 void AddFunctionsToHisto();
92
93 TH1F *fhistoInvMass; // histogram to fit
94 Double_t fminMass; // lower mass limit
95 Double_t fmaxMass; // upper mass limit
96 Int_t fNbin; // number of bins
97 Int_t fParsSize; // size of fFitPars array
98 Float_t *fFitPars; //[fParsSize] array of fit parameters
99 Bool_t fWithBkg; // signal+background (kTRUE) or signal only (kFALSE)
100 Int_t ftypeOfFit4Bkg; // 0 = exponential; 1 = linear; 2 = pol2
101 Int_t ftypeOfFit4Sgn; // 0 = gaus; 1 = gaus+gaus broadened
102 Int_t ffactor; // number to multiply to the sigma of the signal to obtain the reflected gaussian
103 TNtuple *fntuParam; // contains fit parameters
104 Double_t fMass; // signal gaussian mean value
105 Double_t fSigmaSgn; // signal gaussian sigma
106 Bool_t fSideBands; // kTRUE = only side bands considered
107 Int_t fSideBandl; // left side band limit (bin number)
108 Int_t fSideBandr; // right side band limit (bin number)
109 Int_t fcounter; // internal counter
110 TList *fContourGraph; // TList of TGraph containing contour plots
111 ClassDef(AliHFMassFitter,3); // class for invariant mass fit
112};
113
114#endif
115
116