]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/AliBtoJPSItoEle.h
Obsolete CDB entry. The trigger confgiuration is now stored in the GRP/CTP/Config...
[u/mrichter/AliRoot.git] / PWG3 / AliBtoJPSItoEle.h
1 #ifndef AliBtoJPSItoEle_H
2 #define AliBtoJPSItoEle_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 //-------------------------------------------------------------------------
7 //                          Class AliBtoJPSItoEle
8 //                 Reconstructed B -> J\PSI+X --> e+ e- *X  class
9 //      
10 // Note: the two decay tracks are labelled: 0 (positive electron)
11 //                                          1 (negative electron)
12 //
13 //         Origin: G.E. Bruno    giuseppe.bruno@ba.infn.it                  
14 //          based on Class for charm golden channel (D0->Kpi) (thanks to Andrea Dainese!)
15 //-------------------------------------------------------------------------
16
17 #include <TObject.h>
18 #include <TMath.h>
19
20 //----------------------------------------------------------------------------
21 //     Some constants (masses)
22 //
23 // particle masses
24 const Double_t kMJPsi = 3.096916;     // J/Psi mass
25 const Double_t kMe    = 0.000510998902;  // elec  mass
26
27
28
29
30 //-----------------------------------------------------------------------------
31 class AliBtoJPSItoEle : public TObject {
32  public:
33   //
34   AliBtoJPSItoEle();
35   AliBtoJPSItoEle(Int_t ev,Int_t trkNum[2],  
36              Double_t v1[3],Double_t v2[3],Double_t dca,
37              Double_t mom[6],Double_t d0[2]);
38   virtual ~AliBtoJPSItoEle();
39
40   Double_t Alpha() const { return (Ql(0)-Ql(1))/(Ql(0)+Ql(1)); }
41   void     ApplyPID(TString pidScheme="TRDTPCparam");
42   Double_t ChildrenRelAngle() const; 
43   void     ComputeWgts();
44   void     CorrectWgt4BR(Double_t factor);
45   Double_t CosPointing() const;
46   Double_t CosPointingXY() const;
47   void     CosThetaStar(Double_t &ctsJPSI) const;
48   Double_t Ct() const {return Length()*kMJPsi/P();}
49   Double_t Energy() const { return TMath::Sqrt(P()*P()+kMJPsi*kMJPsi); }
50   Double_t Eta() const;
51   Double_t EtaChild(Int_t child) const;
52   Int_t    EventNo() const {return TMath::Abs(fEvent);}
53   Double_t GetDCA() const { return 10000.*fDCA; }
54   Int_t    GetTrkNum(Int_t child) const { return fTrkNum[child]; }
55   Double_t Getd0Child(Int_t child) const { return fd0[child]; }
56   Int_t    GetPdgChild(Int_t child) const { return fPdg[child]; }
57   Int_t    GetPdgMum(Int_t child) const {return fMum[child]; }
58   Int_t    GetPdgGMum(Int_t child) const {return fGMum[child]; }
59   void     GetPIDresponse(Double_t resp0[5],Double_t resp1[5]) const; 
60   void     GetWgts(Double_t &WgtJPsi) const;
61   void     GetPrimaryVtx(Double_t vtx[3]) const 
62     { vtx[0]=fV1x; vtx[1]=fV1y; vtx[2]=fV1z; }
63   void     GetSecondaryVtx(Double_t vtx[3]) const 
64     { vtx[0]=fV2x; vtx[1]=fV2y; vtx[2]=fV2z; }
65
66   Double_t ImpPar() const;
67   void     InvMass(Double_t &mJPsi) const;
68   Bool_t   IsSignal() const { if(fSignal) return kTRUE; return kFALSE; } 
69   Bool_t   IsJpsiPrimary() const { if(fJpsiPrimary) return kTRUE; return kFALSE; } 
70   Double_t Length() const
71     { return TMath::Sqrt((fV1x-fV2x)*(fV1x-fV2x)
72                          +(fV1y-fV2y)*(fV1y-fV2y)+(fV1z-fV2z)*(fV1z-fV2z)); }
73   Double_t P()  const { return TMath::Sqrt(Pt()*Pt()+Pz()*Pz()); } 
74   Double_t PChild(Int_t child)   const { return TMath::Sqrt(fPx[child]*fPx[child]+fPy[child]*fPy[child]+fPz[child]*fPz[child]); } 
75   Double_t ProdImpParams() const { return fd0[0]*fd0[1]; } 
76   Double_t Pt() const { return TMath::Sqrt(Px()*Px()+Py()*Py()); }
77   Double_t PtChild(Int_t child)  const { return TMath::Sqrt(fPx[child]*fPx[child]+fPy[child]*fPy[child]); } 
78   Double_t Px() const { return (fPx[0]+fPx[1]); }
79   Double_t Py() const { return (fPy[0]+fPy[1]); }
80   Double_t Pz() const { return (fPz[0]+fPz[1]); }
81   Double_t Ql(Int_t child) const;
82   Double_t Qt() const;
83   Double_t Rapidity() const { return 0.5*TMath::Log((Energy()+Pz())/(Energy()-Pz()+1.e-13)); }
84   Bool_t   Select(const Double_t* cuts,Int_t& okB) const;
85   void     SetPrimaryVtx(Double_t vtx[3]) 
86     { fV1x=vtx[0]; fV1y=vtx[1]; fV1z=vtx[2]; }
87   void     SetSignal() { fSignal =  kTRUE; }
88   void     SetJpsiPrimary() { fJpsiPrimary =  kTRUE; }
89   void     SetTOFmasses(Double_t mass[2]) 
90     { fTOFmass[0]=mass[0]; fTOFmass[1]=mass[1]; }
91   void     SetPIDresponse(Double_t resp0[5],Double_t resp1[5]); 
92   void     SetPdgCodes(Int_t pdg[2]) {fPdg[0]=pdg[0];fPdg[1]=pdg[1]; }
93   void     SetMumPdgCodes(Int_t mum[2]) {fMum[0]=mum[0];fMum[1]=mum[1]; }
94   void     SetGMumPdgCodes(Int_t gmum[2]) {fGMum[0]=gmum[0];fGMum[1]=gmum[1]; }
95   Double_t TRDTPCCombinedPIDParametrization(Double_t p) const;
96   //
97  private:
98   //
99   Bool_t   fSignal; // TRUE if signal, FALSE if background (for simulation) 
100                     // (background are both combinatorial and primary J/psi)
101   Bool_t   fJpsiPrimary; // TRUE if the current candidate is a primary J/psi, FALSE otherway (for simulation) 
102   Int_t    fEvent;  // number of the event this B comes from
103                     // -1 if the B comes from ev. mixing
104
105   Int_t fTrkNum[2]; // numbers of the two decay tracks
106
107   Double_t fV1x; // X-position of the primary vertex of the event
108   Double_t fV1y; // Y-position of the primary vertex of the event
109   Double_t fV1z; // Z-position of the primary vertex of the event
110   Double_t fV2x; // X-position of the reconstructed secondary vertex
111   Double_t fV2y; // Y-position of the reconstructed secondary vertex
112   Double_t fV2z; // Z-position of the reconstructed secondary vertex
113   Double_t fDCA; // DCA of the two tracks
114
115   Double_t fPx[2];  // X,Y,Z
116   Double_t fPy[2];  // momenta of the two tracks
117   Double_t fPz[2];  // at the reconstructed vertex  
118
119   Double_t fd0[2];  //  impact parameters in the bending plane
120
121   Int_t fPdg[2];  // PDG codes of the two tracks       (for sim.)
122   Int_t fMum[2];  // PDG codes of the mothers          (for sim.)
123   Int_t fGMum[2]; // PDG codes of the grand-mothers    (for sim.)
124
125   Double_t fTagEl[2];  // probability to be tagged as electron 
126   Double_t fTagPi[2];  // probability to be tagged as pion 
127   Double_t fTagKa[2];  // probability to be tagged as kaon 
128   Double_t fTagPr[2];  // probability to be tagged as proton 
129   Double_t fTagNid[2]; // probability to be tagged as "non-identified" 
130
131   Double_t fPIDrespEl[2]; // det. response to be electron
132   Double_t fPIDrespMu[2]; // det. response to be muon
133   Double_t fPIDrespPi[2]; // det. response to be pion
134   Double_t fPIDrespKa[2]; // det. response to be kaon
135   Double_t fPIDrespPr[2]; // det. response to be proton
136   Double_t fTOFmass[2]; // mass estimated by the TOF (-1000. if track not reached TOF)
137   Double_t fWgtJPsi; // weights for the pair
138
139   ClassDef(AliBtoJPSItoEle,1)  // Reconstructed B candidate class
140 };
141
142 #endif
143
144
145
146
147
148
149
150