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