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