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