]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliBtoJPSItoEleCDFfitFCN.h
I improved (by a factor 2.5) the speed of the MatchToMC method
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliBtoJPSItoEleCDFfitFCN.h
1 #ifndef ALIBTOJPSITOELECDFFITFCN_H\r
2 #define ALIBTOJPSITOELECDFFITFCN_H\r
3 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *\r
4  * See cxx source for full Copyright notice                               */\r
5 \r
6 //_________________________________________________________________________\r
7 //                      Class AliBtoJPSItoEleCDFfitFCN\r
8 //                    Definition of main function used in \r
9 //                     unbinned log-likelihood fit for\r
10 //                 the channel B -> JPsi + X -> e+e- + X\r
11 //      \r
12 //                          Origin: C.Di Giglio\r
13 //     Contact: Carmelo.Digiglio@ba.infn.it , Giuseppe.Bruno@ba.infn.it\r
14 //_________________________________________________________________________\r
15 \r
16 #include <TNamed.h>\r
17 #include <TDatabasePDG.h>\r
18 #include "TH1F.h"\r
19 #include "TMath.h"\r
20 \r
21 enum IntegralType {kBkg, \r
22                    kBkgNorm, \r
23                    kSig, \r
24                    kSigNorm};\r
25 \r
26 enum PtBins       {kallpt, \r
27                    kptbin1,kptbin2,kptbin3,\r
28                    kptbin4,kptbin5,kptbin6,\r
29                    kptbin7,kptbin8,kptbin9};\r
30 //_________________________________________________________________________________________________\r
31 class AliBtoJPSItoEleCDFfitFCN : public TNamed {\r
32  public:\r
33   //\r
34   AliBtoJPSItoEleCDFfitFCN();\r
35   AliBtoJPSItoEleCDFfitFCN(const AliBtoJPSItoEleCDFfitFCN& source); \r
36   AliBtoJPSItoEleCDFfitFCN& operator=(const AliBtoJPSItoEleCDFfitFCN& source);  \r
37   virtual ~AliBtoJPSItoEleCDFfitFCN();\r
38 \r
39   Double_t EvaluateLikelihood(const Double_t* pseudoproperdecaytime,\r
40                               const Double_t* invariantmass, const Int_t ncand);\r
41  \r
42   Double_t GetFPlus() const {return fFPlus;}\r
43   Double_t GetFMinus() const {return fFMinus;}\r
44   Double_t GetFSym() const {return fFSym;}\r
45   Double_t GetRadius() const {return fParameters[0];}\r
46   Double_t GetTheta() const {return fParameters[1];}\r
47   Double_t GetPhi() const {return fParameters[2];}\r
48   Double_t GetLamPlus() const {return fParameters[3];}\r
49   Double_t GetLamMinus() const {return fParameters[4];}\r
50   Double_t GetLamSym() const {return fParameters[5];}\r
51   Double_t GetMassSlope() const {return fParameters[6];}\r
52   Double_t GetFractionJpsiFromBeauty() const {return fParameters[7];}\r
53   Double_t GetFsig() const {return fParameters[8];}\r
54   Double_t GetCrystalBallMmean() const {return fParameters[9];}\r
55   Double_t GetCrystalBallNexp() const {return fParameters[10];}\r
56   Double_t GetCrystalBallSigma() const {return fParameters[11];}\r
57   Double_t GetCrystalBallAlpha() const {return fParameters[12];}\r
58   Double_t GetIntegral() const {return fIntegral;}\r
59   Bool_t GetCrystalBallParam() const {return fCrystalBallParam;}\r
60 \r
61   void SetFPlus(Double_t plus) {fFPlus = plus;}\r
62   void SetFMinus(Double_t minus) {fFMinus = minus;}\r
63   void SetFSym(Double_t sym) {fFSym = sym;}\r
64   void SetRadius(Double_t radius) {fParameters[0] = radius;}\r
65   void SetTheta(Double_t theta) {fParameters[1] = theta;}\r
66   void SetPhi(Double_t phi) {fParameters[2] = phi;}\r
67   void SetLamPlus(Double_t lamplus) {fParameters[3] = lamplus;}\r
68   void SetLamMinus(Double_t lamminus) {fParameters[4] = lamminus;}\r
69   void SetLamSym(Double_t lamsym) {fParameters[5] = lamsym;}\r
70   void SetMassSlope(Double_t slope) {fParameters[6] = slope;}\r
71   void SetFractionJpsiFromBeauty(Double_t B) {fParameters[7] = B;}\r
72   void SetFsig(Double_t Fsig) {fParameters[8] = Fsig;}\r
73   void SetCrystalBallMmean(Double_t CrystalBallMmean) {fParameters[9] = CrystalBallMmean;}\r
74   void SetCrystalBallNexp(Double_t CrystalBallNexp) {fParameters[10] = CrystalBallNexp;}\r
75   void SetCrystalBallSigma(Double_t CrystalBallSigma) {fParameters[11] = CrystalBallSigma;}\r
76   void SetCrystalBallAlpha(Double_t CrystalBallAlpha) {fParameters[12] = CrystalBallAlpha;}\r
77 \r
78   void SetAllParameters(const Double_t* parameters);\r
79   void SetIntegral(Double_t integral) {fIntegral = integral;}\r
80   void SetCsiMC(const TH1F* MCtemplate) {fhCsiMC = (TH1F*)MCtemplate->Clone("fhCsiMC");}\r
81   void SetResolutionConstants(Int_t BinNum);\r
82   void SetMassWndHigh(Double_t limit) { fMassWndHigh = TDatabasePDG::Instance()->GetParticle(443)->Mass() + limit ;}//here use pdg code instead\r
83   void SetMassWndLow(Double_t limit) { fMassWndLow = TDatabasePDG::Instance()->GetParticle(443)->Mass() - limit ;}//here use pdg code instead\r
84   void SetCrystalBallParam(Bool_t okCB) {fCrystalBallParam = okCB;}\r
85 \r
86   void ConvertFromSpherical() { fFPlus  = TMath::Power((fParameters[0]*TMath::Cos(fParameters[1])),2.);\r
87                                 fFMinus = TMath::Power((fParameters[0]*TMath::Sin(fParameters[1])*TMath::Sin(fParameters[2])),2.);\r
88                                 fFSym   = TMath::Power((fParameters[0]*TMath::Sin(fParameters[1])*TMath::Cos(fParameters[2])),2.);} \r
89 \r
90   void ComputeIntegral(); \r
91 \r
92   void ReadMCtemplates(Int_t BinNum);\r
93 \r
94   void PrintStatus();\r
95 \r
96  private:  \r
97   //\r
98   Double_t fParameters[13];        /*  par[0]  = fRadius;                \r
99                                        par[1]  = fTheta;\r
100                                        par[2]  = fPhi;\r
101                                        par[3]  = fOneOvLamPlus;\r
102                                        par[4]  = fOneOvLamMinus;\r
103                                        par[5]  = fOneOvLamSym;\r
104                                        par[6]  = fMassBkgSlope;\r
105                                        par[7]  = fFractionJpsiFromBeauty;\r
106                                        par[8]  = fFsig;\r
107                                        par[9]  = fCrystalBallMmean;\r
108                                        par[10] = fCrystalBallNexp;\r
109                                        par[11] = fCrystalBallSigma;\r
110                                        par[12] = fCrystalBallAlpha;*/\r
111 \r
112   Double_t fFPlus;                  // parameters of the log-likelihood function\r
113   Double_t fFMinus;                 // Slopes of the x distributions of the background \r
114   Double_t fFSym;                   // functions \r
115 \r
116   Double_t fIntegral;               // integral values of log-likelihood function terms\r
117   TH1F *fhCsiMC;                    // X distribution used as MC template for JPSI from B\r
118   Double_t fResolutionConstants[7]; // constants for the parametrized resolution function R(X)\r
119   Double_t fMassWndHigh;            // JPSI Mass window higher limit\r
120   Double_t fMassWndLow;             // JPSI Mass window lower limit\r
121   Bool_t fCrystalBallParam;         // Boolean to switch to Crystall Ball parameterisation\r
122 \r
123   ////\r
124 \r
125   Double_t EvaluateCDFfunc(Double_t x, Double_t m) const ;\r
126 \r
127   Double_t EvaluateCDFfuncNorm(Double_t x, Double_t m) const ;\r
128 \r
129   ////\r
130 \r
131   Double_t EvaluateCDFfuncSignalPart(Double_t x, Double_t m) const ;      // Signal part \r
132 \r
133   Double_t EvaluateCDFDecayTimeSigDistr(Double_t x) const ;\r
134 \r
135   Double_t EvaluateCDFInvMassSigDistr(Double_t m) const ;\r
136 \r
137   Double_t EvaluateCDFfuncBkgPart(Double_t x,Double_t m) const ;          // Background part\r
138 \r
139   Double_t EvaluateCDFDecayTimeBkgDistr(Double_t x) const ;\r
140  \r
141   Double_t EvaluateCDFInvMassBkgDistr(Double_t m) const ;\r
142 \r
143   ////\r
144 \r
145   Double_t FunB(Double_t x) const ;\r
146 \r
147   Double_t FunP(Double_t x) const ;\r
148 \r
149   Double_t CsiMC(Double_t x) const ;\r
150 \r
151   Double_t FunBkgPos(Double_t x) const ;\r
152 \r
153   Double_t FunBkgNeg(Double_t x) const ;\r
154 \r
155   Double_t FunBkgSym(Double_t x) const ;\r
156 \r
157   Double_t ResolutionFunc(Double_t x) const ;\r
158   //\r
159   ClassDef (AliBtoJPSItoEleCDFfitFCN,0);         // Unbinned log-likelihood fit \r
160 \r
161 };\r
162 \r
163 #endif\r