]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliV0vertexer.cxx
Initial version, Compiler directive for selection of containers type: either STL...
[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)) continue;
58      if (!iotrack->PropagateTo(2.5,0.,0.)) continue;
59
60      if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
61      else {npos++; postrks.AddLast(iotrack);}
62    }   
63
64
65    for (i=0; i<nneg; i++) {
66       //if (i%10==0) cerr<<nneg-i<<'\r';
67       AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
68
69       if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
70       if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
71
72       for (Int_t k=0; k<npos; k++) {
73          AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
74
75          if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
76          if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
77
78          if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
79          if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
80
81
82          AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
83
84          Double_t dca=PropagateToDCA(pnt,ppt);
85          if (dca > fDCAmax) continue;
86
87          AliV0vertex vertex(*pnt,*ppt);
88          if (vertex.GetChi2() > fChi2max) continue;
89          
90          /*  Think of something better here ! 
91          nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
92          pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
93          */
94
95          Double_t x,y,z; vertex.GetXYZ(x,y,z); 
96          Double_t r2=x*x + y*y; 
97          if (r2 > fRmax*fRmax) continue;
98          if (r2 < fRmin*fRmin) continue;
99
100          Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
101          Double_t p2=px*px+py*py+pz*pz;
102          Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/
103                TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ)));
104
105         //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
106          if (cost < fCPAmax) continue;
107
108          //vertex.ChangeMassHypothesis(); //default is Lambda0 
109
110          event->AddV0(&vertex);
111
112          nvtx++;
113       }
114    }
115
116    Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
117
118    negtrks.Delete();
119    postrks.Delete();
120
121    return 0;
122 }
123
124 Int_t AliV0vertexer::Tracks2V0vertices(TTree *tTree, TTree *vTree) {
125   //--------------------------------------------------------------------
126   //This function reconstructs V0 vertices
127   //--------------------------------------------------------------------
128   Warning("Tracks2V0vertices(TTree*,TTree*)",
129           "Will be removed soon !  Use Tracks2V0vertices(AliESD*) instead");
130
131    TBranch *branch=tTree->GetBranch("tracks");
132    if (!branch) {
133       Error("Tracks2V0vertices","Can't get the branch !");
134       return 1;
135    }
136    Int_t nentr=(Int_t)tTree->GetEntries();
137
138    TObjArray negtrks(nentr/2);
139    TObjArray postrks(nentr/2);
140
141    Int_t nneg=0, npos=0, nvtx=0;
142
143    Int_t i;
144    for (i=0; i<nentr; i++) {
145        AliITStrackV2 *iotrack=new AliITStrackV2;
146        branch->SetAddress(&iotrack);
147        tTree->GetEvent(i);
148
149        if (!iotrack->PropagateTo(3.,0.0023,65.19)) continue; 
150        if (!iotrack->PropagateTo(2.5,0.,0.)) continue;
151
152        if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
153        else {npos++; postrks.AddLast(iotrack);}
154    }   
155
156
157    AliV0vertex *ioVertex=0;
158    branch=vTree->GetBranch("vertices");
159    if (!branch) vTree->Branch("vertices","AliV0vertex",&ioVertex,32000,3);
160    else branch->SetAddress(&ioVertex);
161
162
163    for (i=0; i<nneg; i++) {
164       //if (i%10==0) cerr<<nneg-i<<'\r';
165       AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
166
167       if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
168       if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
169
170       for (Int_t k=0; k<npos; k++) {
171          AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
172
173          if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
174          if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
175
176          if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
177          if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
178
179
180          AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
181
182          Double_t dca=PropagateToDCA(pnt,ppt);
183          if (dca > fDCAmax) continue;
184
185          AliV0vertex vertex(*pnt,*ppt);
186          if (vertex.GetChi2() > fChi2max) continue;
187          
188          /*  Think of something better here ! 
189          nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
190          pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
191          */
192
193          Double_t x,y,z; vertex.GetXYZ(x,y,z); 
194          Double_t r2=x*x + y*y; 
195          if (r2 > fRmax*fRmax) continue;
196          if (r2 < fRmin*fRmin) continue;
197
198          Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
199          Double_t p2=px*px+py*py+pz*pz;
200          Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/
201                TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ)));
202
203         //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
204          if (cost < fCPAmax) continue;
205
206          //vertex.ChangeMassHypothesis(); //default is Lambda0 
207
208          ioVertex=&vertex; vTree->Fill();
209
210          nvtx++;
211       }
212    }
213
214    Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
215
216    negtrks.Delete();
217    postrks.Delete();
218
219    return 0;
220 }
221
222 Double_t 
223 AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) const {
224   //--------------------------------------------------------------------
225   // This function returns the DCA between two tracks
226   // The tracks will be moved to the point of DCA ! 
227   //--------------------------------------------------------------------
228   return n->PropagateToDCA(p);
229 }
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245