Introducing Riostream.h
[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 /*
17 $Log$
18 Revision 1.18  2002/10/14 14:57:44  hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
20
21 Revision 1.14.6.2  2002/07/24 10:09:31  alibrary
22 Updating VirtualMC
23
24 Revision 1.17  2002/06/13 12:09:58  hristov
25 Minor corrections
26
27 Revision 1.16  2002/06/12 09:54:36  cblume
28 Update of tracking code provided by Sergei
29
30 Revision 1.14  2001/11/14 10:50:46  cblume
31 Changes in digits IO. Add merging of summable digits
32
33 Revision 1.13  2001/05/30 12:17:47  hristov
34 Loop variables declared once
35
36 Revision 1.12  2001/05/28 17:07:58  hristov
37 Last minute changes; ExB correction in AliTRDclusterizerV1; taking into account of material in G10 TEC frames and material between TEC planes (C.Blume,S.Sedykh)
38
39 Revision 1.8  2000/12/20 13:00:44  cblume
40 Modifications for the HP-compiler
41
42 Revision 1.7  2000/12/08 16:07:02  cblume
43 Update of the tracking by Sergei
44
45 Revision 1.6  2000/11/30 17:38:08  cblume
46 Changes to get in line with new STEER and EVGEN
47
48 Revision 1.5  2000/11/14 14:40:27  cblume
49 Correction for the Sun compiler (kTRUE and kFALSE)
50
51 Revision 1.4  2000/11/10 14:57:52  cblume
52 Changes in the geometry constants for the DEC compiler
53
54 Revision 1.3  2000/10/15 23:40:01  cblume
55 Remove AliTRDconst
56
57 Revision 1.2  2000/10/06 16:49:46  cblume
58 Made Getters const
59
60 Revision 1.1.2.2  2000/10/04 16:34:58  cblume
61 Replace include files by forward declarations
62
63 Revision 1.1.2.1  2000/09/22 14:47:52  cblume
64 Add the tracking code
65
66 */   
67
68 #include <Riostream.h>
69
70 #include <TFile.h>
71 #include <TBranch.h>
72 #include <TTree.h>  
73 #include <TObjArray.h> 
74
75 #include "AliTRDgeometry.h"
76 #include "AliTRDparameter.h"
77 #include "AliTRDgeometryDetail.h"
78 #include "AliTRDcluster.h" 
79 #include "AliTRDtrack.h"
80 #include "../TPC/AliTPCtrack.h"
81
82 #include "AliTRDtracker.h"
83
84 ClassImp(AliTRDtracker) 
85
86   const  Float_t     AliTRDtracker::fSeedDepth          = 0.5; 
87   const  Float_t     AliTRDtracker::fSeedStep           = 0.10;   
88   const  Float_t     AliTRDtracker::fSeedGap            = 0.25;  
89
90   const  Float_t     AliTRDtracker::fMaxSeedDeltaZ12    = 40.;  
91   const  Float_t     AliTRDtracker::fMaxSeedDeltaZ      = 25.;  
92   const  Float_t     AliTRDtracker::fMaxSeedC           = 0.0052; 
93   const  Float_t     AliTRDtracker::fMaxSeedTan         = 1.2;  
94   const  Float_t     AliTRDtracker::fMaxSeedVertexZ     = 150.; 
95
96   const  Double_t    AliTRDtracker::fSeedErrorSY        = 0.2;
97   const  Double_t    AliTRDtracker::fSeedErrorSY3       = 2.5;
98   const  Double_t    AliTRDtracker::fSeedErrorSZ        = 0.1;
99
100   const  Float_t     AliTRDtracker::fMinClustersInSeed  = 0.7;  
101
102   const  Float_t     AliTRDtracker::fMinClustersInTrack = 0.5;  
103   const  Float_t     AliTRDtracker::fMinFractionOfFoundClusters = 0.8;  
104
105   const  Float_t     AliTRDtracker::fSkipDepth          = 0.05;
106   const  Float_t     AliTRDtracker::fLabelFraction      = 0.8;  
107   const  Float_t     AliTRDtracker::fWideRoad           = 20.;
108
109   const  Double_t    AliTRDtracker::fMaxChi2            = 24.; 
110
111
112 //____________________________________________________________________
113 AliTRDtracker::AliTRDtracker(const TFile *geomfile)
114 {
115   // 
116   //  Main constructor
117   //  
118
119   Float_t fTzero = 0;
120   Int_t   version = 2;
121   
122   fAddTRDseeds = kFALSE;
123
124   fGeom = NULL;
125   
126   TDirectory *savedir=gDirectory; 
127   TFile *in=(TFile*)geomfile;  
128   if (!in->IsOpen()) {
129     printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
130     printf("    DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
131   }
132   else {
133     in->cd();  
134     in->ls();
135     fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
136     fPar  = (AliTRDparameter*) in->Get("TRDparameter");
137     fGeom->Dump();
138   }
139
140   if(fGeom) {
141     //    fTzero = geo->GetT0();
142     fTzero = 0.;
143     version = fGeom->IsVersion();
144     printf("Found geometry version %d on file \n", version);
145   }
146   else { 
147     printf("AliTRDtracker::AliTRDtracker(): cann't find TRD geometry!\n");
148     printf("    DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
149     fGeom = new AliTRDgeometryDetail(); 
150     fPar = new AliTRDparameter();
151   }
152
153   savedir->cd();  
154
155
156   //  fGeom->SetT0(fTzero);
157
158   fEvent     = 0;
159
160   fNclusters = 0;
161   fClusters  = new TObjArray(2000); 
162   fNseeds    = 0;
163   fSeeds     = new TObjArray(2000);
164   fNtracks   = 0;
165   fTracks    = new TObjArray(1000);
166
167   for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
168     Int_t tr_s = CookSectorIndex(geom_s);
169     fTrSec[tr_s] = new AliTRDtrackingSector(fGeom, geom_s, fPar);
170   }
171
172   fSY2corr = 0.025;
173   fSZ2corr = 12.;      
174
175   // calculate max gap on track
176
177   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
178   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
179
180   Double_t dx = (Double_t) fPar->GetTimeBinSize();   
181   Int_t tbAmp = fPar->GetTimeBefore();
182   Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
183   Int_t tbDrift = fPar->GetTimeMax();
184   Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
185
186   tbDrift = TMath::Min(tbDrift,maxDrift);
187   tbAmp = TMath::Min(tbAmp,maxAmp);
188
189   fTimeBinsPerPlane = tbAmp + tbDrift;
190   fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fSkipDepth);
191
192   fVocal = kFALSE;
193
194 }   
195
196 //___________________________________________________________________
197 AliTRDtracker::~AliTRDtracker()
198 {
199   delete fClusters;
200   delete fTracks;
201   delete fSeeds;
202   delete fGeom;  
203   delete fPar;  
204
205   for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
206     delete fTrSec[geom_s];
207   }
208 }   
209
210 //_____________________________________________________________________
211 inline Double_t f1trd(Double_t x1,Double_t y1,
212                       Double_t x2,Double_t y2,
213                       Double_t x3,Double_t y3)
214 {
215   //
216   // Initial approximation of the track curvature
217   //
218   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
219   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
220                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
221   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
222                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
223
224   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
225
226   return -xr*yr/sqrt(xr*xr+yr*yr);
227 }          
228
229 //_____________________________________________________________________
230 inline Double_t f2trd(Double_t x1,Double_t y1,
231                       Double_t x2,Double_t y2,
232                       Double_t x3,Double_t y3)
233 {
234   //
235   // Initial approximation of the track curvature times X coordinate
236   // of the center of curvature
237   //
238
239   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
240   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
241                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
242   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
243                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
244
245   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
246
247   return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
248 }          
249
250 //_____________________________________________________________________
251 inline Double_t f3trd(Double_t x1,Double_t y1,
252                       Double_t x2,Double_t y2,
253                       Double_t z1,Double_t z2)
254 {
255   //
256   // Initial approximation of the tangent of the track dip angle
257   //
258
259   return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
260 }            
261
262 //___________________________________________________________________
263 Int_t AliTRDtracker::Clusters2Tracks(const TFile *inp, TFile *out)
264 {
265   //
266   // Finds tracks within the TRD. File <inp> is expected to contain seeds 
267   // at the outer part of the TRD. If <inp> is NULL, the seeds
268   // are found within the TRD if fAddTRDseeds is TRUE. 
269   // The tracks are propagated to the innermost time bin 
270   // of the TRD and stored in file <out>. 
271   //
272
273   LoadEvent();
274  
275   TDirectory *savedir=gDirectory;
276
277   char   tname[100];
278
279   if (!out->IsOpen()) {
280     cerr<<"AliTRDtracker::Clusters2Tracks(): output file is not open !\n";
281     return 1;
282   }    
283
284   sprintf(tname,"seedTRDtoTPC_%d",fEvent); 
285   TTree tpc_tree(tname,"Tree with seeds from TRD at outer TPC pad row");
286   AliTPCtrack *iotrack=0;
287   tpc_tree.Branch("tracks","AliTPCtrack",&iotrack,32000,0); 
288
289   sprintf(tname,"TreeT%d_TRD",fEvent);
290   TTree trd_tree(tname,"TRD tracks at inner TRD time bin");
291   AliTRDtrack *iotrack_trd=0;
292   trd_tree.Branch("tracks","AliTRDtrack",&iotrack_trd,32000,0);  
293
294   Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
295   Float_t foundMin = fMinClustersInTrack * timeBins; 
296
297   if (inp) {
298      TFile *in=(TFile*)inp;
299      if (!in->IsOpen()) {
300         cerr<<"AliTRDtracker::Clusters2Tracks(): file with seeds is not open !\n";
301         cerr<<" ... going for seeds finding inside the TRD\n";
302      }
303      else {
304        in->cd();
305        sprintf(tname,"TRDb_%d",fEvent);  
306        TTree *seedTree=(TTree*)in->Get(tname);  
307        if (!seedTree) {
308          cerr<<"AliTRDtracker::Clusters2Tracks(): ";
309          cerr<<"can't get a tree with track seeds !\n";
310          return 3;
311        }  
312        AliTRDtrack *seed=new AliTRDtrack;
313        seedTree->SetBranchAddress("tracks",&seed);
314
315        Int_t n=(Int_t)seedTree->GetEntries();
316        for (Int_t i=0; i<n; i++) {
317          seedTree->GetEvent(i);
318          seed->ResetCovariance(); 
319          fSeeds->AddLast(new AliTRDtrack(*seed));
320          fNseeds++;
321        }          
322        delete seed;
323      }
324   }
325
326   out->cd();
327
328   // find tracks from loaded seeds
329
330   Int_t nseed=fSeeds->GetEntriesFast();
331   Int_t i, found = 0;
332   Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
333
334   for (i=0; i<nseed; i++) {   
335     AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt; 
336     FollowProlongation(t, innerTB); 
337     if (t.GetNumberOfClusters() >= foundMin) {
338       UseClusters(&t);
339       CookLabel(pt, 1-fLabelFraction);
340       //      t.CookdEdx();
341     }
342     iotrack_trd = pt;
343     trd_tree.Fill();
344     cerr<<found++<<'\r';     
345
346     if(PropagateToTPC(t)) {
347       AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
348       iotrack = tpc;
349       tpc_tree.Fill();
350       delete tpc;
351     }  
352     delete fSeeds->RemoveAt(i);
353     fNseeds--;
354   }     
355
356   cout<<"Number of loaded seeds: "<<nseed<<endl;  
357   cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
358
359   // after tracks from loaded seeds are found and the corresponding 
360   // clusters are used, look for additional seeds from TRD
361
362   if(fAddTRDseeds) { 
363     // Find tracks for the seeds in the TRD
364     Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
365   
366     Int_t nSteps = (Int_t) (fSeedDepth / fSeedStep);
367     Int_t gap = (Int_t) (timeBins * fSeedGap);
368     Int_t step = (Int_t) (timeBins * fSeedStep);
369   
370     // make a first turn with tight cut on initial curvature
371     for(Int_t turn = 1; turn <= 2; turn++) {
372       if(turn == 2) {
373         nSteps = (Int_t) (fSeedDepth / (3*fSeedStep));
374         step = (Int_t) (timeBins * (3*fSeedStep));
375       }
376       for(Int_t i=0; i<nSteps; i++) {
377         Int_t outer=timeBins-1-i*step; 
378         Int_t inner=outer-gap;
379
380         nseed=fSeeds->GetEntriesFast();
381         printf("\n initial number of seeds %d\n", nseed); 
382       
383         MakeSeeds(inner, outer, turn);
384       
385         nseed=fSeeds->GetEntriesFast();
386         printf("\n number of seeds after MakeSeeds %d\n", nseed); 
387         
388       
389         for (Int_t i=0; i<nseed; i++) {   
390           AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt; 
391           FollowProlongation(t,innerTB); 
392           if (t.GetNumberOfClusters() >= foundMin) {
393             UseClusters(&t);
394             CookLabel(pt, 1-fLabelFraction);
395             t.CookdEdx();
396             cerr<<found++<<'\r';     
397             iotrack_trd = pt;
398             trd_tree.Fill();
399             if(PropagateToTPC(t)) {
400               AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
401               iotrack = tpc;
402               tpc_tree.Fill();
403               delete tpc;
404             }   
405           }
406           delete fSeeds->RemoveAt(i);
407           fNseeds--;
408         }
409       }
410     }
411   }
412   tpc_tree.Write(); 
413   trd_tree.Write(); 
414   
415   cout<<"Total number of found tracks: "<<found<<endl;
416     
417   UnloadEvent();
418     
419   savedir->cd();  
420   
421   return 0;    
422 }     
423      
424   
425
426 //_____________________________________________________________________________
427 Int_t AliTRDtracker::PropagateBack(const TFile *inp, TFile *out) {
428   //
429   // Reads seeds from file <inp>. The seeds are AliTPCtrack's found and
430   // backpropagated by the TPC tracker. Each seed is first propagated 
431   // to the TRD, and then its prolongation is searched in the TRD.
432   // If sufficiently long continuation of the track is found in the TRD
433   // the track is updated, otherwise it's stored as originaly defined 
434   // by the TPC tracker.   
435   //  
436
437
438   LoadEvent();
439
440   TDirectory *savedir=gDirectory;
441
442   TFile *in=(TFile*)inp;
443
444   if (!in->IsOpen()) {
445      cerr<<"AliTRDtracker::PropagateBack(): ";
446      cerr<<"file with back propagated TPC tracks is not open !\n";
447      return 1;
448   }                   
449
450   if (!out->IsOpen()) {
451      cerr<<"AliTRDtracker::PropagateBack(): ";
452      cerr<<"file for back propagated TRD tracks is not open !\n";
453      return 2;
454   }      
455
456   in->cd();
457   char   tname[100];
458   sprintf(tname,"seedsTPCtoTRD_%d",fEvent);       
459   TTree *seedTree=(TTree*)in->Get(tname);
460   if (!seedTree) {
461      cerr<<"AliTRDtracker::PropagateBack(): ";
462      cerr<<"can't get a tree with seeds from TPC !\n";
463      cerr<<"check if your version of TPC tracker creates tree "<<tname<<"\n";
464      return 3;
465   }
466
467   AliTPCtrack *seed=new AliTPCtrack;
468   seedTree->SetBranchAddress("tracks",&seed);
469
470   Int_t n=(Int_t)seedTree->GetEntries();
471   for (Int_t i=0; i<n; i++) {
472      seedTree->GetEvent(i);
473      Int_t lbl = seed->GetLabel();
474      AliTRDtrack *tr = new AliTRDtrack(*seed,seed->GetAlpha());
475      tr->SetSeedLabel(lbl);
476      fSeeds->AddLast(tr);
477      fNseeds++;
478   }
479
480   delete seed;
481   delete seedTree;
482
483   out->cd();
484
485   AliTPCtrack *otrack=0;
486
487   sprintf(tname,"seedsTRDtoTOF1_%d",fEvent);  
488   TTree tofTree1(tname,"Tracks back propagated through TPC and TRD");
489   tofTree1.Branch("tracks","AliTPCtrack",&otrack,32000,0);  
490
491   sprintf(tname,"seedsTRDtoTOF2_%d",fEvent);  
492   TTree tofTree2(tname,"Tracks back propagated through TPC and TRD");
493   tofTree2.Branch("tracks","AliTPCtrack",&otrack,32000,0);  
494
495   sprintf(tname,"seedsTRDtoPHOS_%d",fEvent);  
496   TTree phosTree(tname,"Tracks back propagated through TPC and TRD");
497   phosTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);  
498
499   sprintf(tname,"seedsTRDtoRICH_%d",fEvent);  
500   TTree richTree(tname,"Tracks back propagated through TPC and TRD");
501   richTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);  
502
503   sprintf(tname,"TRDb_%d",fEvent);  
504   TTree trdTree(tname,"Back propagated TRD tracks at outer TRD time bin");
505   AliTRDtrack *otrack_trd=0;
506   trdTree.Branch("tracks","AliTRDtrack",&otrack_trd,32000,0);   
507      
508   Int_t found=0;  
509
510   Int_t nseed=fSeeds->GetEntriesFast();
511
512   Float_t foundMin = fMinClustersInTrack * fTimeBinsPerPlane * fGeom->Nplan(); 
513
514   Int_t outermost_tb  = fTrSec[0]->GetOuterTimeBin();
515
516   for (Int_t i=0; i<nseed; i++) {  
517
518     AliTRDtrack *ps=(AliTRDtrack*)fSeeds->UncheckedAt(i), &s=*ps;
519     Int_t expectedClr = FollowBackProlongation(s);
520     Int_t foundClr = s.GetNumberOfClusters();
521     Int_t last_tb = fTrSec[0]->GetLayerNumber(s.GetX());
522
523     printf("seed %d: found %d out of %d expected clusters, Min is %f\n",
524            i, foundClr, expectedClr, foundMin);
525
526     if (foundClr >= foundMin) {
527       s.CookdEdx(); 
528       CookLabel(ps, 1-fLabelFraction);
529       UseClusters(ps);
530       otrack_trd=ps;
531       trdTree.Fill();
532       cout<<found++<<'\r';
533     }
534
535     if(((expectedClr < 10) && (last_tb == outermost_tb)) ||
536        ((expectedClr >= 10) && 
537         (((Float_t) foundClr) / ((Float_t) expectedClr) >= 
538          fMinFractionOfFoundClusters) && (last_tb == outermost_tb))) {
539
540       Double_t x_tof = 375.5;
541     
542       if(PropagateToOuterPlane(s,x_tof)) {
543         AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
544         otrack = pt;
545         tofTree1.Fill();
546         delete pt;
547
548         x_tof = 381.5;
549     
550         if(PropagateToOuterPlane(s,x_tof)) {
551           AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
552           otrack = pt;
553           tofTree2.Fill();
554           delete pt;
555
556           Double_t x_phos = 460.;
557           
558           if(PropagateToOuterPlane(s,x_phos)) {
559             AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
560             otrack = pt;
561             phosTree.Fill();
562             delete pt;
563             
564             Double_t x_rich = 490+1.267;
565             
566             if(PropagateToOuterPlane(s,x_rich)) {
567               AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
568               otrack = pt;
569               richTree.Fill();
570               delete pt;
571             }   
572           }
573         }
574       }      
575     }
576   }
577   
578   tofTree1.Write(); 
579   tofTree2.Write(); 
580   phosTree.Write(); 
581   richTree.Write(); 
582   trdTree.Write(); 
583
584   savedir->cd();  
585   cerr<<"Number of seeds: "<<nseed<<endl;  
586   cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
587
588   UnloadEvent();
589
590   return 0;
591
592 }
593
594
595 //---------------------------------------------------------------------------
596 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
597 {
598   // Starting from current position on track=t this function tries
599   // to extrapolate the track up to timeBin=0 and to confirm prolongation
600   // if a close cluster is found. Returns the number of clusters
601   // expected to be found in sensitive layers
602
603   Float_t  wIndex, wTB, wChi2;
604   Float_t  wYrt, wYclosest, wYcorrect, wYwindow;
605   Float_t  wZrt, wZclosest, wZcorrect, wZwindow;
606   Float_t  wPx, wPy, wPz, wC;
607   Double_t Px, Py, Pz;
608   Float_t  wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
609
610   Int_t trackIndex = t.GetLabel();  
611
612   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
613
614   Int_t try_again=fMaxGap;
615
616   Double_t alpha=t.GetAlpha();
617
618   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
619   if (alpha < 0.            ) alpha += 2.*TMath::Pi();
620
621   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
622
623   Double_t x0, rho, x, dx, y, ymax, z;
624
625   Int_t expectedNumberOfClusters = 0;
626   Bool_t lookForCluster;
627
628   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
629
630  
631   for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) { 
632
633     y = t.GetY(); z = t.GetZ();
634
635     // first propagate to the inner surface of the current time bin 
636     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
637     x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
638     if(!t.PropagateTo(x,x0,rho,0.139)) break;
639     y = t.GetY();
640     ymax = x*TMath::Tan(0.5*alpha);
641     if (y > ymax) {
642       s = (s+1) % ns;
643       if (!t.Rotate(alpha)) break;
644       if(!t.PropagateTo(x,x0,rho,0.139)) break;
645     } else if (y <-ymax) {
646       s = (s-1+ns) % ns;                           
647       if (!t.Rotate(-alpha)) break;   
648       if(!t.PropagateTo(x,x0,rho,0.139)) break;
649     } 
650
651     y = t.GetY(); z = t.GetZ();
652
653     // now propagate to the middle plane of the next time bin 
654     fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
655     x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
656     if(!t.PropagateTo(x,x0,rho,0.139)) break;
657     y = t.GetY();
658     ymax = x*TMath::Tan(0.5*alpha);
659     if (y > ymax) {
660       s = (s+1) % ns;
661       if (!t.Rotate(alpha)) break;
662       if(!t.PropagateTo(x,x0,rho,0.139)) break;
663     } else if (y <-ymax) {
664       s = (s-1+ns) % ns;                           
665       if (!t.Rotate(-alpha)) break;   
666       if(!t.PropagateTo(x,x0,rho,0.139)) break;
667     } 
668
669
670     if(lookForCluster) {
671
672
673       expectedNumberOfClusters++;       
674
675       wIndex = (Float_t) t.GetLabel();
676       wTB = nr;
677
678
679       AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr-1));
680
681
682       Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
683
684
685       Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
686
687       Double_t road;
688       if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
689       else return expectedNumberOfClusters;
690       
691       wYrt = (Float_t) y;
692       wZrt = (Float_t) z;
693       wYwindow = (Float_t) road;
694       t.GetPxPyPz(Px,Py,Pz);
695       wPx = (Float_t) Px;
696       wPy = (Float_t) Py;
697       wPz = (Float_t) Pz;
698       wC  = (Float_t) t.GetC();
699       wSigmaC2 = (Float_t) t.GetSigmaC2();
700       wSigmaTgl2    = (Float_t) t.GetSigmaTgl2();
701       wSigmaY2 = (Float_t) t.GetSigmaY2();
702       wSigmaZ2 = (Float_t) t.GetSigmaZ2();
703       wChi2 = -1;            
704       
705       if (road>fWideRoad) {
706         if (t.GetNumberOfClusters()>4)
707           cerr<<t.GetNumberOfClusters()
708               <<"FindProlongation warning: Too broad road !\n";
709         return 0;
710       }             
711
712       AliTRDcluster *cl=0;
713       UInt_t index=0;
714
715       Double_t max_chi2=fMaxChi2;
716
717       wYclosest = 12345678;
718       wYcorrect = 12345678;
719       wZclosest = 12345678;
720       wZcorrect = 12345678;
721       wZwindow  = TMath::Sqrt(2.25 * 12 * sz2);   
722
723       // Find the closest correct cluster for debugging purposes
724       if (time_bin) {
725         Float_t minDY = 1000000;
726         for (Int_t i=0; i<time_bin; i++) {
727           AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
728           if((c->GetLabel(0) != trackIndex) &&
729              (c->GetLabel(1) != trackIndex) &&
730              (c->GetLabel(2) != trackIndex)) continue;
731           if(TMath::Abs(c->GetY() - y) > minDY) continue;
732           minDY = TMath::Abs(c->GetY() - y);
733           wYcorrect = c->GetY();
734           wZcorrect = c->GetZ();
735           wChi2 = t.GetPredictedChi2(c);
736         }
737       }                    
738
739       // Now go for the real cluster search
740
741       if (time_bin) {
742
743         for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
744           AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
745           if (c->GetY() > y+road) break;
746           if (c->IsUsed() > 0) continue;
747           if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
748           Double_t chi2=t.GetPredictedChi2(c);
749           
750           if (chi2 > max_chi2) continue;
751           max_chi2=chi2;
752           cl=c;
753           index=time_bin.GetIndex(i);
754         }               
755
756         if(!cl) {
757
758           for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
759             AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
760             
761             if (c->GetY() > y+road) break;
762             if (c->IsUsed() > 0) continue;
763             if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
764             
765             Double_t chi2=t.GetPredictedChi2(c);
766             
767             if (chi2 > max_chi2) continue;
768             max_chi2=chi2;
769             cl=c;
770             index=time_bin.GetIndex(i);
771           }
772         }       
773         
774
775         if (cl) {
776           wYclosest = cl->GetY();
777           wZclosest = cl->GetZ();
778
779           t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters()); 
780           if(!t.Update(cl,max_chi2,index)) {
781             if(!try_again--) return 0;
782           }  
783           else try_again=fMaxGap;
784         }
785         else {
786           if (try_again==0) break; 
787           try_again--;
788         }
789
790         /*
791         if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
792           
793           printf(" %f", wIndex);       //1
794           printf(" %f", wTB);          //2
795           printf(" %f", wYrt);         //3
796           printf(" %f", wYclosest);    //4
797           printf(" %f", wYcorrect);    //5
798           printf(" %f", wYwindow);     //6
799           printf(" %f", wZrt);         //7
800           printf(" %f", wZclosest);    //8
801           printf(" %f", wZcorrect);    //9
802           printf(" %f", wZwindow);     //10
803           printf(" %f", wPx);          //11
804           printf(" %f", wPy);          //12
805           printf(" %f", wPz);          //13
806           printf(" %f", wSigmaC2*1000000);  //14
807           printf(" %f", wSigmaTgl2*1000);   //15
808           printf(" %f", wSigmaY2);     //16
809           //      printf(" %f", wSigmaZ2);     //17
810           printf(" %f", wChi2);     //17
811           printf(" %f", wC);           //18
812           printf("\n");
813         } 
814         */                        
815       }
816     }  
817   }
818   return expectedNumberOfClusters;
819   
820   
821 }                
822
823 //___________________________________________________________________
824
825 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
826 {
827   // Starting from current radial position of track <t> this function
828   // extrapolates the track up to outer timebin and in the sensitive
829   // layers confirms prolongation if a close cluster is found. 
830   // Returns the number of clusters expected to be found in sensitive layers
831
832   Float_t  wIndex, wTB, wChi2;
833   Float_t  wYrt, wYclosest, wYcorrect, wYwindow;
834   Float_t  wZrt, wZclosest, wZcorrect, wZwindow;
835   Float_t  wPx, wPy, wPz, wC;
836   Double_t Px, Py, Pz;
837   Float_t  wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
838
839   Int_t trackIndex = t.GetLabel();  
840
841   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
842
843   Int_t try_again=fMaxGap;
844
845   Double_t alpha=t.GetAlpha();
846
847   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
848   if (alpha < 0.            ) alpha += 2.*TMath::Pi();
849
850   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
851
852   Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
853
854   Double_t x0, rho, x, dx, y, ymax, z;
855
856   Bool_t lookForCluster;
857
858   Int_t expectedNumberOfClusters = 0;
859   x = t.GetX();
860
861   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
862
863
864   for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB; nr++) { 
865
866     y = t.GetY(); z = t.GetZ();
867
868     // first propagate to the outer surface of the current time bin 
869
870     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
871     x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
872
873     if(!t.PropagateTo(x,x0,rho,0.139)) break;
874
875     y = t.GetY();
876     ymax = x*TMath::Tan(0.5*alpha);
877
878     if (y > ymax) {
879       s = (s+1) % ns;
880       if (!t.Rotate(alpha)) break;
881       if(!t.PropagateTo(x,x0,rho,0.139)) break;
882     } else if (y <-ymax) {
883       s = (s-1+ns) % ns;                           
884       if (!t.Rotate(-alpha)) break;   
885       if(!t.PropagateTo(x,x0,rho,0.139)) break;
886     } 
887     y = t.GetY(); z = t.GetZ();
888
889     // now propagate to the middle plane of the next time bin 
890     fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
891
892     x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
893
894     if(!t.PropagateTo(x,x0,rho,0.139)) break;
895
896     y = t.GetY();
897
898     ymax = x*TMath::Tan(0.5*alpha);
899
900     if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
901
902     if (y > ymax) {
903       s = (s+1) % ns;
904       if (!t.Rotate(alpha)) break;
905       if(!t.PropagateTo(x,x0,rho,0.139)) break;
906     } else if (y <-ymax) {
907       s = (s-1+ns) % ns;              
908       if (!t.Rotate(-alpha)) break;   
909       if(!t.PropagateTo(x,x0,rho,0.139)) break;
910     } 
911
912     //    printf("label %d, pl %d, lookForCluster %d \n",
913     //     trackIndex, nr+1, lookForCluster);
914
915     if(lookForCluster) {
916       expectedNumberOfClusters++;       
917
918       wIndex = (Float_t) t.GetLabel();
919       wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
920
921       AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr+1));
922       Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
923       Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
924       if((t.GetSigmaY2() + sy2) < 0) break;
925       Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2); 
926       Double_t y=t.GetY(), z=t.GetZ();
927
928       wYrt = (Float_t) y;
929       wZrt = (Float_t) z;
930       wYwindow = (Float_t) road;
931       t.GetPxPyPz(Px,Py,Pz);
932       wPx = (Float_t) Px;
933       wPy = (Float_t) Py;
934       wPz = (Float_t) Pz;
935       wC  = (Float_t) t.GetC();
936       wSigmaC2 = (Float_t) t.GetSigmaC2();
937       wSigmaTgl2    = (Float_t) t.GetSigmaTgl2();
938       wSigmaY2 = (Float_t) t.GetSigmaY2();
939       wSigmaZ2 = (Float_t) t.GetSigmaZ2();
940       wChi2 = -1;            
941       
942       if (road>fWideRoad) {
943         if (t.GetNumberOfClusters()>4)
944           cerr<<t.GetNumberOfClusters()
945               <<"FindProlongation warning: Too broad road !\n";
946         return 0;
947       }      
948
949       AliTRDcluster *cl=0;
950       UInt_t index=0;
951
952       Double_t max_chi2=fMaxChi2;
953
954       wYclosest = 12345678;
955       wYcorrect = 12345678;
956       wZclosest = 12345678;
957       wZcorrect = 12345678;
958       wZwindow  = TMath::Sqrt(2.25 * 12 * sz2);   
959
960       // Find the closest correct cluster for debugging purposes
961       if (time_bin) {
962         Float_t minDY = 1000000;
963         for (Int_t i=0; i<time_bin; i++) {
964           AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
965           if((c->GetLabel(0) != trackIndex) &&
966              (c->GetLabel(1) != trackIndex) &&
967              (c->GetLabel(2) != trackIndex)) continue;
968           if(TMath::Abs(c->GetY() - y) > minDY) continue;
969           minDY = TMath::Abs(c->GetY() - y);
970           wYcorrect = c->GetY();
971           wZcorrect = c->GetZ();
972           wChi2 = t.GetPredictedChi2(c);
973         }
974       }                    
975
976       // Now go for the real cluster search
977
978       if (time_bin) {
979
980         for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
981           AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
982           if (c->GetY() > y+road) break;
983           if (c->IsUsed() > 0) continue;
984           if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
985           Double_t chi2=t.GetPredictedChi2(c);
986           
987           if (chi2 > max_chi2) continue;
988           max_chi2=chi2;
989           cl=c;
990           index=time_bin.GetIndex(i);
991         }               
992
993         if(!cl) {
994
995           for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
996             AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
997             
998             if (c->GetY() > y+road) break;
999             if (c->IsUsed() > 0) continue;
1000             if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1001             
1002             Double_t chi2=t.GetPredictedChi2(c);
1003             
1004             if (chi2 > max_chi2) continue;
1005             max_chi2=chi2;
1006             cl=c;
1007             index=time_bin.GetIndex(i);
1008           }
1009         }       
1010         
1011         if (cl) {
1012           wYclosest = cl->GetY();
1013           wZclosest = cl->GetZ();
1014
1015           t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters()); 
1016           if(!t.Update(cl,max_chi2,index)) {
1017             if(!try_again--) return 0;
1018           }  
1019           else try_again=fMaxGap;
1020         }
1021         else {
1022           if (try_again==0) break; 
1023           try_again--;
1024         }
1025
1026         /*
1027         if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1028           
1029           printf(" %f", wIndex);       //1
1030           printf(" %f", wTB);          //2
1031           printf(" %f", wYrt);         //3
1032           printf(" %f", wYclosest);    //4
1033           printf(" %f", wYcorrect);    //5
1034           printf(" %f", wYwindow);     //6
1035           printf(" %f", wZrt);         //7
1036           printf(" %f", wZclosest);    //8
1037           printf(" %f", wZcorrect);    //9
1038           printf(" %f", wZwindow);     //10
1039           printf(" %f", wPx);          //11
1040           printf(" %f", wPy);          //12
1041           printf(" %f", wPz);          //13
1042           printf(" %f", wSigmaC2*1000000);  //14
1043           printf(" %f", wSigmaTgl2*1000);   //15
1044           printf(" %f", wSigmaY2);     //16
1045           //      printf(" %f", wSigmaZ2);     //17
1046           printf(" %f", wChi2);     //17
1047           printf(" %f", wC);           //18
1048           printf("\n");
1049         } 
1050         */                        
1051       }
1052     }  
1053   }
1054   return expectedNumberOfClusters;
1055 }         
1056
1057 //___________________________________________________________________
1058
1059 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1060 {
1061   // Starting from current radial position of track <t> this function
1062   // extrapolates the track up to radial position <xToGo>. 
1063   // Returns 1 if track reaches the plane, and 0 otherwise 
1064
1065   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
1066
1067   Double_t alpha=t.GetAlpha();
1068
1069   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1070   if (alpha < 0.            ) alpha += 2.*TMath::Pi();
1071
1072   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
1073
1074   Bool_t lookForCluster;
1075   Double_t x0, rho, x, dx, y, ymax, z;
1076
1077   x = t.GetX();
1078
1079   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1080
1081   Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1082
1083   for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) { 
1084
1085     y = t.GetY(); z = t.GetZ();
1086
1087     // first propagate to the outer surface of the current time bin 
1088     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1089     x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1090     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1091     y = t.GetY();
1092     ymax = x*TMath::Tan(0.5*alpha);
1093     if (y > ymax) {
1094       s = (s+1) % ns;
1095       if (!t.Rotate(alpha)) return 0;
1096     } else if (y <-ymax) {
1097       s = (s-1+ns) % ns;                           
1098       if (!t.Rotate(-alpha)) return 0;   
1099     } 
1100     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1101
1102     y = t.GetY(); z = t.GetZ();
1103
1104     // now propagate to the middle plane of the next time bin 
1105     fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1106     x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1107     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1108     y = t.GetY();
1109     ymax = x*TMath::Tan(0.5*alpha);
1110     if (y > ymax) {
1111       s = (s+1) % ns;
1112       if (!t.Rotate(alpha)) return 0;
1113     } else if (y <-ymax) {
1114       s = (s-1+ns) % ns;                           
1115       if (!t.Rotate(-alpha)) return 0;   
1116     } 
1117     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1118   }
1119   return 1;
1120 }         
1121
1122 //___________________________________________________________________
1123
1124 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1125 {
1126   // Starting from current radial position of track <t> this function
1127   // extrapolates the track up to radial position of the outermost
1128   // padrow of the TPC. 
1129   // Returns 1 if track reaches the TPC, and 0 otherwise 
1130
1131   Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);     
1132
1133   Double_t alpha=t.GetAlpha();
1134
1135   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1136   if (alpha < 0.            ) alpha += 2.*TMath::Pi();
1137
1138   Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;  
1139
1140   Bool_t lookForCluster;
1141   Double_t x0, rho, x, dx, y, ymax, z;
1142
1143   x = t.GetX();
1144
1145   alpha=AliTRDgeometry::GetAlpha();  // note: change in meaning
1146
1147   Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1148
1149   for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) { 
1150
1151     y = t.GetY(); z = t.GetZ();
1152
1153     // first propagate to the outer surface of the current time bin 
1154     fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1155     x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1156     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1157     y = t.GetY();
1158     ymax = x*TMath::Tan(0.5*alpha);
1159     if (y > ymax) {
1160       s = (s+1) % ns;
1161       if (!t.Rotate(alpha)) return 0;
1162     } else if (y <-ymax) {
1163       s = (s-1+ns) % ns;                           
1164       if (!t.Rotate(-alpha)) return 0;   
1165     } 
1166     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1167
1168     y = t.GetY(); z = t.GetZ();
1169
1170     // now propagate to the middle plane of the next time bin 
1171     fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1172     x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1173     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1174     y = t.GetY();
1175     ymax = x*TMath::Tan(0.5*alpha);
1176     if (y > ymax) {
1177       s = (s+1) % ns;
1178       if (!t.Rotate(alpha)) return 0;
1179     } else if (y <-ymax) {
1180       s = (s-1+ns) % ns;                           
1181       if (!t.Rotate(-alpha)) return 0;   
1182     } 
1183     if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1184   }
1185   return 1;
1186 }         
1187
1188
1189 //_____________________________________________________________________________
1190 void AliTRDtracker::LoadEvent()
1191 {
1192   // Fills clusters into TRD tracking_sectors 
1193   // Note that the numbering scheme for the TRD tracking_sectors 
1194   // differs from that of TRD sectors
1195
1196   ReadClusters(fClusters);
1197   Int_t ncl=fClusters->GetEntriesFast();
1198   cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1199               
1200   UInt_t index;
1201   while (ncl--) {
1202     printf("\r %d left  ",ncl); 
1203     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1204     Int_t detector=c->GetDetector(), local_time_bin=c->GetLocalTimeBin();
1205     Int_t sector=fGeom->GetSector(detector);
1206     Int_t plane=fGeom->GetPlane(detector);
1207
1208     Int_t tracking_sector = CookSectorIndex(sector);
1209
1210     Int_t gtb = fTrSec[tracking_sector]->CookTimeBinIndex(plane,local_time_bin);
1211     Int_t layer = fTrSec[tracking_sector]->GetLayerNumber(gtb);
1212
1213     index=ncl;
1214     fTrSec[tracking_sector]->GetLayer(layer)->InsertCluster(c,index);
1215   }    
1216   printf("\r\n");
1217
1218 }
1219
1220 //_____________________________________________________________________________
1221 void AliTRDtracker::UnloadEvent() 
1222
1223   //
1224   // Clears the arrays of clusters and tracks. Resets sectors and timebins 
1225   //
1226
1227   Int_t i, nentr;
1228
1229   nentr = fClusters->GetEntriesFast();
1230   for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1231
1232   nentr = fSeeds->GetEntriesFast();
1233   for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1234
1235   nentr = fTracks->GetEntriesFast();
1236   for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1237
1238   Int_t nsec = AliTRDgeometry::kNsect;
1239
1240   for (i = 0; i < nsec; i++) {    
1241     for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1242       fTrSec[i]->GetLayer(pl)->Clear();
1243     }
1244   }
1245
1246 }
1247
1248 //__________________________________________________________________________
1249 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1250 {
1251   // Creates track seeds using clusters in timeBins=i1,i2
1252
1253   if(turn > 2) {
1254     cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1255     return;
1256   }
1257
1258   Double_t x[5], c[15];
1259   Int_t max_sec=AliTRDgeometry::kNsect;
1260   
1261   Double_t alpha=AliTRDgeometry::GetAlpha();
1262   Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1263   Double_t cs=cos(alpha), sn=sin(alpha);
1264   Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1265     
1266       
1267   Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1268   Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1269       
1270   Double_t x1 =fTrSec[0]->GetX(i1);
1271   Double_t xx2=fTrSec[0]->GetX(i2);
1272       
1273   for (Int_t ns=0; ns<max_sec; ns++) {
1274     
1275     Int_t nl2 = *(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1276     Int_t nl=(*fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1277     Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1278     Int_t nu=(*fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1279     Int_t nu2=(*fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1280     
1281     AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1282     
1283     for (Int_t is=0; is < r1; is++) {
1284       Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1285       
1286       for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1287         
1288         const AliTRDcluster *cl;
1289         Double_t x2,   y2,   z2;
1290         Double_t x3=0., y3=0.;   
1291         
1292         if (js<nl2) {
1293           if(turn != 2) continue;
1294           AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1295           cl=r2[js];
1296           y2=cl->GetY(); z2=cl->GetZ();
1297           
1298           x2= xx2*cs2+y2*sn2;
1299           y2=-xx2*sn2+y2*cs2;
1300         }
1301         else if (js<nl2+nl) {
1302           if(turn != 1) continue;
1303           AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1304           cl=r2[js-nl2];
1305           y2=cl->GetY(); z2=cl->GetZ();
1306           
1307           x2= xx2*cs+y2*sn;
1308           y2=-xx2*sn+y2*cs;
1309         }                                
1310         else if (js<nl2+nl+nm) {
1311           if(turn != 1) continue;
1312           AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1313           cl=r2[js-nl2-nl];
1314           x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1315         }
1316         else if (js<nl2+nl+nm+nu) {
1317           if(turn != 1) continue;
1318           AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1319           cl=r2[js-nl2-nl-nm];
1320           y2=cl->GetY(); z2=cl->GetZ();
1321           
1322           x2=xx2*cs-y2*sn;
1323           y2=xx2*sn+y2*cs;
1324         }              
1325         else {
1326           if(turn != 2) continue;
1327           AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1328           cl=r2[js-nl2-nl-nm-nu];
1329           y2=cl->GetY(); z2=cl->GetZ();
1330           
1331           x2=xx2*cs2-y2*sn2;
1332           y2=xx2*sn2+y2*cs2;
1333         }
1334         
1335         if(TMath::Abs(z1-z2) > fMaxSeedDeltaZ12) continue;
1336         
1337         Double_t zz=z1 - z1/x1*(x1-x2);
1338         
1339         if (TMath::Abs(zz-z2)>fMaxSeedDeltaZ) continue;
1340         
1341         Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1342         if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1343         
1344         x[0]=y1;
1345         x[1]=z1;
1346         x[2]=f1trd(x1,y1,x2,y2,x3,y3);
1347         
1348         if (TMath::Abs(x[2]) > fMaxSeedC) continue;      
1349         
1350         x[3]=f2trd(x1,y1,x2,y2,x3,y3);
1351         
1352         if (TMath::Abs(x[2]*x1-x[3]) >= 0.99999) continue;
1353         
1354         x[4]=f3trd(x1,y1,x2,y2,z1,z2);
1355         
1356         if (TMath::Abs(x[4]) > fMaxSeedTan) continue;
1357         
1358         Double_t a=asin(x[3]);
1359         Double_t zv=z1 - x[4]/x[2]*(a+asin(x[2]*x1-x[3]));
1360         
1361         if (TMath::Abs(zv)>fMaxSeedVertexZ) continue;
1362         
1363         Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1364         Double_t sy2=cl->GetSigmaY2(),     sz2=cl->GetSigmaZ2();
1365         Double_t sy3=fSeedErrorSY3, sy=fSeedErrorSY, sz=fSeedErrorSZ;  
1366         
1367         Double_t f20=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1368         Double_t f22=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1369         Double_t f24=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1370         Double_t f30=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[3])/sy;
1371         Double_t f32=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[3])/sy;
1372         Double_t f34=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[3])/sy;
1373         Double_t f40=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[4])/sy;
1374         Double_t f41=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[4])/sz;
1375         Double_t f42=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[4])/sy;
1376         Double_t f43=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[4])/sz;    
1377         
1378         c[0]=sy1;
1379         c[1]=0.;       c[2]=sz1;
1380         c[3]=f20*sy1;  c[4]=0.;   c[5]=f20*sy1*f20+f22*sy2*f22+f24*sy3*f24;
1381         c[6]=f30*sy1;  c[7]=0.;   c[8]=f30*sy1*f20+f32*sy2*f22+f34*sy3*f24;
1382         c[9]=f30*sy1*f30+f32*sy2*f32+f34*sy3*f34;
1383         c[10]=f40*sy1; c[11]=f41*sz1;   c[12]=f40*sy1*f20+f42*sy2*f22;
1384         c[13]=f40*sy1*f30+f42*sy2*f32;
1385         c[14]=f40*sy1*f40+f41*sz1*f41+f42*sy2*f42+f43*sz2*f43;
1386         
1387         UInt_t index=r1.GetIndex(is);
1388         
1389         AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1390
1391         Int_t rc=FollowProlongation(*track, i2);     
1392         
1393         if ((rc < 1) ||
1394             (track->GetNumberOfClusters() < 
1395              (outer-inner)*fMinClustersInSeed)) delete track;
1396         else {
1397           fSeeds->AddLast(track); fNseeds++;
1398           cerr<<"\r found seed "<<fNseeds;
1399         }
1400       }
1401     }
1402   }
1403 }            
1404
1405 //_____________________________________________________________________________
1406 void AliTRDtracker::ReadClusters(TObjArray *array, const TFile *inp) 
1407 {
1408   //
1409   // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0) 
1410   // from the file. The names of the cluster tree and branches 
1411   // should match the ones used in AliTRDclusterizer::WriteClusters()
1412   //
1413
1414   TDirectory *savedir=gDirectory; 
1415
1416   if (inp) {
1417      TFile *in=(TFile*)inp;
1418      if (!in->IsOpen()) {
1419         cerr<<"AliTRDtracker::ReadClusters(): input file is not open !\n";
1420         return;
1421      }
1422      else{
1423        in->cd();
1424      }
1425   }
1426
1427   Char_t treeName[12];
1428   sprintf(treeName,"TreeR%d_TRD",fEvent);
1429   TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1430   
1431   TObjArray *ClusterArray = new TObjArray(400); 
1432   
1433   ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray); 
1434   
1435   Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1436   printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1437   
1438   // Loop through all entries in the tree
1439   Int_t nbytes;
1440   AliTRDcluster *c = 0;
1441   printf("\n");
1442
1443   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {    
1444     
1445     // Import the tree
1446     nbytes += ClusterTree->GetEvent(iEntry);  
1447     
1448     // Get the number of points in the detector
1449     Int_t nCluster = ClusterArray->GetEntriesFast();  
1450     printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1451     
1452     // Loop through all TRD digits
1453     for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) { 
1454       c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1455       AliTRDcluster *co = new AliTRDcluster(*c);
1456       co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1457       Int_t ltb = co->GetLocalTimeBin();
1458       if(ltb != 0) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1459       
1460       array->AddLast(co);
1461       delete ClusterArray->RemoveAt(iCluster); 
1462     }
1463   }
1464
1465   delete ClusterArray;
1466   savedir->cd();   
1467
1468 }
1469
1470 //______________________________________________________________________
1471 void AliTRDtracker::ReadClusters(TObjArray *array, const Char_t *filename)
1472 {
1473   //
1474   // Reads AliTRDclusters from file <filename>. The names of the cluster
1475   // tree and branches should match the ones used in
1476   // AliTRDclusterizer::WriteClusters()
1477   // if <array> == 0, clusters are added into AliTRDtracker fCluster array
1478   //
1479
1480   TDirectory *savedir=gDirectory;
1481
1482   TFile *file = TFile::Open(filename);
1483   if (!file->IsOpen()) {
1484     cerr<<"Can't open file with TRD clusters"<<endl;
1485     return;
1486   }
1487
1488   Char_t treeName[12];
1489   sprintf(treeName,"TreeR%d_TRD",fEvent);
1490   TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1491
1492   if (!ClusterTree) {
1493      cerr<<"AliTRDtracker::ReadClusters(): ";
1494      cerr<<"can't get a tree with clusters !\n";
1495      return;
1496   }
1497
1498   TObjArray *ClusterArray = new TObjArray(400);
1499
1500   ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray);
1501
1502   Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1503   cout<<"found "<<nEntries<<" in ClusterTree"<<endl;   
1504
1505   // Loop through all entries in the tree
1506   Int_t nbytes;
1507   AliTRDcluster *c = 0;
1508
1509   printf("\n");
1510
1511   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1512
1513     // Import the tree
1514     nbytes += ClusterTree->GetEvent(iEntry);
1515
1516     // Get the number of points in the detector
1517     Int_t nCluster = ClusterArray->GetEntriesFast();
1518     printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1519
1520     // Loop through all TRD digits
1521     for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1522       c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1523       AliTRDcluster *co = new AliTRDcluster(*c);
1524       co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1525       Int_t ltb = co->GetLocalTimeBin();
1526       if(ltb != 0) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1527       array->AddLast(co);
1528       delete ClusterArray->RemoveAt(iCluster);
1529     }
1530   }
1531
1532   file->Close();
1533   delete ClusterArray;
1534   savedir->cd();
1535
1536 }                      
1537
1538
1539 //__________________________________________________________________
1540 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const {
1541
1542   Int_t label=123456789, index, i, j;
1543   Int_t ncl=pt->GetNumberOfClusters();
1544   const Int_t range = fTrSec[0]->GetOuterTimeBin()+1;
1545
1546   Bool_t label_added;
1547
1548   //  Int_t s[range][2];
1549   Int_t **s = new Int_t* [range];
1550   for (i=0; i<range; i++) {
1551     s[i] = new Int_t[2];
1552   }
1553   for (i=0; i<range; i++) {
1554     s[i][0]=-1;
1555     s[i][1]=0;
1556   }
1557
1558   Int_t t0,t1,t2;
1559   for (i=0; i<ncl; i++) {
1560     index=pt->GetClusterIndex(i);
1561     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1562     t0=c->GetLabel(0);
1563     t1=c->GetLabel(1);
1564     t2=c->GetLabel(2);
1565   }
1566
1567   for (i=0; i<ncl; i++) {
1568     index=pt->GetClusterIndex(i);
1569     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1570     for (Int_t k=0; k<3; k++) { 
1571       label=c->GetLabel(k);
1572       label_added=kFALSE; j=0;
1573       if (label >= 0) {
1574         while ( (!label_added) && ( j < range ) ) {
1575           if (s[j][0]==label || s[j][1]==0) {
1576             s[j][0]=label; 
1577             s[j][1]=s[j][1]+1; 
1578             label_added=kTRUE;
1579           }
1580           j++;
1581         }
1582       }
1583     }
1584   }
1585
1586   Int_t max=0;
1587   label = -123456789;
1588
1589   for (i=0; i<range; i++) {
1590     if (s[i][1]>max) {
1591       max=s[i][1]; label=s[i][0];
1592     }
1593   }
1594
1595   for (i=0; i<range; i++) {
1596     delete []s[i];
1597   }        
1598
1599   delete []s;
1600
1601   if ((1.- Float_t(max)/ncl) > wrong) label=-label;   
1602
1603   pt->SetLabel(label); 
1604
1605 }
1606
1607
1608 //__________________________________________________________________
1609 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const {
1610   Int_t ncl=t->GetNumberOfClusters();
1611   for (Int_t i=from; i<ncl; i++) {
1612     Int_t index = t->GetClusterIndex(i);
1613     AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1614     c->Use();
1615   }
1616 }
1617
1618
1619 //_____________________________________________________________________
1620 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t r, Double_t tgl, Double_t pt)
1621 {
1622   // Parametrised "expected" error of the cluster reconstruction in Y 
1623
1624   Double_t s = 0.08 * 0.08;    
1625   return s;
1626 }
1627
1628 //_____________________________________________________________________
1629 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t r, Double_t tgl)
1630 {
1631   // Parametrised "expected" error of the cluster reconstruction in Z 
1632
1633   Double_t s = 6 * 6 /12.;  
1634   return s;
1635 }                  
1636
1637
1638 //_____________________________________________________________________
1639 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t local_tb) const 
1640 {
1641   //
1642   // Returns radial position which corresponds to time bin <local_tb>
1643   // in tracking sector <sector> and plane <plane>
1644   //
1645
1646   Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, local_tb); 
1647   Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1648   return fTrSec[sector]->GetLayer(pl)->GetX();
1649
1650 }
1651
1652
1653 //_______________________________________________________
1654 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x, 
1655                Double_t dx, Double_t rho, Double_t x0, Int_t tb_index)
1656
1657   //
1658   // AliTRDpropagationLayer constructor
1659   //
1660
1661   fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = x0;
1662   fClusters = NULL; fIndex = NULL; fTimeBinIndex = tb_index;
1663
1664
1665   for(Int_t i=0; i < (Int_t) kZONES; i++) {
1666     fZc[i]=0; fZmax[i] = 0;
1667   }
1668
1669   fYmax = 0;
1670
1671   if(fTimeBinIndex >= 0) { 
1672     fClusters = new AliTRDcluster*[kMAX_CLUSTER_PER_TIME_BIN];
1673     fIndex = new UInt_t[kMAX_CLUSTER_PER_TIME_BIN];
1674   }
1675
1676   fHole = kFALSE;
1677   fHoleZc = 0;
1678   fHoleZmax = 0;
1679   fHoleYc = 0;
1680   fHoleYmax = 0;
1681   fHoleRho = 0;
1682   fHoleX0 = 0;
1683
1684 }
1685
1686 //_______________________________________________________
1687 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1688           Double_t Zmax, Double_t Ymax, Double_t rho, 
1689           Double_t x0, Double_t Yc, Double_t Zc) 
1690 {
1691   //
1692   // Sets hole in the layer 
1693   //
1694
1695   fHole = kTRUE;
1696   fHoleZc = Zc;
1697   fHoleZmax = Zmax;
1698   fHoleYc = Yc;
1699   fHoleYmax = Ymax;
1700   fHoleRho = rho;
1701   fHoleX0 = x0;
1702 }
1703   
1704
1705 //_______________________________________________________
1706 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1707 {
1708   //
1709   // AliTRDtrackingSector Constructor
1710   //
1711
1712   fGeom = geo;
1713   fPar = par;
1714   fGeomSector = gs;
1715   fTzeroShift = 0.13;
1716   fN = 0;
1717
1718   for(UInt_t i=0; i < kMAX_TIME_BIN_INDEX; i++) fTimeBinIndex[i] = -1;
1719
1720
1721   AliTRDpropagationLayer* ppl;
1722
1723   Double_t x, xin, xout, dx, rho, x0;
1724   Int_t    steps;
1725
1726   // set time bins in the gas of the TPC
1727
1728   xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
1729   rho = 0.9e-3;  x0 = 28.94;
1730
1731   for(Int_t i=0; i<steps; i++) {
1732     x = xin + i*dx + dx/2;
1733     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1734     InsertLayer(ppl);
1735   }
1736
1737   // set time bins in the outer field cage vessel
1738
1739   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1740   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1741   InsertLayer(ppl);
1742
1743   dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1744   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1745   InsertLayer(ppl);
1746
1747   dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; x0 = 41.28; // Nomex
1748   steps = 5; dx = (xout - xin)/steps;
1749   for(Int_t i=0; i<steps; i++) {
1750     x = xin + i*dx + dx/2;
1751     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1752     InsertLayer(ppl);
1753   }
1754
1755   dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1756   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1757   InsertLayer(ppl);
1758
1759   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1760   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1761   InsertLayer(ppl);
1762
1763
1764   // set time bins in CO2
1765
1766   xin = xout; xout = 275.0; 
1767   steps = 50; dx = (xout - xin)/steps;
1768   rho = 1.977e-3;  x0 = 36.2;
1769   
1770   for(Int_t i=0; i<steps; i++) {
1771     x = xin + i*dx + dx/2;
1772     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1773     InsertLayer(ppl);
1774   }
1775
1776   // set time bins in the outer containment vessel
1777
1778   dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; x0 = 24.01; // Al
1779   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1780   InsertLayer(ppl);
1781
1782   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1783   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1784   InsertLayer(ppl);
1785
1786   dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1787   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1788   InsertLayer(ppl);
1789
1790   dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; x0 = 41.28; // Nomex
1791   steps = 10; dx = (xout - xin)/steps;
1792   for(Int_t i=0; i<steps; i++) {
1793     x = xin + i*dx + dx/2;
1794     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1795     InsertLayer(ppl);
1796   }
1797
1798   dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1799   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1800   InsertLayer(ppl);
1801
1802   dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1803   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1804   InsertLayer(ppl);
1805   
1806   dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; x0 = 24.01; // Al
1807   ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1808   InsertLayer(ppl);
1809
1810   Double_t xtrd = (Double_t) fGeom->Rmin();  
1811
1812   // add layers between TPC and TRD (Air temporarily)
1813   xin = xout; xout = xtrd;
1814   steps = 50; dx = (xout - xin)/steps;
1815   rho = 1.2e-3;  x0 = 36.66;
1816   
1817   for(Int_t i=0; i<steps; i++) {
1818     x = xin + i*dx + dx/2;
1819     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1820     InsertLayer(ppl);
1821   }
1822
1823
1824   Double_t alpha=AliTRDgeometry::GetAlpha();
1825
1826   // add layers for each of the planes
1827
1828   Double_t dxRo = (Double_t) fGeom->CroHght();    // Rohacell 
1829   Double_t dxSpace = (Double_t) fGeom->Cspace();  // Spacing between planes
1830   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
1831   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region  
1832   Double_t dxRad = (Double_t) fGeom->CraHght();   // Radiator
1833   Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo; 
1834   Double_t dxPlane = dxTEC + dxSpace; 
1835
1836   Int_t tbBefore = (Int_t) (dxAmp/fPar->GetTimeBinSize());
1837   if(tbBefore > fPar->GetTimeBefore()) tbBefore = fPar->GetTimeBefore();
1838
1839   Int_t tb, tb_index;
1840   const Int_t  nChambers = AliTRDgeometry::Ncham();
1841   Double_t  Ymax = 0, holeYmax = 0;
1842   Double_t *  Zc  = new Double_t[nChambers];
1843   Double_t * Zmax = new Double_t[nChambers];
1844   Double_t  holeZmax = 1000.;   // the whole sector is missing
1845
1846
1847   for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
1848
1849     // Radiator 
1850     xin = xtrd + plane * dxPlane; xout = xin + dxRad;
1851     steps = 12; dx = (xout - xin)/steps; rho = 0.074; x0 = 40.6; 
1852     for(Int_t i=0; i<steps; i++) {
1853       x = xin + i*dx + dx/2;
1854       ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1855       if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1856         holeYmax = x*TMath::Tan(0.5*alpha);
1857         ppl->SetHole(holeYmax, holeZmax);
1858       }
1859       if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1860         holeYmax = x*TMath::Tan(0.5*alpha);
1861         ppl->SetHole(holeYmax, holeZmax);
1862       }
1863       InsertLayer(ppl);
1864     }
1865
1866     Ymax = fGeom->GetChamberWidth(plane)/2;
1867     for(Int_t ch = 0; ch < nChambers; ch++) {
1868       Zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
1869       Float_t pad = fPar->GetRowPadSize(plane,ch,0);
1870       Float_t row0 = fPar->GetRow0(plane,ch,0);
1871       Int_t nPads = fPar->GetRowMax(plane,ch,0);
1872       Zc[ch] = (pad * nPads)/2 + row0 - pad/2;
1873     }
1874
1875     dx = fPar->GetTimeBinSize(); 
1876     rho = 0.00295 * 0.85; x0 = 11.0;  
1877
1878     Double_t x0 = (Double_t) fPar->GetTime0(plane);
1879     Double_t xbottom = x0 - dxDrift;
1880     Double_t xtop = x0 + dxAmp;
1881
1882     // Amplification region
1883
1884     steps = (Int_t) (dxAmp/dx);
1885
1886     for(tb = 0; tb < steps; tb++) {
1887       x = x0 + tb * dx + dx/2;
1888       tb_index = CookTimeBinIndex(plane, -tb-1);
1889       ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1890       ppl->SetYmax(Ymax);
1891       for(Int_t ch = 0; ch < nChambers; ch++) {
1892         ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1893       }
1894       if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1895         holeYmax = x*TMath::Tan(0.5*alpha);
1896         ppl->SetHole(holeYmax, holeZmax);
1897       }
1898       if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1899         holeYmax = x*TMath::Tan(0.5*alpha);
1900         ppl->SetHole(holeYmax, holeZmax);
1901       }
1902       InsertLayer(ppl);
1903     }
1904     tb_index = CookTimeBinIndex(plane, -steps);
1905     x = (x + dx/2 + xtop)/2;
1906     dx = 2*(xtop-x);
1907     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1908     ppl->SetYmax(Ymax);
1909     for(Int_t ch = 0; ch < nChambers; ch++) {
1910       ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1911     }
1912     if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1913       holeYmax = x*TMath::Tan(0.5*alpha);
1914       ppl->SetHole(holeYmax, holeZmax);
1915     }
1916     if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1917       holeYmax = x*TMath::Tan(0.5*alpha);
1918       ppl->SetHole(holeYmax, holeZmax);
1919     }
1920     InsertLayer(ppl);
1921
1922     // Drift region
1923     dx = fPar->GetTimeBinSize();
1924     steps = (Int_t) (dxDrift/dx);
1925
1926     for(tb = 0; tb < steps; tb++) {
1927       x = x0 - tb * dx - dx/2;
1928       tb_index = CookTimeBinIndex(plane, tb);
1929
1930       ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1931       ppl->SetYmax(Ymax);
1932       for(Int_t ch = 0; ch < nChambers; ch++) {
1933         ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1934       }
1935       if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1936         holeYmax = x*TMath::Tan(0.5*alpha);
1937         ppl->SetHole(holeYmax, holeZmax);
1938       }
1939       if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1940         holeYmax = x*TMath::Tan(0.5*alpha);
1941         ppl->SetHole(holeYmax, holeZmax);
1942       }
1943       InsertLayer(ppl);
1944     }
1945     tb_index = CookTimeBinIndex(plane, steps);
1946     x = (x - dx/2 + xbottom)/2;
1947     dx = 2*(x-xbottom);
1948     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1949     ppl->SetYmax(Ymax);
1950     for(Int_t ch = 0; ch < nChambers; ch++) {
1951       ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1952     }
1953     if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1954       holeYmax = x*TMath::Tan(0.5*alpha);
1955       ppl->SetHole(holeYmax, holeZmax);
1956     }
1957     if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1958       holeYmax = x*TMath::Tan(0.5*alpha);
1959       ppl->SetHole(holeYmax, holeZmax);
1960     }
1961     InsertLayer(ppl);
1962
1963     // Pad Plane
1964     xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; x0 = 33.0;
1965     ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1966     if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1967       holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
1968       ppl->SetHole(holeYmax, holeZmax);
1969     }
1970     if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1971       holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
1972       ppl->SetHole(holeYmax, holeZmax);
1973     }
1974     InsertLayer(ppl);
1975
1976     // Rohacell
1977     xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
1978     steps = 5; dx = (xout - xin)/steps; rho = 0.074; x0 = 40.6; 
1979     for(Int_t i=0; i<steps; i++) {
1980       x = xin + i*dx + dx/2;
1981       ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1982       if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1983         holeYmax = x*TMath::Tan(0.5*alpha);
1984         ppl->SetHole(holeYmax, holeZmax);
1985       }
1986       if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1987         holeYmax = x*TMath::Tan(0.5*alpha);
1988         ppl->SetHole(holeYmax, holeZmax);
1989       }
1990       InsertLayer(ppl);
1991     }
1992
1993     // Space between the chambers, air
1994     xin = xout; xout = xtrd + (plane + 1) * dxPlane;
1995     steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; x0 = 36.66; 
1996     for(Int_t i=0; i<steps; i++) {
1997       x = xin + i*dx + dx/2;
1998       ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1999       if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
2000         holeYmax = x*TMath::Tan(0.5*alpha);
2001         ppl->SetHole(holeYmax, holeZmax);
2002       }
2003       if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2004         holeYmax = x*TMath::Tan(0.5*alpha);
2005         ppl->SetHole(holeYmax, holeZmax);
2006       }
2007       InsertLayer(ppl);
2008     }
2009   }    
2010
2011   // Space between the TRD and RICH
2012   Double_t xRICH = 500.;
2013   xin = xout; xout = xRICH;
2014   steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; x0 = 36.66; 
2015   for(Int_t i=0; i<steps; i++) {
2016     x = xin + i*dx + dx/2;
2017     ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
2018     InsertLayer(ppl);
2019   }
2020
2021   MapTimeBinLayers();
2022   delete [] Zc;
2023   delete [] Zmax;
2024
2025 }
2026
2027 //______________________________________________________
2028
2029 Int_t  AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t local_tb) const
2030 {
2031   //
2032   // depending on the digitization parameters calculates "global"
2033   // time bin index for timebin <local_tb> in plane <plane>
2034   //
2035
2036   Double_t dxAmp = (Double_t) fGeom->CamHght();   // Amplification region
2037   Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region  
2038   Double_t dx = (Double_t) fPar->GetTimeBinSize();  
2039
2040   Int_t tbAmp = fPar->GetTimeBefore();
2041   Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2042   Int_t tbDrift = fPar->GetTimeMax();
2043   Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2044
2045   Int_t tb_per_plane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2046
2047   Int_t gtb = (plane+1) * tb_per_plane - local_tb - 1;
2048
2049   if((local_tb < 0) && 
2050      (TMath::Abs(local_tb) > TMath::Min(tbAmp,maxAmp))) return -1;
2051   if(local_tb >= TMath::Min(tbDrift,maxDrift)) return -1;
2052
2053   return gtb;
2054
2055
2056 }
2057
2058 //______________________________________________________
2059
2060 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers() 
2061 {
2062   //
2063   // For all sensitive time bins sets corresponding layer index
2064   // in the array fTimeBins 
2065   //
2066
2067   Int_t index;
2068
2069   for(Int_t i = 0; i < fN; i++) {
2070     index = fLayers[i]->GetTimeBinIndex();
2071     
2072     //    printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2073
2074     if(index < 0) continue;
2075     if(index >= (Int_t) kMAX_TIME_BIN_INDEX) {
2076       printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2077       printf("    index %d exceeds allowed maximum of %d!\n",
2078              index, kMAX_TIME_BIN_INDEX-1);
2079       continue;
2080     }
2081     fTimeBinIndex[index] = i;
2082   }
2083
2084   Double_t x1, dx1, x2, dx2, gap;
2085
2086   for(Int_t i = 0; i < fN-1; i++) {
2087     x1 = fLayers[i]->GetX();
2088     dx1 = fLayers[i]->GetdX();
2089     x2 = fLayers[i+1]->GetX();
2090     dx2 = fLayers[i+1]->GetdX();
2091     gap = (x2 - dx2/2) - (x1 + dx1/2);
2092     if(gap < -0.01) {
2093       printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2094       printf("             %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2095     }
2096     if(gap > 0.01) { 
2097       printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2098       printf("             (%f - %f) - (%f + %f) = %f\n", 
2099              x2, dx2/2, x1, dx1, gap);
2100     }
2101   }
2102 }
2103   
2104
2105 //______________________________________________________
2106
2107
2108 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2109 {
2110   // 
2111   // Returns the number of time bin which in radial position is closest to <x>
2112   //
2113
2114   if(x >= fLayers[fN-1]->GetX()) return fN-1; 
2115   if(x <= fLayers[0]->GetX()) return 0; 
2116
2117   Int_t b=0, e=fN-1, m=(b+e)/2;
2118   for (; b<e; m=(b+e)/2) {
2119     if (x > fLayers[m]->GetX()) b=m+1;
2120     else e=m;
2121   }
2122   if(TMath::Abs(x - fLayers[m]->GetX()) > 
2123      TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2124   else return m;
2125
2126 }
2127
2128 //______________________________________________________
2129
2130 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const 
2131 {
2132   // 
2133   // Returns number of the innermost SENSITIVE propagation layer
2134   //
2135
2136   return GetLayerNumber(0);
2137 }
2138
2139 //______________________________________________________
2140
2141 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const 
2142 {
2143   // 
2144   // Returns number of the outermost SENSITIVE time bin
2145   //
2146
2147   return GetLayerNumber(GetNumberOfTimeBins() - 1);
2148 }
2149
2150 //______________________________________________________
2151
2152 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const 
2153 {
2154   // 
2155   // Returns number of SENSITIVE time bins
2156   //
2157
2158   Int_t tb, layer;
2159   for(tb = kMAX_TIME_BIN_INDEX-1; tb >=0; tb--) {
2160     layer = GetLayerNumber(tb);
2161     if(layer>=0) break;
2162   }
2163   return tb+1;
2164 }
2165
2166 //______________________________________________________
2167
2168 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2169
2170   //
2171   // Insert layer <pl> in fLayers array.
2172   // Layers are sorted according to X coordinate.
2173
2174   if ( fN == ((Int_t) kMAX_LAYERS_PER_SECTOR)) {
2175     printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2176     return;
2177   }
2178   if (fN==0) {fLayers[fN++] = pl; return;}
2179   Int_t i=Find(pl->GetX());
2180
2181   memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2182   fLayers[i]=pl; fN++;
2183
2184 }              
2185
2186 //______________________________________________________
2187
2188 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const 
2189 {
2190   //
2191   // Returns index of the propagation layer nearest to X 
2192   //
2193
2194   if (x <= fLayers[0]->GetX()) return 0;
2195   if (x > fLayers[fN-1]->GetX()) return fN;
2196   Int_t b=0, e=fN-1, m=(b+e)/2;
2197   for (; b<e; m=(b+e)/2) {
2198     if (x > fLayers[m]->GetX()) b=m+1;
2199     else e=m;
2200   }
2201   return m;
2202 }             
2203
2204 //______________________________________________________
2205
2206 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2207         Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &x0, 
2208         Bool_t &lookForCluster) const
2209 {
2210   //
2211   // Returns radial step <dx>, density <rho>, rad. length <x0>,
2212   // and sensitivity <lookForCluster> in point <y,z>  
2213   //
2214
2215   dx  = fdX;
2216   rho = fRho;
2217   x0  = fX0;
2218   lookForCluster = kFALSE;
2219
2220   // check dead regions
2221   if(fTimeBinIndex >= 0) {
2222     for(Int_t ch = 0; ch < (Int_t) kZONES; ch++) {
2223       if(TMath::Abs(z - fZc[ch]) < fZmax[ch]) 
2224         lookForCluster = kTRUE;
2225     }
2226     if(TMath::Abs(y) > fYmax) lookForCluster = kFALSE;
2227     if(!lookForCluster) { 
2228       //      rho = 1.7; x0 = 33.0; // G10 
2229     }
2230   }
2231
2232   // check hole
2233   if(fHole && (TMath::Abs(y - fHoleYc) < fHoleYmax) && 
2234               (TMath::Abs(z - fHoleZc) < fHoleZmax)) {
2235     lookForCluster = kFALSE;
2236     rho = fHoleRho;
2237     x0  = fHoleX0;
2238   }         
2239
2240   return;
2241 }
2242
2243 //______________________________________________________
2244
2245 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c, 
2246                                                           UInt_t index) {
2247
2248 // Insert cluster in cluster array.
2249 // Clusters are sorted according to Y coordinate.  
2250
2251   if(fTimeBinIndex < 0) { 
2252     printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2253     return;
2254   }
2255
2256   if (fN== (Int_t) kMAX_CLUSTER_PER_TIME_BIN) {
2257     printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n"); 
2258     return;
2259   }
2260   if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2261   Int_t i=Find(c->GetY());
2262   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2263   memmove(fIndex   +i+1 ,fIndex   +i,(fN-i)*sizeof(UInt_t)); 
2264   fIndex[i]=index; fClusters[i]=c; fN++;
2265 }  
2266
2267 //______________________________________________________
2268
2269 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2270
2271 // Returns index of the cluster nearest in Y    
2272
2273   if (y <= fClusters[0]->GetY()) return 0;
2274   if (y > fClusters[fN-1]->GetY()) return fN;
2275   Int_t b=0, e=fN-1, m=(b+e)/2;
2276   for (; b<e; m=(b+e)/2) {
2277     if (y > fClusters[m]->GetY()) b=m+1;
2278     else e=m;
2279   }
2280   return m;
2281 }    
2282
2283
2284
2285
2286
2287
2288
2289
2290