]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCTracklet.cxx
51a2b059ff7f5213cdeb2524ef893d962b27034a
[u/mrichter/AliRoot.git] / TPC / AliTPCTracklet.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 // This class stores a tracklet (a track that lives only in a single TPC
18 // sector). Its objects can be constructed out of TPCseeds, that are
19 // holding the necessary cluster information.
20 ////
21 ////
22 //// 
23
24 #include "AliTPCTracklet.h"
25 #include "TObjArray.h"
26 #include "TLinearFitter.h"
27 #include "AliTPCseed.h"
28 #include "AliESDVertex.h"
29 #include "AliTracker.h"
30 #include "TTreeStream.h"
31 #include "TRandom3.h"
32 #include "TDecompChol.h"
33
34 #include <iostream>
35 using namespace std;
36
37 ClassImp(AliTPCTracklet)
38
39 const Double_t AliTPCTracklet::kB2C=0.299792458e-3;
40 Float_t  AliTPCTracklet::fgEdgeCutY=3;
41 Float_t  AliTPCTracklet::fgEdgeCutX=0;
42
43 AliTPCTracklet::AliTPCTracklet() 
44   : fNClusters(0),fNStoredClusters(0),fClusters(0),fSector(-1),fOuter(0),
45     fInner(0),fPrimary(0) {
46   ////
47   // The default constructor. It is intended to be used for I/O only.
48   ////
49 }
50
51 AliTPCTracklet::AliTPCTracklet(const AliTPCseed *track,Int_t sector,
52                                TrackType type,Bool_t storeClusters)
53   : fNClusters(0),fNStoredClusters(0),fClusters(0),fSector(sector),fOuter(0),
54     fInner(0),fPrimary(0) {
55   ////
56   // Contructor for a tracklet out of a track. Only clusters within a given 
57   // sector are used.
58   ///
59
60   //TODO: only kalman works
61   
62   for (Int_t i=0;i<160;++i) {
63     AliTPCclusterMI *c=track->GetClusterPointer(i);
64     if (c && RejectCluster(c)) continue;
65     if (c&&c->GetDetector()==sector)
66       ++fNClusters;
67   }
68
69   if (storeClusters) {
70     fClusters=new AliTPCclusterMI[fNClusters];
71     for (Int_t i=0;i<160;++i) {
72       AliTPCclusterMI *c=track->GetClusterPointer(i);
73       if (c && RejectCluster(c)) continue;
74       if (c&&c->GetDetector()==sector)
75         fClusters[fNStoredClusters]=*c;
76       ++fNStoredClusters;
77     }
78   }
79
80   switch (type) {
81   case kKalman:
82     FitKalman(track,sector);
83     break;
84   case kLinear:
85   case kQuadratic:
86     FitLinear(track,sector,type);
87     break;
88   case kRiemann:
89     FitRiemann(track,sector);
90     break;
91   }
92
93 }
94
95 AliTPCTracklet::AliTPCTracklet(const TObjArray &/*clusters*/,Int_t sector,
96                                TrackType /*type*/,Bool_t /*storeClusters*/)
97   : fNClusters(0),fNStoredClusters(0),fClusters(0),fSector(sector),fOuter(0),
98     fInner(0),fPrimary(0) {
99   //TODO: write it!
100 }
101
102 AliTPCTracklet::AliTPCTracklet(const AliTPCTracklet &t)
103   : TObject(t),fNClusters(t.fNClusters),fNStoredClusters(t.fNStoredClusters),fClusters(0),
104     fSector(t.fSector),fOuter(0),fInner(0),
105     fPrimary(0) {
106   ////
107   // The copy constructor. You can copy tracklets! 
108   ////
109
110   if (t.fClusters) {
111     fClusters=new AliTPCclusterMI[t.fNStoredClusters];
112     for (int i=0;i<t.fNStoredClusters;++i)
113       fClusters[i]=t.fClusters[i];
114   }
115   if (t.fOuter)
116     fOuter=new AliExternalTrackParam(*t.fOuter);
117   if (t.fInner)
118     fInner=new AliExternalTrackParam(*t.fInner);
119   if (t.fPrimary)
120     fPrimary=new AliExternalTrackParam(*t.fPrimary);
121 }
122
123 AliTPCTracklet& AliTPCTracklet::operator=(const AliTPCTracklet &t) {
124   ////
125   // The assignment constructor. You can assign tracklets!
126   ////
127   if (this!=&t) {
128     fNClusters=t.fNClusters;
129     fNStoredClusters=fNStoredClusters;
130     delete fClusters;
131     if (t.fClusters) {
132       fClusters=new AliTPCclusterMI[t.fNStoredClusters];
133       for (int i=0;i<t.fNStoredClusters;++i)
134         fClusters[i]=t.fClusters[i];
135     }
136     else
137       fClusters=0;
138     fSector=t.fSector;
139     if (t.fOuter) {
140       if (fOuter)
141         *fOuter=*t.fOuter;
142       else
143         fOuter=new AliExternalTrackParam(*t.fOuter);
144     }
145     else {
146       delete fOuter;
147       fOuter=0;
148     }
149
150     if (t.fInner) {
151       if (fInner)
152         *fInner=*t.fInner;
153       else
154         fInner=new AliExternalTrackParam(*t.fInner);
155     }
156     else {
157       delete fInner;
158       fInner=0;
159     }
160
161     if (t.fPrimary) {
162       if (fPrimary)
163         *fPrimary=*t.fPrimary;
164       else
165         fPrimary=new AliExternalTrackParam(*t.fPrimary);
166     }
167     else {
168       delete fPrimary;
169       fPrimary=0;
170     }
171   }
172   return *this;
173 }
174
175 AliTPCTracklet::~AliTPCTracklet() {
176   //
177   // The destructor. Yes, you can even destruct tracklets.
178   //
179   delete fClusters;
180   delete fOuter;
181   delete fInner;
182   delete fPrimary;
183 }
184
185 void AliTPCTracklet::FitKalman(const AliTPCseed *track,Int_t sector) {
186   //
187   // Fit using Kalman filter
188   //
189   AliTPCseed *t=new AliTPCseed(*track);
190   if (!t->Rotate(TMath::DegToRad()*(sector%18*20.+10.)-t->GetAlpha())) {
191     delete t;
192     return;
193   }
194   // fit from inner to outer row
195   AliTPCseed *outerSeed=new AliTPCseed(*t);
196   Int_t n=0;
197   for (Int_t i=0;i<160;++i) {
198     AliTPCclusterMI *c=t->GetClusterPointer(i);
199     if (c && RejectCluster(c,outerSeed)) continue;
200     if (c&&c->GetDetector()==sector) {
201       if (n==1) {
202         outerSeed->ResetCovariance(100.);
203       }
204       ++n;
205       Double_t r[3]={c->GetX(),c->GetY(),c->GetZ()};
206       Double_t cov[3]={0.1,0.,0.1}; //TODO: correct error parametrisation
207       if (!outerSeed->PropagateTo(r[0]) ||
208           !static_cast<AliExternalTrackParam*>(outerSeed)->Update(&r[1],cov)) {
209         delete outerSeed;
210         outerSeed=0;
211         break;
212       }
213     }
214   }
215   if (outerSeed)
216     fOuter=new AliExternalTrackParam(*outerSeed);
217   delete outerSeed;
218   // fit from outer to inner rows
219   AliTPCseed *innerSeed=new AliTPCseed(*t);
220   n=0;
221   for (Int_t i=159;i>=0;--i) {
222     AliTPCclusterMI *c=t->GetClusterPointer(i);
223     if (c && RejectCluster(c, innerSeed)) continue;
224     if (c&&c->GetDetector()==sector) {
225       if (n==1) {
226         innerSeed->ResetCovariance(100.);
227       }
228       ++n;
229       Double_t r[3]={c->GetX(),c->GetY(),c->GetZ()};
230       Double_t cov[3]={0.1,0.,0.1};
231       if (!innerSeed->PropagateTo(r[0]) ||
232           !static_cast<AliExternalTrackParam*>(innerSeed)->Update(&r[1],cov)) {
233         delete innerSeed;
234         innerSeed=0;
235         break;
236       }
237     }
238   }
239   if (innerSeed)
240     fInner=new AliExternalTrackParam(*innerSeed);
241   // propagate to the primary vertex
242   if (innerSeed) {
243     AliTPCseed *primarySeed=new AliTPCseed(*innerSeed);
244     Double_t pos[]={0.,0.,0.};
245     Double_t sigma[]={.1,.1,.1}; //TODO: is this correct?
246     AliESDVertex vertex(pos,sigma);
247     if (primarySeed->PropagateToVertex(&vertex))
248       fPrimary=new AliExternalTrackParam(*primarySeed);
249     delete primarySeed;
250     // for better comparison one does not want to have alpha changed...
251     if (fPrimary) if (!fPrimary->Rotate(fInner->GetAlpha())) {
252       delete fPrimary;
253       fPrimary=0;
254     }
255   }
256   delete innerSeed;
257
258   delete t;
259 }
260
261 void AliTPCTracklet::FitLinear(const AliTPCseed *track,Int_t sector,
262                                TrackType type) {
263   TLinearFitter fy(1);
264   TLinearFitter fz(1);
265   fy.StoreData(kFALSE);
266   fz.StoreData(kFALSE);
267   switch (type) {
268   case kLinear:
269     fy.SetFormula("1 ++ x");
270     fz.SetFormula("1 ++ x");
271     break;
272   case kQuadratic:
273     fy.SetFormula("1 ++ x ++ x*x");
274     fz.SetFormula("1 ++ x");
275     break;
276   case kKalman:
277   case kRiemann:
278     break;
279   }
280   Double_t xmax=-1.;
281   Double_t xmin=1000.;
282   for (Int_t i=0;i<160;++i) {
283     AliTPCclusterMI *c=track->GetClusterPointer(i);
284     if (c && RejectCluster(c)) continue;
285     if (c&&c->GetDetector()==sector) {
286       Double_t x=c->GetX();
287       fy.AddPoint(&x,c->GetY());
288       fz.AddPoint(&x,c->GetZ());
289       xmax=TMath::Max(xmax,x);
290       xmin=TMath::Min(xmin,x);
291     }
292   }
293   fy.Eval();
294   fz.Eval();
295   Double_t a[3]={fy.GetParameter(0),
296                  fy.GetParameter(1),
297                  type==kQuadratic?fy.GetParameter(2):0.};
298   Double_t ca[6]={fy.GetCovarianceMatrixElement(0,0),
299                   fy.GetCovarianceMatrixElement(1,0),
300                   fy.GetCovarianceMatrixElement(1,1),
301                   type==kQuadratic?fy.GetCovarianceMatrixElement(2,0):0.,
302                   type==kQuadratic?fy.GetCovarianceMatrixElement(2,1):0.,
303                   type==kQuadratic?fy.GetCovarianceMatrixElement(2,2):0.};
304   for (int i=0;i<6;++i) ca[i]*=fy.GetChisquare()/fNClusters;
305   Double_t b[2]={fz.GetParameter(0),
306                  fz.GetParameter(1)};
307   Double_t cb[3]={fz.GetCovarianceMatrixElement(0,0),
308                   fz.GetCovarianceMatrixElement(1,0),
309                   fz.GetCovarianceMatrixElement(1,1)};
310   for (int i=0;i<3;++i) cb[i]*=fz.GetChisquare()/fNClusters;
311   Double_t p[5];
312   Double_t c[15];
313   Double_t alpha=track->GetAlpha();
314   Quadratic2Helix(a,ca,b,cb,0.,p,c);
315   fPrimary=new AliExternalTrackParam(0.,alpha,p,c);
316   Quadratic2Helix(a,ca,b,cb,xmin,p,c);
317   fInner=new AliExternalTrackParam(xmin,alpha,p,c);
318   Quadratic2Helix(a,ca,b,cb,xmax,p,c);
319   fOuter=new AliExternalTrackParam(xmax,alpha,p,c);
320 }
321   
322 void AliTPCTracklet::Quadratic2Helix(Double_t *a,Double_t *ca,
323                                      Double_t *b,Double_t *cb,
324                                      Double_t x0,
325                                      Double_t *p,Double_t *c) {
326   // y(x)=a[0]+a[1]*x+a[2]*x^2
327   // z(x)=b[0]+b[1]*x
328   // parametrises the corosponding helix at x0
329
330   // get the polynoms at x0
331   Double_t a0=x0*x0*a[2] + x0*a[1] + a[0];
332   Double_t a1=2.*x0*a[2] +     a[1];
333   Double_t a2=      a[2];
334   Double_t ca00=ca[0]+x0*(2.*ca[1]+x0*(ca[2]+2.*ca[3]+x0*(2.*ca[4]+x0*ca[5])));
335   Double_t ca10=ca[1]+x0*(ca[2]+2.*ca[3]+x0*(3.*ca[4]+x0*2.*ca[5]));
336   Double_t ca11=ca[2]+x0*4.*(ca[4]+x0*ca[5]);
337   Double_t ca20=ca[3]+x0*(ca[4]+x0*ca[5]);
338   Double_t ca21=ca[3]+x0*2.*ca[5];
339   Double_t ca22=ca[5];
340
341   Double_t b0=x0*b[1] + b[0];
342   Double_t b1=   b[1];
343   Double_t cb00=cb[0]+x0*(2.*cb[1]+x0*cb[2]);
344   Double_t cb10=cb[1]+x0*cb[2];
345   Double_t cb11=cb[2];
346
347   // transform to helix parameters
348   Double_t f   =1.+a1*a1;
349   Double_t f2  =f*f;
350   Double_t fi  =1./f; 
351   Double_t fi12=TMath::Sqrt(fi);
352   Double_t fi32=fi*fi12;
353   Double_t fi2 =fi*fi;
354   Double_t fi52=fi2*fi12;
355   Double_t fi3 =fi2*fi;
356   Double_t fi5 =fi2*fi3;
357   
358   Double_t xyz[3]={0.}; // TODO...
359   Double_t fc=1./(GetBz(xyz)*kB2C);
360
361   p[0]=a0;            // y0
362   p[1]=b0;            // z0
363   p[2]=a1*fi12;       // snp
364   p[3]=b1;            // tgl
365   p[4]=2.*a2*fi32*fc; // 1/pt
366
367   c[0] =ca00;      //  y0-y0
368   c[1] =0.;        //  z0-y0
369   c[2] =cb00;      //  z0-z0
370   c[3] =ca10*fi32; // snp-y0
371   c[4] =0.;        // snp-z0
372   c[5] =ca11*fi3;  // snp-snp
373   c[6] =0.;        // tgl-y0
374   c[7] =cb10;      // tgl-z0
375   c[8] =0.;        // tgl-snp
376   c[9] =cb11;      // tgl-tgl
377   c[10]=2.*(-3.*a1*a2*ca10+f*ca20)*fi3*fc;  // 1/pt-y0
378   c[11]=0.;                                 // 1/pt-z0
379   c[12]=2.*(-3.*a1*a2*ca11+f*ca21)*fi52*fc; // 1/pt-snp
380   c[13]=0.;                                 // 1/pt-tgl
381   c[14]=(-12.*a1*a2*(-3.*a1*a2*ca11+2.*f*ca21)+4.*f2*ca22)*fi5
382     *fc*fc;        // 1/pt-1/pt
383 }
384
385
386 void AliTPCTracklet::FitRiemann(const AliTPCseed *track,Int_t sector) {
387   TLinearFitter fy(2);
388   fy.StoreData(kFALSE);
389   fy.SetFormula("hyp2");
390   Double_t xmax=-1.;
391   Double_t xmin=1000.;
392   for (Int_t i=0;i<160;++i) {
393     AliTPCclusterMI *c=track->GetClusterPointer(i);
394     if (c && RejectCluster(c)) continue;
395     if (c&&c->GetDetector()==sector) {
396       Double_t x=c->GetX();
397       Double_t y=c->GetY();
398       Double_t xy[2]={x,y};
399       Double_t r=x*x+y*y;
400       Double_t errx=1.,erry=1.;//TODO!
401       Double_t err=TMath::Sqrt(4.*x*x*errx+4.*y*y*erry);
402       err=1.;
403       fy.AddPoint(xy,r,err);
404       xmax=TMath::Max(xmax,x);
405       xmin=TMath::Min(xmin,x);
406     }
407   }
408   fy.Eval();
409   Double_t a[3]={fy.GetParameter(0),
410                  fy.GetParameter(1),
411                  fy.GetParameter(2)};
412   Double_t ca[6]={fy.GetCovarianceMatrixElement(0,0),
413                   fy.GetCovarianceMatrixElement(1,0),
414                   fy.GetCovarianceMatrixElement(1,1),
415                   fy.GetCovarianceMatrixElement(2,0),
416                   fy.GetCovarianceMatrixElement(2,1),
417                   fy.GetCovarianceMatrixElement(2,2)};
418
419   TLinearFitter fz(1);
420   fz.StoreData(kFALSE);
421   fz.SetFormula("hyp1");
422   Double_t R=.5*TMath::Sqrt(4.*a[0]+a[1]*a[1]+a[2]*a[2]);
423   Double_t oldx=0.;
424   Double_t oldy=R;
425   Double_t phi=0.;
426   for (Int_t i=0;i<160;++i) {
427     AliTPCclusterMI *c=track->GetClusterPointer(i);
428     if (c && RejectCluster(c)) continue;
429     if (c&&c->GetDetector()==sector) {
430       Double_t x=c->GetX();
431       Double_t y=c->GetY();
432       Double_t dx=x-oldx;
433       Double_t dy=y-oldy;
434       phi+=2.*TMath::Abs(TMath::ATan2(.5*TMath::Sqrt(dx*dx+dy*dy),R));
435       Double_t err=1.;
436       fz.AddPoint(&phi,c->GetZ(),err);
437       oldx=x;
438       oldy=y;
439     }
440   }
441   fz.Eval();
442   Double_t b[2]={fz.GetParameter(0),
443                  fz.GetParameter(1)};
444   Double_t cb[3]={fz.GetCovarianceMatrixElement(0,0),
445                   fz.GetCovarianceMatrixElement(1,0),
446                   fz.GetCovarianceMatrixElement(1,1)};
447
448   Double_t p[5];
449   Double_t c[15];
450   Double_t alpha=track->GetAlpha();
451   if (Riemann2Helix(a,ca,b,cb,0.,p,c))
452     fPrimary=new AliExternalTrackParam(0.,alpha,p,c);
453   if (Riemann2Helix(a,ca,b,cb,xmin,p,c))
454     fInner=new AliExternalTrackParam(xmin,alpha,p,c);
455   if (Riemann2Helix(a,ca,b,cb,xmax,p,c))
456     fOuter=new AliExternalTrackParam(xmax,alpha,p,c);
457 }
458
459 Bool_t AliTPCTracklet::Riemann2Helix(Double_t *a,Double_t */*ca*/,
460                                      Double_t *b,Double_t */*cb*/,
461                                      Double_t x0,
462                                      Double_t *p,Double_t *c) {
463   //TODO: signs!
464
465   Double_t xr0=.5*a[1];
466   Double_t yr0=.5*a[2];
467   Double_t R=.5*TMath::Sqrt(4.*a[0]+a[1]*a[1]+a[2]*a[2]);
468   Double_t dx=x0-xr0;
469   if (dx*dx>=R*R) return kFALSE;
470   Double_t dy=TMath::Sqrt(R*R-dx*dx); //sign!!
471   if (TMath::Abs(yr0+dy)>TMath::Abs(yr0-dy))
472     dy=-dy;
473   Double_t y0=yr0+dy; 
474   Double_t tgp=-dx/dy; //TODO: dy!=0
475   Double_t z0=b[0]+TMath::ATan(tgp)*b[1];
476   Double_t xyz[3]={x0,y0,z0};
477   Double_t fc=1./(GetBz(xyz)*kB2C);
478   fc=1;
479   p[0]=y0;  // y0
480   p[1]=z0; // z0
481   p[2]=tgp/TMath::Sqrt(1.+tgp*tgp); // snp
482   p[3]=b[1];       // tgl
483   p[4]=1./R*fc;    // 1/pt
484
485   c[0] =0.;      //  y0-y0
486   c[1] =0.;        //  z0-y0
487   c[2] =0.;      //  z0-z0
488   c[3] =0.; // snp-y0
489   c[4] =0.;        // snp-z0
490   c[5] =0.;  // snp-snp
491   c[6] =0.;        // tgl-y0
492   c[7] =0.;      // tgl-z0
493   c[8] =0.;        // tgl-snp
494   c[9] =0.;      // tgl-tgl
495   c[10]=0.;  // 1/pt-y0
496   c[11]=0.;                                 // 1/pt-z0
497   c[12]=0.; // 1/pt-snp
498   c[13]=0.;                                 // 1/pt-tgl
499   c[14]=0.;        // 1/pt-1/pt
500
501   return kTRUE;
502 }
503
504 TObjArray AliTPCTracklet::CreateTracklets(const AliTPCseed *track,
505                                           TrackType type,
506                                           Bool_t storeClusters,
507                                           Int_t minClusters,
508                                           Int_t maxTracklets) {
509 // The tracklet factory: It creates several tracklets out of a track. They
510 // are created for sectors that fullfill the constraint of having enough
511 // clusters inside. Futhermore you can specify the maximum amount of
512 // tracklets that are to be created.
513 // The tracklets appear in a sorted fashion, beginning with those having the
514 // most clusters.
515
516   Int_t sectors[72]={0};
517   for (Int_t i=0;i<160;++i) {
518     AliTPCclusterMI *c=track->GetClusterPointer(i);
519     if (c && RejectCluster(c)) continue;
520     if (c)
521       ++sectors[c->GetDetector()];
522   }
523   Int_t indices[72];
524   TMath::Sort(72,sectors,indices);
525   TObjArray tracklets;
526   if (maxTracklets>72) maxTracklets=72; // just to protect against "users".
527   for (Int_t i=0;i<maxTracklets&&sectors[indices[i]]>=minClusters;++i) {
528     tracklets.Add(new AliTPCTracklet(track,indices[i],type,storeClusters));
529   }
530   return tracklets;
531 }
532
533 TObjArray AliTPCTracklet::CreateTracklets(const TObjArray &/*clusters*/,
534                                           TrackType /*type*/,
535                                           Bool_t /*storeClusters*/,
536                                           Int_t /*minClusters*/,
537                                           Int_t /*maxTracklets*/) {
538   // TODO!
539
540   TObjArray tracklets;
541   return tracklets;
542 }
543
544 Bool_t AliTPCTracklet::PropagateToMeanX(const AliTPCTracklet &t1,
545                                         const AliTPCTracklet &t2,
546                                         AliExternalTrackParam *&t1m,
547                                         AliExternalTrackParam *&t2m) {
548   // This function propagates two Tracklets to a common x-coordinate. This
549   // x is dermined as the one that is in the middle of the two tracklets (they
550   // are assumed to live on two distinct x-intervalls).
551   // The inner parametrisation of the outer Tracklet and the outer 
552   // parametrisation of the inner Tracklet are used and propagated to this
553   // common x. This result is saved not inside the Tracklets but two new
554   // ExternalTrackParams are created (that means you might want to delete
555   // them).
556   // In the case that the alpha angles of the Tracklets differ both angles
557   // are tried out for this propagation.
558   // In case of any failure kFALSE is returned, no AliExternalTrackParam
559   // is created und the pointers are set to 0.
560
561   if (t1.GetInner() && t1.GetOuter() && 
562       t2.GetInner() && t2.GetOuter()) {
563     if (t1.GetOuter()->GetX()<t2.GetInner()->GetX()) {
564       t1m=new AliExternalTrackParam(*t1.GetOuter());
565       t2m=new AliExternalTrackParam(*t2.GetInner());
566     }
567     else {
568       t1m=new AliExternalTrackParam(*t1.GetInner());
569       t2m=new AliExternalTrackParam(*t2.GetOuter());
570     }
571     Double_t mx=.5*(t1m->GetX()+t2m->GetX());
572     Double_t b1,b2;
573     Double_t xyz[3];
574     t1m->GetXYZ(xyz);
575     b1=GetBz(xyz);
576     t2m->GetXYZ(xyz);
577     b2=GetBz(xyz);
578     if (t1m->Rotate(t2m->GetAlpha()) 
579         && t1m->PropagateTo(mx,b1) 
580         && t2m->PropagateTo(mx,b2));
581     else
582       if (t2m->Rotate(t1m->GetAlpha())
583           && t1m->PropagateTo(mx,b1) 
584           && t2m->PropagateTo(mx,b2));
585       else {
586         delete t1m;
587         delete t2m;
588         t1m=t2m=0;
589       }
590   }
591   else {
592     t1m=t2m=0;
593   }
594   return t1m&&t2m;
595 }
596
597 double AliTPCTracklet::GetBz(Double_t *xyz) {
598   if (AliTracker::UniformField())
599     return AliTracker::GetBz();
600   else
601     return AliTracker::GetBz(xyz);
602 }
603
604 void AliTPCTracklet::RandomND(Int_t ndim,const Double_t *p,const Double_t *c,
605                               Double_t *x) {
606   // This function generates a n-dimensional random variable x with mean
607   // p and covariance c.
608   // That is done using the cholesky decomposition of the covariance matrix,
609   // Begin_Latex C=U^{t} U End_Latex, with Begin_Latex U End_Latex being an
610   // upper triangular matrix. Given a vector v of iid gausian random variables
611   // with variance 1 one obtains the asked result as: Begin_Latex x=U^t v 
612   // End_Latex.
613   // c is expected to be in a lower triangular format:
614   // c[0]
615   // c[1] c[2]
616   // c[3] c[4] c[5]
617   // etc.
618   static TRandom3 random;
619   Double_t *c2= new Double_t[ndim*ndim];
620   Int_t k=0;
621   for (Int_t i=0;i<ndim;++i)
622     for (Int_t j=0;j<=i;++j)
623       c2[i*ndim+j]=c2[j*ndim+i]=c[k++];
624   TMatrixDSym cm(ndim,c2);
625   delete[] c2;
626   TDecompChol chol(cm);
627   chol.Decompose();
628   const TVectorD pv(ndim);
629   const_cast<TVectorD*>(&pv)->Use(ndim,const_cast<Double_t*>(p));
630   TVectorD xv(ndim);
631   xv.Use(ndim,x);
632   for (Int_t i=0;i<ndim;++i)
633     xv[i]=random.Gaus();
634   TMatrixD L=chol.GetU();
635   L.T();
636   xv=L*xv+pv;
637 }
638
639 TEllipse AliTPCTracklet::ErrorEllipse(Double_t x,Double_t y,
640                                       Double_t sx,Double_t sy,Double_t sxy) {
641   /* Begin_Latex
642      r_{1,2}=1/2 (a+c#pm#sqrt{(a-c)^{2}+(2b)^{2}})
643   End_Latex */
644   Double_t det1=1./(sx*sy-sxy*sxy);
645   Double_t a=sy*det1;
646   Double_t b=-sxy*det1;
647   Double_t c=sx*det1;
648   Double_t d=c-a;
649   Double_t s=TMath::Sqrt(d*d+4.*b*b);
650   Double_t r1=TMath::Sqrt(.5*(a+c-s));
651   Double_t r2=TMath::Sqrt(.5*(a+c+s));
652   Double_t alpha=.5*TMath::ATan2(2.*b,d);
653   return TEllipse(x,y,r1,r2,0.,360.,alpha*TMath::RadToDeg());
654 }
655
656 void AliTPCTracklet::Test(const char* filename) {
657   /*
658     aliroot
659     AliTPCTracklet::Test("");
660     TFile f("AliTPCTrackletDebug.root");
661     TTree *t=f.Get("AliTPCTrackletDebug");
662     t->Draw("p0:p4");
663     TEllipse e=AliTPCTracklet::ErrorEllipse(0.,0.,4.,1.,1.8);
664     e.Draw();
665  */
666   TTreeSRedirector ds(filename);
667   Double_t p[5]={0.};
668   Double_t c[15]={4.,
669                   0.,4.,
670                   0.,0.,9.,
671                   0.,0.,0.,16.,
672                   1.8,0.,0.,0.,1.};
673   for (Int_t i=0;i<10000;++i) {
674     Double_t x[5];
675     RandomND(5,p,c,x);
676     ds<<"AliTPCTrackletDebug"
677       <<"p0="<<x[0]
678       <<"p1="<<x[1]
679       <<"p2="<<x[2]
680       <<"p3="<<x[3]
681       <<"p4="<<x[4]
682       <<"\n";
683   }
684
685   /*
686   Double_t b;
687   Double_t x=0.;
688   Double_t alpha=0.;
689   Double_t param[5]={0.};
690   Double_t covar[15]={1.,
691                       0.,1.,
692                       0.,0.,1.,
693                       0.,0.,0.,1.,
694                       0.,0.,0.,0.,1.};
695   AliExternalTrackParam track(x,alpha,param,covar);
696
697   
698
699   for (Int_t i=0;i<points.GetNPoints();++i) {
700     Double_t x=0.;
701     Double_t alpha=0.;
702     Double_t param[5]={0.};
703     Double_t covar[15]={1.,
704                         0.,1.,
705                         0.,0.,1.,
706                         0.,0.,0.,1.,
707                         0.,0.,0.,0.,1.};
708     AliExternalTrackParam track(x,alpha,param,covar);
709     for (x=90.;x<250.;x+=1.) {
710       track.PropagateTo(x,b);
711       AliTPCclusterMI c();
712     }
713   }
714   */
715 }
716
717
718 Bool_t AliTPCTracklet::RejectCluster(AliTPCclusterMI* cl, AliExternalTrackParam * param){
719   //
720   // check the acceptance of cluster
721   // Cut on edge effects
722   //
723   Bool_t isReject = kFALSE;
724   Float_t edgeY = cl->GetX()*TMath::Tan(TMath::Pi()/18);
725   Float_t dist  = edgeY - TMath::Abs(cl->GetY());
726   if (param)  dist  = edgeY - TMath::Abs(param->GetY());
727   if (dist<fgEdgeCutY) isReject=kTRUE;
728   if (cl->GetType()<0) isReject=kTRUE;
729   return isReject;
730 }