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