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)) {
61 if (!iotrack->PropagateTo(2.5,0.,0.)) {
66 if (iotrack->Get1Pt() < 0.) {nneg++; negtrks.AddLast(iotrack);}
67 else {npos++; postrks.AddLast(iotrack);}
71 for (i=0; i<nneg; i++) {
72 //if (i%10==0) cerr<<nneg-i<<'\r';
73 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
75 if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
76 if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
78 for (Int_t k=0; k<npos; k++) {
79 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
81 if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
82 if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
84 if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
85 if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
88 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
90 Double_t dca=PropagateToDCA(pnt,ppt);
91 if (dca > fDCAmax) continue;
93 AliV0vertex vertex(*pnt,*ppt);
94 if (vertex.GetChi2() > fChi2max) continue;
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;
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;
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)));
111 //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
112 if (cost < fCPAmax) continue;
114 //vertex.ChangeMassHypothesis(); //default is Lambda0
116 event->AddV0(&vertex);
122 Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
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");
137 TBranch *branch=tTree->GetBranch("tracks");
139 Error("Tracks2V0vertices","Can't get the branch !");
142 Int_t nentr=(Int_t)tTree->GetEntries();
144 TObjArray negtrks(nentr/2);
145 TObjArray postrks(nentr/2);
147 Int_t nneg=0, npos=0, nvtx=0;
150 for (i=0; i<nentr; i++) {
151 AliITStrackV2 *iotrack=new AliITStrackV2;
152 branch->SetAddress(&iotrack);
155 if (!iotrack->PropagateTo(3.,0.0023,65.19)) {
159 if (!iotrack->PropagateTo(2.5,0.,0.)) {
164 if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
165 else {npos++; postrks.AddLast(iotrack);}
169 AliV0vertex *ioVertex=0;
170 branch=vTree->GetBranch("vertices");
171 if (!branch) vTree->Branch("vertices","AliV0vertex",&ioVertex,32000,3);
172 else branch->SetAddress(&ioVertex);
175 for (i=0; i<nneg; i++) {
176 //if (i%10==0) cerr<<nneg-i<<'\r';
177 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
179 if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
180 if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
182 for (Int_t k=0; k<npos; k++) {
183 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
185 if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
186 if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
188 if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
189 if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
192 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
194 Double_t dca=PropagateToDCA(pnt,ppt);
195 if (dca > fDCAmax) continue;
197 AliV0vertex vertex(*pnt,*ppt);
198 if (vertex.GetChi2() > fChi2max) continue;
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;
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;
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)));
215 //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
216 if (cost < fCPAmax) continue;
218 //vertex.ChangeMassHypothesis(); //default is Lambda0
220 ioVertex=&vertex; vTree->Fill();
226 Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
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);