a9a2d814 |
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 cascade vertexer class |
18 | // |
19 | // Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr |
20 | //------------------------------------------------------------------------- |
a9a2d814 |
21 | #include <TObjArray.h> |
116cbefd |
22 | #include <TTree.h> |
a9a2d814 |
23 | |
e23730c7 |
24 | #include "AliESD.h" |
25 | #include "AliESDv0.h" |
26 | #include "AliESDcascade.h" |
a9a2d814 |
27 | #include "AliCascadeVertex.h" |
28 | #include "AliCascadeVertexer.h" |
a9a2d814 |
29 | #include "AliITStrackV2.h" |
116cbefd |
30 | #include "AliV0vertex.h" |
a9a2d814 |
31 | |
32 | ClassImp(AliCascadeVertexer) |
33 | |
e23730c7 |
34 | Int_t AliCascadeVertexer::V0sTracks2CascadeVertices(AliESD *event) { |
35 | //-------------------------------------------------------------------- |
36 | // This function reconstructs cascade vertices |
37 | // Adapted to the ESD by I.Belikov (Jouri.Belikov@cern.ch) |
38 | //-------------------------------------------------------------------- |
39 | |
40 | Int_t nV0=(Int_t)event->GetNumberOfV0s(); |
41 | TObjArray vtcs(nV0); |
42 | Int_t i; |
43 | for (i=0; i<nV0; i++) { |
44 | const AliESDv0 *esdV0=event->GetV0(i); |
45 | vtcs.AddLast(new AliV0vertex(*esdV0)); |
46 | } |
47 | |
48 | |
49 | Int_t ntr=(Int_t)event->GetNumberOfTracks(); |
50 | TObjArray trks(ntr); |
51 | for (i=0; i<ntr; i++) { |
52 | AliESDtrack *esdtr=event->GetTrack(i); |
53 | AliITStrackV2 *iotrack=new AliITStrackV2(*esdtr); |
54 | iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.); |
55 | trks.AddLast(iotrack); |
56 | } |
57 | |
58 | Double_t massLambda=1.11568; |
59 | Int_t ncasc=0; |
60 | |
61 | // Looking for the cascades... |
62 | for (i=0; i<nV0; i++) { |
63 | AliV0vertex *v=(AliV0vertex*)vtcs.UncheckedAt(i); |
64 | v->ChangeMassHypothesis(kLambda0); // the v0 must be Lambda |
65 | if (TMath::Abs(v->GetEffMass()-massLambda)>fMassWin) continue; |
66 | if (v->GetD(0,0,0)<fDV0min) continue; |
67 | for (Int_t j=0; j<ntr; j++) { |
68 | AliITStrackV2 *b=(AliITStrackV2*)trks.UncheckedAt(j); |
69 | |
70 | if (TMath::Abs(b->GetD())<fDBachMin) continue; |
71 | if (b->Get1Pt()<0.) continue; // bachelor's charge |
72 | |
73 | AliV0vertex v0(*v), *pv0=&v0; |
74 | AliITStrackV2 bt(*b), *pbt=&bt; |
75 | |
76 | Double_t dca=PropagateToDCA(pv0,pbt); |
77 | if (dca > fDCAmax) continue; |
78 | |
79 | AliCascadeVertex cascade(*pv0,*pbt); |
80 | if (cascade.GetChi2() > fChi2max) continue; |
81 | |
82 | Double_t x,y,z; cascade.GetXYZ(x,y,z); |
83 | Double_t r2=x*x + y*y; |
84 | if (r2 > fRmax*fRmax) continue; // condition on fiducial zone |
85 | if (r2 < fRmin*fRmin) continue; |
86 | |
87 | { |
88 | Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1); |
89 | if (r2 > (x1*x1+y1*y1)) continue; |
90 | if (z*z > z1*z1) continue; |
91 | } |
92 | |
93 | Double_t px,py,pz; cascade.GetPxPyPz(px,py,pz); |
94 | Double_t p2=px*px+py*py+pz*pz; |
95 | Double_t cost=(x*px+y*py+z*pz)/TMath::Sqrt(p2*(r2+z*z)); |
96 | |
97 | if (cost<fCPAmax) continue; //condition on the cascade pointing angle |
98 | //cascade.ChangeMassHypothesis(); //default is Xi |
99 | |
100 | event->AddCascade(&cascade); |
101 | |
102 | ncasc++; |
103 | |
104 | } |
105 | } |
106 | |
107 | // Looking for the anti-cascades... |
108 | for (i=0; i<nV0; i++) { |
109 | AliV0vertex *v=(AliV0vertex*)vtcs.UncheckedAt(i); |
110 | v->ChangeMassHypothesis(kLambda0Bar); //the v0 must be anti-Lambda |
111 | if (TMath::Abs(v->GetEffMass()-massLambda)>fMassWin) continue; |
112 | if (v->GetD(0,0,0)<fDV0min) continue; |
113 | for (Int_t j=0; j<ntr; j++) { |
114 | AliITStrackV2 *b=(AliITStrackV2*)trks.UncheckedAt(j); |
115 | |
116 | if (TMath::Abs(b->GetD())<fDBachMin) continue; |
117 | if (b->Get1Pt()>0.) continue; // bachelor's charge |
118 | |
119 | AliV0vertex v0(*v), *pv0=&v0; |
120 | AliITStrackV2 bt(*b), *pbt=&bt; |
121 | |
122 | Double_t dca=PropagateToDCA(pv0,pbt); |
123 | if (dca > fDCAmax) continue; |
124 | |
125 | AliCascadeVertex cascade(*pv0,*pbt); |
126 | if (cascade.GetChi2() > fChi2max) continue; |
127 | |
128 | Double_t x,y,z; cascade.GetXYZ(x,y,z); |
129 | Double_t r2=x*x + y*y; |
130 | if (r2 > fRmax*fRmax) continue; // condition on fiducial zone |
131 | if (r2 < fRmin*fRmin) continue; |
132 | |
133 | { |
134 | Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1); |
135 | if (r2 > (x1*x1+y1*y1)) continue; |
136 | if (z*z > z1*z1) continue; |
137 | } |
138 | |
139 | Double_t px,py,pz; cascade.GetPxPyPz(px,py,pz); |
140 | Double_t p2=px*px+py*py+pz*pz; |
141 | Double_t cost=(x*px+y*py+z*pz)/TMath::Sqrt(p2*(r2+z*z)); |
142 | |
143 | if (cost<fCPAmax) continue; //condition on the cascade pointing angle |
144 | //cascade.ChangeMassHypothesis(); //default is Xi |
145 | |
146 | event->AddCascade(&cascade); |
147 | |
148 | ncasc++; |
149 | |
150 | } |
151 | } |
152 | |
153 | Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc); |
154 | |
155 | trks.Delete(); |
156 | vtcs.Delete(); |
157 | |
158 | return 0; |
159 | } |
160 | |
566abf75 |
161 | Int_t AliCascadeVertexer:: |
162 | V0sTracks2CascadeVertices(TTree *vTree,TTree *tTree, TTree *xTree) { |
a9a2d814 |
163 | //-------------------------------------------------------------------- |
566abf75 |
164 | // This function reconstructs cascade vertices |
a9a2d814 |
165 | //-------------------------------------------------------------------- |
566abf75 |
166 | TBranch *branch=vTree->GetBranch("vertices"); |
167 | if (!branch) { |
168 | Error("V0sTracks2CascadeVertices","Can't get the V0 branch !"); |
169 | return 1; |
170 | } |
171 | Int_t nentrV0=(Int_t)vTree->GetEntries(); |
a9a2d814 |
172 | |
173 | TObjArray vtxV0(nentrV0); |
174 | |
175 | // fill TObjArray vtxV0 with vertices |
176 | |
177 | Int_t i, nV0=0; |
178 | for (i=0; i<nentrV0; i++) { |
179 | |
180 | AliV0vertex *ioVertex=new AliV0vertex; |
181 | branch->SetAddress(&ioVertex); |
566abf75 |
182 | vTree->GetEvent(i); |
a9a2d814 |
183 | nV0++; |
184 | vtxV0.AddLast(ioVertex); |
185 | |
186 | } |
187 | |
566abf75 |
188 | branch=tTree->GetBranch("tracks"); |
189 | if (!branch) { |
190 | Error("V0sTracks2CascadeVertices","Can't get the track branch !"); |
191 | return 2; |
192 | } |
193 | Int_t nentr=(Int_t)tTree->GetEntries(); |
a9a2d814 |
194 | |
195 | TObjArray trks(nentr); |
196 | |
197 | // fill TObjArray trks with tracks |
198 | |
199 | Int_t ntrack=0; |
200 | |
201 | for (i=0; i<nentr; i++) { |
202 | |
203 | AliITStrackV2 *iotrack=new AliITStrackV2; |
204 | branch->SetAddress(&iotrack); |
566abf75 |
205 | tTree->GetEvent(i); |
4ab260d6 |
206 | |
207 | iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.); |
208 | |
a9a2d814 |
209 | ntrack++; trks.AddLast(iotrack); |
210 | |
211 | } |
212 | |
566abf75 |
213 | AliCascadeVertex *ioCascade=0; |
214 | branch=xTree->GetBranch("cascades"); |
215 | if (!branch) xTree->Branch("cascades","AliCascadeVertex",&ioCascade,32000,3); |
216 | else branch->SetAddress(&ioCascade); |
a9a2d814 |
217 | |
218 | // loop on all vertices |
219 | |
220 | Double_t massLambda=1.11568; |
221 | Int_t ncasc=0; |
222 | |
223 | for (i=0; i<nV0; i++) { |
224 | |
93f82b23 |
225 | AliV0vertex *lV0ver=(AliV0vertex *)vtxV0.UncheckedAt(i); |
a9a2d814 |
226 | |
93f82b23 |
227 | lV0ver->ChangeMassHypothesis(kLambda0); //I.B. |
a9a2d814 |
228 | |
93f82b23 |
229 | if (lV0ver->GetEffMass()<massLambda-fMassWin || // condition of the V0 mass window (cut fMassWin) |
230 | lV0ver->GetEffMass()>massLambda+fMassWin) continue; |
a9a2d814 |
231 | |
93f82b23 |
232 | if (lV0ver->GetD(0,0,0)<fDV0min) continue; // condition of minimum impact parameter of the V0 (cut fDV0min) |
a9a2d814 |
233 | // here why not cuting on pointing angle ??? |
234 | |
235 | // for each vertex in the good mass range, loop on all tracks (= bachelor candidates) |
236 | |
237 | for (Int_t k=0; k<ntrack; k++) { |
238 | |
239 | AliITStrackV2 *bachtrk=(AliITStrackV2 *)trks.UncheckedAt(k); |
240 | |
241 | if (TMath::Abs(bachtrk->GetD())<fDBachMin) continue; // eliminate to small impact parameters |
242 | |
93f82b23 |
243 | if (lV0ver->GetPdgCode()==kLambda0 && bachtrk->Get1Pt()<0.) continue; // condition on V0 label |
244 | if (lV0ver->GetPdgCode()==kLambda0Bar && bachtrk->Get1Pt()>0.) continue; // + good sign for bachelor |
a9a2d814 |
245 | |
93f82b23 |
246 | AliV0vertex lV0(*lV0ver), *pV0=&lV0; |
a9a2d814 |
247 | AliITStrackV2 bt(*bachtrk), *pbt=&bt; |
248 | |
249 | // calculation of the distance of closest approach between the V0 and the bachelor |
250 | |
251 | Double_t dca=PropagateToDCA(pV0,pbt); |
252 | if (dca > fDCAmax) continue; // cut on dca |
253 | |
254 | // construction of a cascade object |
255 | |
256 | AliCascadeVertex cascade(*pV0,*pbt); |
257 | if (cascade.GetChi2() > fChi2max) continue; |
258 | |
259 | // get cascade decay position (V0, bachelor) |
260 | |
261 | Double_t x,y,z; cascade.GetXYZ(x,y,z); |
262 | Double_t r2=x*x + y*y; |
263 | if (r2 > fRmax*fRmax) continue; // condition on fiducial zone |
264 | if (r2 < fRmin*fRmin) continue; |
265 | |
4ab260d6 |
266 | { |
267 | //I.B. |
93f82b23 |
268 | Double_t x1,y1,z1; lV0ver->GetXYZ(x1,y1,z1); |
4ab260d6 |
269 | if (r2 > (x1*x1+y1*y1)) continue; |
270 | if (z*z > z1*z1) continue; |
271 | } |
272 | |
a9a2d814 |
273 | // get cascade momentum |
274 | |
275 | Double_t px,py,pz; cascade.GetPxPyPz(px,py,pz); |
276 | Double_t p2=px*px+py*py+pz*pz; |
277 | Double_t cost=(x*px+y*py+z*pz)/TMath::Sqrt(p2*(r2+z*z)); |
278 | |
279 | if (cost<fCPAmax) continue; // condition on the cascade pointing angle |
280 | |
281 | //cascade.ChangeMassHypothesis(); //default is Xi |
282 | |
283 | |
566abf75 |
284 | ioCascade=&cascade; xTree->Fill(); |
a9a2d814 |
285 | |
286 | ncasc++; |
287 | |
288 | } |
289 | } |
290 | |
e23730c7 |
291 | Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc); |
a9a2d814 |
292 | |
a9a2d814 |
293 | trks.Delete(); |
294 | vtxV0.Delete(); |
295 | |
a9a2d814 |
296 | return 0; |
297 | } |
298 | |
299 | inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){ |
300 | // determinant 2x2 |
301 | return a00*a11 - a01*a10; |
302 | } |
303 | |
304 | inline Double_t det (Double_t a00,Double_t a01,Double_t a02, |
305 | Double_t a10,Double_t a11,Double_t a12, |
306 | Double_t a20,Double_t a21,Double_t a22) { |
307 | // determinant 3x3 |
308 | return |
309 | a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21); |
310 | } |
311 | |
312 | |
313 | |
314 | |
315 | Double_t |
316 | AliCascadeVertexer::PropagateToDCA(AliV0vertex *v, AliITStrackV2 *t) { |
317 | //-------------------------------------------------------------------- |
318 | // This function returns the DCA between the V0 and the track |
319 | //-------------------------------------------------------------------- |
320 | |
321 | Double_t phi, x, par[5]; |
322 | Double_t alpha, cs1, sn1; |
323 | |
324 | t->GetExternalParameters(x,par); alpha=t->GetAlpha(); |
325 | phi=TMath::ASin(par[2]) + alpha; |
326 | Double_t px1=TMath::Cos(phi), py1=TMath::Sin(phi), pz1=par[3]; |
327 | |
328 | |
329 | cs1=TMath::Cos(alpha); sn1=TMath::Sin(alpha); |
330 | Double_t x1=x*cs1 - par[0]*sn1; |
331 | Double_t y1=x*sn1 + par[0]*cs1; |
332 | Double_t z1=par[1]; |
333 | |
334 | |
335 | Double_t x2,y2,z2; // position and momentum of V0 |
336 | Double_t px2,py2,pz2; |
337 | |
338 | v->GetXYZ(x2,y2,z2); |
339 | v->GetPxPyPz(px2,py2,pz2); |
340 | |
341 | // calculation dca |
342 | |
343 | Double_t dd= det(x2-x1,y2-y1,z2-z1,px1,py1,pz1,px2,py2,pz2); |
344 | Double_t ax= det(py1,pz1,py2,pz2); |
345 | Double_t ay=-det(px1,pz1,px2,pz2); |
346 | Double_t az= det(px1,py1,px2,py2); |
347 | |
348 | Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az); |
349 | |
350 | //points of the DCA |
351 | Double_t t1 = det(x2-x1,y2-y1,z2-z1,px2,py2,pz2,ax,ay,az)/ |
352 | det(px1,py1,pz1,px2,py2,pz2,ax,ay,az); |
353 | |
354 | x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1; |
355 | |
356 | |
357 | //propagate track to the points of DCA |
358 | |
359 | x1=x1*cs1 + y1*sn1; |
360 | if (!t->PropagateTo(x1,0.,0.)) { |
e23730c7 |
361 | Error("PropagateToDCA","Propagation failed !"); |
a9a2d814 |
362 | return 1.e+33; |
363 | } |
364 | |
365 | return dca; |
366 | } |
367 | |
368 | |
369 | |
370 | |
371 | |
372 | |
373 | |
374 | |
375 | |
376 | |
377 | |