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