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 cascade vertexer class
18 // Reads V0s and tracks, writes out cascade vertices
19 // Fills the ESD with the cascades
20 // Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
21 //-------------------------------------------------------------------------
22 #include <TObjArray.h>
27 #include "AliCascadeVertex.h"
28 #include "AliCascadeVertexer.h"
29 #include "AliITStrackV2.h"
30 #include "AliV0vertex.h"
32 ClassImp(AliCascadeVertexer)
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 //--------------------------------------------------------------------
40 Int_t nV0=(Int_t)event->GetNumberOfV0s();
43 for (i=0; i<nV0; i++) {
44 const AliESDv0 *esdV0=event->GetV0(i);
45 vtcs.AddLast(new AliV0vertex(*esdV0));
49 Int_t ntr=(Int_t)event->GetNumberOfTracks();
51 for (i=0; i<ntr; i++) {
52 AliESDtrack *esdtr=event->GetTrack(i);
53 UInt_t status=esdtr->GetStatus();
54 UInt_t flags=AliESDtrack::kITSin|AliESDtrack::kTPCin;
56 if ((status&AliESDtrack::kITSrefit)==0)
57 if ((status&flags)!=status) continue;
59 AliITStrackV2 *iotrack=new AliITStrackV2(*esdtr);
60 iotrack->SetLabel(i); // now it is the index in array of ESD tracks
61 if ((status&AliESDtrack::kITSrefit)==0) //correction for the beam pipe
62 if (!iotrack->PropagateTo(3.,0.0023,65.19)) continue;
63 if (!iotrack->PropagateTo(2.5,0.,0.)) continue;
64 trks.AddLast(iotrack);
66 ntr=trks.GetEntriesFast();
68 Double_t massLambda=1.11568;
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);
80 if (TMath::Abs(b->GetD())<fDBachMin) continue;
81 if (b->Get1Pt()>0.) continue; // bachelor's charge
83 AliV0vertex v0(*v), *pv0=&v0;
84 AliITStrackV2 bt(*b), *pbt=&bt;
86 Double_t dca=PropagateToDCA(pv0,pbt);
87 if (dca > fDCAmax) continue;
89 AliCascadeVertex cascade(*pv0,*pbt);
90 if (cascade.GetChi2() > fChi2max) continue;
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;
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;
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));
107 if (cost<fCPAmax) continue; //condition on the cascade pointing angle
108 //cascade.ChangeMassHypothesis(); //default is Xi
110 event->AddCascade(&cascade);
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);
126 if (TMath::Abs(b->GetD())<fDBachMin) continue;
127 if (b->Get1Pt()<0.) continue; // bachelor's charge
129 AliV0vertex v0(*v), *pv0=&v0;
130 AliITStrackV2 bt(*b), *pbt=&bt;
132 Double_t dca=PropagateToDCA(pv0,pbt);
133 if (dca > fDCAmax) continue;
135 AliCascadeVertex cascade(*pv0,*pbt);
136 if (cascade.GetChi2() > fChi2max) continue;
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;
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;
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));
153 if (cost<fCPAmax) continue; //condition on the cascade pointing angle
154 //cascade.ChangeMassHypothesis(); //default is Xi
156 event->AddCascade(&cascade);
163 Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);
171 Int_t AliCascadeVertexer::
172 V0sTracks2CascadeVertices(TTree *vTree,TTree *tTree, TTree *xTree) {
173 //--------------------------------------------------------------------
174 // This function reconstructs cascade vertices
175 //--------------------------------------------------------------------
176 Warning("V0sTracks2CascadeVertices(TTree*,TTree*,TTree*)",
177 "Will be removed soon ! Use V0sTracks2CascadeVertices(AliESD*) instead");
179 TBranch *branch=vTree->GetBranch("vertices");
181 Error("V0sTracks2CascadeVertices","Can't get the V0 branch !");
184 Int_t nentrV0=(Int_t)vTree->GetEntries();
186 TObjArray vtxV0(nentrV0);
188 // fill TObjArray vtxV0 with vertices
191 for (i=0; i<nentrV0; i++) {
193 AliV0vertex *ioVertex=new AliV0vertex;
194 branch->SetAddress(&ioVertex);
197 vtxV0.AddLast(ioVertex);
201 branch=tTree->GetBranch("tracks");
203 Error("V0sTracks2CascadeVertices","Can't get the track branch !");
206 Int_t nentr=(Int_t)tTree->GetEntries();
208 TObjArray trks(nentr);
210 // fill TObjArray trks with tracks
214 for (i=0; i<nentr; i++) {
216 AliITStrackV2 *iotrack=new AliITStrackV2;
217 branch->SetAddress(&iotrack);
220 if (!iotrack->PropagateTo(3.,0.0023,65.19)) continue;
221 if (!iotrack->PropagateTo(2.5,0.,0.)) continue;
223 ntrack++; trks.AddLast(iotrack);
227 AliCascadeVertex *ioCascade=0;
228 branch=xTree->GetBranch("cascades");
229 if (!branch) xTree->Branch("cascades","AliCascadeVertex",&ioCascade,32000,3);
230 else branch->SetAddress(&ioCascade);
232 // loop on all vertices
234 Double_t massLambda=1.11568;
237 for (i=0; i<nV0; i++) {
239 AliV0vertex *lV0ver=(AliV0vertex *)vtxV0.UncheckedAt(i);
241 lV0ver->ChangeMassHypothesis(kLambda0); //I.B.
243 if (lV0ver->GetEffMass()<massLambda-fMassWin || // condition of the V0 mass window (cut fMassWin)
244 lV0ver->GetEffMass()>massLambda+fMassWin) continue;
246 if (lV0ver->GetD(0,0,0)<fDV0min) continue; // condition of minimum impact parameter of the V0 (cut fDV0min)
247 // here why not cuting on pointing angle ???
249 // for each vertex in the good mass range, loop on all tracks (= bachelor candidates)
251 for (Int_t k=0; k<ntrack; k++) {
253 AliITStrackV2 *bachtrk=(AliITStrackV2 *)trks.UncheckedAt(k);
255 if (TMath::Abs(bachtrk->GetD())<fDBachMin) continue; // eliminate to small impact parameters
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
260 AliV0vertex lV0(*lV0ver), *pV0=&lV0;
261 AliITStrackV2 bt(*bachtrk), *pbt=&bt;
263 // calculation of the distance of closest approach between the V0 and the bachelor
265 Double_t dca=PropagateToDCA(pV0,pbt);
266 if (dca > fDCAmax) continue; // cut on dca
268 // construction of a cascade object
270 AliCascadeVertex cascade(*pV0,*pbt);
271 if (cascade.GetChi2() > fChi2max) continue;
273 // get cascade decay position (V0, bachelor)
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;
282 Double_t x1,y1,z1; lV0ver->GetXYZ(x1,y1,z1);
283 if (r2 > (x1*x1+y1*y1)) continue;
284 if (z*z > z1*z1) continue;
287 // get cascade momentum
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));
293 if (cost<fCPAmax) continue; // condition on the cascade pointing angle
295 //cascade.ChangeMassHypothesis(); //default is Xi
298 ioCascade=&cascade; xTree->Fill();
305 Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);
313 inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){
315 return a00*a11 - a01*a10;
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) {
323 a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21);
330 AliCascadeVertexer::PropagateToDCA(AliV0vertex *v, AliITStrackV2 *t) {
331 //--------------------------------------------------------------------
332 // This function returns the DCA between the V0 and the track
333 //--------------------------------------------------------------------
335 Double_t phi, x, par[5];
336 Double_t alpha, cs1, sn1;
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];
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;
349 Double_t x2,y2,z2; // position and momentum of V0
350 Double_t px2,py2,pz2;
353 v->GetPxPyPz(px2,py2,pz2);
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);
362 Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);
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);
368 x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
371 //propagate track to the points of DCA
374 if (!t->PropagateTo(x1,0.,0.)) {
375 Error("PropagateToDCA","Propagation failed !");