]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDtracker.cxx
- Return from Gstpar if material is not used.
[u/mrichter/AliRoot.git] / TRD / AliTRDtracker.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 //                                                                           //
20 //  The standard TRD tracker                                                 //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <Riostream.h>
25 #include <TFile.h>
26 #include <TBranch.h>
27 #include <TTree.h>  
28 #include <TObjArray.h> 
29
30 #include "AliTRDgeometry.h"
31 #include "AliTRDparameter.h"
32 #include "AliTRDgeometryHole.h"
33 #include "AliTRDcluster.h" 
34 #include "AliTRDtrack.h"
35 #include "AliBarrelTrack.h"
36 #include "AliESD.h"
37
38 #include "AliTRDtracker.h"
39
40 ClassImp(AliTRDtracker) 
41
42   const  Float_t     AliTRDtracker::fgkSeedDepth          = 0.5; 
43   const  Float_t     AliTRDtracker::fgkSeedStep           = 0.10;   
44   const  Float_t     AliTRDtracker::fgkSeedGap            = 0.25;  
45
46   const  Float_t     AliTRDtracker::fgkMaxSeedDeltaZ12    = 40.;  
47   const  Float_t     AliTRDtracker::fgkMaxSeedDeltaZ      = 25.;  
48   const  Float_t     AliTRDtracker::fgkMaxSeedC           = 0.0052; 
49   const  Float_t     AliTRDtracker::fgkMaxSeedTan         = 1.2;  
50   const  Float_t     AliTRDtracker::fgkMaxSeedVertexZ     = 150.; 
51
52   const  Double_t    AliTRDtracker::fgkSeedErrorSY        = 0.2;
53   const  Double_t    AliTRDtracker::fgkSeedErrorSY3       = 2.5;
54   const  Double_t    AliTRDtracker::fgkSeedErrorSZ        = 0.1;
55
56   const  Float_t     AliTRDtracker::fgkMinClustersInSeed  = 0.7;  
57
58   const  Float_t     AliTRDtracker::fgkMinClustersInTrack = 0.5;  
59   const  Float_t     AliTRDtracker::fgkMinFractionOfFoundClusters = 0.8;  
60
61   const  Float_t     AliTRDtracker::fgkSkipDepth          = 0.3;
62   const  Float_t     AliTRDtracker::fgkLabelFraction      = 0.8;  
63   const  Float_t     AliTRDtracker::fgkWideRoad           = 20.;
64
65   const  Double_t    AliTRDtracker::fgkMaxChi2            = 12.; 
66
67 const Int_t AliTRDtracker::fgkFirstPlane = 5;
68 const Int_t AliTRDtracker::fgkLastPlane = 17;
69
70
71 //____________________________________________________________________
72 AliTRDtracker::AliTRDtracker():AliTracker(),
73                                fGeom(0),
74                                fPar(0),
75                                fNclusters(0),
76                                fClusters(0),
77                                fNseeds(0),
78                                fSeeds(0),
79                                fNtracks(0),
80                                fTracks(0),
81                                fSY2corr(0),
82                                fSZ2corr(0),
83                                fTimeBinsPerPlane(0),
84                                fMaxGap(0),
85                                fVocal(kFALSE),
86                                fAddTRDseeds(kFALSE),
87                                fNoTilt(kFALSE)
88 {
89   // Default constructor
90
91   for(Int_t i=0;i<kTrackingSectors;i++) fTrSec[i]=0;
92   for(Int_t j=0;j<5;j++)
93     for(Int_t k=0;k<18;k++) fHoles[j][k]=kFALSE;
94
95 //____________________________________________________________________
96 AliTRDtracker::AliTRDtracker(const TFile *geomfile):AliTracker()
97 {
98   // 
99   //  Main constructor
100   //  
101
102   //Float_t fTzero = 0;
103    
104   fAddTRDseeds = kFALSE;
105   fGeom = NULL;
106   fNoTilt = kFALSE;
107   
108   TDirectory *savedir=gDirectory; 
109   TFile *in=(TFile*)geomfile;  
110   if (!in->IsOpen()) {
111     printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
112     printf("    DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
113   }
114   else {
115     in->cd();  
116 //    in->ls();
117     fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
118     fPar  = (AliTRDparameter*) in->Get("TRDparameter");
119 //    fGeom->Dump();
120   }
121
122   if(fGeom) {
123     //    fTzero = geo->GetT0();
124     printf("Found geometry version %d on file \n", fGeom->IsVersion());
125   }
126   else { 
127     printf("AliTRDtracker::AliTRDtracker(): can't find TRD geometry!\n");
128     //printf("The DETAIL TRD geometry will be used\n");
129     //fGeom = new AliTRDgeometryDetail();
130     fGeom = new AliTRDgeometryHole();
131     fGeom->SetPHOShole();
132     fGeom->SetRICHhole();    
133   } 
134
135   if (!fPar) {  
136     printf("AliTRDtracker::AliTRDtracker(): can't find TRD parameter!\n");
137     printf("The DEFAULT TRD parameter will be used\n");
138     fPar = new AliTRDparameter();
139   }
140   fPar->ReInit();
141
142   savedir->cd();  
143
144
145   //  fGeom->SetT0(fTzero);
146
147   fNclusters = 0;
148   fClusters  = new TObjArray(2000); 
149   fNseeds    = 0;
150   fSeeds     = new TObjArray(2000);
151   fNtracks   = 0;
152   fTracks    = new TObjArray(1000);
153
154   for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
155     Int_t trS = CookSectorIndex(geomS);
156     fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS, fPar);
157     for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
158       fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
159     }
160   }
161
162   Float_t tiltAngle = TMath::Abs(fPar->GetTiltingAngle()); 
163   if(tiltAngle < 0.1) {
164     fNoTilt = kTRUE;
165   }
166
167   fSY2corr = 0.2;
168   fSZ2corr = 120.;      
169
170   if(fNoTilt && (tiltAngle > 0.1)) fSY2corr = fSY2corr + tiltAngle * 0.05; 
171
172
173   // calculate max gap on track
174
175   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
176   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
177
178   Double_t dx = (Double_t) fPar->GetTimeBinSize();   
179   Int_t tbAmp = fPar->GetTimeBefore();
180   Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
181   if(kTRUE) maxAmp = 0;  // intentional until we change the parameter class 
182   Int_t tbDrift = fPar->GetTimeMax();
183   Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
184
185   tbDrift = TMath::Min(tbDrift,maxDrift);
186   tbAmp = TMath::Min(tbAmp,maxAmp);
187
188   fTimeBinsPerPlane = tbAmp + tbDrift;
189   fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fgkSkipDepth);
190
191   fVocal = kFALSE;
192
193
194   // Barrel Tracks [SR, 03.04.2003]
195
196   fBarrelFile = 0;
197   fBarrelTree = 0;
198   fBarrelArray = 0;
199   fBarrelTrack = 0;
200
201   savedir->cd();
202 }   
203
204 //___________________________________________________________________
205 AliTRDtracker::~AliTRDtracker()
206 {
207   //
208   // Destructor of AliTRDtracker 
209   //
210
211   if (fClusters) {
212     fClusters->Delete();
213     delete fClusters;
214   }
215   if (fTracks) {
216     fTracks->Delete();
217     delete fTracks;
218   }
219   if (fSeeds) {
220     fSeeds->Delete();
221     delete fSeeds;
222   }
223   delete fGeom;  
224   delete fPar;  
225
226   for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
227     delete fTrSec[geomS];
228   }
229 }   
230
231 //_____________________________________________________________________
232
233 void AliTRDtracker::SetBarrelTree(const char *mode) {
234   //
235   //
236   //
237
238   if (!IsStoringBarrel()) return;
239
240   TDirectory *sav = gDirectory;
241   if (!fBarrelFile) fBarrelFile = new TFile("AliBarrelTracks.root", "UPDATE");
242
243   char buff[40];
244   sprintf(buff,  "BarrelTRD_%d_%s", GetEventNumber(), mode);
245
246   fBarrelFile->cd();
247   fBarrelTree = new TTree(buff, "Barrel TPC tracks");
248   
249   Int_t nRefs = fgkLastPlane - fgkFirstPlane + 1;
250
251   if (!fBarrelArray) fBarrelArray = new TClonesArray("AliBarrelTrack", nRefs);
252   for(Int_t i=0; i<nRefs; i++) new((*fBarrelArray)[i]) AliBarrelTrack();
253   
254   fBarrelTree->Branch("tracks", &fBarrelArray);
255   sav->cd();
256 }
257   
258 //_____________________________________________________________________
259
260 void AliTRDtracker::StoreBarrelTrack(AliTRDtrack *ps, Int_t refPlane, Int_t isIn) {
261   //
262   //
263   //
264   
265   if (!IsStoringBarrel()) return;
266   
267   static Int_t nClusters;
268   static Int_t nWrong;
269   static Double_t chi2;
270   static Int_t index;
271   static Bool_t wasLast = kTRUE;
272   
273   Int_t newClusters, newWrong;
274   Double_t newChi2;
275   
276   if (wasLast) {   
277  
278     fBarrelArray->Clear();
279     nClusters = nWrong = 0;
280     chi2 = 0.0;
281     index = 0;
282     wasLast = kFALSE;
283   }
284   
285   fBarrelTrack = (AliBarrelTrack*)(*fBarrelArray)[index++];
286   ps->GetBarrelTrack(fBarrelTrack);
287   
288   newClusters = ps->GetNumberOfClusters() - nClusters; 
289   newWrong = ps->GetNWrong() - nWrong;
290   newChi2 = ps->GetChi2() - chi2;
291   
292   nClusters =  ps->GetNumberOfClusters();
293   nWrong = ps->GetNWrong();
294   chi2 = ps->GetChi2();  
295
296   if (refPlane != fgkLastPlane) {
297     fBarrelTrack->SetNClusters(newClusters, newChi2);
298     fBarrelTrack->SetNWrongClusters(newWrong);
299   } else {
300     wasLast = kTRUE;
301   } 
302
303   fBarrelTrack->SetRefPlane(refPlane, isIn);
304 }
305
306 //_____________________________________________________________________
307
308 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track) {
309   //
310   // Rotates the track when necessary
311   //
312
313   Double_t alpha = AliTRDgeometry::GetAlpha(); 
314   Double_t y = track->GetY();
315   Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
316
317   //Int_t ns = AliTRDgeometry::kNsect;
318   //Int_t s=Int_t(track->GetAlpha()/alpha)%ns; 
319
320   if (y > ymax) {
321     //s = (s+1) % ns;
322     if (!track->Rotate(alpha)) return kFALSE;
323   } else if (y <-ymax) {
324     //s = (s-1+ns) % ns;                           
325     if (!track->Rotate(-alpha)) return kFALSE;   
326   } 
327
328   return kTRUE;
329 }
330
331 //_____________________________________________________________________
332 inline Double_t f1trd(Double_t x1,Double_t y1,
333                       Double_t x2,Double_t y2,
334                       Double_t x3,Double_t y3)
335 {
336   //
337   // Initial approximation of the track curvature
338   //
339   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
340   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
341                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
342   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
343                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
344
345   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
346
347   return -xr*yr/sqrt(xr*xr+yr*yr);
348 }          
349
350 //_____________________________________________________________________
351 inline Double_t f2trd(Double_t x1,Double_t y1,
352                       Double_t x2,Double_t y2,
353                       Double_t x3,Double_t y3)
354 {
355   //
356   // Initial approximation of the track curvature times X coordinate
357   // of the center of curvature
358   //
359
360   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
361   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
362                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
363   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
364                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
365
366   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
367
368   return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
369 }          
370
371 //_____________________________________________________________________
372 inline Double_t f3trd(Double_t x1,Double_t y1,
373                       Double_t x2,Double_t y2,
374                       Double_t z1,Double_t z2)
375 {
376   //
377   // Initial approximation of the tangent of the track dip angle
378   //
379
380   return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
381 }            
382
383
384 AliTRDcluster * AliTRDtracker::GetCluster(AliTRDtrack * track, Int_t plane, Int_t timebin){
385   //
386   //try to find cluster in the backup list
387   //
388   AliTRDcluster * cl =0;
389   UInt_t *indexes = track->GetBackupIndexes();
390   for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
391     if (indexes[i]==0) break;  
392     AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
393     if (!cli) break;
394     if (cli->GetLocalTimeBin()!=timebin) continue;
395     Int_t iplane = fGeom->GetPlane(cli->GetDetector());
396     if (iplane==plane) {
397       cl = cli;
398       break;
399     }
400   }
401   return cl;
402 }
403
404
405 Int_t  AliTRDtracker::GetLastPlane(AliTRDtrack * track){
406   //
407   //return last updated plane
408   Int_t lastplane=0;
409   UInt_t *indexes = track->GetBackupIndexes();
410   for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
411     AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
412     if (!cli) break;
413     Int_t iplane = fGeom->GetPlane(cli->GetDetector());
414     if (iplane>lastplane) {
415       lastplane = iplane;
416     }
417   }
418   return lastplane;
419 }
420 //___________________________________________________________________
421 Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
422 {
423   //
424   // Finds tracks within the TRD. The ESD event is expected to contain seeds 
425   // at the outer part of the TRD. The seeds
426   // are found within the TRD if fAddTRDseeds is TRUE. 
427   // The tracks are propagated to the innermost time bin 
428   // of the TRD and the ESD event is updated
429   //
430
431   Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
432   Float_t foundMin = fgkMinClustersInTrack * timeBins; 
433   Int_t nseed = 0;
434   Int_t found = 0;
435   Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
436
437   Int_t n = event->GetNumberOfTracks();
438   for (Int_t i=0; i<n; i++) {
439     AliESDtrack* seed=event->GetTrack(i);
440     ULong_t status=seed->GetStatus();
441     if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
442     if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
443     nseed++;
444
445     AliTRDtrack* seed2 = new AliTRDtrack(*seed);
446     //seed2->ResetCovariance(); 
447     AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
448     AliTRDtrack &t=*pt; 
449     FollowProlongation(t, innerTB); 
450     if (t.GetNumberOfClusters() >= foundMin) {
451       UseClusters(&t);
452       CookLabel(pt, 1-fgkLabelFraction);
453       //      t.CookdEdx();
454     }
455     found++;
456 //    cout<<found<<'\r';     
457
458     if(PropagateToTPC(t)) {
459       seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
460     }  
461     delete seed2;
462     delete pt;
463   }     
464
465   cout<<"Number of loaded seeds: "<<nseed<<endl;  
466   cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
467
468   // after tracks from loaded seeds are found and the corresponding 
469   // clusters are used, look for additional seeds from TRD
470
471   if(fAddTRDseeds) { 
472     // Find tracks for the seeds in the TRD
473     Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
474   
475     Int_t nSteps = (Int_t) (fgkSeedDepth / fgkSeedStep);
476     Int_t gap = (Int_t) (timeBins * fgkSeedGap);
477     Int_t step = (Int_t) (timeBins * fgkSeedStep);
478   
479     // make a first turn with tight cut on initial curvature
480     for(Int_t turn = 1; turn <= 2; turn++) {
481       if(turn == 2) {
482         nSteps = (Int_t) (fgkSeedDepth / (3*fgkSeedStep));
483         step = (Int_t) (timeBins * (3*fgkSeedStep));
484       }
485       for(Int_t i=0; i<nSteps; i++) {
486         Int_t outer=timeBins-1-i*step; 
487         Int_t inner=outer-gap;
488
489         nseed=fSeeds->GetEntriesFast();
490       
491         MakeSeeds(inner, outer, turn);
492       
493         nseed=fSeeds->GetEntriesFast();
494         //        printf("\n turn %d, step %d: number of seeds for TRD inward %d\n", 
495         //               turn, i, nseed); 
496               
497         for (Int_t i=0; i<nseed; i++) {   
498           AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt; 
499           FollowProlongation(t,innerTB); 
500           if (t.GetNumberOfClusters() >= foundMin) {
501             UseClusters(&t);
502             CookLabel(pt, 1-fgkLabelFraction);
503             t.CookdEdx();
504             found++;
505 //            cout<<found<<'\r';     
506             if(PropagateToTPC(t)) {
507               AliESDtrack track;
508               track.UpdateTrackParams(pt,AliESDtrack::kTRDin);
509               event->AddTrack(&track);
510               //              track.SetTRDtrack(new AliTRDtrack(*pt));
511             }        
512           }
513           delete fSeeds->RemoveAt(i);
514           fNseeds--;
515         }
516       }
517     }
518   }
519   
520   cout<<"Total number of found tracks: "<<found<<endl;
521     
522   return 0;    
523 }     
524      
525   
526
527 //_____________________________________________________________________________
528 Int_t AliTRDtracker::PropagateBack(AliESD* event) {
529   //
530   // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
531   // backpropagated by the TPC tracker. Each seed is first propagated 
532   // to the TRD, and then its prolongation is searched in the TRD.
533   // If sufficiently long continuation of the track is found in the TRD
534   // the track is updated, otherwise it's stored as originaly defined 
535   // by the TPC tracker.   
536   //  
537
538   Int_t found=0;  
539   Float_t foundMin = 20;
540
541   Int_t n = event->GetNumberOfTracks();
542   for (Int_t i=0; i<n; i++) {
543     AliESDtrack* seed=event->GetTrack(i);
544     ULong_t status=seed->GetStatus();
545     if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
546     if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
547
548     Int_t lbl = seed->GetLabel();
549     AliTRDtrack *track = new AliTRDtrack(*seed);
550     track->SetSeedLabel(lbl);
551     fNseeds++;
552     Float_t p4 = track->GetC();
553     //
554     Int_t expectedClr = FollowBackProlongation(*track);
555     if (track->GetNumberOfClusters()<expectedClr/3){
556       AliTRDtrack *track1 = new AliTRDtrack(*seed);
557       track1->SetSeedLabel(lbl);
558       FollowBackProlongation(*track1);
559       AliTRDtrack *track2= new AliTRDtrack(*seed);
560       track->SetSeedLabel(lbl);
561       FollowBackProlongation(*track2);      
562       delete track1;
563       delete track2;
564     }
565     if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)>0.2) {
566       delete track;
567       continue; //too big change of curvature - to be checked
568     }
569     Int_t foundClr = track->GetNumberOfClusters();
570     if (foundClr >= foundMin) {
571       if(foundClr >= 2) {
572         track->CookdEdx(); 
573         CookLabel(track, 1-fgkLabelFraction);
574         UseClusters(track);
575       }
576       
577       // Propagate to outer reference plane [SR, GSI, 18.02.2003]
578       //      track->PropagateTo(364.8);  why?
579       
580       //seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
581       //found++;
582       if (track->GetNCross()==0) seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
583       else{
584         if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
585       }
586     }
587
588     //Propagation to the TOF (I.Belikov)
589     
590     if (track->GetStop()==kFALSE){
591
592       Double_t xtof=371.;
593       Double_t c2=track->GetC()*xtof - track->GetEta();
594       if (TMath::Abs(c2)>=0.85) {
595         delete track;
596         continue;
597       }
598       Double_t xTOF0 = 371. ;          
599       PropagateToOuterPlane(*track,xTOF0); 
600       //      
601       Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
602       Double_t y=track->GetYat(xtof);
603       if (y > ymax) {
604         if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
605           delete track;
606           continue;
607         }
608       } else if (y <-ymax) {
609         if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
610           delete track;
611           continue;
612         }
613       }
614       
615       if (track->PropagateTo(xtof)) {
616         seed->UpdateTrackParams(track, AliESDtrack::kTRDout);    
617         seed->SetTRDtrack(new AliTRDtrack(*track));
618         if (track->GetNumberOfClusters()>foundMin) found++;
619       }
620     }else{
621       if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
622         seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
623         //seed->SetStatus(AliESDtrack::kTRDStop);    
624         seed->SetTRDtrack(new AliTRDtrack(*track));
625         found++;
626       }
627     }
628
629     delete track;
630     
631     //End of propagation to the TOF
632     //if (foundClr>foundMin)
633     //  seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
634     
635
636   }
637   
638   cerr<<"Number of seeds: "<<fNseeds<<endl;  
639   cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
640
641   fSeeds->Clear(); fNseeds=0;
642
643   return 0;
644
645 }
646
647 //_____________________________________________________________________________
648 Int_t AliTRDtracker::RefitInward(AliESD* event)
649 {
650   //
651   // Refits tracks within the TRD. The ESD event is expected to contain seeds 
652   // at the outer part of the TRD. 
653   // The tracks are propagated to the innermost time bin 
654   // of the TRD and the ESD event is updated
655   // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
656   //
657
658   Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
659   Float_t foundMin = fgkMinClustersInTrack * timeBins; 
660   Int_t nseed = 0;
661   Int_t found = 0;
662   Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
663
664   Int_t n = event->GetNumberOfTracks();
665   for (Int_t i=0; i<n; i++) {
666     AliESDtrack* seed=event->GetTrack(i);
667     ULong_t status=seed->GetStatus();
668     if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
669     if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
670     nseed++;
671
672     AliTRDtrack* seed2 = new AliTRDtrack(*seed);
673     seed2->ResetCovariance(5.); 
674     AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
675     UInt_t * indexes2 = seed2->GetIndexes();
676     UInt_t * indexes3 = pt->GetBackupIndexes();
677     for (Int_t i=0;i<200;i++) {
678       if (indexes2[i]==0) break;
679       indexes3[i] = indexes2[i];
680     }          
681     //AliTRDtrack *pt = seed2;
682     AliTRDtrack &t=*pt; 
683     FollowProlongation(t, innerTB); 
684     /*
685     if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
686       // debug  - why we dont go back?
687       AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
688       UInt_t * indexes2 = seed2->GetIndexes();
689       UInt_t * indexes3 = pt2->GetBackupIndexes();
690       for (Int_t i=0;i<200;i++) {
691         if (indexes2[i]==0) break;
692         indexes3[i] = indexes2[i];
693       }  
694       FollowProlongation(*pt2, innerTB);
695       delete pt2;
696     }
697     */
698     if (t.GetNumberOfClusters() >= foundMin) {
699       //      UseClusters(&t);
700       //CookLabel(pt, 1-fgkLabelFraction);
701       //      t.CookdEdx();
702     }
703     found++;
704 //    cout<<found<<'\r';     
705
706     if(PropagateToTPC(t)) {
707       seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
708     }else{
709       //if not prolongation to TPC - propagate without update
710       AliTRDtrack* seed2 = new AliTRDtrack(*seed);
711       seed2->ResetCovariance(5.); 
712       AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
713       delete seed2;
714       if (PropagateToTPC(*pt2)) 
715         seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
716       delete pt2;
717     }  
718
719     delete seed2;
720     delete pt;
721   }     
722
723   cout<<"Number of loaded seeds: "<<nseed<<endl;  
724   cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
725
726   return 0;
727
728 }
729
730
731 //---------------------------------------------------------------------------
732 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
733 {
734   // Starting from current position on track=t this function tries
735   // to extrapolate the track up to timeBin=0 and to confirm prolongation
736   // if a close cluster is found. Returns the number of clusters
737   // expected to be found in sensitive layers
738
739   Float_t  wIndex, wTB, wChi2;
740   Float_t  wYrt, wYclosest, wYcorrect, wYwindow;
741   Float_t  wZrt, wZclosest, wZcorrect, wZwindow;
742   Float_t  wPx, wPy, wPz, wC;
743   Double_t px, py, pz;
744   Float_t  wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
745   Int_t lastplane = GetLastPlane(&t);
746
747   Int_t trackIndex = t.GetLabel();  
748
749   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
750
751   Int_t tryAgain=fMaxGap;
752
753   Double_t alpha=t.GetAlpha();
754   alpha = TVector2::Phi_0_2pi(alpha);
755
756   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
757   Double_t radLength, rho, x, dx, y, ymax, z;
758
759   Int_t expectedNumberOfClusters = 0;
760   Bool_t lookForCluster;
761
762   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
763
764  
765   for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) { 
766
767     y = t.GetY(); z = t.GetZ();
768
769     // first propagate to the inner surface of the current time bin 
770     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
771     x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
772     if(!t.PropagateTo(x,radLength,rho)) break;
773     y = t.GetY();
774     ymax = x*TMath::Tan(0.5*alpha);
775     if (y > ymax) {
776       s = (s+1) % ns;
777       if (!t.Rotate(alpha)) break;
778       if(!t.PropagateTo(x,radLength,rho)) break;
779     } else if (y <-ymax) {
780       s = (s-1+ns) % ns;                           
781       if (!t.Rotate(-alpha)) break;   
782       if(!t.PropagateTo(x,radLength,rho)) break;
783     } 
784
785     y = t.GetY(); z = t.GetZ();
786
787     // now propagate to the middle plane of the next time bin 
788     fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
789     x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
790     if(!t.PropagateTo(x,radLength,rho)) break;
791     y = t.GetY();
792     ymax = x*TMath::Tan(0.5*alpha);
793     if (y > ymax) {
794       s = (s+1) % ns;
795       if (!t.Rotate(alpha)) break;
796       if(!t.PropagateTo(x,radLength,rho)) break;
797     } else if (y <-ymax) {
798       s = (s-1+ns) % ns;                           
799       if (!t.Rotate(-alpha)) break;   
800       if(!t.PropagateTo(x,radLength,rho)) break;
801     } 
802
803
804     if(lookForCluster) {
805
806       expectedNumberOfClusters++;       
807       wIndex = (Float_t) t.GetLabel();
808       wTB = nr;
809
810       AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
811
812       Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
813       Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
814
815       Double_t road;
816       if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
817       else return expectedNumberOfClusters;
818       
819       wYrt = (Float_t) y;
820       wZrt = (Float_t) z;
821       wYwindow = (Float_t) road;
822       t.GetPxPyPz(px,py,pz);
823       wPx = (Float_t) px;
824       wPy = (Float_t) py;
825       wPz = (Float_t) pz;
826       wC  = (Float_t) t.GetC();
827       wSigmaC2 = (Float_t) t.GetSigmaC2();
828       wSigmaTgl2    = (Float_t) t.GetSigmaTgl2();
829       wSigmaY2 = (Float_t) t.GetSigmaY2();
830       wSigmaZ2 = (Float_t) t.GetSigmaZ2();
831       wChi2 = -1;            
832       
833
834       AliTRDcluster *cl=0;
835       UInt_t index=0;
836
837       Double_t maxChi2=fgkMaxChi2;
838
839       wYclosest = 12345678;
840       wYcorrect = 12345678;
841       wZclosest = 12345678;
842       wZcorrect = 12345678;
843       wZwindow  = TMath::Sqrt(2.25 * 12 * sz2);   
844
845       // Find the closest correct cluster for debugging purposes
846       if (timeBin) {
847         Float_t minDY = 1000000;
848         for (Int_t i=0; i<timeBin; i++) {
849           AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
850           if((c->GetLabel(0) != trackIndex) &&
851              (c->GetLabel(1) != trackIndex) &&
852              (c->GetLabel(2) != trackIndex)) continue;
853           if(TMath::Abs(c->GetY() - y) > minDY) continue;
854           minDY = TMath::Abs(c->GetY() - y);
855           wYcorrect = c->GetY();
856           wZcorrect = c->GetZ();
857
858           Double_t h01 = GetTiltFactor(c);
859           wChi2 = t.GetPredictedChi2(c, h01);
860         }
861       }                    
862
863       // Now go for the real cluster search
864
865       if (timeBin) {
866         //
867         //find cluster in history
868         cl =0;
869         
870         AliTRDcluster * cl0 = timeBin[0];
871         if (!cl0) {
872           continue;
873         }
874         Int_t plane = fGeom->GetPlane(cl0->GetDetector());
875         if (plane>lastplane) continue;
876         Int_t timebin = cl0->GetLocalTimeBin();
877         AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
878         if (cl2) {
879           cl =cl2;      
880           Double_t h01 = GetTiltFactor(cl);
881           maxChi2=t.GetPredictedChi2(cl,h01);
882         }
883         if ((!cl) && road>fgkWideRoad) {
884           //if (t.GetNumberOfClusters()>4)
885           //  cerr<<t.GetNumberOfClusters()
886           //    <<"FindProlongation warning: Too broad road !\n";
887           continue;
888         }             
889
890
891         if(!cl){
892
893           for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
894             AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
895             if (c->GetY() > y+road) break;
896             if (c->IsUsed() > 0) continue;
897             if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
898             
899             Double_t h01 = GetTiltFactor(c);
900             Double_t chi2=t.GetPredictedChi2(c,h01);
901             
902             if (chi2 > maxChi2) continue;
903             maxChi2=chi2;
904             cl=c;
905             index=timeBin.GetIndex(i);
906           }               
907         }
908
909         if(!cl) {
910
911           for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
912             AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
913             
914             if (c->GetY() > y+road) break;
915             if (c->IsUsed() > 0) continue;
916             if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
917             
918             Double_t h01 = GetTiltFactor(c);
919             Double_t chi2=t.GetPredictedChi2(c, h01);
920             
921             if (chi2 > maxChi2) continue;
922             maxChi2=chi2;
923             cl=c;
924             index=timeBin.GetIndex(i);
925           }
926         }        
927         if (cl) {
928           wYclosest = cl->GetY();
929           wZclosest = cl->GetZ();
930           Double_t h01 = GetTiltFactor(cl);
931
932           t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters()); 
933           //printf("Track   position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
934           //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
935           Int_t det = cl->GetDetector();    
936           Int_t plane = fGeom->GetPlane(det);
937
938           if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
939             //if(!t.Update(cl,maxChi2,index,h01)) {
940             //if(!tryAgain--) return 0;
941           }  
942           else tryAgain=fMaxGap;
943         }
944         else {
945           //if (tryAgain==0) break; 
946           tryAgain--;
947         }
948
949         /*
950         if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
951           
952           printf(" %f", wIndex);       //1
953           printf(" %f", wTB);          //2
954           printf(" %f", wYrt);         //3
955           printf(" %f", wYclosest);    //4
956           printf(" %f", wYcorrect);    //5
957           printf(" %f", wYwindow);     //6
958           printf(" %f", wZrt);         //7
959           printf(" %f", wZclosest);    //8
960           printf(" %f", wZcorrect);    //9
961           printf(" %f", wZwindow);     //10
962           printf(" %f", wPx);          //11
963           printf(" %f", wPy);          //12
964           printf(" %f", wPz);          //13
965           printf(" %f", wSigmaC2*1000000);  //14
966           printf(" %f", wSigmaTgl2*1000);   //15
967           printf(" %f", wSigmaY2);     //16
968           //      printf(" %f", wSigmaZ2);     //17
969           printf(" %f", wChi2);     //17
970           printf(" %f", wC);           //18
971           printf("\n");
972         } 
973         */                        
974       }
975     }  
976   }
977   return expectedNumberOfClusters;
978   
979   
980 }                
981
982 //___________________________________________________________________
983
984 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
985 {
986   // Starting from current radial position of track <t> this function
987   // extrapolates the track up to outer timebin and in the sensitive
988   // layers confirms prolongation if a close cluster is found. 
989   // Returns the number of clusters expected to be found in sensitive layers
990
991
992   Float_t  wIndex, wTB, wChi2;
993   Float_t  wYrt, wYclosest, wYcorrect, wYwindow;
994   Float_t  wZrt, wZclosest, wZcorrect, wZwindow;
995   Float_t  wPx, wPy, wPz, wC;
996   Double_t px, py, pz;
997   Float_t  wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
998
999   Int_t trackIndex = t.GetLabel();  
1000   Int_t tryAgain=fMaxGap;
1001
1002   Double_t alpha=t.GetAlpha();
1003   TVector2::Phi_0_2pi(alpha);
1004
1005   Int_t s;
1006
1007   Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
1008   Double_t radLength, rho, x, dx, y, ymax = 0, z;
1009   Bool_t lookForCluster;
1010
1011   Int_t expectedNumberOfClusters = 0;
1012   x = t.GetX();
1013
1014   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1015
1016   Int_t nRefPlane = fgkFirstPlane;
1017   Bool_t isNewLayer = kFALSE; 
1018
1019   Double_t chi2;
1020   Double_t minDY;
1021   Int_t zone =-10;
1022   Int_t nr;
1023   for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) { 
1024     
1025     y = t.GetY(); 
1026     z = t.GetZ();
1027
1028     // first propagate to the outer surface of the current time bin 
1029
1030     s = t.GetSector();
1031     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1032     x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; 
1033     y = t.GetY(); 
1034     z = t.GetZ();
1035
1036     if(!t.PropagateTo(x,radLength,rho)) break;
1037     //    if (!AdjustSector(&t)) break;
1038     //
1039     // MI -fix untill correct material desription will be implemented
1040     //
1041     Float_t angle =  t.GetAlpha();  // MI - if rotation - we go through the material 
1042     if (!AdjustSector(&t)) break;
1043     Int_t cross = kFALSE;
1044     
1045     if (TMath::Abs(angle -  t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
1046     Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1047     if (currentzone==-10) cross = kTRUE;  // we are in the frame
1048     if (currentzone>-10){   // layer knows where we are
1049       if (zone==-10) zone = currentzone;
1050       if (zone!=currentzone) cross=kTRUE;  
1051     }
1052     if (cross) {
1053       t.IncCross();
1054       if (t.GetNCross()==1) t.MakeBackupTrack();
1055       if (t.GetNCross()>2) break;
1056     }
1057     
1058     //
1059     //
1060     s = t.GetSector();
1061     if (!t.PropagateTo(x,radLength,rho)) break;
1062
1063     y = t.GetY();
1064     z = t.GetZ();
1065
1066     // Barrel Tracks [SR, 04.04.2003]
1067
1068     s = t.GetSector();
1069     if (fTrSec[s]->GetLayer(nr)->IsSensitive() != 
1070         fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1071
1072 //      if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1073     }
1074
1075     if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() && 
1076           ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1077       isNewLayer = kTRUE;
1078     } else {isNewLayer = kFALSE;}
1079
1080     y = t.GetY();
1081     z = t.GetZ();
1082
1083     // now propagate to the middle plane of the next time bin 
1084     fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1085
1086     x = fTrSec[s]->GetLayer(nr+1)->GetX(); 
1087       if(!t.PropagateTo(x,radLength,rho)) break;
1088     if (!AdjustSector(&t)) break;
1089     s = t.GetSector();
1090       if(!t.PropagateTo(x,radLength,rho)) break;
1091
1092     y = t.GetY();
1093     z = t.GetZ();
1094
1095     if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1096     //    printf("label %d, pl %d, lookForCluster %d \n",
1097     //     trackIndex, nr+1, lookForCluster);
1098
1099     if(lookForCluster) {
1100       expectedNumberOfClusters++;       
1101
1102       wIndex = (Float_t) t.GetLabel();
1103       wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1104
1105       AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1106       Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1107       Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1108       if((t.GetSigmaY2() + sy2) < 0) break;
1109       Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2); 
1110       Double_t y=t.GetY(), z=t.GetZ();
1111
1112       wYrt = (Float_t) y;
1113       wZrt = (Float_t) z;
1114       wYwindow = (Float_t) road;
1115       t.GetPxPyPz(px,py,pz);
1116       wPx = (Float_t) px;
1117       wPy = (Float_t) py;
1118       wPz = (Float_t) pz;
1119       wC  = (Float_t) t.GetC();
1120       wSigmaC2 = (Float_t) t.GetSigmaC2();
1121       wSigmaTgl2    = (Float_t) t.GetSigmaTgl2();
1122       wSigmaY2 = (Float_t) t.GetSigmaY2();
1123       wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1124       wChi2 = -1;            
1125       
1126       if (road>fgkWideRoad) {
1127         if (t.GetNumberOfClusters()>4)
1128           cerr<<t.GetNumberOfClusters()
1129               <<"FindProlongation warning: Too broad road !\n";
1130         return 0;
1131       }      
1132
1133       AliTRDcluster *cl=0;
1134       UInt_t index=0;
1135
1136       Double_t maxChi2=fgkMaxChi2;
1137
1138       if (isNewLayer) { 
1139         road = 3 * road;
1140         //sz2 = 3 * sz2;
1141         maxChi2 = 10 * fgkMaxChi2;
1142       }
1143       
1144       if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2; 
1145       if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1146       if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1147       
1148
1149       wYclosest = 12345678;
1150       wYcorrect = 12345678;
1151       wZclosest = 12345678;
1152       wZcorrect = 12345678;
1153       wZwindow  = TMath::Sqrt(2.25 * 12 * sz2);   
1154
1155       // Find the closest correct cluster for debugging purposes
1156       if (timeBin) {
1157         minDY = 1000000;
1158         for (Int_t i=0; i<timeBin; i++) {
1159           AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1160           if((c->GetLabel(0) != trackIndex) &&
1161              (c->GetLabel(1) != trackIndex) &&
1162              (c->GetLabel(2) != trackIndex)) continue;
1163           if(TMath::Abs(c->GetY() - y) > minDY) continue;
1164           //minDY = TMath::Abs(c->GetY() - y);
1165           minDY = c->GetY() - y;
1166           wYcorrect = c->GetY();
1167           wZcorrect = c->GetZ();
1168
1169           Double_t h01 = GetTiltFactor(c);
1170           wChi2 = t.GetPredictedChi2(c, h01);
1171         }
1172       }                    
1173
1174       // Now go for the real cluster search
1175
1176       if (timeBin) {
1177
1178         for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1179           AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1180           if (c->GetY() > y+road) break;
1181           if (c->IsUsed() > 0) continue;
1182           if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1183
1184           Double_t h01 = GetTiltFactor(c);
1185           chi2=t.GetPredictedChi2(c,h01);
1186           
1187           if (chi2 > maxChi2) continue;
1188           maxChi2=chi2;
1189           cl=c;
1190           index=timeBin.GetIndex(i);
1191
1192           //check is correct
1193           if((c->GetLabel(0) != trackIndex) &&
1194              (c->GetLabel(1) != trackIndex) &&
1195              (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1196         }               
1197         
1198         if(!cl) {
1199
1200           for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1201             AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1202             
1203             if (c->GetY() > y+road) break;
1204             if (c->IsUsed() > 0) continue;
1205             if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1206             
1207             Double_t h01 = GetTiltFactor(c);
1208             chi2=t.GetPredictedChi2(c,h01);
1209             
1210             if (chi2 > maxChi2) continue;
1211             maxChi2=chi2;
1212             cl=c;
1213             index=timeBin.GetIndex(i);
1214           }
1215         }        
1216         
1217         if (cl) {
1218           wYclosest = cl->GetY();
1219           wZclosest = cl->GetZ();
1220
1221           t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters()); 
1222           Double_t h01 = GetTiltFactor(cl);
1223           Int_t det = cl->GetDetector();    
1224           Int_t plane = fGeom->GetPlane(det);
1225
1226           if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1227           //if(!t.Update(cl,maxChi2,index,h01)) {
1228             if(!tryAgain--) return 0;
1229           }  
1230           else tryAgain=fMaxGap;
1231         }
1232         else {
1233           if (tryAgain==0) break; 
1234           tryAgain--;
1235           
1236           //if (minDY < 1000000 && isNewLayer) 
1237             //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() <<  "\t" << 
1238             //  road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
1239                                                                      
1240         }
1241
1242         isNewLayer = kFALSE;
1243
1244         /*
1245         if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1246           
1247           printf(" %f", wIndex);       //1
1248           printf(" %f", wTB);          //2
1249           printf(" %f", wYrt);         //3
1250           printf(" %f", wYclosest);    //4
1251           printf(" %f", wYcorrect);    //5
1252           printf(" %f", wYwindow);     //6
1253           printf(" %f", wZrt);         //7
1254           printf(" %f", wZclosest);    //8
1255           printf(" %f", wZcorrect);    //9
1256           printf(" %f", wZwindow);     //10
1257           printf(" %f", wPx);          //11
1258           printf(" %f", wPy);          //12
1259           printf(" %f", wPz);          //13
1260           printf(" %f", wSigmaC2*1000000);  //14
1261           printf(" %f", wSigmaTgl2*1000);   //15
1262           printf(" %f", wSigmaY2);     //16
1263           //      printf(" %f", wSigmaZ2);     //17
1264           printf(" %f", wChi2);     //17
1265           printf(" %f", wC);           //18
1266           printf("\n");
1267         } 
1268         */                        
1269       }
1270     }  
1271   }
1272   if (nr<outerTB) 
1273     t.SetStop(kTRUE);
1274   else
1275     t.SetStop(kFALSE);
1276   return expectedNumberOfClusters;
1277
1278
1279 }         
1280
1281 //---------------------------------------------------------------------------
1282 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1283 {
1284   // Starting from current position on track=t this function tries
1285   // to extrapolate the track up to timeBin=0 and to reuse already
1286   // assigned clusters. Returns the number of clusters
1287   // expected to be found in sensitive layers
1288   // get indices of assigned clusters for each layer
1289   // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1290
1291   Int_t iCluster[90];
1292   for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1293   for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1294     Int_t index = t.GetClusterIndex(i);
1295     AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1296     if (!cl) continue;
1297     Int_t detector=cl->GetDetector();
1298     Int_t localTimeBin=cl->GetLocalTimeBin();
1299     Int_t sector=fGeom->GetSector(detector);
1300     Int_t plane=fGeom->GetPlane(detector);
1301
1302     Int_t trackingSector = CookSectorIndex(sector);
1303
1304     Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1305     if(gtb < 0) continue; 
1306     Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1307     iCluster[layer] = index;
1308   }
1309   t.ResetClusters();
1310
1311   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
1312
1313   Double_t alpha=t.GetAlpha();
1314   alpha = TVector2::Phi_0_2pi(alpha);
1315
1316   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
1317   Double_t radLength, rho, x, dx, y, ymax, z;
1318
1319   Int_t expectedNumberOfClusters = 0;
1320   Bool_t lookForCluster;
1321
1322   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1323
1324  
1325   for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) { 
1326
1327     y = t.GetY(); z = t.GetZ();
1328
1329     // first propagate to the inner surface of the current time bin 
1330     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1331     x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1332     if(!t.PropagateTo(x,radLength,rho)) break;
1333     y = t.GetY();
1334     ymax = x*TMath::Tan(0.5*alpha);
1335     if (y > ymax) {
1336       s = (s+1) % ns;
1337       if (!t.Rotate(alpha)) break;
1338       if(!t.PropagateTo(x,radLength,rho)) break;
1339     } else if (y <-ymax) {
1340       s = (s-1+ns) % ns;                           
1341       if (!t.Rotate(-alpha)) break;   
1342       if(!t.PropagateTo(x,radLength,rho)) break;
1343     } 
1344
1345     y = t.GetY(); z = t.GetZ();
1346
1347     // now propagate to the middle plane of the next time bin 
1348     fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1349     x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1350     if(!t.PropagateTo(x,radLength,rho)) break;
1351     y = t.GetY();
1352     ymax = x*TMath::Tan(0.5*alpha);
1353     if (y > ymax) {
1354       s = (s+1) % ns;
1355       if (!t.Rotate(alpha)) break;
1356       if(!t.PropagateTo(x,radLength,rho)) break;
1357     } else if (y <-ymax) {
1358       s = (s-1+ns) % ns;                           
1359       if (!t.Rotate(-alpha)) break;   
1360       if(!t.PropagateTo(x,radLength,rho)) break;
1361     } 
1362
1363     if(lookForCluster) expectedNumberOfClusters++;       
1364
1365     // use assigned cluster
1366     if (!iCluster[nr-1]) continue;
1367     AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1368     Double_t h01 = GetTiltFactor(cl);
1369     Double_t chi2=t.GetPredictedChi2(cl, h01);
1370     t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters()); 
1371     t.Update(cl,chi2,iCluster[nr-1],h01);
1372   }
1373
1374   return expectedNumberOfClusters;
1375 }                
1376
1377 //___________________________________________________________________
1378
1379 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1380 {
1381   // Starting from current radial position of track <t> this function
1382   // extrapolates the track up to radial position <xToGo>. 
1383   // Returns 1 if track reaches the plane, and 0 otherwise 
1384
1385   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
1386
1387   Double_t alpha=t.GetAlpha();
1388
1389   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1390   if (alpha < 0.            ) alpha += 2.*TMath::Pi();
1391
1392   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
1393
1394   Bool_t lookForCluster;
1395   Double_t radLength, rho, x, dx, y, ymax, z;
1396
1397   x = t.GetX();
1398
1399   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1400
1401   Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1402
1403   for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) { 
1404
1405     y = t.GetY(); z = t.GetZ();
1406
1407     // first propagate to the outer surface of the current time bin 
1408     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1409     x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1410     if(!t.PropagateTo(x,radLength,rho)) return 0;
1411     y = t.GetY();
1412     ymax = x*TMath::Tan(0.5*alpha);
1413     if (y > ymax) {
1414       s = (s+1) % ns;
1415       if (!t.Rotate(alpha)) return 0;
1416     } else if (y <-ymax) {
1417       s = (s-1+ns) % ns;                           
1418       if (!t.Rotate(-alpha)) return 0;   
1419     } 
1420     if(!t.PropagateTo(x,radLength,rho)) return 0;
1421
1422     y = t.GetY(); z = t.GetZ();
1423
1424     // now propagate to the middle plane of the next time bin 
1425     fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1426     x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1427     if(!t.PropagateTo(x,radLength,rho)) return 0;
1428     y = t.GetY();
1429     ymax = x*TMath::Tan(0.5*alpha);
1430     if (y > ymax) {
1431       s = (s+1) % ns;
1432       if (!t.Rotate(alpha)) return 0;
1433     } else if (y <-ymax) {
1434       s = (s-1+ns) % ns;                           
1435       if (!t.Rotate(-alpha)) return 0;   
1436     } 
1437     if(!t.PropagateTo(x,radLength,rho)) return 0;
1438   }
1439   return 1;
1440 }         
1441
1442 //___________________________________________________________________
1443
1444 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1445 {
1446   // Starting from current radial position of track <t> this function
1447   // extrapolates the track up to radial position of the outermost
1448   // padrow of the TPC. 
1449   // Returns 1 if track reaches the TPC, and 0 otherwise 
1450
1451   //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
1452
1453   Double_t alpha=t.GetAlpha();
1454   alpha = TVector2::Phi_0_2pi(alpha);
1455
1456   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
1457
1458   Bool_t lookForCluster;
1459   Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1460
1461   x = t.GetX();
1462
1463   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1464   Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1465
1466   for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) { 
1467
1468     y = t.GetY(); 
1469     z = t.GetZ();
1470
1471     // first propagate to the outer surface of the current time bin 
1472     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1473     x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; 
1474     
1475     if(!t.PropagateTo(x,radLength,rho)) return 0;
1476     AdjustSector(&t);
1477     if(!t.PropagateTo(x,radLength,rho)) return 0;
1478
1479     y = t.GetY(); 
1480     z = t.GetZ();
1481
1482     // now propagate to the middle plane of the next time bin 
1483     fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1484     x = fTrSec[s]->GetLayer(nr-1)->GetX(); 
1485     
1486     if(!t.PropagateTo(x,radLength,rho)) return 0;
1487     AdjustSector(&t);
1488     if(!t.PropagateTo(x,radLength,rho)) return 0;
1489   } 
1490   return 1;
1491 }         
1492
1493 //_____________________________________________________________________________
1494 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1495 {
1496   // Fills clusters into TRD tracking_sectors 
1497   // Note that the numbering scheme for the TRD tracking_sectors 
1498   // differs from that of TRD sectors
1499
1500   if (ReadClusters(fClusters,cTree)) {
1501      Error("LoadClusters","Problem with reading the clusters !");
1502      return 1;
1503   }
1504   Int_t ncl=fClusters->GetEntriesFast();
1505   fNclusters=ncl;
1506   cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1507               
1508   UInt_t index;
1509   for (Int_t ichamber=0;ichamber<5;ichamber++)
1510     for (Int_t isector=0;isector<18;isector++){
1511       fHoles[ichamber][isector]=kTRUE;
1512     }
1513
1514
1515   while (ncl--) {
1516 //    printf("\r %d left  ",ncl); 
1517     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1518     Int_t detector=c->GetDetector();
1519     Int_t localTimeBin=c->GetLocalTimeBin();
1520     Int_t sector=fGeom->GetSector(detector);
1521     Int_t plane=fGeom->GetPlane(detector);
1522       
1523     Int_t trackingSector = CookSectorIndex(sector);
1524     if (c->GetLabel(0)>0){
1525       Int_t chamber = fGeom->GetChamber(detector);
1526       fHoles[chamber][trackingSector]=kFALSE;
1527     }
1528
1529     Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1530     if(gtb < 0) continue; 
1531     Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1532
1533     index=ncl;
1534     fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1535   }    
1536   //  printf("\r\n");
1537   //
1538   //
1539   /*
1540   for (Int_t isector=0;isector<18;isector++){
1541     for (Int_t ichamber=0;ichamber<5;ichamber++)      
1542       if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector)) 
1543         printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1544              fGeom->IsHole(0,ichamber,17-isector));
1545   }
1546   */
1547   return 0;
1548 }
1549
1550 //_____________________________________________________________________________
1551 void AliTRDtracker::UnloadClusters() 
1552
1553   //
1554   // Clears the arrays of clusters and tracks. Resets sectors and timebins 
1555   //
1556
1557   Int_t i, nentr;
1558
1559   nentr = fClusters->GetEntriesFast();
1560   for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1561   fNclusters = 0;
1562
1563   nentr = fSeeds->GetEntriesFast();
1564   for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1565
1566   nentr = fTracks->GetEntriesFast();
1567   for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1568
1569   Int_t nsec = AliTRDgeometry::kNsect;
1570
1571   for (i = 0; i < nsec; i++) {    
1572     for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1573       fTrSec[i]->GetLayer(pl)->Clear();
1574     }
1575   }
1576
1577 }
1578
1579 //__________________________________________________________________________
1580 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1581 {
1582   // Creates track seeds using clusters in timeBins=i1,i2
1583
1584   if(turn > 2) {
1585     cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1586     return;
1587   }
1588
1589   Double_t x[5], c[15];
1590   Int_t maxSec=AliTRDgeometry::kNsect;
1591   
1592   Double_t alpha=AliTRDgeometry::GetAlpha();
1593   Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1594   Double_t cs=cos(alpha), sn=sin(alpha);
1595   Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1596     
1597       
1598   Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1599   Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1600       
1601   Double_t x1 =fTrSec[0]->GetX(i1);
1602   Double_t xx2=fTrSec[0]->GetX(i2);
1603       
1604   for (Int_t ns=0; ns<maxSec; ns++) {
1605     
1606     Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1607     Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1608     Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1609     Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1610     Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1611     
1612     AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1613     
1614     for (Int_t is=0; is < r1; is++) {
1615       Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1616       
1617       for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1618         
1619         const AliTRDcluster *cl;
1620         Double_t x2,   y2,   z2;
1621         Double_t x3=0., y3=0.;   
1622         
1623         if (js<nl2) {
1624           if(turn != 2) continue;
1625           AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1626           cl=r2[js];
1627           y2=cl->GetY(); z2=cl->GetZ();
1628           
1629           x2= xx2*cs2+y2*sn2;
1630           y2=-xx2*sn2+y2*cs2;
1631         }
1632         else if (js<nl2+nl) {
1633           if(turn != 1) continue;
1634           AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1635           cl=r2[js-nl2];
1636           y2=cl->GetY(); z2=cl->GetZ();
1637           
1638           x2= xx2*cs+y2*sn;
1639           y2=-xx2*sn+y2*cs;
1640         }                                
1641         else if (js<nl2+nl+nm) {
1642           if(turn != 1) continue;
1643           AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1644           cl=r2[js-nl2-nl];
1645           x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1646         }
1647         else if (js<nl2+nl+nm+nu) {
1648           if(turn != 1) continue;
1649           AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1650           cl=r2[js-nl2-nl-nm];
1651           y2=cl->GetY(); z2=cl->GetZ();
1652           
1653           x2=xx2*cs-y2*sn;
1654           y2=xx2*sn+y2*cs;
1655         }              
1656         else {
1657           if(turn != 2) continue;
1658           AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1659           cl=r2[js-nl2-nl-nm-nu];
1660           y2=cl->GetY(); z2=cl->GetZ();
1661           
1662           x2=xx2*cs2-y2*sn2;
1663           y2=xx2*sn2+y2*cs2;
1664         }
1665         
1666         if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1667         
1668         Double_t zz=z1 - z1/x1*(x1-x2);
1669         
1670         if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1671         
1672         Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1673         if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1674         
1675         x[0]=y1;
1676         x[1]=z1;
1677         x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1678         
1679         if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;      
1680         
1681         x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1682         
1683         if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1684         
1685         x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1686         
1687         if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1688         
1689         Double_t a=asin(x[2]);
1690         Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1691         
1692         if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1693         
1694         Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1695         Double_t sy2=cl->GetSigmaY2(),     sz2=cl->GetSigmaZ2();
1696         Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;  
1697
1698         // Tilt changes
1699         Double_t h01 = GetTiltFactor(r1[is]);
1700         Double_t xuFactor = 100.;
1701         if(fNoTilt) { 
1702           h01 = 0;
1703           xuFactor = 1;
1704         }
1705
1706         sy1=sy1+sz1*h01*h01;
1707         Double_t syz=sz1*(-h01);
1708         // end of tilt changes
1709         
1710         Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1711         Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1712         Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1713         Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1714         Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1715         Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1716         Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1717         Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1718         Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1719         Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;    
1720
1721         
1722         c[0]=sy1;
1723         //        c[1]=0.;       c[2]=sz1;
1724         c[1]=syz;       c[2]=sz1*xuFactor;
1725         c[3]=f20*sy1;  c[4]=0.;       c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1726         c[6]=f30*sy1;  c[7]=f31*sz1;  c[8]=f30*sy1*f20+f32*sy2*f22;
1727                        c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1728         c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1729         c[13]=f30*sy1*f40+f32*sy2*f42;
1730         c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;      
1731         
1732         UInt_t index=r1.GetIndex(is);
1733         
1734         AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1735
1736         Int_t rc=FollowProlongation(*track, i2);     
1737         
1738         if ((rc < 1) ||
1739             (track->GetNumberOfClusters() < 
1740              (outer-inner)*fgkMinClustersInSeed)) delete track;
1741         else {
1742           fSeeds->AddLast(track); fNseeds++;
1743 //          cerr<<"\r found seed "<<fNseeds;
1744         }
1745       }
1746     }
1747   }
1748 }            
1749
1750 //_____________________________________________________________________________
1751 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1752 {
1753   //
1754   // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0) 
1755   // from the file. The names of the cluster tree and branches 
1756   // should match the ones used in AliTRDclusterizer::WriteClusters()
1757   //
1758   TObjArray *clusterArray = new TObjArray(400); 
1759   
1760   TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1761   if (!branch) {
1762     Error("ReadClusters","Can't get the branch !");
1763     return 1;
1764   }
1765   branch->SetAddress(&clusterArray); 
1766   
1767   Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1768   //  printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1769   
1770   // Loop through all entries in the tree
1771   Int_t nbytes = 0;
1772   AliTRDcluster *c = 0;
1773   //  printf("\n");
1774
1775   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {    
1776     
1777     // Import the tree
1778     nbytes += ClusterTree->GetEvent(iEntry);  
1779     
1780     // Get the number of points in the detector
1781     Int_t nCluster = clusterArray->GetEntriesFast();  
1782 //    printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1783     
1784     // Loop through all TRD digits
1785     for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) { 
1786       c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1787       AliTRDcluster *co = new AliTRDcluster(*c);
1788       co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1789       Int_t ltb = co->GetLocalTimeBin();
1790       if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1791       else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1792       array->AddLast(co);
1793       delete clusterArray->RemoveAt(iCluster); 
1794     }
1795   }
1796
1797   delete clusterArray;
1798
1799   return 0;
1800 }
1801
1802 //__________________________________________________________________
1803 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const 
1804 {
1805   //
1806   // This cooks a label. Mmmmh, smells good...
1807   //
1808
1809   Int_t label=123456789, index, i, j;
1810   Int_t ncl=pt->GetNumberOfClusters();
1811   const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1812
1813   Bool_t labelAdded;
1814
1815   //  Int_t s[kRange][2];
1816   Int_t **s = new Int_t* [kRange];
1817   for (i=0; i<kRange; i++) {
1818     s[i] = new Int_t[2];
1819   }
1820   for (i=0; i<kRange; i++) {
1821     s[i][0]=-1;
1822     s[i][1]=0;
1823   }
1824
1825   Int_t t0,t1,t2;
1826   for (i=0; i<ncl; i++) {
1827     index=pt->GetClusterIndex(i);
1828     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1829     t0=c->GetLabel(0);
1830     t1=c->GetLabel(1);
1831     t2=c->GetLabel(2);
1832   }
1833
1834   for (i=0; i<ncl; i++) {
1835     index=pt->GetClusterIndex(i);
1836     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1837     for (Int_t k=0; k<3; k++) { 
1838       label=c->GetLabel(k);
1839       labelAdded=kFALSE; j=0;
1840       if (label >= 0) {
1841         while ( (!labelAdded) && ( j < kRange ) ) {
1842           if (s[j][0]==label || s[j][1]==0) {
1843             s[j][0]=label; 
1844             s[j][1]=s[j][1]+1; 
1845             labelAdded=kTRUE;
1846           }
1847           j++;
1848         }
1849       }
1850     }
1851   }
1852
1853   Int_t max=0;
1854   label = -123456789;
1855
1856   for (i=0; i<kRange; i++) {
1857     if (s[i][1]>max) {
1858       max=s[i][1]; label=s[i][0];
1859     }
1860   }
1861
1862   for (i=0; i<kRange; i++) {
1863     delete []s[i];
1864   }        
1865
1866   delete []s;
1867
1868   if ((1.- Float_t(max)/ncl) > wrong) label=-label;   
1869
1870   pt->SetLabel(label); 
1871
1872 }
1873
1874
1875 //__________________________________________________________________
1876 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const 
1877 {
1878   //
1879   // Use clusters, but don't abuse them!
1880   //
1881
1882   Int_t ncl=t->GetNumberOfClusters();
1883   for (Int_t i=from; i<ncl; i++) {
1884     Int_t index = t->GetClusterIndex(i);
1885     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1886     c->Use();
1887   }
1888 }
1889
1890
1891 //_____________________________________________________________________
1892 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1893 {
1894   // Parametrised "expected" error of the cluster reconstruction in Y 
1895
1896   Double_t s = 0.08 * 0.08;    
1897   return s;
1898 }
1899
1900 //_____________________________________________________________________
1901 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1902 {
1903   // Parametrised "expected" error of the cluster reconstruction in Z 
1904
1905   Double_t s = 9 * 9 /12.;  
1906   return s;
1907 }                  
1908
1909 //_____________________________________________________________________
1910 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const 
1911 {
1912   //
1913   // Returns radial position which corresponds to time bin <localTB>
1914   // in tracking sector <sector> and plane <plane>
1915   //
1916
1917   Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB); 
1918   Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1919   return fTrSec[sector]->GetLayer(pl)->GetX();
1920
1921 }
1922
1923
1924 //_______________________________________________________
1925 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x, 
1926                Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1927
1928   //
1929   // AliTRDpropagationLayer constructor
1930   //
1931
1932   fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1933   fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1934
1935
1936   for(Int_t i=0; i < (Int_t) kZones; i++) {
1937     fZc[i]=0; fZmax[i] = 0;
1938   }
1939
1940   fYmax = 0;
1941
1942   if(fTimeBinIndex >= 0) { 
1943     fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1944     fIndex = new UInt_t[kMaxClusterPerTimeBin];
1945   }
1946
1947   for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
1948   fHole = kFALSE;
1949   fHoleZc = 0;
1950   fHoleZmax = 0;
1951   fHoleYc = 0;
1952   fHoleYmax = 0;
1953   fHoleRho = 0;
1954   fHoleX0 = 0;
1955
1956 }
1957
1958 //_______________________________________________________
1959 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1960           Double_t Zmax, Double_t Ymax, Double_t rho, 
1961           Double_t radLength, Double_t Yc, Double_t Zc) 
1962 {
1963   //
1964   // Sets hole in the layer 
1965   //
1966   fHole = kTRUE;
1967   fHoleZc = Zc;
1968   fHoleZmax = Zmax;
1969   fHoleYc = Yc;
1970   fHoleYmax = Ymax;
1971   fHoleRho = rho;
1972   fHoleX0 = radLength;
1973 }
1974   
1975
1976 //_______________________________________________________
1977 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1978 {
1979   //
1980   // AliTRDtrackingSector Constructor
1981   //
1982
1983   fGeom = geo;
1984   fPar = par;
1985   fGeomSector = gs;
1986   fTzeroShift = 0.13;
1987   fN = 0;
1988   //
1989   // get holes description from geometry
1990   Bool_t holes[AliTRDgeometry::kNcham];
1991   //printf("sector\t%d\t",gs);
1992   for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
1993     holes[icham] = fGeom->IsHole(0,icham,gs);
1994     //printf("%d",holes[icham]);
1995   } 
1996   //printf("\n");
1997   
1998   for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
1999
2000
2001   AliTRDpropagationLayer* ppl;
2002
2003   Double_t x, xin, xout, dx, rho, radLength;
2004   Int_t    steps;
2005
2006   // set time bins in the gas of the TPC
2007
2008   xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
2009   rho = 0.9e-3;  radLength = 28.94;
2010
2011   for(Int_t i=0; i<steps; i++) {
2012     x = xin + i*dx + dx/2;
2013     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2014     InsertLayer(ppl);
2015   }
2016
2017   // set time bins in the outer field cage vessel
2018
2019   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2020   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2021   InsertLayer(ppl);
2022
2023   dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2024   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2025   InsertLayer(ppl);
2026
2027   dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2028   steps = 5; dx = (xout - xin)/steps;
2029   for(Int_t i=0; i<steps; i++) {
2030     x = xin + i*dx + dx/2;
2031     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2032     InsertLayer(ppl);
2033   }
2034
2035   dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2036   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2037   InsertLayer(ppl);
2038
2039   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2040   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2041   InsertLayer(ppl);
2042
2043
2044   // set time bins in CO2
2045
2046   xin = xout; xout = 275.0; 
2047   steps = 50; dx = (xout - xin)/steps;
2048   rho = 1.977e-3;  radLength = 36.2;
2049   
2050   for(Int_t i=0; i<steps; i++) {
2051     x = xin + i*dx + dx/2;
2052     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2053     InsertLayer(ppl);
2054   }
2055
2056   // set time bins in the outer containment vessel
2057
2058   dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2059   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2060   InsertLayer(ppl);
2061
2062   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2063   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2064   InsertLayer(ppl);
2065
2066   dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2067   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2068   InsertLayer(ppl);
2069
2070   dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2071   steps = 10; dx = (xout - xin)/steps;
2072   for(Int_t i=0; i<steps; i++) {
2073     x = xin + i*dx + dx/2;
2074     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2075     InsertLayer(ppl);
2076   }
2077
2078   dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2079   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2080   InsertLayer(ppl);
2081
2082   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2083   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2084   InsertLayer(ppl);
2085   
2086   dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2087   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2088   InsertLayer(ppl);
2089
2090   Double_t xtrd = (Double_t) fGeom->Rmin();  
2091
2092   // add layers between TPC and TRD (Air temporarily)
2093   xin = xout; xout = xtrd;
2094   steps = 50; dx = (xout - xin)/steps;
2095   rho = 1.2e-3;  radLength = 36.66;
2096   
2097   for(Int_t i=0; i<steps; i++) {
2098     x = xin + i*dx + dx/2;
2099     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2100     InsertLayer(ppl);
2101   }
2102
2103
2104   //  Double_t alpha=AliTRDgeometry::GetAlpha();
2105
2106   // add layers for each of the planes
2107
2108   Double_t dxRo = (Double_t) fGeom->CroHght();    // Rohacell 
2109   Double_t dxSpace = (Double_t) fGeom->Cspace();  // Spacing between planes
2110   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
2111   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region  
2112   Double_t dxRad = (Double_t) fGeom->CraHght();   // Radiator
2113   Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo; 
2114   Double_t dxPlane = dxTEC + dxSpace; 
2115
2116   Int_t tb, tbIndex;
2117   const Int_t  kNchambers = AliTRDgeometry::Ncham();
2118   Double_t  ymax = 0;
2119   //, holeYmax = 0;
2120   Double_t ymaxsensitive=0;
2121   Double_t *zc = new Double_t[kNchambers];
2122   Double_t *zmax = new Double_t[kNchambers];
2123   Double_t *zmaxsensitive = new Double_t[kNchambers];  
2124   //  Double_t  holeZmax = 1000.;   // the whole sector is missing
2125
2126   for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2127     //
2128     // Radiator 
2129     xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2130     steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6; 
2131     for(Int_t i=0; i<steps; i++) {
2132       x = xin + i*dx + dx/2;
2133       ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);      
2134       InsertLayer(ppl);
2135     }
2136
2137     ymax          = fGeom->GetChamberWidth(plane)/2.;
2138     ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2139     
2140     for(Int_t ch = 0; ch < kNchambers; ch++) {
2141       zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2142       Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2143       Float_t row0 = fPar->GetRow0(plane,ch,0);
2144       Int_t nPads = fPar->GetRowMax(plane,ch,0);
2145       zmaxsensitive[ch] = Float_t(nPads)*pad/2.;      
2146       //      zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2147       zc[ch] = (pad * nPads)/2 + row0;
2148       //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2149
2150     }
2151
2152     dx = fPar->GetTimeBinSize(); 
2153     rho = 0.00295 * 0.85; radLength = 11.0;  
2154
2155     Double_t x0 = (Double_t) fPar->GetTime0(plane);
2156     Double_t xbottom = x0 - dxDrift;
2157     Double_t xtop = x0 + dxAmp;
2158     //
2159     // Amplification region
2160     steps = (Int_t) (dxAmp/dx);
2161
2162     for(tb = 0; tb < steps; tb++) {
2163       x = x0 + tb * dx + dx/2;
2164       tbIndex = CookTimeBinIndex(plane, -tb-1);
2165       ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2166       ppl->SetYmax(ymax,ymaxsensitive);
2167       ppl->SetZ(zc, zmax, zmaxsensitive);
2168       ppl->SetHoles(holes);
2169       InsertLayer(ppl);
2170     }
2171     tbIndex = CookTimeBinIndex(plane, -steps);
2172     x = (x + dx/2 + xtop)/2;
2173     dx = 2*(xtop-x);
2174     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2175     ppl->SetYmax(ymax,ymaxsensitive);
2176     ppl->SetZ(zc, zmax,zmaxsensitive);
2177     ppl->SetHoles(holes);
2178     InsertLayer(ppl);
2179
2180     // Drift region
2181     dx = fPar->GetTimeBinSize();
2182     steps = (Int_t) (dxDrift/dx);
2183
2184     for(tb = 0; tb < steps; tb++) {
2185       x = x0 - tb * dx - dx/2;
2186       tbIndex = CookTimeBinIndex(plane, tb);
2187
2188       ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2189       ppl->SetYmax(ymax,ymaxsensitive);
2190       ppl->SetZ(zc, zmax, zmaxsensitive);
2191       ppl->SetHoles(holes);
2192       InsertLayer(ppl);
2193     }
2194     tbIndex = CookTimeBinIndex(plane, steps);
2195     x = (x - dx/2 + xbottom)/2;
2196     dx = 2*(x-xbottom);
2197     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2198     ppl->SetYmax(ymax,ymaxsensitive);
2199     ppl->SetZ(zc, zmax, zmaxsensitive);
2200     ppl->SetHoles(holes);    
2201     InsertLayer(ppl);
2202
2203     // Pad Plane
2204     xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2205     ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2206     ppl->SetYmax(ymax,ymaxsensitive);
2207     ppl->SetZ(zc, zmax,zmax);
2208     ppl->SetHoles(holes);         
2209     InsertLayer(ppl);
2210
2211     // Rohacell
2212     xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2213     steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6; 
2214     for(Int_t i=0; i<steps; i++) {
2215       x = xin + i*dx + dx/2;
2216       ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2217       ppl->SetYmax(ymax,ymaxsensitive);
2218       ppl->SetZ(zc, zmax,zmax);
2219       ppl->SetHoles(holes);
2220       InsertLayer(ppl);
2221     }
2222
2223     // Space between the chambers, air
2224     xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2225     steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66; 
2226     for(Int_t i=0; i<steps; i++) {
2227       x = xin + i*dx + dx/2;
2228       ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2229       InsertLayer(ppl);
2230     }
2231   }    
2232
2233   // Space between the TRD and RICH
2234   Double_t xRICH = 500.;
2235   xin = xout; xout = xRICH;
2236   steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66; 
2237   for(Int_t i=0; i<steps; i++) {
2238     x = xin + i*dx + dx/2;
2239     ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2240     InsertLayer(ppl);
2241   }
2242
2243   MapTimeBinLayers();
2244   delete [] zc;
2245   delete [] zmax;
2246
2247 }
2248
2249 //______________________________________________________
2250
2251 Int_t  AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2252 {
2253   //
2254   // depending on the digitization parameters calculates "global"
2255   // time bin index for timebin <localTB> in plane <plane>
2256   //
2257
2258   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
2259   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region  
2260   Double_t dx = (Double_t) fPar->GetTimeBinSize();  
2261
2262   Int_t tbAmp = fPar->GetTimeBefore();
2263   Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2264   if(kTRUE) maxAmp = 0;   // intentional until we change parameter class 
2265   Int_t tbDrift = fPar->GetTimeMax();
2266   Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2267
2268   Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2269
2270   Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2271
2272   if((localTB < 0) && 
2273      (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2274   if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2275
2276   return gtb;
2277
2278
2279 }
2280
2281 //______________________________________________________
2282
2283 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers() 
2284 {
2285   //
2286   // For all sensitive time bins sets corresponding layer index
2287   // in the array fTimeBins 
2288   //
2289
2290   Int_t index;
2291
2292   for(Int_t i = 0; i < fN; i++) {
2293     index = fLayers[i]->GetTimeBinIndex();
2294     
2295     //    printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2296
2297     if(index < 0) continue;
2298     if(index >= (Int_t) kMaxTimeBinIndex) {
2299       printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2300       printf("    index %d exceeds allowed maximum of %d!\n",
2301              index, kMaxTimeBinIndex-1);
2302       continue;
2303     }
2304     fTimeBinIndex[index] = i;
2305   }
2306
2307   Double_t x1, dx1, x2, dx2, gap;
2308
2309   for(Int_t i = 0; i < fN-1; i++) {
2310     x1 = fLayers[i]->GetX();
2311     dx1 = fLayers[i]->GetdX();
2312     x2 = fLayers[i+1]->GetX();
2313     dx2 = fLayers[i+1]->GetdX();
2314     gap = (x2 - dx2/2) - (x1 + dx1/2);
2315     if(gap < -0.01) {
2316       printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2317       printf("             %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2318     }
2319     if(gap > 0.01) { 
2320       printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2321       printf("             (%f - %f) - (%f + %f) = %f\n", 
2322              x2, dx2/2, x1, dx1, gap);
2323     }
2324   }
2325 }
2326   
2327
2328 //______________________________________________________
2329
2330
2331 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2332 {
2333   // 
2334   // Returns the number of time bin which in radial position is closest to <x>
2335   //
2336
2337   if(x >= fLayers[fN-1]->GetX()) return fN-1; 
2338   if(x <= fLayers[0]->GetX()) return 0; 
2339
2340   Int_t b=0, e=fN-1, m=(b+e)/2;
2341   for (; b<e; m=(b+e)/2) {
2342     if (x > fLayers[m]->GetX()) b=m+1;
2343     else e=m;
2344   }
2345   if(TMath::Abs(x - fLayers[m]->GetX()) > 
2346      TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2347   else return m;
2348
2349 }
2350
2351 //______________________________________________________
2352
2353 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const 
2354 {
2355   // 
2356   // Returns number of the innermost SENSITIVE propagation layer
2357   //
2358
2359   return GetLayerNumber(0);
2360 }
2361
2362 //______________________________________________________
2363
2364 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const 
2365 {
2366   // 
2367   // Returns number of the outermost SENSITIVE time bin
2368   //
2369
2370   return GetLayerNumber(GetNumberOfTimeBins() - 1);
2371 }
2372
2373 //______________________________________________________
2374
2375 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const 
2376 {
2377   // 
2378   // Returns number of SENSITIVE time bins
2379   //
2380
2381   Int_t tb, layer;
2382   for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2383     layer = GetLayerNumber(tb);
2384     if(layer>=0) break;
2385   }
2386   return tb+1;
2387 }
2388
2389 //______________________________________________________
2390
2391 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2392
2393   //
2394   // Insert layer <pl> in fLayers array.
2395   // Layers are sorted according to X coordinate.
2396
2397   if ( fN == ((Int_t) kMaxLayersPerSector)) {
2398     printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2399     return;
2400   }
2401   if (fN==0) {fLayers[fN++] = pl; return;}
2402   Int_t i=Find(pl->GetX());
2403
2404   memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2405   fLayers[i]=pl; fN++;
2406
2407 }              
2408
2409 //______________________________________________________
2410
2411 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const 
2412 {
2413   //
2414   // Returns index of the propagation layer nearest to X 
2415   //
2416
2417   if (x <= fLayers[0]->GetX()) return 0;
2418   if (x > fLayers[fN-1]->GetX()) return fN;
2419   Int_t b=0, e=fN-1, m=(b+e)/2;
2420   for (; b<e; m=(b+e)/2) {
2421     if (x > fLayers[m]->GetX()) b=m+1;
2422     else e=m;
2423   }
2424   return m;
2425 }             
2426
2427 //______________________________________________________
2428 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2429 {
2430   //
2431   // set centers and the width of sectors
2432   for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2433     fZc[icham] = center[icham];  
2434     fZmax[icham] = w[icham];
2435     fZmaxSensitive[icham] = wsensitive[icham];
2436     //   printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2437   }  
2438 }
2439 //______________________________________________________
2440
2441 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2442 {
2443   //
2444   // set centers and the width of sectors
2445   fHole = kFALSE;
2446   for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2447     fIsHole[icham] = holes[icham]; 
2448     if (holes[icham]) fHole = kTRUE;
2449   }  
2450 }
2451
2452
2453
2454 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2455         Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength, 
2456         Bool_t &lookForCluster) const
2457 {
2458   //
2459   // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2460   // and sensitivity <lookForCluster> in point <y,z>  
2461   //
2462
2463   dx  = fdX;
2464   rho = fRho;
2465   radLength  = fX0;
2466   lookForCluster = kFALSE;
2467   //
2468   // check dead regions in sensitive volume 
2469   if(fTimeBinIndex >= 0) {
2470     //
2471     Int_t zone=-1;
2472     for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2473       if  (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){ 
2474         zone = ch;
2475         lookForCluster = !(fIsHole[zone]);
2476         if(TMath::Abs(y) > fYmaxSensitive){  
2477           lookForCluster = kFALSE;
2478         }
2479         if (fIsHole[zone]) {
2480           //if hole
2481           rho = 1.29e-3;
2482           radLength = 36.66;
2483         }
2484       }    
2485     }
2486     return;
2487   }
2488   //
2489   //
2490   // check hole
2491   if (fHole==kFALSE) return;
2492   //
2493   for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2494     if  (TMath::Abs(z - fZc[ch]) < fZmax[ch]){ 
2495       if (fIsHole[ch]) {
2496         //if hole
2497         rho = 1.29e-3;
2498         radLength = 36.66;
2499       }
2500     }
2501   }
2502   return;
2503 }
2504
2505 Int_t  AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2506 {
2507   //
2508   //
2509   if (fTimeBinIndex < 0) return -20;  //unknown 
2510   Int_t zone=-10;   // dead zone
2511   for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2512     if(TMath::Abs(z - fZc[ch]) < fZmax[ch]) 
2513       zone = ch;
2514   }
2515   return zone;
2516 }
2517
2518
2519 //______________________________________________________
2520
2521 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c, 
2522                                                           UInt_t index) {
2523
2524 // Insert cluster in cluster array.
2525 // Clusters are sorted according to Y coordinate.  
2526
2527   if(fTimeBinIndex < 0) { 
2528     printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2529     return;
2530   }
2531
2532   if (fN== (Int_t) kMaxClusterPerTimeBin) {
2533     printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n"); 
2534     return;
2535   }
2536   if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2537   Int_t i=Find(c->GetY());
2538   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2539   memmove(fIndex   +i+1 ,fIndex   +i,(fN-i)*sizeof(UInt_t)); 
2540   fIndex[i]=index; fClusters[i]=c; fN++;
2541 }  
2542
2543 //______________________________________________________
2544
2545 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2546
2547 // Returns index of the cluster nearest in Y    
2548
2549   if (y <= fClusters[0]->GetY()) return 0;
2550   if (y > fClusters[fN-1]->GetY()) return fN;
2551   Int_t b=0, e=fN-1, m=(b+e)/2;
2552   for (; b<e; m=(b+e)/2) {
2553     if (y > fClusters[m]->GetY()) b=m+1;
2554     else e=m;
2555   }
2556   return m;
2557 }    
2558
2559 //---------------------------------------------------------
2560
2561 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2562 //
2563 //  Returns correction factor for tilted pads geometry 
2564 //
2565
2566   Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2567   Int_t det = c->GetDetector();    
2568   Int_t plane = fGeom->GetPlane(det);
2569
2570   //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2571   if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
2572
2573   if(fNoTilt) h01 = 0;
2574   
2575   return h01;
2576 }
2577
2578
2579
2580
2581
2582