]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliV0vertexer.cxx
Fixing memory leaks
[u/mrichter/AliRoot.git] / ITS / AliV0vertexer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //               Implementation of the V0 vertexer class
18 //                  reads tracks writes out V0 vertices
19 //                      fills the ESD with the V0s       
20 //     Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch
21 //-------------------------------------------------------------------------
22 #include <TObjArray.h>
23 #include <TTree.h>
24
25 #include "AliESD.h"
26 #include "AliESDtrack.h"
27 #include "AliITStrackV2.h"
28 #include "AliV0vertex.h"
29 #include "AliV0vertexer.h"
30
31 ClassImp(AliV0vertexer)
32
33 Int_t AliV0vertexer::Tracks2V0vertices(AliESD *event) {
34   //--------------------------------------------------------------------
35   //This function reconstructs V0 vertices
36   //--------------------------------------------------------------------
37
38    Int_t nentr=event->GetNumberOfTracks();
39
40    TObjArray negtrks(nentr/2);
41    TObjArray postrks(nentr/2);
42
43    Int_t nneg=0, npos=0, nvtx=0;
44
45    Int_t i;
46    for (i=0; i<nentr; i++) {
47      AliESDtrack *esd=event->GetTrack(i);
48      UInt_t status=esd->GetStatus();
49      UInt_t flags=AliESDtrack::kITSin|AliESDtrack::kTPCin;
50
51      if ((status&AliESDtrack::kITSrefit)==0)
52         if ((status&flags)!=status) continue;
53
54      AliITStrackV2 *iotrack=new AliITStrackV2(*esd);
55      iotrack->SetLabel(i);  // now it is the index in array of ESD tracks
56      if ((status&AliESDtrack::kITSrefit)==0)   //correction for the beam pipe
57         if (!iotrack->PropagateTo(3.,0.0023,65.19)) {
58           delete iotrack;
59           continue;
60         }
61      if (!iotrack->PropagateTo(2.5,0.,0.)) {
62        delete iotrack;
63        continue;
64      }
65
66      if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
67      else {npos++; postrks.AddLast(iotrack);}
68    }   
69
70
71    for (i=0; i<nneg; i++) {
72       //if (i%10==0) cerr<<nneg-i<<'\r';
73       AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
74
75       if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
76       if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
77
78       for (Int_t k=0; k<npos; k++) {
79          AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
80
81          if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
82          if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
83
84          if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
85          if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
86
87
88          AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
89
90          Double_t dca=PropagateToDCA(pnt,ppt);
91          if (dca > fDCAmax) continue;
92
93          AliV0vertex vertex(*pnt,*ppt);
94          if (vertex.GetChi2() > fChi2max) continue;
95          
96          /*  Think of something better here ! 
97          nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
98          pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
99          */
100
101          Double_t x,y,z; vertex.GetXYZ(x,y,z); 
102          Double_t r2=x*x + y*y; 
103          if (r2 > fRmax*fRmax) continue;
104          if (r2 < fRmin*fRmin) continue;
105
106          Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
107          Double_t p2=px*px+py*py+pz*pz;
108          Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/
109                TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ)));
110
111         //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
112          if (cost < fCPAmax) continue;
113
114          //vertex.ChangeMassHypothesis(); //default is Lambda0 
115
116          event->AddV0(&vertex);
117
118          nvtx++;
119       }
120    }
121
122    Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
123
124    negtrks.Delete();
125    postrks.Delete();
126
127    return 0;
128 }
129
130 Int_t AliV0vertexer::Tracks2V0vertices(TTree *tTree, TTree *vTree) {
131   //--------------------------------------------------------------------
132   //This function reconstructs V0 vertices
133   //--------------------------------------------------------------------
134   Warning("Tracks2V0vertices(TTree*,TTree*)",
135           "Will be removed soon !  Use Tracks2V0vertices(AliESD*) instead");
136
137    TBranch *branch=tTree->GetBranch("tracks");
138    if (!branch) {
139       Error("Tracks2V0vertices","Can't get the branch !");
140       return 1;
141    }
142    Int_t nentr=(Int_t)tTree->GetEntries();
143
144    TObjArray negtrks(nentr/2);
145    TObjArray postrks(nentr/2);
146
147    Int_t nneg=0, npos=0, nvtx=0;
148
149    Int_t i;
150    for (i=0; i<nentr; i++) {
151        AliITStrackV2 *iotrack=new AliITStrackV2;
152        branch->SetAddress(&iotrack);
153        tTree->GetEvent(i);
154
155        if (!iotrack->PropagateTo(3.,0.0023,65.19)) {
156          delete iotrack;
157          continue; 
158        }
159        if (!iotrack->PropagateTo(2.5,0.,0.)) {
160          delete iotrack;
161          continue;
162        }
163
164        if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
165        else {npos++; postrks.AddLast(iotrack);}
166    }   
167
168
169    AliV0vertex *ioVertex=0;
170    branch=vTree->GetBranch("vertices");
171    if (!branch) vTree->Branch("vertices","AliV0vertex",&ioVertex,32000,3);
172    else branch->SetAddress(&ioVertex);
173
174
175    for (i=0; i<nneg; i++) {
176       //if (i%10==0) cerr<<nneg-i<<'\r';
177       AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
178
179       if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
180       if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
181
182       for (Int_t k=0; k<npos; k++) {
183          AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
184
185          if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
186          if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
187
188          if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
189          if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
190
191
192          AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
193
194          Double_t dca=PropagateToDCA(pnt,ppt);
195          if (dca > fDCAmax) continue;
196
197          AliV0vertex vertex(*pnt,*ppt);
198          if (vertex.GetChi2() > fChi2max) continue;
199          
200          /*  Think of something better here ! 
201          nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
202          pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
203          */
204
205          Double_t x,y,z; vertex.GetXYZ(x,y,z); 
206          Double_t r2=x*x + y*y; 
207          if (r2 > fRmax*fRmax) continue;
208          if (r2 < fRmin*fRmin) continue;
209
210          Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
211          Double_t p2=px*px+py*py+pz*pz;
212          Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/
213                TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ)));
214
215         //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
216          if (cost < fCPAmax) continue;
217
218          //vertex.ChangeMassHypothesis(); //default is Lambda0 
219
220          ioVertex=&vertex; vTree->Fill();
221
222          nvtx++;
223       }
224    }
225
226    Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
227
228    negtrks.Delete();
229    postrks.Delete();
230
231    return 0;
232 }
233
234 Double_t 
235 AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) const {
236   //--------------------------------------------------------------------
237   // This function returns the DCA between two tracks
238   // The tracks will be moved to the point of DCA ! 
239   //--------------------------------------------------------------------
240   return n->PropagateToDCA(p);
241 }
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257