1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //-------------------------------------------------------------------------
17 // Implementation of the V0 vertexer class
19 // Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch
20 //-------------------------------------------------------------------------
23 #include <TObjArray.h>
26 #include "AliV0vertex.h"
27 #include "AliV0vertexer.h"
28 #include "AliITStrackV2.h"
30 ClassImp(AliV0vertexer)
32 Int_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;
40 cerr<<"AliV0vertexer::Tracks2V0vertices(): ";
41 cerr<<"file with ITS tracks has not been open !\n";
46 cerr<<"AliV0vertexer::Tracks2V0vertices(): ";
47 cerr<<"file for V0 vertices has not been open !\n";
53 TTree *trkTree=(TTree*)in->Get("TreeT_ITS_0");
54 TBranch *branch=trkTree->GetBranch("tracks");
55 Int_t nentr=(Int_t)trkTree->GetEntries();
57 TObjArray negtrks(nentr/2);
58 TObjArray postrks(nentr/2);
60 Int_t nneg=0, npos=0, nvtx=0;
63 for (i=0; i<nentr; i++) {
64 AliITStrackV2 *iotrack=new AliITStrackV2;
65 branch->SetAddress(&iotrack);
67 if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);}
68 else {npos++; postrks.AddLast(iotrack);}
73 TTree vtxTree("TreeV","Tree with V0 vertices");
74 AliV0vertex *ioVertex=0;
75 vtxTree.Branch("vertices","AliV0vertex",&ioVertex,32000,0);
78 for (i=0; i<nneg; i++) {
79 if (i%10==0) cerr<<nneg-i<<'\r';
80 AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i);
81 for (Int_t k=0; k<npos; k++) {
82 AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k);
84 if (TMath::Abs(ntrk->GetD())<fDNmin) continue;
85 if (TMath::Abs(ptrk->GetD())<fDPmin) continue;
87 AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt;
89 Double_t dca=PropagateToDCA(pnt,ppt);
90 if (dca > fDCAmax) continue;
92 AliV0vertex vertex(*pnt,*ppt);
93 if (vertex.GetChi2() > fChi2max) continue;
95 /* Think of something better here !
96 nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue;
97 pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue;
100 Double_t x,y,z; vertex.GetXYZ(x,y,z);
101 Double_t r2=x*x + y*y;
102 if (r2 > fRmax*fRmax) continue;
103 if (r2 < fRmin*fRmin) continue;
105 Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz);
106 Double_t p2=px*px+py*py+pz*pz;
107 Double_t cost=(x*px+y*py+z*pz)/TMath::Sqrt(p2*(r2+z*z));
108 if (cost<fCPAmax) continue;
110 //vertex.ChangeMassHypothesis(); //default is Lambda0
112 ioVertex=&vertex; vtxTree.Fill();
118 cerr<<"Number of reconstructed V0 vertices: "<<nvtx<<endl;
133 static void External2Helix(const AliITStrackV2 *t, Double_t helix[6]) {
134 //--------------------------------------------------------------------
135 // External track parameters -> helix parameters
136 //--------------------------------------------------------------------
137 Double_t alpha,x,cs,sn;
138 t->GetExternalParameters(x,helix); alpha=t->GetAlpha();
140 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
141 helix[5]=x*cs - helix[0]*sn; // x0
142 helix[0]=x*sn + helix[0]*cs; // y0
144 helix[2]=TMath::ASin(helix[2]) + alpha; // phi0
146 helix[4]=helix[4]/t->GetConvConst(); // C
149 static void Evaluate(const Double_t *h, Double_t t,
150 Double_t r[3], //radius vector
151 Double_t g[3], //first defivatives
152 Double_t gg[3]) //second derivatives
154 //--------------------------------------------------------------------
155 // Calculate position of a point on a track and some derivatives
156 //--------------------------------------------------------------------
157 Double_t phase=h[4]*t+h[2];
158 Double_t sn=TMath::Sin(phase), cs=TMath::Cos(phase);
160 r[0] = h[5] + (sn - h[6])/h[4];
161 r[1] = h[0] - (cs - h[7])/h[4];
162 r[2] = h[1] + h[3]*t;
164 g[0] = cs; g[1]=sn; g[2]=h[3];
166 gg[0]=-h[4]*sn; gg[1]=h[4]*cs; gg[2]=0.;
169 Double_t AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) {
170 //--------------------------------------------------------------------
171 // This function returns the DCA between two tracks
172 // The tracks will be moved to the point of DCA !
173 //--------------------------------------------------------------------
174 Double_t p1[8]; External2Helix(n,p1);
175 p1[6]=TMath::Sin(p1[2]); p1[7]=TMath::Cos(p1[2]);
176 Double_t p2[8]; External2Helix(p,p2);
177 p2[6]=TMath::Sin(p2[2]); p2[7]=TMath::Cos(p2[2]);
180 Double_t r1[3],g1[3],gg1[3]; Double_t t1=0.;
181 Evaluate(p1,t1,r1,g1,gg1);
182 Double_t r2[3],g2[3],gg2[3]; Double_t t2=0.;
183 Evaluate(p2,t2,r2,g2,gg2);
185 Double_t dx=r2[0]-r1[0], dy=r2[1]-r1[1], dz=r2[2]-r1[2];
186 Double_t dm=dx*dx + dy*dy + dz*dz;
190 Double_t gt1=-(dx*g1[0] + dy*g1[1] + dz*g1[2]);
191 Double_t gt2=+(dx*g2[0] + dy*g2[1] + dz*g2[2]);
192 Double_t h11=g1[0]*g1[0] - dx*gg1[0] +
193 g1[1]*g1[1] - dy*gg1[1] +
194 g1[2]*g1[2] - dz*gg1[2];
195 Double_t h22=g2[0]*g2[0] + dx*gg2[0] +
196 g2[1]*g2[1] + dy*gg2[1] +
197 g2[2]*g2[2] + dz*gg2[2];
198 Double_t h12=-(g1[0]*g2[0] + g1[1]*g2[1] + g1[2]*g2[2]);
200 Double_t det=h11*h22-h12*h12;
203 if (TMath::Abs(det)<1.e-33) {
204 //(quasi)singular Hessian
207 dt1=-(gt1*h22 - gt2*h12)/det;
208 dt2=-(h11*gt2 - h12*gt1)/det;
211 if ((dt1*gt1+dt2*gt2)>0) {dt1=-dt1; dt2=-dt2;}
213 //check delta(phase1) ?
214 //check delta(phase2) ?
216 if (TMath::Abs(dt1)/(TMath::Abs(t1)+1.e-3) < 1.e-4)
217 if (TMath::Abs(dt2)/(TMath::Abs(t2)+1.e-3) < 1.e-4) {
218 if ((gt1*gt1+gt2*gt2) > 1.e-4)
219 cerr<<"AliV0vertexer::PropagateToDCA:"
220 " stopped at not a stationary point !\n";
221 Double_t lmb=h11+h22; lmb=lmb-TMath::Sqrt(lmb*lmb-4*det);
223 cerr<<"AliV0vertexer::PropagateToDCA:"
224 " stopped at not a minimum !\n";
229 for (Int_t div=1 ; ; div*=2) {
230 Evaluate(p1,t1+dt1,r1,g1,gg1);
231 Evaluate(p2,t2+dt2,r2,g2,gg2);
232 dx=r2[0]-r1[0]; dy=r2[1]-r1[1]; dz=r2[2]-r1[2];
233 dd=dx*dx + dy*dy + dz*dz;
237 cerr<<"AliV0vertexer::PropagateToDCA: overshoot !\n"; break;
247 if (max<=0) cerr<<"AliV0vertexer::PropagateToDCA: too many iterations !\n";
249 //propagate tracks to the points of DCA
250 Double_t cs=TMath::Cos(n->GetAlpha());
251 Double_t sn=TMath::Sin(n->GetAlpha());
252 Double_t x=r1[0]*cs + r1[1]*sn;
253 if (!n->PropagateTo(x,0.,0.)) {
254 //cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
258 cs=TMath::Cos(p->GetAlpha());
259 sn=TMath::Sin(p->GetAlpha());
260 x=r2[0]*cs + r2[1]*sn;
261 if (!p->PropagateTo(x,0.,0.)) {
262 //cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
266 return TMath::Sqrt(dm);