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