]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliHFMassFitter.h
Added like-sign 3Prong
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliHFMassFitter.h
CommitLineData
fabf4d8e 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
24class 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(){if(fhistoInvMass) delete fhistoInvMass; if(fntuParam) delete fntuParam;}
30
31
32 //setters
33 void SetHisto(TH1F *histoToFit){fhistoInvMass=histoToFit;};
34 void SetRangeFit(Double_t minvalue, Double_t maxvalue){fminMass=minvalue; fmaxMass=maxvalue;};
35 void SetMinRangeFit(Double_t minvalue){fminMass=minvalue;};
36 void SetMaxRangeFit(Double_t maxvalue){fmaxMass=maxvalue;};
37 void SetBinN(Int_t newbinN){fNbin=newbinN;};
38 void SetType(Int_t fittypeb, Int_t fittypes) {ftypeOfFit4Bkg = fittypeb; ftypeOfFit4Sgn = fittypes; };
39 void SetReflectionSigmaFactor(Int_t constant) {ffactor=constant;}
40 void SetInitialGaussianMean(Double_t mean) {fMass=mean;}
41 void SetInitialGaussianSigma(Double_t sigma) {fSigmaSgn=sigma;}
42 void SetSideBands(Bool_t onlysidebands=kTRUE) {fSideBands=onlysidebands;}
43
44 //getters
45 TH1F* GetHisto() const {return fhistoInvMass;};
46 Double_t GetMinRangeFit()const {return fminMass;};
47 Double_t GetMaxRangeFit()const {return fmaxMass;};
48 Int_t GetBinN() const {return fNbin;};
49 void GetFitPars(Float_t*) const;
50 void GetTypeOfFit(Bool_t background, Int_t typeb) const {background = fWithBkg; typeb = ftypeOfFit4Bkg;}
51 Int_t GetReflectionSigmaFactor() const {return ffactor;}
52
53 void InitNtuParam(char *ntuname="ntupar");
54 void FillNtuParam();
55 TNtuple* GetNtuParam() const {return fntuParam;}
56 TNtuple* NtuParamOneShot(char *ntuname="ntupar");
57
58 void Signal(Double_t nOfSigma,Double_t &signal,Double_t &errsignal) const;
59 void Background(Double_t nOfSigma,Double_t &background,Double_t &errbackground) const;
60 void Significance(Double_t nOfSigma,Double_t &significance,Double_t &errsignificance) const;
61
62 Double_t FitFunction4MassDistr (Double_t*, Double_t*);
63 Double_t FitFunction4Sgn (Double_t*, Double_t*);
64 Double_t FitFunction4Bkg (Double_t*, Double_t*);
65 void MassFitter(Bool_t draw=kTRUE);
66 void RebinMass(Int_t binground=1);
67
68
69 private:
70
71 TH1F* fhistoInvMass; // histogram to fit
72 Double_t fminMass; // lower mass limit
73 Double_t fmaxMass; // upper mass limit
74 Int_t fNbin; // number of bins
75 Float_t fFitPars[30]; // array of fit parameters
76 Bool_t fWithBkg; // signal+background (kTRUE) or signal only (kFALSE)
77 Int_t ftypeOfFit4Bkg; // 0 = exponential; 1 = linear; 2 = pol2
78 Int_t ftypeOfFit4Sgn; // 0 = gaus; 1 = gaus+gaus broadened
79 Int_t ffactor; // number to multiply to the sigma of the signal to obtain the reflected gaussian
80 TNtuple *fntuParam; // contains fit parameters
81 Double_t fMass; // signal gaussian mean value
82 Double_t fSigmaSgn; // signal gaussian sigma
83 Bool_t fSideBands; // kTRUE = only side bands considered
84
85 ClassDef(AliHFMassFitter,0); // class for invariant mass fit
86};
87
88#endif
89