]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliHFMassFitter.h
Fix in copy constructor, added draw functionality (Chiara)
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliHFMassFitter.h
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 <Riostream.h>
16 #include <TH1F.h>
17 #include <TF1.h>
18 #include <TMath.h>
19 #include <TNtuple.h>
20 #include <TFile.h>
21 #include <TString.h>
22
23
24 class AliHFMassFitter : public TNamed {
25
26  public:
27   AliHFMassFitter();
28   AliHFMassFitter(TH1F* histoToFit, Double_t minvalue, Double_t maxvalue, Int_t rebin=1,Int_t fittypeb=0,Int_t fittypes=0);
29   virtual ~AliHFMassFitter();
30
31   AliHFMassFitter(const AliHFMassFitter &mfit);
32   AliHFMassFitter& operator=(const AliHFMassFitter &mfit);
33
34   //setters
35   void     SetHisto(TH1F *histoToFit) {fhistoInvMass=histoToFit;}
36   void     SetRangeFit(Double_t minvalue, Double_t maxvalue){fminMass=minvalue; fmaxMass=maxvalue;}
37   void     SetMinRangeFit(Double_t minvalue){fminMass=minvalue;}
38   void     SetMaxRangeFit(Double_t maxvalue){fmaxMass=maxvalue;}
39   void     SetBinN(Int_t newbinN){fNbin=newbinN;}
40   void     SetType(Int_t fittypeb, Int_t fittypes);
41   void     SetReflectionSigmaFactor(Int_t constant) {ffactor=constant;}
42   void     SetInitialGaussianMean(Double_t mean) {fMass=mean;}
43   void     SetInitialGaussianSigma(Double_t sigma) {fSigmaSgn=sigma;}
44   void     SetSideBands(Bool_t onlysidebands=kTRUE) {fSideBands=onlysidebands;}
45
46   //getters
47   TH1F*    GetHisto()      const {return fhistoInvMass;}
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*) const;
52   void     GetTypeOfFit(Bool_t background, Int_t typeb) const {background = fWithBkg; typeb = ftypeOfFit4Bkg;}
53   Int_t    GetReflectionSigmaFactor() const {return ffactor;} 
54
55   void     InitNtuParam(char *ntuname="ntupar");
56   void     FillNtuParam();
57   TNtuple* GetNtuParam() const {return fntuParam;}
58   TNtuple* NtuParamOneShot(char *ntuname="ntupar");
59
60   void     Signal(Double_t nOfSigma,Double_t &signal,Double_t &errsignal) const; 
61   void     Background(Double_t nOfSigma,Double_t &background,Double_t &errbackground) const; 
62   void     Significance(Double_t nOfSigma,Double_t &significance,Double_t &errsignificance) const;
63
64   Double_t FitFunction4MassDistr (Double_t*, Double_t*);
65   Double_t FitFunction4Sgn (Double_t*, Double_t*);
66   Double_t FitFunction4Bkg (Double_t*, Double_t*);
67   void     MassFitter(Bool_t draw=kTRUE);
68   void     RebinMass(Int_t binground=1);
69   
70
71  private:
72
73   void     ComputeParSize();
74
75   TH1F    *fhistoInvMass;  // histogram to fit
76   Double_t fminMass;       // lower mass limit
77   Double_t fmaxMass;       // upper mass limit
78   Int_t    fNbin;          // number of bins
79   Int_t    fParsSize;      // size of fFitPars array
80   Float_t *fFitPars;       //[fParsSize] array of fit parameters
81   Bool_t   fWithBkg;       // signal+background (kTRUE) or signal only (kFALSE)
82   Int_t    ftypeOfFit4Bkg; // 0 = exponential; 1 = linear; 2 = pol2
83   Int_t    ftypeOfFit4Sgn; // 0 = gaus; 1 = gaus+gaus broadened
84   Int_t    ffactor;         // number to multiply to the sigma of the signal to obtain the reflected gaussian
85   TNtuple *fntuParam;      // contains fit parameters
86   Double_t fMass;          // signal gaussian mean value
87   Double_t fSigmaSgn;      // signal gaussian sigma
88   Bool_t   fSideBands;     // kTRUE = only side bands considered
89
90   ClassDef(AliHFMassFitter,1); // class for invariant mass fit
91 };
92
93 #endif
94
95