]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackV2.cxx
Change in the initialisation of two data members of AliITSTrackleterSPDEff (G. E...
[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 xOverX0,xTimesRho; 
186   xOverX0 = d; xTimesRho = d*x0;
187   if (!CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kFALSE;
188
189   Double_t x=GetX(), y=GetY(), z=GetZ();
190   if (IsStartedTimeIntegral() && x>oldX) {
191     Double_t l2 = (x-oldX)*(x-oldX) + (y-oldY)*(y-oldY) + (z-oldZ)*(z-oldZ);
192     AddTimeStep(TMath::Sqrt(l2));
193   }
194
195   return kTRUE;
196 }
197
198 //____________________________________________________________________________
199 Bool_t AliITStrackV2::PropagateToTGeo(Double_t xToGo, Int_t nstep, Double_t &xOverX0, Double_t &xTimesRho, Bool_t addTime) {
200   //-------------------------------------------------------------------
201   //  Propagates the track to a reference plane x=xToGo in n steps.
202   //  These n steps are only used to take into account the curvature.
203   //  The material is calculated with TGeo. (L.Gaudichet)
204   //-------------------------------------------------------------------
205   
206   Double_t startx = GetX(), starty = GetY(), startz = GetZ();
207   Double_t sign = (startx<xToGo) ? -1.:1.;
208   Double_t step = (xToGo-startx)/TMath::Abs(nstep);
209
210   Double_t start[3], end[3], mparam[7], bz = GetBz();
211   Double_t x = startx;
212   
213   for (Int_t i=0; i<nstep; i++) {
214     
215     GetXYZ(start);   //starting global position
216     x += step;
217     if (!GetXYZAt(x, bz, end)) return kFALSE;
218     if (!AliExternalTrackParam::PropagateTo(x, bz)) return kFALSE;
219     AliTracker::MeanMaterialBudget(start, end, mparam);
220     xTimesRho = sign*mparam[4]*mparam[0];
221     xOverX0   = mparam[1];
222     if (mparam[1]<900000) {
223       if (!AliExternalTrackParam::CorrectForMeanMaterial(xOverX0,
224                            xTimesRho,GetMass())) return kFALSE;
225     } else { // this happens when MeanMaterialBudget cannot cross a boundary
226       return kFALSE;
227     }
228   }
229
230   if (addTime && IsStartedTimeIntegral() && GetX()>startx) {
231     Double_t l2 = ( (GetX()-startx)*(GetX()-startx) +
232                     (GetY()-starty)*(GetY()-starty) +
233                     (GetZ()-startz)*(GetZ()-startz) );
234     AddTimeStep(TMath::Sqrt(l2));
235   }
236
237   return kTRUE;
238 }
239
240 //____________________________________________________________________________
241 Bool_t AliITStrackV2::Update(const AliCluster* c, Double_t chi2, Int_t index) 
242 {
243   //------------------------------------------------------------------
244   //This function updates track parameters
245   //------------------------------------------------------------------
246   Double_t p[2]={c->GetY(), c->GetZ()};
247   Double_t cov[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};
248
249   if (!AliExternalTrackParam::Update(p,cov)) return kFALSE;
250
251   Int_t n=GetNumberOfClusters();
252   if (!Invariant()) {
253      if (n>fgkWARN) AliWarning("Wrong invariant !");
254      return kFALSE;
255   }
256
257   if (chi2<0) return kTRUE;
258
259   // fill residuals for ITS+TPC tracks 
260   if (fESDtrack) {
261     if (fESDtrack->GetStatus()&AliESDtrack::kTPCin) {
262       AliTracker::FillResiduals(this,p,cov,c->GetVolumeId());
263     }
264   }
265
266   fIndex[n]=index;
267   SetNumberOfClusters(n+1);
268   SetChi2(GetChi2()+chi2);
269
270   return kTRUE;
271 }
272
273 Bool_t AliITStrackV2::Invariant() const {
274   //------------------------------------------------------------------
275   // This function is for debugging purpose only
276   //------------------------------------------------------------------
277   Int_t n=GetNumberOfClusters();
278
279   // take into account the misalignment error
280   Float_t maxMisalErrY2=0,maxMisalErrZ2=0;
281   for (Int_t lay=0; lay<AliITSgeomTGeo::kNLayers; lay++) {
282     maxMisalErrY2 = TMath::Max(maxMisalErrY2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorY(lay,GetBz()));
283     maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorZ(lay,GetBz()));
284   }
285   maxMisalErrY2 *= maxMisalErrY2;
286   maxMisalErrZ2 *= maxMisalErrZ2;
287   // this is because when we reset before refitting, we multiply the
288   // matrix by 10
289   maxMisalErrY2 *= 10.; 
290   maxMisalErrZ2 *= 10.;
291
292   Double_t sP2=GetParameter()[2];
293   if (TMath::Abs(sP2) >= kAlmost1){
294      if (n>fgkWARN) Warning("Invariant","fP2=%f\n",sP2);
295      return kFALSE;
296   }
297   Double_t sC00=GetCovariance()[0];
298   if (sC00<=0 || sC00>(9.+maxMisalErrY2)) {
299      if (n>fgkWARN) Warning("Invariant","fC00=%f\n",sC00); 
300      return kFALSE;
301   }
302   Double_t sC11=GetCovariance()[2];
303   if (sC11<=0 || sC11>(9.+maxMisalErrZ2)) {
304      if (n>fgkWARN) Warning("Invariant","fC11=%f\n",sC11); 
305      return kFALSE;
306   }
307   Double_t sC22=GetCovariance()[5];
308   if (sC22<=0 || sC22>1.) {
309      if (n>fgkWARN) Warning("Invariant","fC22=%f\n",sC22); 
310      return kFALSE;
311   }
312   Double_t sC33=GetCovariance()[9];
313   if (sC33<=0 || sC33>1.) {
314      if (n>fgkWARN) Warning("Invariant","fC33=%f\n",sC33); 
315      return kFALSE;
316   }
317   Double_t sC44=GetCovariance()[14];
318   if (sC44<=0 /*|| sC44>6e-5*/) {
319      if (n>fgkWARN) Warning("Invariant","fC44=%f\n",sC44);
320      return kFALSE;
321   }
322
323   return kTRUE;
324 }
325
326 //____________________________________________________________________________
327 Bool_t AliITStrackV2::Propagate(Double_t alp,Double_t xk) {
328   //------------------------------------------------------------------
329   //This function propagates a track
330   //------------------------------------------------------------------
331   Double_t bz=GetBz();
332   if (!AliExternalTrackParam::Propagate(alp,xk,bz)) return kFALSE;
333
334   if (!Invariant()) {
335     Int_t n=GetNumberOfClusters();
336     if (n>fgkWARN) AliWarning("Wrong invariant !");
337     return kFALSE;
338   }
339
340   return kTRUE;
341 }
342
343 Bool_t AliITStrackV2::MeanBudgetToPrimVertex(Double_t xyz[3], Double_t step, Double_t &d) const {
344
345   //-------------------------------------------------------------------
346   //  Get the mean material budget between the actual point and the
347   //  primary vertex. (L.Gaudichet)
348   //-------------------------------------------------------------------
349
350   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
351   Double_t vertexX = xyz[0]*cs + xyz[1]*sn;
352
353   Int_t nstep = Int_t((GetX()-vertexX)/step);
354   if (nstep<1) nstep = 1;
355   step = (GetX()-vertexX)/nstep;
356
357   //  Double_t mparam[7], densMean=0, radLength=0, length=0;
358   Double_t mparam[7];
359   Double_t p1[3], p2[3], x = GetX(), bz = GetBz();
360   GetXYZ(p1);
361
362   d=0.;
363
364   for (Int_t i=0; i<nstep; i++) {
365     x  += step;
366     if (!GetXYZAt(x, bz, p2)) return kFALSE;
367     AliTracker::MeanMaterialBudget(p1, p2, mparam);
368     if (mparam[1]>900000) return kFALSE;
369     d  += mparam[1];
370
371     p1[0] = p2[0];
372     p1[1] = p2[1];
373     p1[2] = p2[2];
374   }
375
376   return kTRUE;
377 }
378
379 Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) {
380   //------------------------------------------------------------------
381   //This function improves angular track parameters
382   //------------------------------------------------------------------
383   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
384 //Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex
385   Double_t yv =-xyz[0]*sn + xyz[1]*cs; // in the
386   Double_t zv = xyz[2];                // local frame
387
388   Double_t dy = Par(0) - yv, dz = Par(1) - zv;
389   Double_t r2=GetX()*GetX() + dy*dy;
390   Double_t p2=(1.+ GetTgl()*GetTgl())/(GetSigned1Pt()*GetSigned1Pt());
391   Double_t beta2=p2/(p2 + GetMass()*GetMass());
392   x0*=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(1.- GetSnp()*GetSnp()));
393   Double_t theta2=14.1*14.1/(beta2*p2*1e6)*x0;
394   //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*x0*9.36*2.33;
395
396   Double_t cnv=GetBz()*kB2C;
397   {
398     Double_t dummy = 4/r2 - GetC()*GetC();
399     if (dummy < 0) return kFALSE;
400     Double_t parp = 0.5*(GetC()*GetX() + dy*TMath::Sqrt(dummy));
401     Double_t sigma2p = theta2*(1.- GetSnp()*GetSnp())*(1. + GetTgl()*GetTgl());
402     sigma2p += Cov(0)/r2*(1.- dy*dy/r2)*(1.- dy*dy/r2);
403     sigma2p += ers[1]*ers[1]/r2;
404     sigma2p += 0.25*Cov(14)*cnv*cnv*GetX()*GetX();
405     Double_t eps2p=sigma2p/(Cov(5) + sigma2p);
406     Par(0) += Cov(3)/(Cov(5) + sigma2p)*(parp - GetSnp());
407     Par(2)  = eps2p*GetSnp() + (1 - eps2p)*parp;
408     Cov(5) *= eps2p;
409     Cov(3) *= eps2p;
410   }
411   {
412     Double_t parl=0.5*GetC()*dz/TMath::ASin(0.5*GetC()*TMath::Sqrt(r2));
413     Double_t sigma2l=theta2;
414     sigma2l += Cov(2)/r2 + Cov(0)*dy*dy*dz*dz/(r2*r2*r2);
415     sigma2l += ers[2]*ers[2]/r2;
416     Double_t eps2l = sigma2l/(Cov(9) + sigma2l);
417     Par(1) += Cov(7 )/(Cov(9) + sigma2l)*(parl - Par(3));
418     Par(4) += Cov(13)/(Cov(9) + sigma2l)*(parl - Par(3));
419     Par(3)  = eps2l*Par(3) + (1-eps2l)*parl;
420     Cov(9) *= eps2l; 
421     Cov(13)*= eps2l; 
422     Cov(7) *= eps2l; 
423   }
424   if (!Invariant()) return kFALSE;
425
426   return kTRUE;
427 }
428
429 void AliITStrackV2::CookdEdx(Double_t low, Double_t up) {
430   //-----------------------------------------------------------------
431   // This function calculates dE/dX within the "low" and "up" cuts.
432   // Origin: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
433   // Updated: F. Prino 8-June-2009
434   //-----------------------------------------------------------------
435   // The cluster order is: SDD-1, SDD-2, SSD-1, SSD-2
436
437   Int_t nc=0;
438   Float_t dedx[4];
439   for (Int_t il=0; il<4; il++) { // count good (>0) dE/dx values
440     if(fdEdxSample[il]>0.){
441       dedx[nc]= fdEdxSample[il];
442       nc++;
443     }
444   }
445   if(nc<1){
446     SetdEdx(0.);
447     return;
448   }
449
450   Int_t swap; // sort in ascending order
451   do {
452     swap=0;
453     for (Int_t i=0; i<nc-1; i++) {
454       if (dedx[i]<=dedx[i+1]) continue;
455       Float_t tmp=dedx[i];
456       dedx[i]=dedx[i+1]; 
457       dedx[i+1]=tmp;
458       swap++;
459     }
460   } while (swap);
461
462
463   Double_t nl=low*nc, nu =up*nc;
464   Float_t sumamp = 0;
465   Float_t sumweight =0;
466   for (Int_t i=0; i<nc; i++) {
467     Float_t weight =1;
468     if (i<nl+0.1)     weight = TMath::Max(1.-(nl-i),0.);
469     if (i>nu-1)     weight = TMath::Max(nu-i,0.);
470     sumamp+= dedx[i]*weight;
471     sumweight+=weight;
472   }
473   SetdEdx(sumamp/sumweight);
474 }
475
476 //____________________________________________________________________________
477 Bool_t AliITStrackV2::
478 GetPhiZat(Double_t r, Double_t &phi, Double_t &z) const {
479   //------------------------------------------------------------------
480   // This function returns the global cylindrical (phi,z) of the track 
481   // position estimated at the radius r. 
482   // The track curvature is neglected.
483   //------------------------------------------------------------------
484   Double_t d=GetD(0.,0.);
485   if (TMath::Abs(d) > r) {
486     if (r>1e-1) return kFALSE;
487     r = TMath::Abs(d);
488   }
489
490   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
491   if (TMath::Abs(d) > rcurr) return kFALSE;
492   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
493   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
494
495   if (GetX()>=0.) {
496     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
497   } else {
498     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
499   }
500
501   // return a phi in [0,2pi[ 
502   if (phi<0.) phi+=2.*TMath::Pi();
503   else if (phi>=2.*TMath::Pi()) phi-=2.*TMath::Pi();
504   z=GetZ()+GetTgl()*(TMath::Sqrt((r-d)*(r+d))-TMath::Sqrt((rcurr-d)*(rcurr+d)));
505   return kTRUE;
506 }
507 //____________________________________________________________________________
508 Bool_t AliITStrackV2::
509 GetLocalXat(Double_t r,Double_t &xloc) const {
510   //------------------------------------------------------------------
511   // This function returns the local x of the track 
512   // position estimated at the radius r. 
513   // The track curvature is neglected.
514   //------------------------------------------------------------------
515   Double_t d=GetD(0.,0.);
516   if (TMath::Abs(d) > r) { 
517     if (r>1e-1) return kFALSE; 
518     r = TMath::Abs(d); 
519   } 
520
521   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
522   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
523   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
524   Double_t phi;
525   if (GetX()>=0.) {
526     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
527   } else {
528     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
529   }
530
531   xloc=r*(TMath::Cos(phi)*TMath::Cos(GetAlpha())
532          +TMath::Sin(phi)*TMath::Sin(GetAlpha())); 
533
534   return kTRUE;
535 }