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