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