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