]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDtrackV1.cxx
protect against overstepping ADC array boundary
[u/mrichter/AliRoot.git] / TRD / AliTRDtrackV1.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19 #include "TVectorT.h"
20 #include "AliLog.h"
21 #include "AliESDtrack.h"
22 #include "AliTracker.h"
23
24 #include "AliTRDtrackV1.h"
25 #include "AliTRDcluster.h"
26 #include "AliTRDcalibDB.h"
27 #include "AliTRDReconstructor.h"
28 #include "AliTRDPIDResponse.h"
29 #include "AliTRDrecoParam.h"
30 #include "AliTRDdEdxUtils.h"
31
32 ClassImp(AliTRDtrackV1)
33
34 ///////////////////////////////////////////////////////////////////////////////
35 //                                                                           //
36 //  Represents a reconstructed TRD track                                     //
37 //  Local TRD Kalman track                                                   //
38 //                                                                           //
39 //  Authors:                                                                 //
40 //    Alex Bercuci <A.Bercuci@gsi.de>                                        //
41 //    Markus Fasel <M.Fasel@gsi.de>                                          //
42 //                                                                           //
43 ///////////////////////////////////////////////////////////////////////////////
44
45 //_______________________________________________________________
46 AliTRDtrackV1::AliTRDtrackV1() : AliKalmanTrack()
47   ,fStatus(0)
48   ,fESDid(0)
49   ,fDE(0.)
50   ,fTruncatedMean(0)
51   ,fkReconstructor(NULL)
52   ,fBackupTrack(NULL)
53   ,fTrackLow(NULL)
54   ,fTrackHigh(NULL)
55 {
56   //
57   // Default constructor
58   //
59   //printf("AliTRDtrackV1::AliTRDtrackV1()\n");
60
61   for(int i =0; i<3; i++) fBudget[i] = 0.;
62   
63   Float_t pid = 1./AliPID::kSPECIES;
64   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
65
66   for(int ip=0; ip<kNplane; ip++){
67     fTrackletIndex[ip] = -1;
68     fTracklet[ip]      = NULL;
69   }
70 }
71
72 //_______________________________________________________________
73 AliTRDtrackV1::AliTRDtrackV1(const AliTRDtrackV1 &ref) : AliKalmanTrack(ref)
74   ,fStatus(ref.fStatus)
75   ,fESDid(ref.fESDid)
76   ,fDE(ref.fDE)
77   ,fTruncatedMean(ref.fTruncatedMean)
78   ,fkReconstructor(ref.fkReconstructor)
79   ,fBackupTrack(NULL)
80   ,fTrackLow(NULL)
81   ,fTrackHigh(NULL)
82 {
83   //
84   // Copy constructor
85   //
86
87   //printf("AliTRDtrackV1::AliTRDtrackV1(const AliTRDtrackV1 &)\n");
88   SetBit(kOwner, kFALSE);
89   for(int ip=0; ip<kNplane; ip++){ 
90     fTrackletIndex[ip] = ref.fTrackletIndex[ip];
91     fTracklet[ip]      = ref.fTracklet[ip];
92   }
93   if(ref.fTrackLow) fTrackLow = new AliExternalTrackParam(*ref.fTrackLow);
94   if(ref.fTrackHigh) fTrackHigh = new AliExternalTrackParam(*ref.fTrackHigh);
95  
96   for (Int_t i = 0; i < 3;i++) fBudget[i]      = ref.fBudget[i];
97
98         for(Int_t is = 0; is<AliPID::kSPECIES; is++) fPID[is] = ref.fPID[is];
99   
100   AliKalmanTrack::SetNumberOfClusters(ref.GetNumberOfClusters());
101 }
102
103 //_______________________________________________________________
104 AliTRDtrackV1::AliTRDtrackV1(const AliESDtrack &t) : AliKalmanTrack()
105   ,fStatus(0)
106   ,fESDid(0)
107   ,fDE(0.)
108   ,fTruncatedMean(0)
109   ,fkReconstructor(NULL)
110   ,fBackupTrack(NULL)
111   ,fTrackLow(NULL)
112   ,fTrackHigh(NULL)
113 {
114   //
115   // Constructor from AliESDtrack
116   //
117
118   SetESDid(t.GetID());
119   SetLabel(t.GetLabel());
120   SetChi2(0.0);
121
122   SetMass(t.GetMass(kTRUE));
123   AliKalmanTrack::SetNumberOfClusters(t.GetTRDncls()); 
124   Int_t ti[]={-1, -1, -1, -1, -1, -1}; t.GetTRDtracklets(&ti[0]);
125   for(int ip=0; ip<kNplane; ip++){ 
126     fTrackletIndex[ip] = ti[ip];
127     fTracklet[ip]      = NULL;
128   }
129   for(int i =0; i<3; i++) fBudget[i] = 0.;
130   
131   Float_t pid = 1./AliPID::kSPECIES;
132   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
133
134   const AliExternalTrackParam *par = &t;
135   if (t.GetStatus() & AliESDtrack::kTRDbackup) { 
136     par = t.GetOuterParam();
137     if (!par) {
138       AliError("No backup info!"); 
139       par = &t;
140     }
141   }
142   Set(par->GetX() 
143      ,par->GetAlpha()
144      ,par->GetParameter()
145      ,par->GetCovariance());
146
147   if(t.GetStatus() & AliESDtrack::kTIME) {
148     StartTimeIntegral();
149     Double_t times[10]; 
150     t.GetIntegratedTimes(times); 
151     SetIntegratedTimes(times);
152     SetIntegratedLength(t.GetIntegratedLength());
153   }
154 }
155
156 //_______________________________________________________________
157 AliTRDtrackV1::AliTRDtrackV1(AliTRDseedV1 * const trklts, const Double_t p[5], const Double_t cov[15]
158              , Double_t x, Double_t alpha) : AliKalmanTrack()
159   ,fStatus(0)
160   ,fESDid(0)
161   ,fDE(0.)
162   ,fTruncatedMean(0)
163   ,fkReconstructor(NULL)
164   ,fBackupTrack(NULL)
165   ,fTrackLow(NULL)
166   ,fTrackHigh(NULL)
167 {
168   //
169   // The stand alone tracking constructor
170   // TEMPORARY !!!!!!!!!!!
171   // to check :
172   // 1. covariance matrix
173   // 2. dQdl calculation
174   //
175
176   Double_t b(GetBz());
177   Double_t cnv = (TMath::Abs(b) < 1.e-5) ? 1.e5 : 1./GetBz()/kB2C;
178   
179   Double_t pp[5] = { p[0]    
180                     , p[1]
181                     , p[2]
182                     , p[3]
183                     , p[4]*cnv      };
184   
185   Double_t c22 = x*x*cov[14] - 2*x*cov[12] + cov[ 5];
186   Double_t c32 =   x*cov[13] -     cov[ 8];
187   Double_t c20 =   x*cov[10] -     cov[ 3];
188   Double_t c21 =   x*cov[11] -     cov[ 4];
189   Double_t c42 =   x*cov[14] -     cov[12];
190   
191   Double_t cc[15] = { cov[ 0]
192                     , cov[ 1],     cov[ 2]
193                     , c20,         c21,         c22
194                     , cov[ 6],     cov[ 7],     c32,     cov[ 9]
195                     , cov[10]*cnv, cov[11]*cnv, c42*cnv, cov[13]*cnv, cov[14]*cnv*cnv };
196   
197   Double_t mostProbablePt=AliExternalTrackParam::GetMostProbablePt();
198   Double_t p0=TMath::Sign(1/mostProbablePt,pp[4]);
199   Double_t w0=cc[14]/(cc[14] + p0*p0), w1=p0*p0/(cc[14] + p0*p0);
200   AliDebug(3, Form("Pt mixing : w0[%4.2f] pt0[%5.3f] w1[%4.2f] pt[%5.3f]", w0, 1./p0, w1, 1./pp[4]));
201   
202   pp[4] = w0*p0 + w1*pp[4];
203   cc[10]*=w1; cc[11]*=w1; cc[12]*=w1; cc[13]*=w1; cc[14]*=w1;
204   Set(x,alpha,pp,cc);
205   AliDebug(2, Form("Init @ x[%6.2f] pt[%5.3f]", x, 1./pp[4]));
206   Int_t ncls = 0;
207   for(int iplane=0; iplane<kNplane; iplane++){
208     fTrackletIndex[iplane] = -1;
209           if(!trklts[iplane].IsOK()) fTracklet[iplane] = NULL;
210     else{ 
211       fTracklet[iplane] = &trklts[iplane];
212       ncls += fTracklet[iplane]->GetN();
213     }
214   }
215   AliKalmanTrack::SetNumberOfClusters(ncls);            
216   for(int i =0; i<3; i++) fBudget[i] = 0.;
217   
218   Float_t pid = 1./AliPID::kSPECIES;
219   for(int is =0; is<AliPID::kSPECIES; is++) fPID[is] = pid;
220 }
221
222 //_______________________________________________________________
223 AliTRDtrackV1::~AliTRDtrackV1()
224 {
225   // Clean up all objects allocated by the track during its lifetime.
226   AliDebug(2, Form("Deleting track[%d]\n   fBackupTrack[%p] fTrackLow[%p] fTrackHigh[%p] Owner[%c].", fESDid, (void*)fBackupTrack, (void*)fTrackLow, (void*)fTrackHigh, TestBit(kOwner)?'y':'n'));
227
228   if(fBackupTrack) delete fBackupTrack; fBackupTrack = NULL;
229
230   if(fTrackLow) delete fTrackLow; fTrackLow = NULL;
231   if(fTrackHigh) delete fTrackHigh; fTrackHigh = NULL;
232
233   for(Int_t ip=0; ip<kNplane; ip++){
234     if(TestBit(kOwner) && fTracklet[ip]) delete fTracklet[ip];
235     fTracklet[ip] = NULL;
236     fTrackletIndex[ip] = -1;
237   }
238 }
239         
240 //_______________________________________________________________
241 AliTRDtrackV1 &AliTRDtrackV1::operator=(const AliTRDtrackV1 &t)
242 {
243   //
244   // Assignment operator
245   //
246
247   if (this != &t) {
248     AliKalmanTrack::operator=(t);
249     ((AliTRDtrackV1 &) t).Copy(*this);
250   }
251
252   return *this;
253
254 }
255
256 //_____________________________________________________________________________
257 void AliTRDtrackV1::Copy(TObject &t) const
258 {
259   //
260   // Copy function
261   //
262
263   ((AliTRDtrackV1 &) t).fStatus         = fStatus;
264   ((AliTRDtrackV1 &) t).fESDid          = fESDid;
265   ((AliTRDtrackV1 &) t).fDE             = fDE;
266   ((AliTRDtrackV1 &) t).fkReconstructor = fkReconstructor;
267   ((AliTRDtrackV1 &) t).fBackupTrack    = 0x0;
268   ((AliTRDtrackV1 &) t).fTrackLow       = 0x0;
269   ((AliTRDtrackV1 &) t).fTrackHigh      = 0x0;
270
271   for(Int_t ip = 0; ip < kNplane; ip++) { 
272     ((AliTRDtrackV1 &) t).fTrackletIndex[ip] = fTrackletIndex[ip];
273     ((AliTRDtrackV1 &) t).fTracklet[ip]      = fTracklet[ip];
274   }
275   if (fTrackLow) {
276     ((AliTRDtrackV1 &) t).fTrackLow  = new AliExternalTrackParam(*fTrackLow);
277   }
278   if (fTrackHigh){
279     ((AliTRDtrackV1 &) t).fTrackHigh = new AliExternalTrackParam(*fTrackHigh);
280   }
281  
282   for (Int_t i = 0; i < 3; i++) {
283     ((AliTRDtrackV1 &) t).fBudget[i] = fBudget[i];
284   }
285   for (Int_t is = 0; is < AliPID::kSPECIES; is++) {
286     ((AliTRDtrackV1 &) t).fPID[is] = fPID[is];
287   }  
288
289 }
290
291 //_______________________________________________________________
292 Int_t AliTRDtrackV1::CookLabel(Float_t wrong, Int_t *labs, Float_t *freq)
293 {
294 // Set MC label for this track
295 // On demand i.e. if arrays "labs" and "freq" are allocated by user returns :
296 //   nlabs = the no of distinct labels
297 //   labs  = array of distinct labels in decreasing order of frequency
298 //   freq  = frequency of each label  in decreasing order
299
300   Int_t ncl(0);
301   if(!(ncl = GetNumberOfClusters())) return 0;
302
303   Int_t s[2][kMAXCLUSTERSPERTRACK];
304   for (Int_t i = 0; i < kMAXCLUSTERSPERTRACK; i++) {
305     s[0][i] = -1;
306     s[1][i] =  0;
307   }
308
309   Int_t label(-123456789), nlabels(0);
310   AliTRDcluster *c(NULL);
311   for (Int_t ip(0); ip < AliTRDgeometry::kNlayer; ip++) {
312     if(fTrackletIndex[ip]<0 || !fTracklet[ip]) continue;
313     for (Int_t ic(0); ic < AliTRDseedV1::kNclusters; ic++) {
314       if(!(c = fTracklet[ip]->GetClusters(ic))) continue;
315       for (Int_t k(0); k < 3; k++) {
316         if ((label = c->GetLabel(k)) < 0) continue;
317         Int_t j(0);
318         while(j < kMAXCLUSTERSPERTRACK){
319           if(s[0][j]!=label && s[1][j]!=0){j++; continue;}
320           if(!s[1][j]) nlabels++;
321           s[0][j] = label; s[1][j]++;
322           break;
323         }
324       }
325     }
326   }
327   //printf("  Found %4d labels\n", nlabels);
328   Float_t prob(1.);
329   if(!nlabels){
330     AliError(Form("No MC labels found for track %d.", fESDid));
331     return 0;
332   } else if(nlabels==1) {
333     label = s[0][0];
334     if(labs && freq){labs[0]=label; freq[0]=1.;}
335   } else {
336     Int_t idx[kMAXCLUSTERSPERTRACK];
337     TMath::Sort(nlabels, s[1], idx);
338     label = s[0][idx[0]]; prob = s[1][idx[0]]/Float_t(ncl);
339     if(labs && freq){
340       for (Int_t i(0); i<nlabels; i++){
341         labs[i] = s[0][idx[i]];
342         freq[i] = s[1][idx[i]]/Float_t(ncl);
343       }
344     }
345   }
346   SetLabel((1.-prob > wrong)?-label:label);
347   return nlabels;
348 }
349
350 //_______________________________________________________________
351 Bool_t AliTRDtrackV1::CookPID()
352 {
353 //
354 // Cook the PID information for the track by delegating the omonim function of the tracklets. 
355 // Computes the number of tracklets used. The tracklet information are considered independent. 
356 // For the moment no global track measurement of PID is performed as for example to estimate 
357 // bremsstrahlung probability based on global chi2 of the track.
358 //
359 // The status bit AliESDtrack::kTRDpid is set during the call of AliTRDtrackV1::UpdateESDtrack().The PID performance of the 
360 //TRD for tracks with 6 tacklets is displayed below.
361 //Begin_Html
362 //<img src="TRD/trackPID.gif">
363 //End_Html
364 //
365   const AliTRDPIDResponse *pidResponse = AliTRDcalibDB::Instance()->GetPIDResponse(fkReconstructor->GetRecoParam()->GetPIDmethod());
366   if(!pidResponse){
367     AliError("PID Response not available");
368     return kFALSE;
369   }
370   Int_t nslices = pidResponse->GetNumberOfSlices();
371   Double_t dEdx[kNplane * (Int_t)AliTRDPIDResponse::kNslicesNN];
372   Float_t trackletP[kNplane];
373   memset(dEdx, 0, sizeof(Double_t) * kNplane * (Int_t)AliTRDPIDResponse::kNslicesNN);
374   memset(trackletP, 0, sizeof(Float_t)*kNplane);
375   for(Int_t iseed = 0; iseed < kNplane; iseed++){
376     if(!fTracklet[iseed]) continue;
377     trackletP[iseed] = fTracklet[iseed]->GetMomentum();
378     fTracklet[iseed]->SetPID();
379     if(pidResponse->GetPIDmethod() == AliTRDPIDResponse::kLQ1D){
380       dEdx[iseed] = fTracklet[iseed]->GetdQdl();
381     } else {
382       fTracklet[iseed]->CookdEdx(nslices);
383       const Float_t *trackletdEdx = fTracklet[iseed]->GetdEdx();
384       for(Int_t islice = 0; islice < nslices; islice++){
385         dEdx[iseed*nslices + islice] = trackletdEdx[islice];
386       }
387     }
388   }
389   pidResponse->GetResponse(nslices, dEdx, trackletP, fPID);
390
391   //do truncated mean
392   //ncls needs to be included!! todo!!
393   AliTRDdEdxUtils::SetObjPHQ(AliTRDcalibDB::Instance()->GetPHQ());
394   const Double_t mag    = AliTRDdEdxUtils::IsExBOn() ? GetBz()  : -1;
395   const Double_t charge = AliTRDdEdxUtils::IsExBOn() ? Charge() : -1;
396   fTruncatedMean = CookTruncatedMean(0, mag, charge, kTRUE);
397
398   return kTRUE;
399 }
400
401 //___________________________________________________________
402 UChar_t AliTRDtrackV1::GetNumberOfTrackletsPID() const
403 {
404 // Retrieve number of tracklets used for PID calculation. 
405
406   UChar_t nPID = 0;
407   for(int ip=0; ip<kNplane; ip++){
408     if(fTrackletIndex[ip]<0 || !fTracklet[ip]) continue;
409     if(!fTracklet[ip]->IsOK()) continue;
410     
411     nPID++;
412   }
413   return nPID;
414 }
415
416 //_______________________________________________________________
417 AliTRDcluster* AliTRDtrackV1::GetCluster(Int_t id)
418 {
419   // Get the cluster at a certain position in the track
420   Int_t n = 0;
421   for(Int_t ip=0; ip<kNplane; ip++){
422     if(!fTracklet[ip]) continue;
423     if(n+fTracklet[ip]->GetN() <= id){ 
424       n+=fTracklet[ip]->GetN();
425       continue;
426     }
427     AliTRDcluster *c = NULL;
428     for(Int_t ic=AliTRDseedV1::kNclusters; ic--;){
429       if(!(c = fTracklet[ip]->GetClusters(ic))) continue;
430
431       if(n<id){n++; continue;}
432       return c;
433     }
434   }
435   return NULL;
436 }
437
438 //_______________________________________________________________
439 Int_t  AliTRDtrackV1::GetClusterIndex(Int_t id) const
440 {
441   // Get the cluster index at a certain position in the track
442   Int_t n = 0;
443   for(Int_t ip=0; ip<kNplane; ip++){
444     if(!fTracklet[ip]) continue;
445     if(n+fTracklet[ip]->GetN() <= id){ 
446       n+=fTracklet[ip]->GetN();
447       continue;
448     }
449     for(Int_t ic=AliTRDseedV1::kNclusters; ic--;){
450       if(!(fTracklet[ip]->GetClusters(ic))) continue;
451       if(n<id){n++; continue;}
452       return fTracklet[ip]->GetIndexes(ic);
453     }
454   }
455   return -1;
456 }
457
458 //_______________________________________________________________
459 Double_t AliTRDtrackV1::GetPredictedChi2(const AliTRDseedV1 *trklt, Double_t *cov) const
460 {
461 // Compute chi2 between tracklet and track. The value is calculated at the radial position of the track
462 // equal to the reference radial position of the tracklet (see AliTRDseedV1)
463 // 
464 // The chi2 estimator is computed according to the following formula
465 // BEGIN_LATEX
466 // #chi^{2}=(X_{trklt}-X_{track})(C_{trklt}+C_{track})^{-1}(X_{trklt}-X_{track})^{T}
467 // END_LATEX
468 // where X=(y z), the position of the track/tracklet in the yz plane
469 // 
470
471   Double_t p[2]   = { trklt->GetY(), trklt->GetZ()};
472   trklt->GetCovAt(trklt->GetX(), cov);
473   return AliExternalTrackParam::GetPredictedChi2(p, cov);
474 }
475
476 //_______________________________________________________________
477 Int_t AliTRDtrackV1::GetSector() const
478 {
479   return Int_t(GetAlpha()/AliTRDgeometry::GetAlpha() + (GetAlpha()>0. ? 0 : AliTRDgeometry::kNsector));
480 }
481
482 //_______________________________________________________________
483 Bool_t AliTRDtrackV1::IsEqual(const TObject *o) const
484 {
485   // Checks whether two tracks are equal
486   if (!o) return kFALSE;
487   const AliTRDtrackV1 *inTrack = dynamic_cast<const AliTRDtrackV1*>(o);
488   if (!inTrack) return kFALSE;
489   
490   //if ( fPIDquality != inTrack->GetPIDquality() ) return kFALSE;
491   
492   if(memcmp(fPID, inTrack->fPID, AliPID::kSPECIES*sizeof(Double32_t))) return kFALSE;
493   if(memcmp(fBudget, inTrack->fBudget, 3*sizeof(Double32_t))) return kFALSE;
494   if(memcmp(&fDE, &inTrack->fDE, sizeof(Double32_t))) return kFALSE;
495   if(memcmp(&fFakeRatio, &inTrack->fFakeRatio, sizeof(Double32_t))) return kFALSE;
496   if(memcmp(&fChi2, &inTrack->fChi2, sizeof(Double32_t))) return kFALSE;
497   if(memcmp(&fMass, &inTrack->fMass, sizeof(Double32_t))) return kFALSE;
498   if( fLab != inTrack->fLab ) return kFALSE;
499   if( fN != inTrack->fN ) return kFALSE;
500   Double32_t l(0.), in(0.);
501   l = GetIntegratedLength(); in = inTrack->GetIntegratedLength();
502   if(memcmp(&l, &in, sizeof(Double32_t))) return kFALSE;
503   l=GetX(); in=inTrack->GetX();
504   if(memcmp(&l, &in, sizeof(Double32_t))) return kFALSE;
505   l = GetAlpha(); in = inTrack->GetAlpha();
506   if(memcmp(&l, &in, sizeof(Double32_t))) return kFALSE;
507   if(memcmp(GetParameter(), inTrack->GetParameter(), 5*sizeof(Double32_t))) return kFALSE;
508   if(memcmp(GetCovariance(), inTrack->GetCovariance(), 15*sizeof(Double32_t))) return kFALSE;
509   
510   for (Int_t iTracklet = 0; iTracklet < kNplane; iTracklet++){
511     AliTRDseedV1 *curTracklet = fTracklet[iTracklet];
512     AliTRDseedV1 *inTracklet = inTrack->GetTracklet(iTracklet);
513     if (curTracklet && inTracklet){
514       if (! curTracklet->IsEqual(inTracklet) ) {
515         curTracklet->Print();
516         inTracklet->Print();
517         return kFALSE;
518       }
519     } else {
520       // if one tracklet exists, and corresponding 
521       // in other track doesn't - return kFALSE
522       if(inTracklet || curTracklet) return kFALSE;
523     }
524   }
525
526   return kTRUE;
527 }
528
529 //_______________________________________________________________
530 Bool_t AliTRDtrackV1::IsElectron() const
531 {
532   if(GetPID(0) > fkReconstructor->GetRecoParam()->GetPIDThreshold(GetP())) return kTRUE;
533   return kFALSE;
534 }
535
536         
537 //_____________________________________________________________________________
538 Int_t AliTRDtrackV1::MakeBackupTrack()
539 {
540 //
541 // Creates a backup track
542 // TO DO update quality check of the track.
543 //
544
545   Float_t occupancy(0.); Int_t n(0), ncls(0);
546   for(Int_t il(AliTRDgeometry::kNlayer); il--;){ 
547     if(!fTracklet[il]) continue;
548     n++; 
549     occupancy+=fTracklet[il]->GetTBoccupancy()/AliTRDseedV1::kNtb;
550     ncls += fTracklet[il]->GetN();
551   }
552   if(!n) return -1;
553   occupancy/=n;
554
555   //Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);  
556   
557   Int_t failedCutId(0);
558   if(GetChi2()/n > 5.0) failedCutId=1; 
559   if(occupancy < 0.7) failedCutId=2;
560   //if(ratio1 >   0.6) && 
561   //if(ratio0+ratio1  >   1.5) && 
562   if(GetNCross() != 0)  failedCutId=3;
563   if(TMath::Abs(GetSnp()) > 0.85) failedCutId=4;
564   if(ncls < 20) failedCutId=5;
565
566   if(failedCutId){ 
567     AliDebug(2, Form("\n"
568       "chi2/tracklet < 5.0   [%c] %5.2f\n"
569       "occupancy > 0.7       [%c] %4.2f\n"
570       "NCross == 0           [%c] %d\n"
571       "Abs(snp) < 0.85       [%c] %4.2f\n"
572       "NClusters > 20        [%c] %d"
573       ,(GetChi2()/n<5.0)?'y':'n', GetChi2()/n
574       ,(occupancy>0.7)?'y':'n', occupancy
575       ,(GetNCross()==0)?'y':'n', GetNCross()
576       ,(TMath::Abs(GetSnp())<0.85)?'y':'n', TMath::Abs(GetSnp())
577       ,(ncls>20)?'y':'n', ncls
578     ));
579     return failedCutId;
580   }
581
582   if(fBackupTrack) {
583     fBackupTrack->~AliTRDtrackV1();
584     new(fBackupTrack) AliTRDtrackV1((AliTRDtrackV1&)(*this));
585     return 0;
586   }
587   fBackupTrack = new AliTRDtrackV1((AliTRDtrackV1&)(*this));
588   return 0;
589 }
590
591 //_____________________________________________________________________________
592 Int_t AliTRDtrackV1::GetProlongation(Double_t xk, Double_t &y, Double_t &z) const
593 {
594   //
595   // Find a prolongation at given x
596   // Return -1 if it does not exist
597   //  
598
599   Double_t bz = GetBz();
600   if (!AliExternalTrackParam::GetYAt(xk,bz,y)) return -1;
601   if (!AliExternalTrackParam::GetZAt(xk,bz,z)) return -1;
602
603   return 1;  
604
605 }
606
607 //_____________________________________________________________________________
608 Bool_t AliTRDtrackV1::PropagateTo(Double_t xk, Double_t xx0, Double_t xrho)
609 {
610   //
611   // Propagates this track to a reference plane defined by "xk" [cm] 
612   // correcting for the mean crossed material.
613   //
614   // "xx0"  - thickness/rad.length [units of the radiation length] 
615   // "xrho" - thickness*density    [g/cm^2] 
616   // 
617
618   if (TMath::Abs(xk - GetX())<AliTRDReconstructor::GetEpsilon()*0.1) return kTRUE; // 10% of the tracker precision
619
620   Double_t xyz0[3] = {GetX(), GetY(), GetZ()}, // track position BEFORE propagation 
621            b[3];    // magnetic field 
622   GetBxByBz(b);
623   if(!AliExternalTrackParam::PropagateToBxByBz(xk,b)) return kFALSE;
624  
625   // local track position AFTER propagation 
626   Double_t xyz1[3] = {GetX(), GetY(), GetZ()};
627 //  printf("x0[%6.2f] -> x1[%6.2f] dx[%6.2f] rho[%f]\n", xyz0[0], xyz1[0], xyz0[0]-xk, xrho/TMath::Abs(xyz0[0]-xk));
628   if(xyz0[0] < xk) {
629     xrho = -xrho;
630     if (IsStartedTimeIntegral()) {
631       Double_t l2  = TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) 
632                                + (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]) 
633                                + (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
634       Double_t crv = AliExternalTrackParam::GetC(b[2]);
635       if (TMath::Abs(l2*crv) > 0.0001) {
636         // Make correction for curvature if neccesary
637         l2 = 0.5 * TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) 
638                              + (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]));
639         l2 = 2.0 * TMath::ASin(l2 * crv) / crv;
640         l2 = TMath::Sqrt(l2*l2 + (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
641       }
642       AddTimeStep(l2);
643     }
644   }
645   if (!AliExternalTrackParam::CorrectForMeanMaterial(xx0, xrho, fMass)) return kFALSE;
646
647
648   {
649
650     // Energy losses
651     Double_t p2    = (1.0 + GetTgl()*GetTgl()) / (GetSigned1Pt()*GetSigned1Pt());
652     Double_t beta2 = p2 / (p2 + fMass*fMass);
653     if ((beta2 < 1.0e-10) || 
654         ((5940.0 * beta2/(1.0 - beta2 + 1.0e-10) - beta2) < 0.0)) {
655       return kFALSE;
656     }
657
658     Double_t dE    = 0.153e-3 / beta2 
659                    * (TMath::Log(5940.0 * beta2/(1.0 - beta2 + 1.0e-10)) - beta2)
660                    * xrho;
661     fBudget[0] += xrho;
662
663     /*
664     // Suspicious part - think about it ?
665     Double_t kinE =  TMath::Sqrt(p2);
666     if (dE > 0.8*kinE) dE = 0.8 * kinE;  //      
667     if (dE < 0)        dE = 0.0;         // Not valid region for Bethe bloch 
668     */
669  
670     fDE += dE;
671
672     /*
673     // Suspicious ! I.B.
674     Double_t sigmade = 0.07 * TMath::Sqrt(TMath::Abs(dE));   // Energy loss fluctuation 
675     Double_t sigmac2 = sigmade*sigmade*fC*fC*(p2+fMass*fMass)/(p2*p2);
676     fCcc += sigmac2;
677     fCee += fX*fX * sigmac2;  
678     */
679
680   }
681
682   return kTRUE;
683 }
684
685 //_____________________________________________________________________________
686 Int_t   AliTRDtrackV1::PropagateToR(Double_t r,Double_t step)
687 {
688   //
689   // Propagate track to the radial position
690   // Rotation always connected to the last track position
691   //
692
693   Double_t xyz0[3];
694   Double_t xyz1[3];
695   Double_t y;
696   Double_t z; 
697
698   Double_t radius = TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
699   // Direction +-
700   Double_t dir    = (radius > r) ? -1.0 : 1.0;   
701
702   for (Double_t x = radius+dir*step; dir*x < dir*r; x += dir*step) {
703
704     GetXYZ(xyz0);       
705     Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
706     Rotate(alpha,kTRUE);
707     GetXYZ(xyz0);       
708     if(GetProlongation(x,y,z)<0) return -1;
709     xyz1[0] = x * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
710     xyz1[1] = x * TMath::Sin(alpha) - y * TMath::Cos(alpha);
711     xyz1[2] = z;
712     Double_t param[7];
713     if(AliTracker::MeanMaterialBudget(xyz0,xyz1,param)<=0.) return -1;
714     if (param[1] <= 0) {
715       param[1] = 100000000;
716     }
717     PropagateTo(x,param[1],param[0]*param[4]);
718
719   } 
720
721   GetXYZ(xyz0); 
722   Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
723   Rotate(alpha,kTRUE);
724   GetXYZ(xyz0); 
725   if(GetProlongation(r,y,z)<0) return -1;
726   xyz1[0] = r * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
727   xyz1[1] = r * TMath::Sin(alpha) - y * TMath::Cos(alpha);
728   xyz1[2] = z;
729   Double_t param[7];
730   if(AliTracker::MeanMaterialBudget(xyz0,xyz1,param) <= 0.) return -1;
731
732   if (param[1] <= 0) {
733     param[1] = 100000000;
734   }
735   PropagateTo(r,param[1],param[0]*param[4]);
736
737   return 0;
738
739 }
740
741 //_____________________________________________________________________________
742 void AliTRDtrackV1::Print(Option_t *o) const
743 {
744   // Print track status
745   AliInfo(Form("PID [%4.1f %4.1f %4.1f %4.1f %4.1f]", 1.E2*fPID[0], 1.E2*fPID[1], 1.E2*fPID[2], 1.E2*fPID[3], 1.E2*fPID[4]));
746   AliInfo(Form("Material[%5.2f %5.2f %5.2f]", fBudget[0], fBudget[1], fBudget[2]));
747
748   AliInfo(Form("x[%7.2f] t[%7.4f] alpha[%f] mass[%f]", GetX(), GetIntegratedLength(), GetAlpha(), fMass));
749   AliInfo(Form("Ntr[%1d] NtrPID[%1d] Ncl[%3d] lab[%3d]", GetNumberOfTracklets(), GetNumberOfTrackletsPID(), fN, fLab));
750
751   printf("|X| = (");
752   const Double_t *curP = GetParameter();
753   for (Int_t i = 0; i < 5; i++) printf("%7.2f ", curP[i]);
754   printf(")\n");
755
756   printf("|V| = \n");
757   const Double_t *curC = GetCovariance();
758   for (Int_t i = 0, j=4, k=0; i<15; i++, k++){
759     printf("%7.2f ", curC[i]);
760     if(k==j){ 
761       printf("\n");
762       k=-1; j--;
763     }
764   }
765   if(strcmp(o, "a")!=0) return;
766
767   for(Int_t ip=0; ip<kNplane; ip++){
768     if(!fTracklet[ip]) continue;
769     fTracklet[ip]->Print(o);
770   }
771 }
772
773
774 //_____________________________________________________________________________
775 Bool_t AliTRDtrackV1::Rotate(Double_t alpha, Bool_t absolute)
776 {
777   //
778   // Rotates track parameters in R*phi plane
779   // if absolute rotation alpha is in global system
780   // otherwise alpha rotation is relative to the current rotation angle
781   //  
782
783   if (absolute) alpha -= GetAlpha();
784   //else fNRotate++;
785
786   return AliExternalTrackParam::Rotate(GetAlpha()+alpha);
787 }
788
789 //___________________________________________________________
790 void AliTRDtrackV1::SetNumberOfClusters() 
791 {
792 // Calculate the number of clusters attached to this track
793         
794   Int_t ncls = 0;
795   for(int ip=0; ip<kNplane; ip++){
796     if(fTracklet[ip] && fTrackletIndex[ip] >= 0) ncls += fTracklet[ip]->GetN();
797   }
798   AliKalmanTrack::SetNumberOfClusters(ncls);    
799 }
800
801         
802 //_______________________________________________________________
803 void AliTRDtrackV1::SetOwner()
804 {
805   //
806   // Toggle ownership of tracklets
807   //
808
809   if(TestBit(kOwner)) return;
810   for (Int_t ip = 0; ip < kNplane; ip++) {
811     if(fTrackletIndex[ip]<0 || !fTracklet[ip]) continue;
812     fTracklet[ip] = new AliTRDseedV1(*fTracklet[ip]);
813     fTracklet[ip]->SetOwner();
814   }
815   SetBit(kOwner);
816 }
817
818 //_______________________________________________________________
819 void AliTRDtrackV1::SetTracklet(AliTRDseedV1 *const trklt, Int_t index)
820 {
821   //
822   // Set the tracklets
823   //
824   Int_t plane = trklt->GetPlane();
825
826   fTracklet[plane]      = trklt;
827   fTrackletIndex[plane] = index;
828 }
829
830 //_______________________________________________________________
831 void AliTRDtrackV1::SetTrackIn()
832 {
833 //  Save location of birth for the TRD track
834 //  If the pointer is not valid allocate memory
835 //
836   const AliExternalTrackParam *op = dynamic_cast<const AliExternalTrackParam*>(this);
837
838   //printf("SetTrackIn() : fTrackLow[%p]\n", (void*)fTrackLow);
839   if(fTrackLow){
840     fTrackLow->~AliExternalTrackParam();
841     new(fTrackLow) AliExternalTrackParam(*op);
842   } else fTrackLow = new AliExternalTrackParam(*op);
843 }
844
845 //_______________________________________________________________
846 void AliTRDtrackV1::SetTrackOut(const AliExternalTrackParam *op)
847 {
848 //  Save location of death for the TRD track
849 //  If the pointer is not valid allocate memory
850 //
851   if(!op) op = dynamic_cast<const AliExternalTrackParam*>(this);
852   if(fTrackHigh){
853     fTrackHigh->~AliExternalTrackParam();
854     new(fTrackHigh) AliExternalTrackParam(*op);
855   } else fTrackHigh = new AliExternalTrackParam(*op);
856 }
857
858 //_______________________________________________________________
859 void AliTRDtrackV1::UnsetTracklet(Int_t plane)
860 {
861   if(plane<0) return;
862   fTrackletIndex[plane] = -1;
863   fTracklet[plane] = NULL;
864 }
865
866
867 //_______________________________________________________________
868 void AliTRDtrackV1::UpdateChi2(Float_t chi2)
869 {
870 // Update chi2/track with one tracklet contribution
871   SetChi2(GetChi2() + chi2);
872 }
873
874 //_______________________________________________________________
875 void AliTRDtrackV1::UpdateESDtrack(AliESDtrack *track)
876 {
877   //
878   // Update the TRD PID information in the ESD track
879   //
880
881 //   Int_t nslices = AliTRDcalibDB::Instance()->GetPIDResponse(fkReconstructor->GetRecoParam()->GetPIDmethod())->GetNumberOfSlices();
882   // number of tracklets used for PID calculation
883   UChar_t nPID = GetNumberOfTrackletsPID();
884   // number of tracklets attached to the track
885   UChar_t nTrk = GetNumberOfTracklets();
886   // pack the two numbers together and store them in the ESD
887   track->SetTRDntracklets(nPID | (nTrk<<3));
888   // allocate space to store raw PID signals dEdx & momentum
889   track->SetNumberOfTRDslices((AliTRDPIDResponse::kNslicesNN+3)*AliTRDgeometry::kNlayer);
890   // store raw signals
891   Float_t p, sp; Double_t spd;
892   for (Int_t ip = 0; ip < kNplane; ip++) {
893     if(fTrackletIndex[ip]<0 || !fTracklet[ip]) continue;
894     if(!fTracklet[ip]->HasPID()) continue;
895     fTracklet[ip]->CookdEdx(AliTRDPIDResponse::kNslicesNN);
896     const Float_t *dedx = fTracklet[ip]->GetdEdx();
897     for (Int_t js = 0; js < AliTRDPIDResponse::kNslicesNN; js++, dedx++){
898       track->SetTRDslice(*dedx, ip, js+1);
899     }
900     p = fTracklet[ip]->GetMomentum(&sp);
901     spd = sp; track->SetTRDmomentum(p, ip, &spd);
902     // store global quality per tracklet instead of momentum error
903     // 26.09.11 A.Bercuci
904     // first implementation store no. of time bins filled in tracklet (5bits  see "y" bits) and
905     // no. of double clusters in case of pad row cross (4bits see "x" bits)
906     // bit map for tracklet quality xxxxyyyyy
907     // 27.10.11 A.Bercuci
908     // add chamber status bit "z" bit
909     // bit map for tracklet quality zxxxxyyyyy
910     // 12.11.11 A.Bercuci
911     // fit tracklet quality into the field fTRDTimeBin [Char_t]
912     // bit map for tracklet quality zxxyyyyy
913     // The information should be retrieved by the following functions of AliESDtrack for each TRD layer
914     // GetTRDtrkltOccupancy(layer) -> no of TB filled in tracklet
915     // GetTRDtrkltClCross(layer)   -> no of TB filled in crossing pad rows
916     // IsTRDtrkltChmbGood(layer)   -> status of the chamber from which the tracklet is found
917     Int_t nCross(fTracklet[ip]->IsRowCross()?fTracklet[ip]->GetTBcross():0); if(nCross>3) nCross = 3;
918     Char_t trackletQ = Char_t(fTracklet[ip]->GetTBoccupancy() | (nCross<<5) | (fTracklet[ip]->IsChmbGood()<<7));
919     track->SetTRDTimBin(trackletQ, ip);
920     track->SetTRDslice(fTracklet[ip]->GetdQdl(), ip, 0); // Set Summed dEdx into the first slice
921   }
922   // store PID probabilities
923   track->SetTRDpid(fPID);
924
925   //store truncated mean
926   track->SetTRDsignal(fTruncatedMean);
927 }
928
929 //_______________________________________________________________
930 Double_t  AliTRDtrackV1::CookTruncatedMean(const Bool_t kinvq, const Double_t mag, const Int_t charge, const Int_t kcalib, TVectorD *Qs, TVectorD *Xs, Int_t timeBin0, Int_t timeBin1, Int_t tstep) const
931 {
932   //
933   //Origin: Xianguo Lu <xianguo.lu@cern.ch>, Marian Ivanov <marian.ivanov@cern.ch>
934   //
935
936   TVectorD arrayQ(200), arrayX(200);
937   Int_t ncls = AliTRDdEdxUtils::GetArrayClusterQ(kinvq, &arrayQ, &arrayX, this, timeBin0, timeBin1, tstep);
938
939   const TObjArray *cobj = kcalib ? AliTRDdEdxUtils::GetObjPHQ(kinvq, mag, charge) : NULL;
940   
941   const Double_t tmean = AliTRDdEdxUtils::ToyCook(kinvq, ncls, &arrayQ, &arrayX, cobj);
942
943   const Int_t nch = AliTRDdEdxUtils::UpdateArrayX(ncls, &arrayX);
944
945   if(Qs && Xs){
946     (*Qs)=arrayQ;
947     (*Xs)=arrayX;
948   }
949
950   return AliTRDdEdxUtils::GetSignal(nch, ncls, tmean);
951 }
952