1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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>
26 #include "AliESDtrack.h"
27 #include "AliITStrackV2.h"
28 #include "AliV0vertex.h"
29 #include "AliV0vertexer.h"
31 ClassImp(AliV0vertexer)
33 Int_t AliV0vertexer::Tracks2V0vertices(AliESD *event) {
34 //--------------------------------------------------------------------
35 //This function reconstructs V0 vertices
36 //--------------------------------------------------------------------
38 Int_t nentr=event->GetNumberOfTracks();
40 TObjArray negtrks(nentr/2);
41 TObjArray postrks(nentr/2);
43 Int_t nneg=0, npos=0, nvtx=0;
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;
51 if ((status&AliESDtrack::kITSrefit)==0)
52 if ((status&flags)!=status) continue;
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;
60 if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
61 else {npos++; postrks.AddLast(iotrack);}
65 for (i=0; i<nneg; i++) {
66 //if (i%10==0) cerr<<nneg-i<<'\r';
67 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
69 if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
70 if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
72 for (Int_t k=0; k<npos; k++) {
73 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
75 if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
76 if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
78 if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
79 if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
82 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
84 Double_t dca=PropagateToDCA(pnt,ppt);
85 if (dca > fDCAmax) continue;
87 AliV0vertex vertex(*pnt,*ppt);
88 if (vertex.GetChi2() > fChi2max) continue;
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;
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;
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)));
105 //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
106 if (cost < fCPAmax) continue;
108 //vertex.ChangeMassHypothesis(); //default is Lambda0
110 event->AddV0(&vertex);
116 Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
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");
131 TBranch *branch=tTree->GetBranch("tracks");
133 Error("Tracks2V0vertices","Can't get the branch !");
136 Int_t nentr=(Int_t)tTree->GetEntries();
138 TObjArray negtrks(nentr/2);
139 TObjArray postrks(nentr/2);
141 Int_t nneg=0, npos=0, nvtx=0;
144 for (i=0; i<nentr; i++) {
145 AliITStrackV2 *iotrack=new AliITStrackV2;
146 branch->SetAddress(&iotrack);
149 if (!iotrack->PropagateTo(3.,0.0023,65.19)) continue;
150 if (!iotrack->PropagateTo(2.5,0.,0.)) continue;
152 if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
153 else {npos++; postrks.AddLast(iotrack);}
157 AliV0vertex *ioVertex=0;
158 branch=vTree->GetBranch("vertices");
159 if (!branch) vTree->Branch("vertices","AliV0vertex",&ioVertex,32000,3);
160 else branch->SetAddress(&ioVertex);
163 for (i=0; i<nneg; i++) {
164 //if (i%10==0) cerr<<nneg-i<<'\r';
165 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
167 if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
168 if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
170 for (Int_t k=0; k<npos; k++) {
171 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
173 if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
174 if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
176 if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
177 if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
180 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
182 Double_t dca=PropagateToDCA(pnt,ppt);
183 if (dca > fDCAmax) continue;
185 AliV0vertex vertex(*pnt,*ppt);
186 if (vertex.GetChi2() > fChi2max) continue;
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;
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;
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)));
203 //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
204 if (cost < fCPAmax) continue;
206 //vertex.ChangeMassHypothesis(); //default is Lambda0
208 ioVertex=&vertex; vTree->Fill();
214 Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
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);