]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackV2.cxx
impoved num precision
[u/mrichter/AliRoot.git] / ITS / AliITStrackV2.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////
19 //                Implementation of the ITS track class
20 //
21 //          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
22 //     dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
23 ///////////////////////////////////////////////////////////////////////////
24 #include <TMath.h>
25
26 #include "AliCluster.h"
27 #include "AliESDtrack.h"
28 #include "AliESDVertex.h"
29 #include "AliITSReconstructor.h"
30 #include "AliITStrackV2.h"
31 #include "AliTracker.h"
32
33 const Int_t AliITStrackV2::fgkWARN = 5;
34
35 ClassImp(AliITStrackV2)
36
37
38 //____________________________________________________________________________
39 AliITStrackV2::AliITStrackV2() : AliKalmanTrack(),
40   fdEdx(0),
41   fESDtrack(0)
42 {
43   for(Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) {fIndex[i]=-1; fModule[i]=-1;}
44   for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
45 }
46
47
48 //____________________________________________________________________________
49 AliITStrackV2::AliITStrackV2(AliESDtrack& t,Bool_t c) throw (const Char_t *) :
50   AliKalmanTrack(),
51   fdEdx(t.GetITSsignal()),
52   fESDtrack(&t)
53 {
54   //------------------------------------------------------------------
55   // Conversion ESD track -> ITS track.
56   // If c==kTRUE, create the ITS track out of the constrained params.
57   //------------------------------------------------------------------
58   const AliExternalTrackParam *par=&t;
59   if (c) {
60     par=t.GetConstrainedParam();
61     if (!par) throw "AliITStrackV2: conversion failed !\n";
62   }
63   Set(par->GetX(),par->GetAlpha(),par->GetParameter(),par->GetCovariance());
64
65   //if (!Invariant()) throw "AliITStrackV2: conversion failed !\n";
66
67   SetLabel(t.GetLabel());
68   SetMass(t.GetMass());
69   SetNumberOfClusters(t.GetITSclusters(fIndex));
70
71   if (t.GetStatus()&AliESDtrack::kTIME) {
72     StartTimeIntegral();
73     Double_t times[10]; t.GetIntegratedTimes(times); SetIntegratedTimes(times);
74     SetIntegratedLength(t.GetIntegratedLength());
75   }
76
77   for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
78 }
79
80 //____________________________________________________________________________
81 void AliITStrackV2::ResetClusters() {
82   //------------------------------------------------------------------
83   // Reset the array of attached clusters.
84   //------------------------------------------------------------------
85   for (Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) fIndex[i]=-1;
86   SetChi2(0.); 
87   SetNumberOfClusters(0);
88
89
90 //____________________________________________________________________________
91 void AliITStrackV2::UpdateESDtrack(ULong_t flags) const {
92   // Update track params
93   fESDtrack->UpdateTrackParams(this,flags);
94   // copy the module indices
95   for(Int_t i=0;i<12;i++) {
96     //   printf("     %d\n",GetModuleIndex(i));
97     fESDtrack->SetITSModuleIndex(i,GetModuleIndex(i));
98   }
99   // copy the 4 dedx samples
100   Double_t sdedx[4]={0.,0.,0.,0.};
101   for(Int_t i=0; i<4; i++) sdedx[i]=fdEdxSample[i];
102   fESDtrack->SetITSdEdxSamples(sdedx);
103 }
104
105 //____________________________________________________________________________
106 AliITStrackV2::AliITStrackV2(const AliITStrackV2& t) : 
107   AliKalmanTrack(t),
108   fdEdx(t.fdEdx),
109   fESDtrack(t.fESDtrack) 
110 {
111   //------------------------------------------------------------------
112   //Copy constructor
113   //------------------------------------------------------------------
114   Int_t i;
115   for (i=0; i<4; i++) fdEdxSample[i]=t.fdEdxSample[i];
116   for (i=0; i<2*AliITSgeomTGeo::GetNLayers(); i++) {
117     fIndex[i]=t.fIndex[i];
118     fModule[i]=t.fModule[i];
119   }
120 }
121
122 //_____________________________________________________________________________
123 Int_t AliITStrackV2::Compare(const TObject *o) const {
124   //-----------------------------------------------------------------
125   // This function compares tracks according to the their curvature
126   //-----------------------------------------------------------------
127   AliITStrackV2 *t=(AliITStrackV2*)o;
128   //Double_t co=OneOverPt();
129   //Double_t c =OneOverPt();
130   Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
131   Double_t c =GetSigmaY2()*GetSigmaZ2();
132   if (c>co) return 1;
133   else if (c<co) return -1;
134   return 0;
135 }
136
137 //____________________________________________________________________________
138 Bool_t 
139 AliITStrackV2::PropagateToVertex(const AliESDVertex *v,Double_t d,Double_t x0) 
140 {
141   //------------------------------------------------------------------
142   //This function propagates a track to the minimal distance from the origin
143   //------------------------------------------------------------------  
144   Double_t bz=GetBz();
145   if (PropagateToDCA(v,bz,kVeryBig)) {
146     Double_t xOverX0,xTimesRho; 
147     xOverX0 = d; xTimesRho = d*x0;
148     if (CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kTRUE;
149   }
150   return kFALSE;
151 }
152
153 //____________________________________________________________________________
154 Bool_t AliITStrackV2::
155 GetGlobalXYZat(Double_t xloc, Double_t &x, Double_t &y, Double_t &z) const {
156   //------------------------------------------------------------------
157   //This function returns a track position in the global system
158   //------------------------------------------------------------------
159   Double_t r[3];
160   Bool_t rc=GetXYZAt(xloc, GetBz(), r);
161   x=r[0]; y=r[1]; z=r[2]; 
162   return rc;
163 }
164
165 //_____________________________________________________________________________
166 Double_t AliITStrackV2::GetPredictedChi2(const AliCluster *c) const {
167   //-----------------------------------------------------------------
168   // This function calculates a predicted chi2 increment.
169   //-----------------------------------------------------------------
170   Double_t p[2]={c->GetY(), c->GetZ()};
171   Double_t cov[3]={c->GetSigmaY2(), 0., c->GetSigmaZ2()};
172   return AliExternalTrackParam::GetPredictedChi2(p,cov);
173 }
174
175 //____________________________________________________________________________
176 Bool_t AliITStrackV2::PropagateTo(Double_t xk, Double_t d, Double_t x0) {
177   //------------------------------------------------------------------
178   //This function propagates a track
179   //------------------------------------------------------------------
180
181   Double_t oldX=GetX(), oldY=GetY(), oldZ=GetZ();
182   
183   //Double_t bz=GetBz();
184   //if (!AliExternalTrackParam::PropagateTo(xk,bz)) return kFALSE;
185   Double_t b[3]; GetBxByBz(b);
186   if (!AliExternalTrackParam::PropagateToBxByBz(xk,b)) return kFALSE;
187   Double_t xOverX0,xTimesRho; 
188   xOverX0 = d; xTimesRho = d*x0;
189   if (!CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kFALSE;
190
191   Double_t x=GetX(), y=GetY(), z=GetZ();
192   if (IsStartedTimeIntegral() && x>oldX) {
193     Double_t l2 = (x-oldX)*(x-oldX) + (y-oldY)*(y-oldY) + (z-oldZ)*(z-oldZ);
194     AddTimeStep(TMath::Sqrt(l2));
195   }
196
197   return kTRUE;
198 }
199
200 //____________________________________________________________________________
201 Bool_t AliITStrackV2::PropagateToTGeo(Double_t xToGo, Int_t nstep, Double_t &xOverX0, Double_t &xTimesRho, Bool_t addTime) {
202   //-------------------------------------------------------------------
203   //  Propagates the track to a reference plane x=xToGo in n steps.
204   //  These n steps are only used to take into account the curvature.
205   //  The material is calculated with TGeo. (L.Gaudichet)
206   //-------------------------------------------------------------------
207   
208   Double_t startx = GetX(), starty = GetY(), startz = GetZ();
209   Double_t sign = (startx<xToGo) ? -1.:1.;
210   Double_t step = (xToGo-startx)/TMath::Abs(nstep);
211
212   Double_t start[3], end[3], mparam[7];
213   //Double_t bz = GetBz();
214   Double_t b[3]; GetBxByBz(b);
215   Double_t bz = b[2];
216
217   Double_t x = startx;
218   
219   for (Int_t i=0; i<nstep; i++) {
220     
221     GetXYZ(start);   //starting global position
222     x += step;
223     if (!GetXYZAt(x, bz, end)) return kFALSE;
224     //if (!AliExternalTrackParam::PropagateTo(x, bz)) return kFALSE;
225     if (!AliExternalTrackParam::PropagateToBxByBz(x, b)) return kFALSE;
226     AliTracker::MeanMaterialBudget(start, end, mparam);
227     xTimesRho = sign*mparam[4]*mparam[0];
228     xOverX0   = mparam[1];
229     if (mparam[1]<900000) {
230       if (!AliExternalTrackParam::CorrectForMeanMaterial(xOverX0,
231                            xTimesRho,GetMass())) return kFALSE;
232     } else { // this happens when MeanMaterialBudget cannot cross a boundary
233       return kFALSE;
234     }
235   }
236
237   if (addTime && IsStartedTimeIntegral() && GetX()>startx) {
238     Double_t l2 = ( (GetX()-startx)*(GetX()-startx) +
239                     (GetY()-starty)*(GetY()-starty) +
240                     (GetZ()-startz)*(GetZ()-startz) );
241     AddTimeStep(TMath::Sqrt(l2));
242   }
243
244   return kTRUE;
245 }
246
247 //____________________________________________________________________________
248 Bool_t AliITStrackV2::Update(const AliCluster* c, Double_t chi2, Int_t index) 
249 {
250   //------------------------------------------------------------------
251   //This function updates track parameters
252   //------------------------------------------------------------------
253   Double_t p[2]={c->GetY(), c->GetZ()};
254   Double_t cov[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};
255
256   if (!AliExternalTrackParam::Update(p,cov)) return kFALSE;
257
258   Int_t n=GetNumberOfClusters();
259   if (!Invariant()) {
260      if (n>fgkWARN) AliWarning("Wrong invariant !");
261      return kFALSE;
262   }
263
264   if (chi2<0) return kTRUE;
265
266   // fill residuals for ITS+TPC tracks 
267   if (fESDtrack) {
268     if (fESDtrack->GetStatus()&AliESDtrack::kTPCin) {
269       AliTracker::FillResiduals(this,p,cov,c->GetVolumeId());
270     }
271   }
272
273   fIndex[n]=index;
274   SetNumberOfClusters(n+1);
275   SetChi2(GetChi2()+chi2);
276
277   return kTRUE;
278 }
279
280 Bool_t AliITStrackV2::Invariant() const {
281   //------------------------------------------------------------------
282   // This function is for debugging purpose only
283   //------------------------------------------------------------------
284   Int_t n=GetNumberOfClusters();
285
286   // take into account the misalignment error
287   Float_t maxMisalErrY2=0,maxMisalErrZ2=0;
288   for (Int_t lay=0; lay<AliITSgeomTGeo::kNLayers; lay++) {
289     maxMisalErrY2 = TMath::Max(maxMisalErrY2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorY(lay,GetBz()));
290     maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorZ(lay,GetBz()));
291   }
292   maxMisalErrY2 *= maxMisalErrY2;
293   maxMisalErrZ2 *= maxMisalErrZ2;
294   // this is because when we reset before refitting, we multiply the
295   // matrix by 10
296   maxMisalErrY2 *= 10.; 
297   maxMisalErrZ2 *= 10.;
298
299   Double_t sP2=GetParameter()[2];
300   if (TMath::Abs(sP2) >= kAlmost1){
301      if (n>fgkWARN) Warning("Invariant","fP2=%f\n",sP2);
302      return kFALSE;
303   }
304   Double_t sC00=GetCovariance()[0];
305   if (sC00<=0 || sC00>(9.+maxMisalErrY2)) {
306      if (n>fgkWARN) Warning("Invariant","fC00=%f\n",sC00); 
307      return kFALSE;
308   }
309   Double_t sC11=GetCovariance()[2];
310   if (sC11<=0 || sC11>(9.+maxMisalErrZ2)) {
311      if (n>fgkWARN) Warning("Invariant","fC11=%f\n",sC11); 
312      return kFALSE;
313   }
314   Double_t sC22=GetCovariance()[5];
315   if (sC22<=0 || sC22>1.) {
316      if (n>fgkWARN) Warning("Invariant","fC22=%f\n",sC22); 
317      return kFALSE;
318   }
319   Double_t sC33=GetCovariance()[9];
320   if (sC33<=0 || sC33>1.) {
321      if (n>fgkWARN) Warning("Invariant","fC33=%f\n",sC33); 
322      return kFALSE;
323   }
324   Double_t sC44=GetCovariance()[14];
325   if (sC44<=0 /*|| sC44>6e-5*/) {
326      if (n>fgkWARN) Warning("Invariant","fC44=%f\n",sC44);
327      return kFALSE;
328   }
329
330   return kTRUE;
331 }
332
333 //____________________________________________________________________________
334 Bool_t AliITStrackV2::Propagate(Double_t alp,Double_t xk) {
335   //------------------------------------------------------------------
336   //This function propagates a track
337   //------------------------------------------------------------------
338   //Double_t bz=GetBz();
339   //if (!AliExternalTrackParam::Propagate(alp,xk,bz)) return kFALSE;
340   Double_t b[3]; GetBxByBz(b);
341   if (!AliExternalTrackParam::PropagateBxByBz(alp,xk,b)) return kFALSE;
342
343   if (!Invariant()) {
344     Int_t n=GetNumberOfClusters();
345     if (n>fgkWARN) AliWarning("Wrong invariant !");
346     return kFALSE;
347   }
348
349   return kTRUE;
350 }
351
352 Bool_t AliITStrackV2::MeanBudgetToPrimVertex(Double_t xyz[3], Double_t step, Double_t &d) const {
353
354   //-------------------------------------------------------------------
355   //  Get the mean material budget between the actual point and the
356   //  primary vertex. (L.Gaudichet)
357   //-------------------------------------------------------------------
358
359   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
360   Double_t vertexX = xyz[0]*cs + xyz[1]*sn;
361
362   Int_t nstep = Int_t((GetX()-vertexX)/step);
363   if (nstep<1) nstep = 1;
364   step = (GetX()-vertexX)/nstep;
365
366   //  Double_t mparam[7], densMean=0, radLength=0, length=0;
367   Double_t mparam[7];
368   Double_t p1[3], p2[3], x = GetX(), bz = GetBz();
369   GetXYZ(p1);
370
371   d=0.;
372
373   for (Int_t i=0; i<nstep; i++) {
374     x  += step;
375     if (!GetXYZAt(x, bz, p2)) return kFALSE;
376     AliTracker::MeanMaterialBudget(p1, p2, mparam);
377     if (mparam[1]>900000) return kFALSE;
378     d  += mparam[1];
379
380     p1[0] = p2[0];
381     p1[1] = p2[1];
382     p1[2] = p2[2];
383   }
384
385   return kTRUE;
386 }
387
388 Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) {
389   //------------------------------------------------------------------
390   //This function improves angular track parameters
391   //------------------------------------------------------------------
392   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
393 //Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex
394   Double_t yv =-xyz[0]*sn + xyz[1]*cs; // in the
395   Double_t zv = xyz[2];                // local frame
396
397   Double_t dy = Par(0) - yv, dz = Par(1) - zv;
398   Double_t r2=GetX()*GetX() + dy*dy;
399   Double_t p2=(1.+ GetTgl()*GetTgl())/(GetSigned1Pt()*GetSigned1Pt());
400   Double_t beta2=p2/(p2 + GetMass()*GetMass());
401   x0*=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(1.- GetSnp()*GetSnp()));
402   Double_t theta2=14.1*14.1/(beta2*p2*1e6)*x0;
403   //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*x0*9.36*2.33;
404
405   Double_t cnv=GetBz()*kB2C;
406   {
407     Double_t dummy = 4/r2 - GetC()*GetC();
408     if (dummy < 0) return kFALSE;
409     Double_t parp = 0.5*(GetC()*GetX() + dy*TMath::Sqrt(dummy));
410     Double_t sigma2p = theta2*(1.-GetSnp())*(1.+GetSnp())*(1. + GetTgl()*GetTgl());
411     Double_t ovSqr2 = 1./TMath::Sqrt(r2);
412     Double_t tfact = ovSqr2*(1.-dy*ovSqr2)*(1.+dy*ovSqr2);
413     sigma2p += Cov(0)*tfact*tfact;
414     sigma2p += ers[1]*ers[1]/r2;
415     sigma2p += 0.25*Cov(14)*cnv*cnv*GetX()*GetX();
416     Double_t eps2p=sigma2p/(Cov(5) + sigma2p);
417     Par(0) += Cov(3)/(Cov(5) + sigma2p)*(parp - GetSnp());
418     Par(2)  = eps2p*GetSnp() + (1 - eps2p)*parp;
419     Cov(5) *= eps2p;
420     Cov(3) *= eps2p;
421   }
422   {
423     Double_t parl=0.5*GetC()*dz/TMath::ASin(0.5*GetC()*TMath::Sqrt(r2));
424     Double_t sigma2l=theta2;
425     sigma2l += Cov(2)/r2 + Cov(0)*dy*dy*dz*dz/(r2*r2*r2);
426     sigma2l += ers[2]*ers[2]/r2;
427     Double_t eps2l = sigma2l/(Cov(9) + sigma2l);
428     Par(1) += Cov(7 )/(Cov(9) + sigma2l)*(parl - Par(3));
429     Par(4) += Cov(13)/(Cov(9) + sigma2l)*(parl - Par(3));
430     Par(3)  = eps2l*Par(3) + (1-eps2l)*parl;
431     Cov(9) *= eps2l; 
432     Cov(13)*= eps2l; 
433     Cov(7) *= eps2l; 
434   }
435   if (!Invariant()) return kFALSE;
436
437   return kTRUE;
438 }
439
440 void AliITStrackV2::CookdEdx(Double_t low, Double_t up) {
441   //-----------------------------------------------------------------
442   // This function calculates dE/dX within the "low" and "up" cuts.
443   // Origin: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
444   // Updated: F. Prino 8-June-2009
445   //-----------------------------------------------------------------
446   // The cluster order is: SDD-1, SDD-2, SSD-1, SSD-2
447
448   Int_t nc=0;
449   Float_t dedx[4];
450   for (Int_t il=0; il<4; il++) { // count good (>0) dE/dx values
451     if(fdEdxSample[il]>0.){
452       dedx[nc]= fdEdxSample[il];
453       nc++;
454     }
455   }
456   if(nc<1){
457     SetdEdx(0.);
458     return;
459   }
460
461   Int_t swap; // sort in ascending order
462   do {
463     swap=0;
464     for (Int_t i=0; i<nc-1; i++) {
465       if (dedx[i]<=dedx[i+1]) continue;
466       Float_t tmp=dedx[i];
467       dedx[i]=dedx[i+1]; 
468       dedx[i+1]=tmp;
469       swap++;
470     }
471   } while (swap);
472
473
474   Double_t nl=low*nc, nu =up*nc;
475   Float_t sumamp = 0;
476   Float_t sumweight =0;
477   for (Int_t i=0; i<nc; i++) {
478     Float_t weight =1;
479     if (i<nl+0.1)     weight = TMath::Max(1.-(nl-i),0.);
480     if (i>nu-1)     weight = TMath::Max(nu-i,0.);
481     sumamp+= dedx[i]*weight;
482     sumweight+=weight;
483   }
484   SetdEdx(sumamp/sumweight);
485 }
486
487 //____________________________________________________________________________
488 Bool_t AliITStrackV2::
489 GetPhiZat(Double_t r, Double_t &phi, Double_t &z) const {
490   //------------------------------------------------------------------
491   // This function returns the global cylindrical (phi,z) of the track 
492   // position estimated at the radius r. 
493   // The track curvature is neglected.
494   //------------------------------------------------------------------
495   Double_t d=GetD(0.,0.);
496   if (TMath::Abs(d) > r) {
497     if (r>1e-1) return kFALSE;
498     r = TMath::Abs(d);
499   }
500
501   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
502   if (TMath::Abs(d) > rcurr) return kFALSE;
503   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
504   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
505
506   if (GetX()>=0.) {
507     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
508   } else {
509     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
510   }
511
512   // return a phi in [0,2pi[ 
513   if (phi<0.) phi+=2.*TMath::Pi();
514   else if (phi>=2.*TMath::Pi()) phi-=2.*TMath::Pi();
515   z=GetZ()+GetTgl()*(TMath::Sqrt((r-d)*(r+d))-TMath::Sqrt((rcurr-d)*(rcurr+d)));
516   return kTRUE;
517 }
518 //____________________________________________________________________________
519 Bool_t AliITStrackV2::
520 GetLocalXat(Double_t r,Double_t &xloc) const {
521   //------------------------------------------------------------------
522   // This function returns the local x of the track 
523   // position estimated at the radius r. 
524   // The track curvature is neglected.
525   //------------------------------------------------------------------
526   Double_t d=GetD(0.,0.);
527   if (TMath::Abs(d) > r) { 
528     if (r>1e-1) return kFALSE; 
529     r = TMath::Abs(d); 
530   } 
531
532   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
533   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
534   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
535   Double_t phi;
536   if (GetX()>=0.) {
537     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
538   } else {
539     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
540   }
541
542   xloc=r*(TMath::Cos(phi)*TMath::Cos(GetAlpha())
543          +TMath::Sin(phi)*TMath::Sin(GetAlpha())); 
544
545   return kTRUE;
546 }