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