]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODv0.cxx
Handler for AOD input for analysis.
[u/mrichter/AliRoot.git] / STEER / AliAODv0.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //-------------------------------------------------------------------------
17 //     Analysis Oriented Data (AOD) V0 vertex class
18 //     Authors: B.Hippolyte, IReS, hippolyt@in2p3.fr 
19 //              G.Van Buren, BNL,  gene@bnl.gov      (original STAR MuDsts)
20 //-------------------------------------------------------------------------
21
22 //#include "AliESDEvent.h"
23 //#include "AliESDv0.h"
24
25 #include "AliAODv0.h"
26
27 //#include "AliAODTrack.h"
28
29 ClassImp(AliAODv0)
30
31   AliAODv0::AliAODv0() : 
32     AliAODRecoDecay(),
33     fDcaV0ToPrimVertex(999)
34 {
35   //--------------------------------------------------------------------
36   // Default constructor
37   //--------------------------------------------------------------------
38   fCharge  = 0;
39   fNProngs = 2;
40   fNDCA    = 1;
41   fNPID    = 2;
42
43   fDCA = new Double_t[fNDCA];
44   fDCA[0] = 999;
45
46   fPx = new Double_t[GetNProngs()];
47   fPy = new Double_t[GetNProngs()];
48   fPz = new Double_t[GetNProngs()];
49   fPx[0] = 999;
50   fPy[0] = 999;
51   fPz[0] = 999;
52
53   fPx[1] = 999;
54   fPy[1] = 999;
55   fPz[1] = 999;
56
57   fd0 = new Double_t[GetNProngs()];
58   fd0[0] = 999;
59   fd0[1] = 999;
60 }
61
62 //AliAODv0::AliAODv0(AliESDv0* rV0Vertex ,AliESDEvent* rEvent) :
63 //  AliAODRecoDecay(),
64 //  fDcaV0ToPrimVertex(999)
65 //{
66   //--------------------------------------------------------------------
67   // Constructor via fill to be removed eventually
68   //--------------------------------------------------------------------
69 //  fCharge  = 0;
70 //  fNProngs = 2;
71 //  fNDCA    = 1;
72 //  fNPID    = 2;
73
74 //  fDCA = new Double_t[fNDCA];
75
76 //   fPx = new Double_t[GetNProngs()];
77 //   fPy = new Double_t[GetNProngs()];
78 //   fPz = new Double_t[GetNProngs()];
79
80 //   fd0 = new Double_t[GetNProngs()];
81
82 //   this->Fill(rV0Vertex,rEvent);
83 // }
84
85 AliAODv0::AliAODv0(AliAODVertex* rAODVertex, Double_t rDcaV0Daughters, Double_t rDcaV0ToPrimVertex,
86            Double_t *rMomPos, Double_t *rMomNeg, Double_t *rDcaDaughterToPrimVertex) :
87   AliAODRecoDecay(rAODVertex,2,0,rDcaDaughterToPrimVertex),
88   fDcaV0ToPrimVertex(rDcaV0ToPrimVertex)
89 {
90   //--------------------------------------------------------------------
91   // Constructor via setting each data member
92   //--------------------------------------------------------------------
93   fCharge  = 0;
94   fNProngs = 2;
95   fNDCA    = 1;
96   fNPID    = 2;
97
98   fDCA = new Double_t[fNDCA];
99
100   fDCA[0] = rDcaV0Daughters;
101   fDcaV0ToPrimVertex = rDcaV0ToPrimVertex;
102
103   fPx = new Double_t[GetNProngs()];
104   fPy = new Double_t[GetNProngs()];
105   fPz = new Double_t[GetNProngs()];
106
107   fPx[0] = rMomPos[0] ;
108   fPy[0] = rMomPos[1];
109   fPz[0] = rMomPos[2];
110
111   fPx[1] = rMomNeg[0];
112   fPy[1] = rMomNeg[1];
113   fPz[1] = rMomNeg[2];
114 }
115
116 AliAODv0::AliAODv0(const AliAODv0& rAliAODv0) :
117   AliAODRecoDecay(rAliAODv0),
118   fDcaV0ToPrimVertex(rAliAODv0.fDcaV0ToPrimVertex)
119  {
120   //--------------------------------------------------------------------
121   // Copy constructor
122   //--------------------------------------------------------------------
123 }
124
125 AliAODv0& AliAODv0::operator=(const AliAODv0& rAliAODv0){
126   //--------------------------------------------------------------------
127   // Assignment overload
128   //--------------------------------------------------------------------
129   this->fDcaV0ToPrimVertex  = rAliAODv0.fDcaV0ToPrimVertex ;
130   return *this;
131 }
132
133 AliAODv0::~AliAODv0(){
134   //--------------------------------------------------------------------
135   // Empty destructor
136   //--------------------------------------------------------------------
137 }
138
139
140 // void AliAODv0::Fill(AliESDv0* rV0Vertex ,AliESDEvent* rEvent){
141
142 //   Double_t tDecayVertexV0[3]; rV0Vertex->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]); 
143 //   GetSecondaryVtx()->SetX(tDecayVertexV0[0]);
144 //   GetSecondaryVtx()->SetY(tDecayVertexV0[1]);
145 //   GetSecondaryVtx()->SetZ(tDecayVertexV0[2]);
146
147 //   Double_t lCovVtx[6];
148 //   rV0Vertex->GetPosCov(lCovVtx);
149 //   GetSecondaryVtx()->SetCovMatrix(lCovVtx);
150
151 //   GetSecondaryVtx()->SetChi2perNDF(rV0Vertex->GetChi2V0());
152 //   GetSecondaryVtx()->SetType(AliAODVertex::kV0);
153
154 //   UInt_t lKeyPos = (UInt_t)TMath::Abs(rV0Vertex->GetPindex());// need to ask why Abs
155 //   UInt_t lKeyNeg = (UInt_t)TMath::Abs(rV0Vertex->GetNindex());
156 //   GetSecondaryVtx()->AddDaughter(rEvent->GetTrack(lKeyPos));
157 //   GetSecondaryVtx()->AddDaughter(rEvent->GetTrack(lKeyNeg));
158
159 //   fDCA[0] = rV0Vertex->GetDcaV0Daughters();
160 //   fDcaV0ToPrimVertex = rV0Vertex->GetD();
161
162 //   Double_t tMomPos[3]; rV0Vertex->GetPPxPyPz(tMomPos[0],tMomPos[1],tMomPos[2]); 
163 //   fPx[0] = tMomPos[0];
164 //   fPy[0] = tMomPos[1];
165 //   fPz[0] = tMomPos[2];
166
167 //   Double_t tMomNeg[3]; rV0Vertex->GetNPxPyPz(tMomNeg[0],tMomNeg[1],tMomNeg[2]); 
168 //   fPx[1] = tMomNeg[0];
169 //   fPy[1] = tMomNeg[1];
170 //   fPz[1] = tMomNeg[2];
171
172 //   AliESDtrack *pTrack=rEvent->GetTrack(lKeyPos);
173 //   AliESDtrack *nTrack=rEvent->GetTrack(lKeyNeg);
174
175 //   Float_t tDcaPosToPrimVertex[2];
176 //   if(pTrack) pTrack->GetImpactParameters(tDcaPosToPrimVertex[0],tDcaPosToPrimVertex[1]);
177 //   else { tDcaPosToPrimVertex[0]=999.;  tDcaPosToPrimVertex[1]=999.;}
178 //   fd0[0] = TMath::Sqrt(tDcaPosToPrimVertex[0]*tDcaPosToPrimVertex[0]+tDcaPosToPrimVertex[1]*tDcaPosToPrimVertex[1]);
179
180 //   Float_t tDcaNegToPrimVertex[2];
181 //   if(nTrack) nTrack->GetImpactParameters(tDcaNegToPrimVertex[0],tDcaNegToPrimVertex[1]);
182 //   else { tDcaNegToPrimVertex[0]=999.;  tDcaNegToPrimVertex[1]=999.;}
183
184 //   fd0[1] = TMath::Sqrt(tDcaNegToPrimVertex[0]*tDcaNegToPrimVertex[0]+tDcaNegToPrimVertex[1]*tDcaNegToPrimVertex[1]);
185 // }
186
187 void AliAODv0::Fill(AliAODVertex *rAODVertex, Double_t rDcaV0Daughters, Double_t rDcaV0ToPrimVertex,
188                     Double_t *rMomPos, Double_t *rMomNeg, Double_t *rDcaDaughterToPrimVertex){
189
190   this->SetSecondaryVtx(rAODVertex);
191
192   fDCA[0] = rDcaV0Daughters;
193   fDcaV0ToPrimVertex = rDcaV0ToPrimVertex;
194
195   fPx[0] = rMomPos[0] ;
196   fPy[0] = rMomPos[1];
197   fPz[0] = rMomPos[2];
198
199   fPx[1] = rMomNeg[0];
200   fPy[1] = rMomNeg[1];
201   fPz[1] = rMomNeg[2];
202
203   fd0[0] = rDcaDaughterToPrimVertex[0];
204   fd0[1] = rDcaDaughterToPrimVertex[1];
205 }
206
207 void AliAODv0::ResetV0(){
208
209   GetSecondaryVtx()->SetX(999);
210   GetSecondaryVtx()->SetY(999);
211   GetSecondaryVtx()->SetZ(999);
212   GetSecondaryVtx()->SetChi2perNDF(999);
213   GetSecondaryVtx()->SetType(AliAODVertex::kUndef);
214
215   Int_t lNumDaughters = GetSecondaryVtx()->GetNDaughters();
216   for(Int_t iDaughterIndex = 0; iDaughterIndex<lNumDaughters;iDaughterIndex++)
217     GetSecondaryVtx()->RemoveDaughter(GetSecondaryVtx()->GetDaughter(iDaughterIndex));
218
219   fDCA[0] = 999;
220   fDcaV0ToPrimVertex  = 999;
221
222   fPx[0] = 999;
223   fPy[0] = 999;
224   fPz[0] = 999;
225
226   fPx[1] = 999;
227   fPy[1] = 999;
228   fPz[1] = 999;
229
230   fd0[0] = 999;
231   fd0[1] = 999;
232 }
233
234 void AliAODv0::Print(Option_t* /*option*/) const {
235   //
236   // Print some information
237   //
238   AliAODRecoDecay::Print();
239   printf("AliAODv0: invariant mass (k0s %.6f, lambda %.6f, anti-lambda %.6f) \n",MassK0Short(),MassLambda(),MassAntiLambda());
240   printf("AliAODv0: dca (v0d %.6f, v0tpv %.6f, postpv %.6f, negtpv %.6f ) \n",DcaV0Daughters(),DcaV0ToPrimVertex(),DcaPosToPrimVertex(),DcaNegToPrimVertex());
241   printf("AliAODv0: mom (ptot2 %.6f, pt2 %.6f, rapk0 %.6f, rapla %.6f ) \n",Ptot2V0(),Pt2V0(),RapK0Short(),RapLambda());
242   printf("AliAODv0: cin (mpav0 %.6f, mnav0 %.6f, alpha %.6f, ptarm %.6f ) \n",MomPosAlongV0(),MomNegAlongV0(),AlphaV0(),PtArmV0());
243   printf("AliAODv0: nrg (eppro %.6f, enpro %.6f, eppio %.6f, enpio %.6f ) \n",EPosProton(),ENegProton(),EPosPion(),ENegPion());
244
245   return;
246 }