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