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