]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliCascadeVertexer.cxx
Hyperon reconstruction ported to newIO (Y. 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 //
19 //    Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
20 //-------------------------------------------------------------------------
21 #include <Riostream.h>
22 #include <TObjArray.h>
23 #include <TPDGCode.h>
24 #include <TTree.h>
25
26 #include "AliCascadeVertex.h"
27 #include "AliCascadeVertexer.h"
28 #include "AliITStrackV2.h"
29 #include "AliV0vertex.h"
30
31 ClassImp(AliCascadeVertexer)
32
33 Int_t AliCascadeVertexer::
34 V0sTracks2CascadeVertices(TTree *vTree,TTree *tTree, TTree *xTree) {
35   //--------------------------------------------------------------------
36   // This function reconstructs cascade vertices
37   //--------------------------------------------------------------------
38    TBranch *branch=vTree->GetBranch("vertices");
39    if (!branch) {
40       Error("V0sTracks2CascadeVertices","Can't get the V0 branch !");
41       return 1;
42    }   
43    Int_t nentrV0=(Int_t)vTree->GetEntries();
44    
45    TObjArray vtxV0(nentrV0);
46
47    // fill TObjArray vtxV0 with vertices
48
49    Int_t i, nV0=0;
50    for (i=0; i<nentrV0; i++) {
51
52        AliV0vertex *ioVertex=new AliV0vertex;
53        branch->SetAddress(&ioVertex);
54        vTree->GetEvent(i);
55        nV0++; 
56        vtxV0.AddLast(ioVertex);
57        
58    }
59
60    branch=tTree->GetBranch("tracks");
61    if (!branch) {
62       Error("V0sTracks2CascadeVertices","Can't get the track branch !");
63       return 2;
64    }
65    Int_t nentr=(Int_t)tTree->GetEntries();
66
67    TObjArray trks(nentr);
68
69    // fill TObjArray trks with tracks
70
71    Int_t ntrack=0;
72
73    for (i=0; i<nentr; i++) {
74
75        AliITStrackV2 *iotrack=new AliITStrackV2;
76        branch->SetAddress(&iotrack);
77        tTree->GetEvent(i);
78
79        iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.);
80
81        ntrack++; trks.AddLast(iotrack);
82        
83    }   
84
85   AliCascadeVertex *ioCascade=0;
86   branch=xTree->GetBranch("cascades");
87   if (!branch) xTree->Branch("cascades","AliCascadeVertex",&ioCascade,32000,3);
88   else branch->SetAddress(&ioCascade); 
89
90    // loop on all vertices
91
92    Double_t massLambda=1.11568;
93    Int_t ncasc=0;
94
95    for (i=0; i<nV0; i++) {
96
97        AliV0vertex *V0ver=(AliV0vertex *)vtxV0.UncheckedAt(i);
98
99        V0ver->ChangeMassHypothesis(kLambda0); //I.B.
100
101        if (V0ver->GetEffMass()<massLambda-fMassWin ||       // condition of the V0 mass window (cut fMassWin)
102            V0ver->GetEffMass()>massLambda+fMassWin) continue; 
103
104        if (V0ver->GetD(0,0,0)<fDV0min) continue;          // condition of minimum impact parameter of the V0 (cut fDV0min) 
105                                                           // here why not cuting on pointing angle ???
106
107    // for each vertex in the good mass range, loop on all tracks (= bachelor candidates)
108
109        for (Int_t k=0; k<ntrack; k++) {
110
111           AliITStrackV2 *bachtrk=(AliITStrackV2 *)trks.UncheckedAt(k);
112
113           if (TMath::Abs(bachtrk->GetD())<fDBachMin) continue;        // eliminate to small impact parameters
114
115           if (V0ver->GetPdgCode()==kLambda0 && bachtrk->Get1Pt()<0.) continue;     // condition on V0 label 
116           if (V0ver->GetPdgCode()==kLambda0Bar && bachtrk->Get1Pt()>0.) continue;  // + good sign for bachelor
117           
118           AliV0vertex V0(*V0ver), *pV0=&V0;
119           AliITStrackV2 bt(*bachtrk), *pbt=&bt;
120
121    // calculation of the distance of closest approach between the V0 and the bachelor
122
123           Double_t dca=PropagateToDCA(pV0,pbt);
124           if (dca > fDCAmax) continue;                         // cut on dca
125
126    // construction of a cascade object
127
128           AliCascadeVertex cascade(*pV0,*pbt);
129           if (cascade.GetChi2() > fChi2max) continue;
130
131    // get cascade decay position (V0, bachelor)
132           
133          Double_t x,y,z; cascade.GetXYZ(x,y,z); 
134          Double_t r2=x*x + y*y; 
135          if (r2 > fRmax*fRmax) continue;   // condition on fiducial zone
136          if (r2 < fRmin*fRmin) continue;
137
138          {
139    //I.B.
140          Double_t x1,y1,z1; V0ver->GetXYZ(x1,y1,z1);
141          if (r2 > (x1*x1+y1*y1)) continue;
142          if (z*z > z1*z1) continue;
143          }
144
145    // get cascade momentum
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
153          //cascade.ChangeMassHypothesis(); //default is Xi
154
155
156          ioCascade=&cascade; xTree->Fill();
157
158          ncasc++;
159
160       }
161    }
162
163    cerr<<"Number of reconstructed cascades: "<<ncasc<<endl;
164
165    trks.Delete();
166    vtxV0.Delete();
167
168    return 0;
169 }
170
171 inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){
172   // determinant 2x2
173   return a00*a11 - a01*a10;
174 }
175
176 inline Double_t det (Double_t a00,Double_t a01,Double_t a02,
177                      Double_t a10,Double_t a11,Double_t a12,
178                      Double_t a20,Double_t a21,Double_t a22) {
179   // determinant 3x3
180   return 
181   a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21);
182 }
183
184
185
186
187 Double_t 
188 AliCascadeVertexer::PropagateToDCA(AliV0vertex *v, AliITStrackV2 *t) {
189   //--------------------------------------------------------------------
190   // This function returns the DCA between the V0 and the track
191   //--------------------------------------------------------------------
192
193   Double_t phi, x, par[5];
194   Double_t alpha, cs1, sn1;
195
196   t->GetExternalParameters(x,par); alpha=t->GetAlpha();
197   phi=TMath::ASin(par[2]) + alpha;  
198   Double_t px1=TMath::Cos(phi), py1=TMath::Sin(phi), pz1=par[3];
199
200  
201   cs1=TMath::Cos(alpha); sn1=TMath::Sin(alpha);
202   Double_t x1=x*cs1 - par[0]*sn1;
203   Double_t y1=x*sn1 + par[0]*cs1;
204   Double_t z1=par[1];
205
206   
207   Double_t x2,y2,z2;     // position and momentum of V0
208   Double_t px2,py2,pz2;
209   
210   v->GetXYZ(x2,y2,z2);
211   v->GetPxPyPz(px2,py2,pz2);
212  
213 // calculation dca
214    
215   Double_t dd= det(x2-x1,y2-y1,z2-z1,px1,py1,pz1,px2,py2,pz2);
216   Double_t ax= det(py1,pz1,py2,pz2);
217   Double_t ay=-det(px1,pz1,px2,pz2);
218   Double_t az= det(px1,py1,px2,py2);
219
220   Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);
221
222 //points of the DCA
223   Double_t t1 = det(x2-x1,y2-y1,z2-z1,px2,py2,pz2,ax,ay,az)/
224                 det(px1,py1,pz1,px2,py2,pz2,ax,ay,az);
225   
226   x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
227   
228
229   //propagate track to the points of DCA
230
231   x1=x1*cs1 + y1*sn1;
232   if (!t->PropagateTo(x1,0.,0.)) {
233     cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
234     return 1.e+33;
235   }  
236
237   return dca;
238 }
239
240
241
242
243
244
245
246
247
248
249