]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDtrackV1.cxx
overwrite a bug fix by mistake. Sorry. Thank you Markus
[u/mrichter/AliRoot.git] / TRD / AliTRDtrackV1.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 #include "AliESDtrack.h"
19 #include "AliTracker.h"
20
21 #include "AliTRDtrackV1.h"
22 #include "AliTRDcluster.h"
23 #include "AliTRDcalibDB.h"
24 #include "AliTRDReconstructor.h"
25 #include "AliTRDrecoParam.h"
26
27 ClassImp(AliTRDtrackV1)
28
29 ///////////////////////////////////////////////////////////////////////////////
30 //                                                                           //
31 //  Represents a reconstructed TRD track                                     //
32 //  Local TRD Kalman track                                                   //
33 //                                                                           //
34 //  Authors:                                                                 //
35 //    Alex Bercuci <A.Bercuci@gsi.de>                                        //
36 //    Markus Fasel <M.Fasel@gsi.de>                                          //
37 //                                                                           //
38 ///////////////////////////////////////////////////////////////////////////////
39
40 //_______________________________________________________________
41 AliTRDtrackV1::AliTRDtrackV1() : AliKalmanTrack()
42   ,fPIDquality(0)
43   ,fDE(0.)
44   ,fBackupTrack(0x0)
45 {
46   //
47   // Default constructor
48   //
49   //printf("AliTRDtrackV1::AliTRDtrackV1()\n");
50
51   for(int i =0; i<3; i++) fBudget[i] = 0.;
52   
53   Float_t pid = 1./AliPID::kSPECIES;
54   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
55
56   for(int ip=0; ip<kNplane; ip++){
57     fTrackletIndex[ip] = 0xffff;
58     fTracklet[ip]      = 0x0;
59   }
60 }
61
62 //_______________________________________________________________
63 AliTRDtrackV1::AliTRDtrackV1(const AliTRDtrackV1 &ref) : AliKalmanTrack(ref)
64   ,fPIDquality(ref.fPIDquality)
65   ,fDE(ref.fDE)
66   ,fBackupTrack(0x0)
67 {
68   //
69   // Copy constructor
70   //
71
72   //printf("AliTRDtrackV1::AliTRDtrackV1(const AliTRDtrackV1 &)\n");
73   SetBit(kOwner, kFALSE);
74   for(int ip=0; ip<kNplane; ip++){ 
75     fTrackletIndex[ip] = ref.fTrackletIndex[ip];
76     fTracklet[ip]      = ref.fTracklet[ip];
77   }
78
79   for (Int_t i = 0; i < 3;i++) fBudget[i]      = ref.fBudget[i];
80
81         for(Int_t is = 0; is<AliPID::kSPECIES; is++) fPID[is] = ref.fPID[is];
82   
83   AliKalmanTrack::SetNumberOfClusters(ref.GetNumberOfClusters());
84 }
85
86 //_______________________________________________________________
87 AliTRDtrackV1::AliTRDtrackV1(const AliESDtrack &t) : AliKalmanTrack()
88   ,fPIDquality(0)
89   ,fDE(0.)
90   ,fBackupTrack(0x0)
91 {
92   //
93   // Constructor from AliESDtrack
94   //
95
96   SetLabel(t.GetLabel());
97   SetChi2(0.0);
98   SetMass(t.GetMass());
99   AliKalmanTrack::SetNumberOfClusters(t.GetTRDncls()); 
100   Int_t ti[kNplane]; t.GetTRDtracklets(&ti[0]);
101   for(int ip=0; ip<kNplane; ip++){ 
102     fTrackletIndex[ip] = ti[ip] < 0 ? 0xffff : ti[ip];
103     fTracklet[ip]      = 0x0;
104   }
105   for(int i =0; i<3; i++) fBudget[i] = 0.;
106   
107   Float_t pid = 1./AliPID::kSPECIES;
108   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
109
110   const AliExternalTrackParam *par = &t;
111   if (t.GetStatus() & AliESDtrack::kTRDbackup) { 
112     par = t.GetOuterParam();
113     if (!par) {
114       AliError("No backup info!"); 
115       par = &t;
116     }
117   }
118   Set(par->GetX() 
119      ,par->GetAlpha()
120      ,par->GetParameter()
121      ,par->GetCovariance());
122
123   if(t.GetStatus() & AliESDtrack::kTIME) {
124     StartTimeIntegral();
125     Double_t times[10]; 
126     t.GetIntegratedTimes(times); 
127     SetIntegratedTimes(times);
128     SetIntegratedLength(t.GetIntegratedLength());
129   }
130 }
131
132 //_______________________________________________________________
133 AliTRDtrackV1::AliTRDtrackV1(AliTRDseedV1 *trklts, const Double_t p[5], const Double_t cov[15]
134              , Double_t x, Double_t alpha) : AliKalmanTrack()
135   ,fPIDquality(0)
136   ,fDE(0.)
137   ,fBackupTrack(0x0)
138 {
139   //
140   // The stand alone tracking constructor
141   // TEMPORARY !!!!!!!!!!!
142   // to check :
143   // 1. covariance matrix
144   // 2. dQdl calculation
145   //
146
147   Double_t cnv = GetBz() < 1.e-5 ? 1.e5 : 1.0 / (GetBz() * kB2C);
148   //  Double_t cnv   = 1.0 / (GetBz() * kB2C);
149
150   Double_t pp[5] = { p[0]    
151                    , p[1]
152                    , p[2]
153                    , p[3]
154                    , p[4]*cnv      };
155
156   Double_t c22 = x*x*cov[14] - 2*x*cov[12] + cov[ 5];
157   Double_t c32 =   x*cov[13] -     cov[ 8];
158   Double_t c20 =   x*cov[10] -     cov[ 3];
159   Double_t c21 =   x*cov[11] -     cov[ 4];
160   Double_t c42 =   x*cov[14] -     cov[12];
161
162   Double_t cc[15] = { cov[ 0]
163                     , cov[ 1],     cov[ 2]
164                     , c20,         c21,         c22
165                     , cov[ 6],     cov[ 7],     c32,     cov[ 9]
166                     , cov[10]*cnv, cov[11]*cnv, c42*cnv, cov[13]*cnv, cov[14]*cnv*cnv };
167
168   Double_t mostProbablePt=AliExternalTrackParam::GetMostProbablePt();
169   Double_t p0=TMath::Sign(1/mostProbablePt,pp[4]);
170   Double_t w0=cc[14]/(cc[14] + p0*p0), w1=p0*p0/(cc[14] + p0*p0);
171   pp[4] = w0*p0 + w1*pp[4];
172   cc[10]*=w1; cc[11]*=w1; cc[12]*=w1; cc[13]*=w1; cc[14]*=w1;
173
174         Set(x,alpha,pp,cc);
175   Int_t ncls = 0;
176         for(int iplane=0; iplane<kNplane; iplane++){
177     fTrackletIndex[iplane] = 0xffff;
178                 if(!trklts[iplane].IsOK()) fTracklet[iplane] = 0x0;
179     else{ 
180       fTracklet[iplane] = &trklts[iplane];
181       ncls += fTracklet[iplane]->GetN();
182     }
183         }
184   AliKalmanTrack::SetNumberOfClusters(ncls);            
185   for(int i =0; i<3; i++) fBudget[i] = 0.;
186   
187   Float_t pid = 1./AliPID::kSPECIES;
188   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
189
190 }
191
192 //_______________________________________________________________
193 AliTRDtrackV1::~AliTRDtrackV1()
194 {
195   //AliInfo("");
196   //printf("I-AliTRDtrackV1::~AliTRDtrackV1() : Owner[%s]\n", TestBit(kOwner)?"YES":"NO");
197
198   if(fBackupTrack) delete fBackupTrack;
199   fBackupTrack = 0x0;
200
201   for(Int_t ip=0; ip<kNplane; ip++){
202     if(TestBit(kOwner) && fTracklet[ip]) delete fTracklet[ip];
203     fTracklet[ip] = 0x0;
204     fTrackletIndex[ip] = 0xffff;
205   }
206 }
207         
208 //_______________________________________________________________
209 Bool_t AliTRDtrackV1::CookLabel(Float_t wrong)
210 {
211   // set MC label for this track
212   
213   Int_t s[kMAXCLUSTERSPERTRACK][2];
214   for (Int_t i = 0; i < kMAXCLUSTERSPERTRACK; i++) {
215     s[i][0] = -1;
216     s[i][1] =  0;
217   }
218   
219   Bool_t labelAdded;
220   Int_t label;
221   AliTRDcluster *c    = 0x0;
222   for (Int_t ip = 0; ip < kNplane; ip++) {
223     if(fTrackletIndex[ip] == 0xffff) continue;
224     for (Int_t ic = 0; ic < AliTRDseed::knTimebins; ic++) {
225       if(!(c = fTracklet[ip]->GetClusters(ic))) continue;
226       for (Int_t k = 0; k < 3; k++) { 
227         label      = c->GetLabel(k);
228         labelAdded = kFALSE; 
229         Int_t j = 0;
230         if (label >= 0) {
231           while ((!labelAdded) && (j < kMAXCLUSTERSPERTRACK)) {
232             if ((s[j][0] == label) || 
233                 (s[j][1] ==     0)) {
234               s[j][0] = label; 
235               s[j][1]++; 
236               labelAdded = kTRUE;
237             }
238             j++;
239           }
240         }
241       }
242     }
243   }
244   
245   Int_t max = 0;
246   label = -123456789;
247   for (Int_t i = 0; i < kMAXCLUSTERSPERTRACK; i++) {
248     if (s[i][1] <= max) continue;
249     max   = s[i][1]; 
250     label = s[i][0];
251   }
252   
253   if ((1. - Float_t(max)/GetNumberOfClusters()) > wrong) label = -label;
254   
255   SetLabel(label); 
256   
257   return kTRUE;
258 }
259
260 //_______________________________________________________________
261 Bool_t AliTRDtrackV1::CookPID()
262 {
263   //
264   // Cook the PID information
265   //
266   
267   // Reset the a priori probabilities
268   Double_t pid = 1. / AliPID::kSPECIES;
269   for(int ispec=0; ispec<AliPID::kSPECIES; ispec++) {
270     fPID[ispec] = pid;  
271   }
272   fPIDquality = 0;
273   
274   // steer PID calculation @ tracklet level
275   Double_t *prob = 0x0;
276   for(int ip=0; ip<kNplane; ip++){
277     if(fTrackletIndex[ip] == 0xffff) continue;
278     if(!fTracklet[ip]->IsOK()) continue;
279     if(!(prob = fTracklet[ip]->GetProbability())) return kFALSE;
280     
281     Int_t nspec = 0; // quality check of tracklet dEdx
282     for(int ispec=0; ispec<AliPID::kSPECIES; ispec++){
283       if(prob[ispec] < 0.) continue;
284       fPID[ispec] *= prob[ispec];
285       nspec++;
286     }
287     if(!nspec) continue;
288     
289     fPIDquality++;
290   }
291   
292   // no tracklet found for PID calculations
293   if(!fPIDquality) return kTRUE;
294   
295   // slot for PID calculation @ track level
296   
297   
298   // normalize probabilities
299   Double_t probTotal = 0.0;
300   for (Int_t is = 0; is < AliPID::kSPECIES; is++) probTotal += fPID[is];
301   
302   
303   if (probTotal <= 0.0) {
304     AliWarning("The total probability over all species <= 0. This may be caused by some error in the reference data.");
305     return kFALSE;
306   }
307   
308   for (Int_t iSpecies = 0; iSpecies < AliPID::kSPECIES; iSpecies++) fPID[iSpecies] /= probTotal;
309   
310   return kTRUE;
311 }
312
313 //_____________________________________________________________________________
314 Double_t AliTRDtrackV1::GetBz() const 
315 {
316   //
317   // Returns Bz component of the magnetic field (kG)
318   //
319
320   if (AliTracker::UniformField()) return AliTracker::GetBz();
321
322   Double_t r[3]; 
323   GetXYZ(r);
324   return AliTracker::GetBz(r);
325 }
326
327 //_______________________________________________________________
328 Int_t  AliTRDtrackV1::GetClusterIndex(Int_t id) const
329 {
330   Int_t n = 0;
331   for(Int_t ip=0; ip<kNplane; ip++){
332     if(!fTracklet[ip]) continue;
333     if(n+fTracklet[ip]->GetN() <= id){ 
334       n+=fTracklet[ip]->GetN();
335       continue;
336     }
337     AliTRDcluster *c = 0x0;
338     for(Int_t ic=AliTRDseed::knTimebins-1; ic>=0; ic--){
339       if(!(c = fTracklet[ip]->GetClusters(ic))) continue;
340
341       if(n<id){n++; continue;}
342       return fTracklet[ip]->GetIndexes(ic);
343     }
344   }
345   return -1;
346 }
347
348 //_______________________________________________________________
349 Double_t AliTRDtrackV1::GetPredictedChi2(const AliTRDseedV1 *trklt) const
350 {
351   //
352   // Get the predicted chi2
353   //
354   
355   Double_t x      = trklt->GetX0();
356   Double_t p[2]   = { trklt->GetYat(x)
357                     , trklt->GetZat(x) };
358   Double_t cov[3];
359   trklt->GetCovAt(x, cov);
360   
361   return AliExternalTrackParam::GetPredictedChi2(p, cov);
362 }
363
364         
365 //_____________________________________________________________________________
366 void AliTRDtrackV1::MakeBackupTrack()
367 {
368   //
369   // Creates a backup track
370   //
371
372   if(fBackupTrack) {
373     fBackupTrack->~AliTRDtrackV1();
374     new(fBackupTrack) AliTRDtrackV1((AliTRDtrackV1&)(*this));
375     return;
376   }
377   fBackupTrack = new AliTRDtrackV1((AliTRDtrackV1&)(*this));
378 }
379
380 //_____________________________________________________________________________
381 Int_t AliTRDtrackV1::GetProlongation(Double_t xk, Double_t &y, Double_t &z)
382 {
383   //
384   // Find a prolongation at given x
385   // Return 0 if it does not exist
386   //  
387
388   Double_t bz = GetBz();
389   if (!AliExternalTrackParam::GetYAt(xk,bz,y)) return 0;
390   if (!AliExternalTrackParam::GetZAt(xk,bz,z)) return 0;
391
392   return 1;  
393
394 }
395
396 //_____________________________________________________________________________
397 Bool_t AliTRDtrackV1::PropagateTo(Double_t xk, Double_t xx0, Double_t xrho)
398 {
399   //
400   // Propagates this track to a reference plane defined by "xk" [cm] 
401   // correcting for the mean crossed material.
402   //
403   // "xx0"  - thickness/rad.length [units of the radiation length] 
404   // "xrho" - thickness*density    [g/cm^2] 
405   // 
406
407   if (xk == GetX()) {
408     return kTRUE;
409   }
410
411   Double_t oldX = GetX();
412   Double_t oldY = GetY();
413   Double_t oldZ = GetZ();
414
415   Double_t bz   = GetBz();
416
417   if (!AliExternalTrackParam::PropagateTo(xk,bz)) {
418     return kFALSE;
419   }
420
421   Double_t x = GetX();
422   Double_t y = GetY();
423   Double_t z = GetZ();
424
425   if (oldX < xk) {
426     xrho = -xrho;
427     if (IsStartedTimeIntegral()) {
428       Double_t l2  = TMath::Sqrt((x-oldX)*(x-oldX) 
429                                + (y-oldY)*(y-oldY) 
430                                + (z-oldZ)*(z-oldZ));
431       Double_t crv = AliExternalTrackParam::GetC(bz);
432       if (TMath::Abs(l2*crv) > 0.0001) {
433         // Make correction for curvature if neccesary
434         l2 = 0.5 * TMath::Sqrt((x-oldX)*(x-oldX) 
435                              + (y-oldY)*(y-oldY));
436         l2 = 2.0 * TMath::ASin(l2 * crv) / crv;
437         l2 = TMath::Sqrt(l2*l2 + (z-oldZ)*(z-oldZ));
438       }
439       AddTimeStep(l2);
440     }
441   }
442
443   if (!AliExternalTrackParam::CorrectForMeanMaterial(xx0,xrho,GetMass())) { 
444     return kFALSE;
445   }
446
447   {
448
449     // Energy losses
450     Double_t p2    = (1.0 + GetTgl()*GetTgl()) / (GetSigned1Pt()*GetSigned1Pt());
451     Double_t beta2 = p2 / (p2 + GetMass()*GetMass());
452     if ((beta2 < 1.0e-10) || 
453         ((5940.0 * beta2/(1.0 - beta2 + 1.0e-10) - beta2) < 0.0)) {
454       return kFALSE;
455     }
456
457     Double_t dE    = 0.153e-3 / beta2 
458                    * (TMath::Log(5940.0 * beta2/(1.0 - beta2 + 1.0e-10)) - beta2)
459                    * xrho;
460     fBudget[0] += xrho;
461
462     /*
463     // Suspicious part - think about it ?
464     Double_t kinE =  TMath::Sqrt(p2);
465     if (dE > 0.8*kinE) dE = 0.8 * kinE;  //      
466     if (dE < 0)        dE = 0.0;         // Not valid region for Bethe bloch 
467     */
468  
469     fDE += dE;
470
471     /*
472     // Suspicious ! I.B.
473     Double_t sigmade = 0.07 * TMath::Sqrt(TMath::Abs(dE));   // Energy loss fluctuation 
474     Double_t sigmac2 = sigmade*sigmade*fC*fC*(p2+GetMass()*GetMass())/(p2*p2);
475     fCcc += sigmac2;
476     fCee += fX*fX * sigmac2;  
477     */
478
479   }
480
481   return kTRUE;
482 }
483
484 //_____________________________________________________________________________
485 Int_t   AliTRDtrackV1::PropagateToR(Double_t r,Double_t step)
486 {
487   //
488   // Propagate track to the radial position
489   // Rotation always connected to the last track position
490   //
491
492   Double_t xyz0[3];
493   Double_t xyz1[3];
494   Double_t y;
495   Double_t z; 
496
497   Double_t radius = TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
498   // Direction +-
499   Double_t dir    = (radius > r) ? -1.0 : 1.0;   
500
501   for (Double_t x = radius+dir*step; dir*x < dir*r; x += dir*step) {
502
503     GetXYZ(xyz0);       
504     Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
505     Rotate(alpha,kTRUE);
506     GetXYZ(xyz0);       
507     GetProlongation(x,y,z);
508     xyz1[0] = x * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
509     xyz1[1] = x * TMath::Sin(alpha) - y * TMath::Cos(alpha);
510     xyz1[2] = z;
511     Double_t param[7];
512     AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
513     if (param[1] <= 0) {
514       param[1] = 100000000;
515     }
516     PropagateTo(x,param[1],param[0]*param[4]);
517
518   } 
519
520   GetXYZ(xyz0); 
521   Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
522   Rotate(alpha,kTRUE);
523   GetXYZ(xyz0); 
524   GetProlongation(r,y,z);
525   xyz1[0] = r * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
526   xyz1[1] = r * TMath::Sin(alpha) - y * TMath::Cos(alpha);
527   xyz1[2] = z;
528   Double_t param[7];
529   AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
530
531   if (param[1] <= 0) {
532     param[1] = 100000000;
533   }
534   PropagateTo(r,param[1],param[0]*param[4]);
535
536   return 0;
537
538 }
539
540 //_____________________________________________________________________________
541 Bool_t AliTRDtrackV1::Rotate(Double_t alpha, Bool_t absolute)
542 {
543   //
544   // Rotates track parameters in R*phi plane
545   // if absolute rotation alpha is in global system
546   // otherwise alpha rotation is relative to the current rotation angle
547   //  
548
549   if (absolute) alpha -= GetAlpha();
550   //else fNRotate++;
551
552   return AliExternalTrackParam::Rotate(GetAlpha()+alpha);
553 }
554
555 //___________________________________________________________
556 void AliTRDtrackV1::SetNumberOfClusters() 
557 {
558 // Calculate the number of clusters attached to this track
559         
560   Int_t ncls = 0;
561   for(int ip=0; ip<kNplane; ip++){
562     if(fTracklet[ip] && fTrackletIndex[ip] != 0xffff) ncls += fTracklet[ip]->GetN();
563   }
564   AliKalmanTrack::SetNumberOfClusters(ncls);    
565 }
566
567         
568 //_______________________________________________________________
569 void AliTRDtrackV1::SetOwner()
570 {
571   //
572   // Toggle ownership of tracklets
573   //
574
575   if(TestBit(kOwner)) return;
576   for (Int_t ip = 0; ip < kNplane; ip++) {
577     if(fTrackletIndex[ip] == 0xffff) continue;
578     fTracklet[ip] = new AliTRDseedV1(*fTracklet[ip]);
579     fTracklet[ip]->SetOwner();
580   }
581   SetBit(kOwner);
582 }
583
584 //_______________________________________________________________
585 void AliTRDtrackV1::SetTracklet(AliTRDseedV1 *trklt, Int_t index)
586 {
587   //
588   // Set the tracklets
589   //
590   Int_t plane = trklt->GetPlane();
591
592   fTracklet[plane]      = trklt;
593   fTrackletIndex[plane] = index;
594 }
595
596 //_______________________________________________________________
597 void AliTRDtrackV1::UnsetTracklet(Int_t plane)
598 {
599   if(plane<0 && plane >= kNplane) return;
600   fTrackletIndex[plane] = 0xffff;
601   fTracklet[plane] = 0x0;
602 }
603
604
605 //_______________________________________________________________
606 Bool_t  AliTRDtrackV1::Update(AliTRDseedV1 *trklt, Double_t chisq)
607 {
608   //
609   // Update track and tracklet parameters 
610   //
611   
612   Double_t x      = trklt->GetX0();
613   Double_t p[2]   = { trklt->GetYat(x)
614                     , trklt->GetZat(x) };
615   Double_t cov[3];
616   trklt->GetCovAt(x, cov);
617   
618   if(!AliExternalTrackParam::Update(p, cov)) return kFALSE;
619   
620   AliTRDcluster *c = 0x0;
621   Int_t ic = 0; while(!(c = trklt->GetClusters(ic))) ic++;
622   AliTracker::FillResiduals(this, p, cov, c->GetVolumeId());
623   
624   
625   // Register info to track
626   SetNumberOfClusters();
627   SetChi2(GetChi2() + chisq);
628   
629   // update tracklet
630   trklt->SetMomentum(GetP());
631   trklt->SetSnp(GetSnp());
632   trklt->SetTgl(GetTgl());
633   return kTRUE;
634 }
635
636 //_______________________________________________________________
637 void AliTRDtrackV1::UpdateESDtrack(AliESDtrack *track)
638 {
639   //
640   // Update the TRD PID information in the ESD track
641   //
642
643   track->SetNumberOfTRDslices(kNslice);
644
645   for (Int_t ip = 0; ip < kNplane; ip++) {
646     if(fTrackletIndex[ip] == 0xffff) continue;
647     fTracklet[ip]->CookdEdx(kNslice);
648     Float_t *dedx = fTracklet[ip]->GetdEdx();
649     for (Int_t js = 0; js < kNslice; js++) track->SetTRDslice(dedx[js], ip, js);
650   }
651
652   // copy PID to ESD
653   if(!fPIDquality) return;
654   track->SetTRDpid(fPID);
655   track->SetTRDpidQuality(fPIDquality);
656 }