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
19 // Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
20 //-------------------------------------------------------------------------
21 #include <TObjArray.h>
26 #include "AliESDcascade.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 AliITStrackV2 *iotrack=new AliITStrackV2(*esdtr);
54 iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.);
55 trks.AddLast(iotrack);
58 Double_t massLambda=1.11568;
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);
70 if (TMath::Abs(b->GetD())<fDBachMin) continue;
71 if (b->Get1Pt()<0.) continue; // bachelor's charge
73 AliV0vertex v0(*v), *pv0=&v0;
74 AliITStrackV2 bt(*b), *pbt=&bt;
76 Double_t dca=PropagateToDCA(pv0,pbt);
77 if (dca > fDCAmax) continue;
79 AliCascadeVertex cascade(*pv0,*pbt);
80 if (cascade.GetChi2() > fChi2max) continue;
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;
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;
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));
97 if (cost<fCPAmax) continue; //condition on the cascade pointing angle
98 //cascade.ChangeMassHypothesis(); //default is Xi
100 event->AddCascade(&cascade);
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);
116 if (TMath::Abs(b->GetD())<fDBachMin) continue;
117 if (b->Get1Pt()>0.) continue; // bachelor's charge
119 AliV0vertex v0(*v), *pv0=&v0;
120 AliITStrackV2 bt(*b), *pbt=&bt;
122 Double_t dca=PropagateToDCA(pv0,pbt);
123 if (dca > fDCAmax) continue;
125 AliCascadeVertex cascade(*pv0,*pbt);
126 if (cascade.GetChi2() > fChi2max) continue;
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;
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;
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));
143 if (cost<fCPAmax) continue; //condition on the cascade pointing angle
144 //cascade.ChangeMassHypothesis(); //default is Xi
146 event->AddCascade(&cascade);
153 Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);
161 Int_t AliCascadeVertexer::
162 V0sTracks2CascadeVertices(TTree *vTree,TTree *tTree, TTree *xTree) {
163 //--------------------------------------------------------------------
164 // This function reconstructs cascade vertices
165 //--------------------------------------------------------------------
166 TBranch *branch=vTree->GetBranch("vertices");
168 Error("V0sTracks2CascadeVertices","Can't get the V0 branch !");
171 Int_t nentrV0=(Int_t)vTree->GetEntries();
173 TObjArray vtxV0(nentrV0);
175 // fill TObjArray vtxV0 with vertices
178 for (i=0; i<nentrV0; i++) {
180 AliV0vertex *ioVertex=new AliV0vertex;
181 branch->SetAddress(&ioVertex);
184 vtxV0.AddLast(ioVertex);
188 branch=tTree->GetBranch("tracks");
190 Error("V0sTracks2CascadeVertices","Can't get the track branch !");
193 Int_t nentr=(Int_t)tTree->GetEntries();
195 TObjArray trks(nentr);
197 // fill TObjArray trks with tracks
201 for (i=0; i<nentr; i++) {
203 AliITStrackV2 *iotrack=new AliITStrackV2;
204 branch->SetAddress(&iotrack);
207 iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.);
209 ntrack++; trks.AddLast(iotrack);
213 AliCascadeVertex *ioCascade=0;
214 branch=xTree->GetBranch("cascades");
215 if (!branch) xTree->Branch("cascades","AliCascadeVertex",&ioCascade,32000,3);
216 else branch->SetAddress(&ioCascade);
218 // loop on all vertices
220 Double_t massLambda=1.11568;
223 for (i=0; i<nV0; i++) {
225 AliV0vertex *lV0ver=(AliV0vertex *)vtxV0.UncheckedAt(i);
227 lV0ver->ChangeMassHypothesis(kLambda0); //I.B.
229 if (lV0ver->GetEffMass()<massLambda-fMassWin || // condition of the V0 mass window (cut fMassWin)
230 lV0ver->GetEffMass()>massLambda+fMassWin) continue;
232 if (lV0ver->GetD(0,0,0)<fDV0min) continue; // condition of minimum impact parameter of the V0 (cut fDV0min)
233 // here why not cuting on pointing angle ???
235 // for each vertex in the good mass range, loop on all tracks (= bachelor candidates)
237 for (Int_t k=0; k<ntrack; k++) {
239 AliITStrackV2 *bachtrk=(AliITStrackV2 *)trks.UncheckedAt(k);
241 if (TMath::Abs(bachtrk->GetD())<fDBachMin) continue; // eliminate to small impact parameters
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
246 AliV0vertex lV0(*lV0ver), *pV0=&lV0;
247 AliITStrackV2 bt(*bachtrk), *pbt=&bt;
249 // calculation of the distance of closest approach between the V0 and the bachelor
251 Double_t dca=PropagateToDCA(pV0,pbt);
252 if (dca > fDCAmax) continue; // cut on dca
254 // construction of a cascade object
256 AliCascadeVertex cascade(*pV0,*pbt);
257 if (cascade.GetChi2() > fChi2max) continue;
259 // get cascade decay position (V0, bachelor)
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;
268 Double_t x1,y1,z1; lV0ver->GetXYZ(x1,y1,z1);
269 if (r2 > (x1*x1+y1*y1)) continue;
270 if (z*z > z1*z1) continue;
273 // get cascade momentum
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));
279 if (cost<fCPAmax) continue; // condition on the cascade pointing angle
281 //cascade.ChangeMassHypothesis(); //default is Xi
284 ioCascade=&cascade; xTree->Fill();
291 Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);
299 inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){
301 return a00*a11 - a01*a10;
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) {
309 a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21);
316 AliCascadeVertexer::PropagateToDCA(AliV0vertex *v, AliITStrackV2 *t) {
317 //--------------------------------------------------------------------
318 // This function returns the DCA between the V0 and the track
319 //--------------------------------------------------------------------
321 Double_t phi, x, par[5];
322 Double_t alpha, cs1, sn1;
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];
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;
335 Double_t x2,y2,z2; // position and momentum of V0
336 Double_t px2,py2,pz2;
339 v->GetPxPyPz(px2,py2,pz2);
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);
348 Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);
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);
354 x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
357 //propagate track to the points of DCA
360 if (!t->PropagateTo(x1,0.,0.)) {
361 Error("PropagateToDCA","Propagation failed !");