4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
9 //-------------------------------------------------------------------------
10 // ESD V0 Vertex Class
11 // This class is part of the Event Summary Data set of classes
12 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
13 // Modified by: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
14 // and Boris Hippolyte,IPHC, hippolyt@in2p3.fr
15 //-------------------------------------------------------------------------
19 #include "AliESDV0Params.h"
20 #include "AliExternalTrackParam.h"
21 #include "AliVParticle.h"
23 class AliESDv0 : public AliVParticle {
26 AliESDv0(const AliExternalTrackParam &t1, Int_t i1,
27 const AliExternalTrackParam &t2, Int_t i2);
29 AliESDv0(const AliESDv0&);
31 AliESDv0& operator=(const AliESDv0&);
32 virtual void Copy(TObject &obj) const;
34 // Start with AliVParticle functions
35 virtual Double_t Px() const { return fNmom[0]+fPmom[0]; }
36 virtual Double_t Py() const { return fNmom[1]+fPmom[1]; }
37 virtual Double_t Pz() const { return fNmom[2]+fPmom[2]; }
38 virtual Double_t Pt() const { return TMath::Sqrt(Px()*Px()+Py()*Py()); }
39 virtual Double_t P() const {
40 return TMath::Sqrt(Px()*Px()+Py()*Py()+Pz()*Pz());
42 virtual Bool_t PxPyPz(Double_t p[3]) const { p[0] = Px(); p[1] = Py(); p[2] = Pz(); return kTRUE; }
43 virtual Double_t Xv() const { return fPos[0]; }
44 virtual Double_t Yv() const { return fPos[1]; }
45 virtual Double_t Zv() const { return fPos[2]; }
46 virtual Bool_t XvYvZv(Double_t x[3]) const { x[0] = Xv(); x[1] = Yv(); x[2] = Zv(); return kTRUE; }
47 virtual Double_t OneOverPt() const { return (Pt() != 0.) ? 1./Pt() : -999.; }
48 virtual Double_t Phi() const {return TMath::Pi()+TMath::ATan2(-Py(),-Px()); }
49 virtual Double_t Theta() const {return 0.5*TMath::Pi()-TMath::ATan(Pz()/(Pt()+1.e-13)); }
50 virtual Double_t E() const; // default is KOs but can be changed via ChangeMassHypothesis (defined in the .cxx)
51 virtual Double_t M() const { return GetEffMass(); }
52 virtual Double_t Eta() const { return 0.5*TMath::Log((P()+Pz())/(P()-Pz()+1.e-13)); }
53 virtual Double_t Y() const { return 0.5*TMath::Log((E()+Pz())/(E()-Pz()+1.e-13)); }
54 virtual Short_t Charge() const { return 0; }
55 virtual Int_t GetLabel() const { return -1; } // temporary
56 virtual const Double_t *PID() const { return 0; } // return PID object ? (to be discussed!)
58 // Then the older functions
59 Double_t ChangeMassHypothesis(Int_t code=kK0Short);
61 Int_t GetPdgCode() const {return fPdgCode;}
62 Float_t GetEffMass(UInt_t p1, UInt_t p2);
63 Double_t GetEffMass() const {return fEffMass;}
64 Double_t GetChi2V0() const {return fChi2V0;}
65 void GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const;
66 void GetNPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const;
67 void GetPPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const;
68 void GetXYZ(Double_t &x, Double_t &y, Double_t &z) const;
69 Float_t GetD(Double_t x0=0.,Double_t y0=0.,Double_t z0=0.) const;
70 Int_t GetNindex() const {return fNidx;}
71 Int_t GetPindex() const {return fPidx;}
72 void SetDcaV0Daughters(Double_t rDcaV0Daughters=0.);
73 Double_t GetDcaV0Daughters() {return fDcaV0Daughters;}
74 Float_t GetV0CosineOfPointingAngle(Double_t&, Double_t&, Double_t&) const;
75 Double_t GetV0CosineOfPointingAngle() const {return fPointAngle;}
76 void SetV0CosineOfPointingAngle(Double_t cpa) {fPointAngle=cpa;}
77 void SetOnFlyStatus(Bool_t status){fOnFlyStatus=status;}
78 Bool_t GetOnFlyStatus() const {return fOnFlyStatus;}
79 const AliExternalTrackParam *GetParamP() const {return &fParamP;}
80 const AliExternalTrackParam *GetParamN() const {return &fParamN;}
84 // **** The following member functions need to be revised ***
86 void GetPosCov(Double_t cov[6])const ; // getter for the covariance matrix of the V0 position
87 Double_t GetSigmaY(); // sigma of y coordinate at vertex posistion
88 Double_t GetSigmaZ(); // sigma of z coordinate at vertex posistion
89 Double_t GetSigmaAP0(); // calculate sigma of Point angle resolution at vertex pos.
90 Double_t GetSigmaD0(); // calculate sigma of position resolution at vertex pos.
91 Double_t GetEffectiveSigmaAP0(); // calculate sigma of point angle resolution at vertex pos. effecive parameterization
92 Double_t GetEffectiveSigmaD0(); // calculate sigma of position resolution at vertex pos.
93 Double_t GetMinimaxSigmaAP0(); // calculate mini-max sigma of point angle resolution
94 Double_t GetMinimaxSigmaD0(); // calculate mini-max sigma of dca resolution
95 Double_t GetLikelihoodAP(Int_t mode0, Int_t mode1); // get likelihood for point angle
96 Double_t GetLikelihoodD(Int_t mode0, Int_t mode1); // get likelihood for DCA
97 Double_t GetLikelihoodC(Int_t mode0, Int_t mode1); // get likelihood for Causality
100 static const AliESDV0Params & GetParameterization(){return fgkParams;}
101 void SetParamP(const AliExternalTrackParam & paramP) {fParamP = paramP;}
102 void SetParamN(const AliExternalTrackParam & paramN) {fParamN = paramN;}
103 void SetStatus(Int_t status){fStatus=status;}
104 Int_t GetStatus() const {return fStatus;}
105 Int_t GetIndex(Int_t i) const {return (i==0) ? fNidx : fPidx;}
106 void SetIndex(Int_t i, Int_t ind) {(i==0) ? (fNidx=ind) : (fPidx=ind);}
107 Double_t *GetAnglep() {return fAngle;}
108 Double_t GetRr() const {return fRr;}
109 Double_t GetDistSigma() const {return fDistSigma;}
110 void SetDistSigma(Double_t ds) {fDistSigma=ds;}
111 Float_t GetChi2Before() const {return fChi2Before;}
112 void SetChi2Before(Float_t cb) {fChi2Before=cb;}
113 Float_t GetChi2After() const {return fChi2After;}
114 void SetChi2After(Float_t ca) {fChi2After=ca;}
115 Float_t GetNAfter() const {return fNAfter;}
116 void SetNAfter(Short_t na) {fNAfter=na;}
117 Short_t GetNBefore() const {return fNBefore;}
118 void SetNBefore(Short_t nb) {fNBefore=nb;}
119 void SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1);
120 const Double_t * GetCausalityP() const {return fCausality;}
121 void SetClusters(Int_t *clp, Int_t *clm);
122 const Int_t * GetClusters(Int_t i) const {return fClusters[i];}
123 void SetNormDCAPrim(Float_t nd0, Float_t nd1){fNormDCAPrim[0] = nd0; fNormDCAPrim[1]=nd1;}
124 const Double_t *GetNormDCAPrimP() const {return fNormDCAPrim;}
127 AliExternalTrackParam fParamN; // external parameters of negative particle
128 AliExternalTrackParam fParamP; // external parameters of positive particle
130 // CKBrev: tkink about revision
132 Double32_t fEffMass; // reconstructed V0's effective mass
133 Double32_t fDcaV0Daughters; // dca between V0's daughters
134 Double32_t fChi2V0; // V0's chi2 value
135 Double32_t fPos[3]; // V0's position (global)
136 Double32_t fPosCov[6]; // covariance matrix of the vertex position
137 Double32_t fNmom[3]; // momentum of the negative daughter (global)
138 Double32_t fPmom[3]; // momentum of the positive daughter (global)
139 Double32_t fNormDCAPrim[2]; // normalize distance to the primary vertex CKBrev
140 Double32_t fRr; //rec position of the vertex CKBrev
141 Double32_t fDistSigma; //sigma of distance CKBrev
142 Double32_t fChi2Before; //chi2 of the tracks before V0 CKBrev
143 Double32_t fChi2After; // chi2 of the tracks after V0 CKBrev
146 Double32_t fCausality[4]; //[0,1,8] causality information - see comments in SetCausality CKBrev
147 Double32_t fAngle[3]; //[-2*pi,2*pi,16]three angles CKBrev
148 Double32_t fPointAngleFi; //[-1,1,16]point angle fi CKBrev
149 Double32_t fPointAngleTh; //[-1,1,16]point angle theta CKBrev
150 Double32_t fPointAngle; //[-1,1,16] cosine of the pointing angle
153 Int_t fPdgCode; // reconstructed V0's type (PDG code)
154 Int_t fClusters[2][6]; //! its clusters CKBrev
155 Int_t fNidx; // index of the negative daughter
156 Int_t fPidx; // index of the positive daughter
160 Short_t fStatus; //status CKBrev
161 Short_t fNBefore; // number of possible points before V0 CKBrev
162 Short_t fNAfter; // number of possible points after V0 CKBrev
164 Bool_t fOnFlyStatus; // if kTRUE, then this V0 is recontructed
165 // "on fly" during the tracking
168 // parameterization coefficients
169 static const AliESDV0Params fgkParams; //! resolution and likelihood parameterization
173 ClassDef(AliESDv0,5) // ESD V0 vertex
177 void AliESDv0::GetNPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
178 px=fNmom[0]; py=fNmom[1]; pz=fNmom[2];
182 void AliESDv0::GetPPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
183 px=fPmom[0]; py=fPmom[1]; pz=fPmom[2];
187 void AliESDv0::SetDcaV0Daughters(Double_t rDcaV0Daughters){
188 fDcaV0Daughters=rDcaV0Daughters;