]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliV0vertexer.cxx
Track time measurement (S.Radomski)
[u/mrichter/AliRoot.git] / ITS / AliV0vertexer.cxx
CommitLineData
7f6ddf58 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 V0 vertexer class
18//
19// Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch
20//-------------------------------------------------------------------------
116cbefd 21#include <Riostream.h>
7f6ddf58 22#include <TFile.h>
116cbefd 23#include <TPDGCode.h>
7f6ddf58 24#include <TObjArray.h>
116cbefd 25#include <TTree.h>
7f6ddf58 26
116cbefd 27#include "AliITStrackV2.h"
7f6ddf58 28#include "AliV0vertex.h"
29#include "AliV0vertexer.h"
7f6ddf58 30
31ClassImp(AliV0vertexer)
32
33Int_t AliV0vertexer::Tracks2V0vertices(const TFile *inp, TFile *out) {
34 //--------------------------------------------------------------------
35 //This function reconstructs V0 vertices
36 //--------------------------------------------------------------------
37 TFile *in=(TFile*)inp;
38 TDirectory *savedir=gDirectory;
39
40 if (!in->IsOpen()) {
41 cerr<<"AliV0vertexer::Tracks2V0vertices(): ";
42 cerr<<"file with ITS tracks has not been open !\n";
43 return 1;
44 }
45
46 if (!out->IsOpen()) {
47 cerr<<"AliV0vertexer::Tracks2V0vertices(): ";
48 cerr<<"file for V0 vertices has not been open !\n";
49 return 2;
50 }
51
52 in->cd();
53
5d390eda 54 Char_t name[100];
55 sprintf(name,"TreeT_ITS_%d",fEventN);
56 TTree *trkTree=(TTree*)in->Get(name);
7f6ddf58 57 TBranch *branch=trkTree->GetBranch("tracks");
58 Int_t nentr=(Int_t)trkTree->GetEntries();
59
60 TObjArray negtrks(nentr/2);
61 TObjArray postrks(nentr/2);
62
63 Int_t nneg=0, npos=0, nvtx=0;
64
65 Int_t i;
66 for (i=0; i<nentr; i++) {
67 AliITStrackV2 *iotrack=new AliITStrackV2;
68 branch->SetAddress(&iotrack);
69 trkTree->GetEvent(i);
4ab260d6 70
71 iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.);
72
7f6ddf58 73 if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
74 else {npos++; postrks.AddLast(iotrack);}
75 }
76
77
78 out->cd();
5d390eda 79 sprintf(name,"TreeV%d",fEventN);
80 TTree vtxTree(name,"Tree with V0 vertices");
7f6ddf58 81 AliV0vertex *ioVertex=0;
82 vtxTree.Branch("vertices","AliV0vertex",&ioVertex,32000,0);
83
84
85 for (i=0; i<nneg; i++) {
86 if (i%10==0) cerr<<nneg-i<<'\r';
87 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
a9a2d814 88
5d390eda 89 if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue;
90 if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue;
4ab260d6 91
7f6ddf58 92 for (Int_t k=0; k<npos; k++) {
93 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
94
5d390eda 95 if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue;
96 if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue;
4ab260d6 97
5d390eda 98 if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin)
99 if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue;
4ab260d6 100
7f6ddf58 101
102 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
103
104 Double_t dca=PropagateToDCA(pnt,ppt);
105 if (dca > fDCAmax) continue;
106
107 AliV0vertex vertex(*pnt,*ppt);
108 if (vertex.GetChi2() > fChi2max) continue;
109
110 /* Think of something better here !
111 nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
112 pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
113 */
114
115 Double_t x,y,z; vertex.GetXYZ(x,y,z);
116 Double_t r2=x*x + y*y;
117 if (r2 > fRmax*fRmax) continue;
118 if (r2 < fRmin*fRmin) continue;
119
120 Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
121 Double_t p2=px*px+py*py+pz*pz;
5d390eda 122 Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/
123 TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ)));
4ab260d6 124
5d390eda 125 //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue;
126 if (cost < fCPAmax) continue;
7f6ddf58 127
128 //vertex.ChangeMassHypothesis(); //default is Lambda0
129
130 ioVertex=&vertex; vtxTree.Fill();
131
132 nvtx++;
133 }
134 }
135
136 cerr<<"Number of reconstructed V0 vertices: "<<nvtx<<endl;
137
138 vtxTree.Write();
139
140 negtrks.Delete();
141 postrks.Delete();
142
143 savedir->cd();
144
145 delete trkTree;
146
147 return 0;
148}
149
150
151static void External2Helix(const AliITStrackV2 *t, Double_t helix[6]) {
152 //--------------------------------------------------------------------
153 // External track parameters -> helix parameters
154 //--------------------------------------------------------------------
155 Double_t alpha,x,cs,sn;
156 t->GetExternalParameters(x,helix); alpha=t->GetAlpha();
157
158 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
159 helix[5]=x*cs - helix[0]*sn; // x0
160 helix[0]=x*sn + helix[0]*cs; // y0
161//helix[1]= // z0
162 helix[2]=TMath::ASin(helix[2]) + alpha; // phi0
163//helix[3]= // tgl
164 helix[4]=helix[4]/t->GetConvConst(); // C
165}
166
167static void Evaluate(const Double_t *h, Double_t t,
168 Double_t r[3], //radius vector
169 Double_t g[3], //first defivatives
170 Double_t gg[3]) //second derivatives
171{
172 //--------------------------------------------------------------------
173 // Calculate position of a point on a track and some derivatives
174 //--------------------------------------------------------------------
175 Double_t phase=h[4]*t+h[2];
176 Double_t sn=TMath::Sin(phase), cs=TMath::Cos(phase);
177
178 r[0] = h[5] + (sn - h[6])/h[4];
179 r[1] = h[0] - (cs - h[7])/h[4];
180 r[2] = h[1] + h[3]*t;
181
182 g[0] = cs; g[1]=sn; g[2]=h[3];
183
184 gg[0]=-h[4]*sn; gg[1]=h[4]*cs; gg[2]=0.;
185}
186
187Double_t AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) {
188 //--------------------------------------------------------------------
189 // This function returns the DCA between two tracks
190 // The tracks will be moved to the point of DCA !
191 //--------------------------------------------------------------------
a9a2d814 192 Double_t dy2=n->GetSigmaY2() + p->GetSigmaY2();
193 Double_t dz2=n->GetSigmaZ2() + p->GetSigmaZ2();
194 Double_t dx2=dy2;
195
196 //dx2=dy2=dz2=1.;
197
7f6ddf58 198 Double_t p1[8]; External2Helix(n,p1);
199 p1[6]=TMath::Sin(p1[2]); p1[7]=TMath::Cos(p1[2]);
200 Double_t p2[8]; External2Helix(p,p2);
201 p2[6]=TMath::Sin(p2[2]); p2[7]=TMath::Cos(p2[2]);
202
203
204 Double_t r1[3],g1[3],gg1[3]; Double_t t1=0.;
205 Evaluate(p1,t1,r1,g1,gg1);
206 Double_t r2[3],g2[3],gg2[3]; Double_t t2=0.;
207 Evaluate(p2,t2,r2,g2,gg2);
208
209 Double_t dx=r2[0]-r1[0], dy=r2[1]-r1[1], dz=r2[2]-r1[2];
a9a2d814 210 Double_t dm=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
7f6ddf58 211
212 Int_t max=27;
213 while (max--) {
a9a2d814 214 Double_t gt1=-(dx*g1[0]/dx2 + dy*g1[1]/dy2 + dz*g1[2]/dz2);
215 Double_t gt2=+(dx*g2[0]/dx2 + dy*g2[1]/dy2 + dz*g2[2]/dz2);
216 Double_t h11=(g1[0]*g1[0] - dx*gg1[0])/dx2 +
217 (g1[1]*g1[1] - dy*gg1[1])/dy2 +
218 (g1[2]*g1[2] - dz*gg1[2])/dz2;
219 Double_t h22=(g2[0]*g2[0] + dx*gg2[0])/dx2 +
220 (g2[1]*g2[1] + dy*gg2[1])/dy2 +
221 (g2[2]*g2[2] + dz*gg2[2])/dz2;
222 Double_t h12=-(g1[0]*g2[0]/dx2 + g1[1]*g2[1]/dy2 + g1[2]*g2[2]/dz2);
7f6ddf58 223
224 Double_t det=h11*h22-h12*h12;
225
226 Double_t dt1,dt2;
227 if (TMath::Abs(det)<1.e-33) {
228 //(quasi)singular Hessian
229 dt1=-gt1; dt2=-gt2;
230 } else {
231 dt1=-(gt1*h22 - gt2*h12)/det;
232 dt2=-(h11*gt2 - h12*gt1)/det;
233 }
234
235 if ((dt1*gt1+dt2*gt2)>0) {dt1=-dt1; dt2=-dt2;}
236
237 //check delta(phase1) ?
238 //check delta(phase2) ?
239
240 if (TMath::Abs(dt1)/(TMath::Abs(t1)+1.e-3) < 1.e-4)
241 if (TMath::Abs(dt2)/(TMath::Abs(t2)+1.e-3) < 1.e-4) {
a9a2d814 242 if ((gt1*gt1+gt2*gt2) > 1.e-4/dy2/dy2)
7f6ddf58 243 cerr<<"AliV0vertexer::PropagateToDCA:"
244 " stopped at not a stationary point !\n";
245 Double_t lmb=h11+h22; lmb=lmb-TMath::Sqrt(lmb*lmb-4*det);
246 if (lmb < 0.)
247 cerr<<"AliV0vertexer::PropagateToDCA:"
248 " stopped at not a minimum !\n";
249 break;
250 }
251
252 Double_t dd=dm;
253 for (Int_t div=1 ; ; div*=2) {
254 Evaluate(p1,t1+dt1,r1,g1,gg1);
255 Evaluate(p2,t2+dt2,r2,g2,gg2);
256 dx=r2[0]-r1[0]; dy=r2[1]-r1[1]; dz=r2[2]-r1[2];
a9a2d814 257 dd=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
7f6ddf58 258 if (dd<dm) break;
259 dt1*=0.5; dt2*=0.5;
260 if (div>512) {
261 cerr<<"AliV0vertexer::PropagateToDCA: overshoot !\n"; break;
262 }
263 }
264 dm=dd;
265
266 t1+=dt1;
267 t2+=dt2;
268
269 }
270
271 if (max<=0) cerr<<"AliV0vertexer::PropagateToDCA: too many iterations !\n";
272
273 //propagate tracks to the points of DCA
274 Double_t cs=TMath::Cos(n->GetAlpha());
275 Double_t sn=TMath::Sin(n->GetAlpha());
276 Double_t x=r1[0]*cs + r1[1]*sn;
277 if (!n->PropagateTo(x,0.,0.)) {
278 //cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
279 return 1.e+33;
280 }
281
282 cs=TMath::Cos(p->GetAlpha());
283 sn=TMath::Sin(p->GetAlpha());
284 x=r2[0]*cs + r2[1]*sn;
285 if (!p->PropagateTo(x,0.,0.)) {
286 //cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
287 return 1.e+33;
288 }
289
4ab260d6 290 return TMath::Sqrt(dm*TMath::Sqrt(dy2*dz2));
291 //return TMath::Sqrt(dx*dx + dy*dy + dz*dz);
7f6ddf58 292}
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308