1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.14.6.2 2002/07/24 10:09:31 alibrary
21 Revision 1.17 2002/06/13 12:09:58 hristov
24 Revision 1.16 2002/06/12 09:54:36 cblume
25 Update of tracking code provided by Sergei
27 Revision 1.14 2001/11/14 10:50:46 cblume
28 Changes in digits IO. Add merging of summable digits
30 Revision 1.13 2001/05/30 12:17:47 hristov
31 Loop variables declared once
33 Revision 1.12 2001/05/28 17:07:58 hristov
34 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)
36 Revision 1.8 2000/12/20 13:00:44 cblume
37 Modifications for the HP-compiler
39 Revision 1.7 2000/12/08 16:07:02 cblume
40 Update of the tracking by Sergei
42 Revision 1.6 2000/11/30 17:38:08 cblume
43 Changes to get in line with new STEER and EVGEN
45 Revision 1.5 2000/11/14 14:40:27 cblume
46 Correction for the Sun compiler (kTRUE and kFALSE)
48 Revision 1.4 2000/11/10 14:57:52 cblume
49 Changes in the geometry constants for the DEC compiler
51 Revision 1.3 2000/10/15 23:40:01 cblume
54 Revision 1.2 2000/10/06 16:49:46 cblume
57 Revision 1.1.2.2 2000/10/04 16:34:58 cblume
58 Replace include files by forward declarations
60 Revision 1.1.2.1 2000/09/22 14:47:52 cblume
70 #include <TObjArray.h>
72 #include "AliTRDgeometry.h"
73 #include "AliTRDparameter.h"
74 #include "AliTRDgeometryDetail.h"
75 #include "AliTRDcluster.h"
76 #include "AliTRDtrack.h"
77 #include "../TPC/AliTPCtrack.h"
79 #include "AliTRDtracker.h"
81 ClassImp(AliTRDtracker)
83 const Float_t AliTRDtracker::fSeedDepth = 0.5;
84 const Float_t AliTRDtracker::fSeedStep = 0.10;
85 const Float_t AliTRDtracker::fSeedGap = 0.25;
87 const Float_t AliTRDtracker::fMaxSeedDeltaZ12 = 40.;
88 const Float_t AliTRDtracker::fMaxSeedDeltaZ = 25.;
89 const Float_t AliTRDtracker::fMaxSeedC = 0.0052;
90 const Float_t AliTRDtracker::fMaxSeedTan = 1.2;
91 const Float_t AliTRDtracker::fMaxSeedVertexZ = 150.;
93 const Double_t AliTRDtracker::fSeedErrorSY = 0.2;
94 const Double_t AliTRDtracker::fSeedErrorSY3 = 2.5;
95 const Double_t AliTRDtracker::fSeedErrorSZ = 0.1;
97 const Float_t AliTRDtracker::fMinClustersInSeed = 0.7;
99 const Float_t AliTRDtracker::fMinClustersInTrack = 0.5;
100 const Float_t AliTRDtracker::fMinFractionOfFoundClusters = 0.8;
102 const Float_t AliTRDtracker::fSkipDepth = 0.05;
103 const Float_t AliTRDtracker::fLabelFraction = 0.8;
104 const Float_t AliTRDtracker::fWideRoad = 20.;
106 const Double_t AliTRDtracker::fMaxChi2 = 24.;
109 //____________________________________________________________________
110 AliTRDtracker::AliTRDtracker(const TFile *geomfile)
119 fAddTRDseeds = kFALSE;
123 TDirectory *savedir=gDirectory;
124 TFile *in=(TFile*)geomfile;
126 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
127 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
132 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
133 fPar = (AliTRDparameter*) in->Get("TRDparameter");
138 // fTzero = geo->GetT0();
140 version = fGeom->IsVersion();
141 printf("Found geometry version %d on file \n", version);
144 printf("AliTRDtracker::AliTRDtracker(): cann't find TRD geometry!\n");
145 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
146 fGeom = new AliTRDgeometryDetail();
147 fPar = new AliTRDparameter();
153 // fGeom->SetT0(fTzero);
158 fClusters = new TObjArray(2000);
160 fSeeds = new TObjArray(2000);
162 fTracks = new TObjArray(1000);
164 for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
165 Int_t tr_s = CookSectorIndex(geom_s);
166 fTrSec[tr_s] = new AliTRDtrackingSector(fGeom, geom_s, fPar);
172 // calculate max gap on track
174 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
175 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
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 Int_t tbDrift = fPar->GetTimeMax();
181 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
183 tbDrift = TMath::Min(tbDrift,maxDrift);
184 tbAmp = TMath::Min(tbAmp,maxAmp);
186 fTimeBinsPerPlane = tbAmp + tbDrift;
187 fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fSkipDepth);
193 //___________________________________________________________________
194 AliTRDtracker::~AliTRDtracker()
202 for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
203 delete fTrSec[geom_s];
207 //_____________________________________________________________________
208 inline Double_t f1trd(Double_t x1,Double_t y1,
209 Double_t x2,Double_t y2,
210 Double_t x3,Double_t y3)
213 // Initial approximation of the track curvature
215 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
216 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
217 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
218 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
219 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
221 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
223 return -xr*yr/sqrt(xr*xr+yr*yr);
226 //_____________________________________________________________________
227 inline Double_t f2trd(Double_t x1,Double_t y1,
228 Double_t x2,Double_t y2,
229 Double_t x3,Double_t y3)
232 // Initial approximation of the track curvature times X coordinate
233 // of the center of curvature
236 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
237 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
238 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
239 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
240 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
242 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
244 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
247 //_____________________________________________________________________
248 inline Double_t f3trd(Double_t x1,Double_t y1,
249 Double_t x2,Double_t y2,
250 Double_t z1,Double_t z2)
253 // Initial approximation of the tangent of the track dip angle
256 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
259 //___________________________________________________________________
260 Int_t AliTRDtracker::Clusters2Tracks(const TFile *inp, TFile *out)
263 // Finds tracks within the TRD. File <inp> is expected to contain seeds
264 // at the outer part of the TRD. If <inp> is NULL, the seeds
265 // are found within the TRD if fAddTRDseeds is TRUE.
266 // The tracks are propagated to the innermost time bin
267 // of the TRD and stored in file <out>.
272 TDirectory *savedir=gDirectory;
276 if (!out->IsOpen()) {
277 cerr<<"AliTRDtracker::Clusters2Tracks(): output file is not open !\n";
281 sprintf(tname,"seedTRDtoTPC_%d",fEvent);
282 TTree tpc_tree(tname,"Tree with seeds from TRD at outer TPC pad row");
283 AliTPCtrack *iotrack=0;
284 tpc_tree.Branch("tracks","AliTPCtrack",&iotrack,32000,0);
286 sprintf(tname,"TreeT%d_TRD",fEvent);
287 TTree trd_tree(tname,"TRD tracks at inner TRD time bin");
288 AliTRDtrack *iotrack_trd=0;
289 trd_tree.Branch("tracks","AliTRDtrack",&iotrack_trd,32000,0);
291 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
292 Float_t foundMin = fMinClustersInTrack * timeBins;
295 TFile *in=(TFile*)inp;
297 cerr<<"AliTRDtracker::Clusters2Tracks(): file with seeds is not open !\n";
298 cerr<<" ... going for seeds finding inside the TRD\n";
302 sprintf(tname,"TRDb_%d",fEvent);
303 TTree *seedTree=(TTree*)in->Get(tname);
305 cerr<<"AliTRDtracker::Clusters2Tracks(): ";
306 cerr<<"can't get a tree with track seeds !\n";
309 AliTRDtrack *seed=new AliTRDtrack;
310 seedTree->SetBranchAddress("tracks",&seed);
312 Int_t n=(Int_t)seedTree->GetEntries();
313 for (Int_t i=0; i<n; i++) {
314 seedTree->GetEvent(i);
315 seed->ResetCovariance();
316 fSeeds->AddLast(new AliTRDtrack(*seed));
325 // find tracks from loaded seeds
327 Int_t nseed=fSeeds->GetEntriesFast();
329 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
331 for (i=0; i<nseed; i++) {
332 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
333 FollowProlongation(t, innerTB);
334 if (t.GetNumberOfClusters() >= foundMin) {
336 CookLabel(pt, 1-fLabelFraction);
343 if(PropagateToTPC(t)) {
344 AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
349 delete fSeeds->RemoveAt(i);
353 cout<<"Number of loaded seeds: "<<nseed<<endl;
354 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
356 // after tracks from loaded seeds are found and the corresponding
357 // clusters are used, look for additional seeds from TRD
360 // Find tracks for the seeds in the TRD
361 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
363 Int_t nSteps = (Int_t) (fSeedDepth / fSeedStep);
364 Int_t gap = (Int_t) (timeBins * fSeedGap);
365 Int_t step = (Int_t) (timeBins * fSeedStep);
367 // make a first turn with tight cut on initial curvature
368 for(Int_t turn = 1; turn <= 2; turn++) {
370 nSteps = (Int_t) (fSeedDepth / (3*fSeedStep));
371 step = (Int_t) (timeBins * (3*fSeedStep));
373 for(Int_t i=0; i<nSteps; i++) {
374 Int_t outer=timeBins-1-i*step;
375 Int_t inner=outer-gap;
377 nseed=fSeeds->GetEntriesFast();
378 printf("\n initial number of seeds %d\n", nseed);
380 MakeSeeds(inner, outer, turn);
382 nseed=fSeeds->GetEntriesFast();
383 printf("\n number of seeds after MakeSeeds %d\n", nseed);
386 for (Int_t i=0; i<nseed; i++) {
387 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
388 FollowProlongation(t,innerTB);
389 if (t.GetNumberOfClusters() >= foundMin) {
391 CookLabel(pt, 1-fLabelFraction);
396 if(PropagateToTPC(t)) {
397 AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
403 delete fSeeds->RemoveAt(i);
412 cout<<"Total number of found tracks: "<<found<<endl;
423 //_____________________________________________________________________________
424 Int_t AliTRDtracker::PropagateBack(const TFile *inp, TFile *out) {
426 // Reads seeds from file <inp>. The seeds are AliTPCtrack's found and
427 // backpropagated by the TPC tracker. Each seed is first propagated
428 // to the TRD, and then its prolongation is searched in the TRD.
429 // If sufficiently long continuation of the track is found in the TRD
430 // the track is updated, otherwise it's stored as originaly defined
431 // by the TPC tracker.
437 TDirectory *savedir=gDirectory;
439 TFile *in=(TFile*)inp;
442 cerr<<"AliTRDtracker::PropagateBack(): ";
443 cerr<<"file with back propagated TPC tracks is not open !\n";
447 if (!out->IsOpen()) {
448 cerr<<"AliTRDtracker::PropagateBack(): ";
449 cerr<<"file for back propagated TRD tracks is not open !\n";
455 sprintf(tname,"seedsTPCtoTRD_%d",fEvent);
456 TTree *seedTree=(TTree*)in->Get(tname);
458 cerr<<"AliTRDtracker::PropagateBack(): ";
459 cerr<<"can't get a tree with seeds from TPC !\n";
460 cerr<<"check if your version of TPC tracker creates tree "<<tname<<"\n";
464 AliTPCtrack *seed=new AliTPCtrack;
465 seedTree->SetBranchAddress("tracks",&seed);
467 Int_t n=(Int_t)seedTree->GetEntries();
468 for (Int_t i=0; i<n; i++) {
469 seedTree->GetEvent(i);
470 Int_t lbl = seed->GetLabel();
471 AliTRDtrack *tr = new AliTRDtrack(*seed,seed->GetAlpha());
472 tr->SetSeedLabel(lbl);
482 AliTPCtrack *otrack=0;
484 sprintf(tname,"seedsTRDtoTOF1_%d",fEvent);
485 TTree tofTree1(tname,"Tracks back propagated through TPC and TRD");
486 tofTree1.Branch("tracks","AliTPCtrack",&otrack,32000,0);
488 sprintf(tname,"seedsTRDtoTOF2_%d",fEvent);
489 TTree tofTree2(tname,"Tracks back propagated through TPC and TRD");
490 tofTree2.Branch("tracks","AliTPCtrack",&otrack,32000,0);
492 sprintf(tname,"seedsTRDtoPHOS_%d",fEvent);
493 TTree phosTree(tname,"Tracks back propagated through TPC and TRD");
494 phosTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
496 sprintf(tname,"seedsTRDtoRICH_%d",fEvent);
497 TTree richTree(tname,"Tracks back propagated through TPC and TRD");
498 richTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
500 sprintf(tname,"TRDb_%d",fEvent);
501 TTree trdTree(tname,"Back propagated TRD tracks at outer TRD time bin");
502 AliTRDtrack *otrack_trd=0;
503 trdTree.Branch("tracks","AliTRDtrack",&otrack_trd,32000,0);
507 Int_t nseed=fSeeds->GetEntriesFast();
509 Float_t foundMin = fMinClustersInTrack * fTimeBinsPerPlane * fGeom->Nplan();
511 Int_t outermost_tb = fTrSec[0]->GetOuterTimeBin();
513 for (Int_t i=0; i<nseed; i++) {
515 AliTRDtrack *ps=(AliTRDtrack*)fSeeds->UncheckedAt(i), &s=*ps;
516 Int_t expectedClr = FollowBackProlongation(s);
517 Int_t foundClr = s.GetNumberOfClusters();
518 Int_t last_tb = fTrSec[0]->GetLayerNumber(s.GetX());
520 printf("seed %d: found %d out of %d expected clusters, Min is %f\n",
521 i, foundClr, expectedClr, foundMin);
523 if (foundClr >= foundMin) {
525 CookLabel(ps, 1-fLabelFraction);
532 if(((expectedClr < 10) && (last_tb == outermost_tb)) ||
533 ((expectedClr >= 10) &&
534 (((Float_t) foundClr) / ((Float_t) expectedClr) >=
535 fMinFractionOfFoundClusters) && (last_tb == outermost_tb))) {
537 Double_t x_tof = 375.5;
539 if(PropagateToOuterPlane(s,x_tof)) {
540 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
547 if(PropagateToOuterPlane(s,x_tof)) {
548 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
553 Double_t x_phos = 460.;
555 if(PropagateToOuterPlane(s,x_phos)) {
556 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
561 Double_t x_rich = 490+1.267;
563 if(PropagateToOuterPlane(s,x_rich)) {
564 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
582 cerr<<"Number of seeds: "<<nseed<<endl;
583 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
592 //---------------------------------------------------------------------------
593 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
595 // Starting from current position on track=t this function tries
596 // to extrapolate the track up to timeBin=0 and to confirm prolongation
597 // if a close cluster is found. Returns the number of clusters
598 // expected to be found in sensitive layers
600 Float_t wIndex, wTB, wChi2;
601 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
602 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
603 Float_t wPx, wPy, wPz, wC;
605 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
607 Int_t trackIndex = t.GetLabel();
609 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
611 Int_t try_again=fMaxGap;
613 Double_t alpha=t.GetAlpha();
615 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
616 if (alpha < 0. ) alpha += 2.*TMath::Pi();
618 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
620 Double_t x0, rho, x, dx, y, ymax, z;
622 Int_t expectedNumberOfClusters = 0;
623 Bool_t lookForCluster;
625 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
628 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
630 y = t.GetY(); z = t.GetZ();
632 // first propagate to the inner surface of the current time bin
633 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
634 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
635 if(!t.PropagateTo(x,x0,rho,0.139)) break;
637 ymax = x*TMath::Tan(0.5*alpha);
640 if (!t.Rotate(alpha)) break;
641 if(!t.PropagateTo(x,x0,rho,0.139)) break;
642 } else if (y <-ymax) {
644 if (!t.Rotate(-alpha)) break;
645 if(!t.PropagateTo(x,x0,rho,0.139)) break;
648 y = t.GetY(); z = t.GetZ();
650 // now propagate to the middle plane of the next time bin
651 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
652 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
653 if(!t.PropagateTo(x,x0,rho,0.139)) break;
655 ymax = x*TMath::Tan(0.5*alpha);
658 if (!t.Rotate(alpha)) break;
659 if(!t.PropagateTo(x,x0,rho,0.139)) break;
660 } else if (y <-ymax) {
662 if (!t.Rotate(-alpha)) break;
663 if(!t.PropagateTo(x,x0,rho,0.139)) break;
670 expectedNumberOfClusters++;
672 wIndex = (Float_t) t.GetLabel();
676 AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr-1));
679 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
682 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
685 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
686 else return expectedNumberOfClusters;
690 wYwindow = (Float_t) road;
691 t.GetPxPyPz(Px,Py,Pz);
695 wC = (Float_t) t.GetC();
696 wSigmaC2 = (Float_t) t.GetSigmaC2();
697 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
698 wSigmaY2 = (Float_t) t.GetSigmaY2();
699 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
702 if (road>fWideRoad) {
703 if (t.GetNumberOfClusters()>4)
704 cerr<<t.GetNumberOfClusters()
705 <<"FindProlongation warning: Too broad road !\n";
712 Double_t max_chi2=fMaxChi2;
714 wYclosest = 12345678;
715 wYcorrect = 12345678;
716 wZclosest = 12345678;
717 wZcorrect = 12345678;
718 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
720 // Find the closest correct cluster for debugging purposes
722 Float_t minDY = 1000000;
723 for (Int_t i=0; i<time_bin; i++) {
724 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
725 if((c->GetLabel(0) != trackIndex) &&
726 (c->GetLabel(1) != trackIndex) &&
727 (c->GetLabel(2) != trackIndex)) continue;
728 if(TMath::Abs(c->GetY() - y) > minDY) continue;
729 minDY = TMath::Abs(c->GetY() - y);
730 wYcorrect = c->GetY();
731 wZcorrect = c->GetZ();
732 wChi2 = t.GetPredictedChi2(c);
736 // Now go for the real cluster search
740 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
741 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
742 if (c->GetY() > y+road) break;
743 if (c->IsUsed() > 0) continue;
744 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
745 Double_t chi2=t.GetPredictedChi2(c);
747 if (chi2 > max_chi2) continue;
750 index=time_bin.GetIndex(i);
755 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
756 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
758 if (c->GetY() > y+road) break;
759 if (c->IsUsed() > 0) continue;
760 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
762 Double_t chi2=t.GetPredictedChi2(c);
764 if (chi2 > max_chi2) continue;
767 index=time_bin.GetIndex(i);
773 wYclosest = cl->GetY();
774 wZclosest = cl->GetZ();
776 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
777 if(!t.Update(cl,max_chi2,index)) {
778 if(!try_again--) return 0;
780 else try_again=fMaxGap;
783 if (try_again==0) break;
788 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
790 printf(" %f", wIndex); //1
791 printf(" %f", wTB); //2
792 printf(" %f", wYrt); //3
793 printf(" %f", wYclosest); //4
794 printf(" %f", wYcorrect); //5
795 printf(" %f", wYwindow); //6
796 printf(" %f", wZrt); //7
797 printf(" %f", wZclosest); //8
798 printf(" %f", wZcorrect); //9
799 printf(" %f", wZwindow); //10
800 printf(" %f", wPx); //11
801 printf(" %f", wPy); //12
802 printf(" %f", wPz); //13
803 printf(" %f", wSigmaC2*1000000); //14
804 printf(" %f", wSigmaTgl2*1000); //15
805 printf(" %f", wSigmaY2); //16
806 // printf(" %f", wSigmaZ2); //17
807 printf(" %f", wChi2); //17
808 printf(" %f", wC); //18
815 return expectedNumberOfClusters;
820 //___________________________________________________________________
822 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
824 // Starting from current radial position of track <t> this function
825 // extrapolates the track up to outer timebin and in the sensitive
826 // layers confirms prolongation if a close cluster is found.
827 // Returns the number of clusters expected to be found in sensitive layers
829 Float_t wIndex, wTB, wChi2;
830 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
831 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
832 Float_t wPx, wPy, wPz, wC;
834 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
836 Int_t trackIndex = t.GetLabel();
838 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
840 Int_t try_again=fMaxGap;
842 Double_t alpha=t.GetAlpha();
844 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
845 if (alpha < 0. ) alpha += 2.*TMath::Pi();
847 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
849 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
851 Double_t x0, rho, x, dx, y, ymax, z;
853 Bool_t lookForCluster;
855 Int_t expectedNumberOfClusters = 0;
858 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
861 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB; nr++) {
863 y = t.GetY(); z = t.GetZ();
865 // first propagate to the outer surface of the current time bin
867 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
868 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
870 if(!t.PropagateTo(x,x0,rho,0.139)) break;
873 ymax = x*TMath::Tan(0.5*alpha);
877 if (!t.Rotate(alpha)) break;
878 if(!t.PropagateTo(x,x0,rho,0.139)) break;
879 } else if (y <-ymax) {
881 if (!t.Rotate(-alpha)) break;
882 if(!t.PropagateTo(x,x0,rho,0.139)) break;
884 y = t.GetY(); z = t.GetZ();
886 // now propagate to the middle plane of the next time bin
887 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
889 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
891 if(!t.PropagateTo(x,x0,rho,0.139)) break;
895 ymax = x*TMath::Tan(0.5*alpha);
897 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
901 if (!t.Rotate(alpha)) break;
902 if(!t.PropagateTo(x,x0,rho,0.139)) break;
903 } else if (y <-ymax) {
905 if (!t.Rotate(-alpha)) break;
906 if(!t.PropagateTo(x,x0,rho,0.139)) break;
909 // printf("label %d, pl %d, lookForCluster %d \n",
910 // trackIndex, nr+1, lookForCluster);
913 expectedNumberOfClusters++;
915 wIndex = (Float_t) t.GetLabel();
916 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
918 AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr+1));
919 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
920 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
921 if((t.GetSigmaY2() + sy2) < 0) break;
922 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
923 Double_t y=t.GetY(), z=t.GetZ();
927 wYwindow = (Float_t) road;
928 t.GetPxPyPz(Px,Py,Pz);
932 wC = (Float_t) t.GetC();
933 wSigmaC2 = (Float_t) t.GetSigmaC2();
934 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
935 wSigmaY2 = (Float_t) t.GetSigmaY2();
936 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
939 if (road>fWideRoad) {
940 if (t.GetNumberOfClusters()>4)
941 cerr<<t.GetNumberOfClusters()
942 <<"FindProlongation warning: Too broad road !\n";
949 Double_t max_chi2=fMaxChi2;
951 wYclosest = 12345678;
952 wYcorrect = 12345678;
953 wZclosest = 12345678;
954 wZcorrect = 12345678;
955 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
957 // Find the closest correct cluster for debugging purposes
959 Float_t minDY = 1000000;
960 for (Int_t i=0; i<time_bin; i++) {
961 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
962 if((c->GetLabel(0) != trackIndex) &&
963 (c->GetLabel(1) != trackIndex) &&
964 (c->GetLabel(2) != trackIndex)) continue;
965 if(TMath::Abs(c->GetY() - y) > minDY) continue;
966 minDY = TMath::Abs(c->GetY() - y);
967 wYcorrect = c->GetY();
968 wZcorrect = c->GetZ();
969 wChi2 = t.GetPredictedChi2(c);
973 // Now go for the real cluster search
977 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
978 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
979 if (c->GetY() > y+road) break;
980 if (c->IsUsed() > 0) continue;
981 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
982 Double_t chi2=t.GetPredictedChi2(c);
984 if (chi2 > max_chi2) continue;
987 index=time_bin.GetIndex(i);
992 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
993 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
995 if (c->GetY() > y+road) break;
996 if (c->IsUsed() > 0) continue;
997 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
999 Double_t chi2=t.GetPredictedChi2(c);
1001 if (chi2 > max_chi2) continue;
1004 index=time_bin.GetIndex(i);
1009 wYclosest = cl->GetY();
1010 wZclosest = cl->GetZ();
1012 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1013 if(!t.Update(cl,max_chi2,index)) {
1014 if(!try_again--) return 0;
1016 else try_again=fMaxGap;
1019 if (try_again==0) break;
1024 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1026 printf(" %f", wIndex); //1
1027 printf(" %f", wTB); //2
1028 printf(" %f", wYrt); //3
1029 printf(" %f", wYclosest); //4
1030 printf(" %f", wYcorrect); //5
1031 printf(" %f", wYwindow); //6
1032 printf(" %f", wZrt); //7
1033 printf(" %f", wZclosest); //8
1034 printf(" %f", wZcorrect); //9
1035 printf(" %f", wZwindow); //10
1036 printf(" %f", wPx); //11
1037 printf(" %f", wPy); //12
1038 printf(" %f", wPz); //13
1039 printf(" %f", wSigmaC2*1000000); //14
1040 printf(" %f", wSigmaTgl2*1000); //15
1041 printf(" %f", wSigmaY2); //16
1042 // printf(" %f", wSigmaZ2); //17
1043 printf(" %f", wChi2); //17
1044 printf(" %f", wC); //18
1051 return expectedNumberOfClusters;
1054 //___________________________________________________________________
1056 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1058 // Starting from current radial position of track <t> this function
1059 // extrapolates the track up to radial position <xToGo>.
1060 // Returns 1 if track reaches the plane, and 0 otherwise
1062 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1064 Double_t alpha=t.GetAlpha();
1066 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1067 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1069 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1071 Bool_t lookForCluster;
1072 Double_t x0, rho, x, dx, y, ymax, z;
1076 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1078 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1080 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1082 y = t.GetY(); z = t.GetZ();
1084 // first propagate to the outer surface of the current time bin
1085 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1086 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1087 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1089 ymax = x*TMath::Tan(0.5*alpha);
1092 if (!t.Rotate(alpha)) return 0;
1093 } else if (y <-ymax) {
1095 if (!t.Rotate(-alpha)) return 0;
1097 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1099 y = t.GetY(); z = t.GetZ();
1101 // now propagate to the middle plane of the next time bin
1102 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1103 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1104 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1106 ymax = x*TMath::Tan(0.5*alpha);
1109 if (!t.Rotate(alpha)) return 0;
1110 } else if (y <-ymax) {
1112 if (!t.Rotate(-alpha)) return 0;
1114 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1119 //___________________________________________________________________
1121 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1123 // Starting from current radial position of track <t> this function
1124 // extrapolates the track up to radial position of the outermost
1125 // padrow of the TPC.
1126 // Returns 1 if track reaches the TPC, and 0 otherwise
1128 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1130 Double_t alpha=t.GetAlpha();
1132 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1133 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1135 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1137 Bool_t lookForCluster;
1138 Double_t x0, rho, x, dx, y, ymax, z;
1142 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1144 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1146 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1148 y = t.GetY(); z = t.GetZ();
1150 // first propagate to the outer surface of the current time bin
1151 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1152 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1153 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1155 ymax = x*TMath::Tan(0.5*alpha);
1158 if (!t.Rotate(alpha)) return 0;
1159 } else if (y <-ymax) {
1161 if (!t.Rotate(-alpha)) return 0;
1163 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1165 y = t.GetY(); z = t.GetZ();
1167 // now propagate to the middle plane of the next time bin
1168 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,x0,lookForCluster);
1169 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1170 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1172 ymax = x*TMath::Tan(0.5*alpha);
1175 if (!t.Rotate(alpha)) return 0;
1176 } else if (y <-ymax) {
1178 if (!t.Rotate(-alpha)) return 0;
1180 if(!t.PropagateTo(x,x0,rho,0.139)) return 0;
1186 //_____________________________________________________________________________
1187 void AliTRDtracker::LoadEvent()
1189 // Fills clusters into TRD tracking_sectors
1190 // Note that the numbering scheme for the TRD tracking_sectors
1191 // differs from that of TRD sectors
1193 ReadClusters(fClusters);
1194 Int_t ncl=fClusters->GetEntriesFast();
1195 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1199 printf("\r %d left ",ncl);
1200 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1201 Int_t detector=c->GetDetector(), local_time_bin=c->GetLocalTimeBin();
1202 Int_t sector=fGeom->GetSector(detector);
1203 Int_t plane=fGeom->GetPlane(detector);
1205 Int_t tracking_sector = CookSectorIndex(sector);
1207 Int_t gtb = fTrSec[tracking_sector]->CookTimeBinIndex(plane,local_time_bin);
1208 Int_t layer = fTrSec[tracking_sector]->GetLayerNumber(gtb);
1211 fTrSec[tracking_sector]->GetLayer(layer)->InsertCluster(c,index);
1217 //_____________________________________________________________________________
1218 void AliTRDtracker::UnloadEvent()
1221 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1226 nentr = fClusters->GetEntriesFast();
1227 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1229 nentr = fSeeds->GetEntriesFast();
1230 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1232 nentr = fTracks->GetEntriesFast();
1233 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1235 Int_t nsec = AliTRDgeometry::kNsect;
1237 for (i = 0; i < nsec; i++) {
1238 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1239 fTrSec[i]->GetLayer(pl)->Clear();
1245 //__________________________________________________________________________
1246 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1248 // Creates track seeds using clusters in timeBins=i1,i2
1251 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1255 Double_t x[5], c[15];
1256 Int_t max_sec=AliTRDgeometry::kNsect;
1258 Double_t alpha=AliTRDgeometry::GetAlpha();
1259 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1260 Double_t cs=cos(alpha), sn=sin(alpha);
1261 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1264 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1265 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1267 Double_t x1 =fTrSec[0]->GetX(i1);
1268 Double_t xx2=fTrSec[0]->GetX(i2);
1270 for (Int_t ns=0; ns<max_sec; ns++) {
1272 Int_t nl2 = *(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1273 Int_t nl=(*fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1274 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1275 Int_t nu=(*fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1276 Int_t nu2=(*fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1278 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1280 for (Int_t is=0; is < r1; is++) {
1281 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1283 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1285 const AliTRDcluster *cl;
1286 Double_t x2, y2, z2;
1287 Double_t x3=0., y3=0.;
1290 if(turn != 2) continue;
1291 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1293 y2=cl->GetY(); z2=cl->GetZ();
1298 else if (js<nl2+nl) {
1299 if(turn != 1) continue;
1300 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1302 y2=cl->GetY(); z2=cl->GetZ();
1307 else if (js<nl2+nl+nm) {
1308 if(turn != 1) continue;
1309 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1311 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1313 else if (js<nl2+nl+nm+nu) {
1314 if(turn != 1) continue;
1315 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1316 cl=r2[js-nl2-nl-nm];
1317 y2=cl->GetY(); z2=cl->GetZ();
1323 if(turn != 2) continue;
1324 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1325 cl=r2[js-nl2-nl-nm-nu];
1326 y2=cl->GetY(); z2=cl->GetZ();
1332 if(TMath::Abs(z1-z2) > fMaxSeedDeltaZ12) continue;
1334 Double_t zz=z1 - z1/x1*(x1-x2);
1336 if (TMath::Abs(zz-z2)>fMaxSeedDeltaZ) continue;
1338 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1339 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1343 x[2]=f1trd(x1,y1,x2,y2,x3,y3);
1345 if (TMath::Abs(x[2]) > fMaxSeedC) continue;
1347 x[3]=f2trd(x1,y1,x2,y2,x3,y3);
1349 if (TMath::Abs(x[2]*x1-x[3]) >= 0.99999) continue;
1351 x[4]=f3trd(x1,y1,x2,y2,z1,z2);
1353 if (TMath::Abs(x[4]) > fMaxSeedTan) continue;
1355 Double_t a=asin(x[3]);
1356 Double_t zv=z1 - x[4]/x[2]*(a+asin(x[2]*x1-x[3]));
1358 if (TMath::Abs(zv)>fMaxSeedVertexZ) continue;
1360 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1361 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1362 Double_t sy3=fSeedErrorSY3, sy=fSeedErrorSY, sz=fSeedErrorSZ;
1364 Double_t f20=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1365 Double_t f22=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1366 Double_t f24=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1367 Double_t f30=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[3])/sy;
1368 Double_t f32=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[3])/sy;
1369 Double_t f34=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[3])/sy;
1370 Double_t f40=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[4])/sy;
1371 Double_t f41=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[4])/sz;
1372 Double_t f42=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[4])/sy;
1373 Double_t f43=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[4])/sz;
1377 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f24*sy3*f24;
1378 c[6]=f30*sy1; c[7]=0.; c[8]=f30*sy1*f20+f32*sy2*f22+f34*sy3*f24;
1379 c[9]=f30*sy1*f30+f32*sy2*f32+f34*sy3*f34;
1380 c[10]=f40*sy1; c[11]=f41*sz1; c[12]=f40*sy1*f20+f42*sy2*f22;
1381 c[13]=f40*sy1*f30+f42*sy2*f32;
1382 c[14]=f40*sy1*f40+f41*sz1*f41+f42*sy2*f42+f43*sz2*f43;
1384 UInt_t index=r1.GetIndex(is);
1386 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1388 Int_t rc=FollowProlongation(*track, i2);
1391 (track->GetNumberOfClusters() <
1392 (outer-inner)*fMinClustersInSeed)) delete track;
1394 fSeeds->AddLast(track); fNseeds++;
1395 cerr<<"\r found seed "<<fNseeds;
1402 //_____________________________________________________________________________
1403 void AliTRDtracker::ReadClusters(TObjArray *array, const TFile *inp)
1406 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1407 // from the file. The names of the cluster tree and branches
1408 // should match the ones used in AliTRDclusterizer::WriteClusters()
1411 TDirectory *savedir=gDirectory;
1414 TFile *in=(TFile*)inp;
1415 if (!in->IsOpen()) {
1416 cerr<<"AliTRDtracker::ReadClusters(): input file is not open !\n";
1424 Char_t treeName[12];
1425 sprintf(treeName,"TreeR%d_TRD",fEvent);
1426 TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1428 TObjArray *ClusterArray = new TObjArray(400);
1430 ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray);
1432 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1433 printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1435 // Loop through all entries in the tree
1437 AliTRDcluster *c = 0;
1440 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1443 nbytes += ClusterTree->GetEvent(iEntry);
1445 // Get the number of points in the detector
1446 Int_t nCluster = ClusterArray->GetEntriesFast();
1447 printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1449 // Loop through all TRD digits
1450 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1451 c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1452 AliTRDcluster *co = new AliTRDcluster(*c);
1453 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1454 Int_t ltb = co->GetLocalTimeBin();
1455 if(ltb != 0) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1458 delete ClusterArray->RemoveAt(iCluster);
1462 delete ClusterArray;
1467 //______________________________________________________________________
1468 void AliTRDtracker::ReadClusters(TObjArray *array, const Char_t *filename)
1471 // Reads AliTRDclusters from file <filename>. The names of the cluster
1472 // tree and branches should match the ones used in
1473 // AliTRDclusterizer::WriteClusters()
1474 // if <array> == 0, clusters are added into AliTRDtracker fCluster array
1477 TDirectory *savedir=gDirectory;
1479 TFile *file = TFile::Open(filename);
1480 if (!file->IsOpen()) {
1481 cerr<<"Can't open file with TRD clusters"<<endl;
1485 Char_t treeName[12];
1486 sprintf(treeName,"TreeR%d_TRD",fEvent);
1487 TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1490 cerr<<"AliTRDtracker::ReadClusters(): ";
1491 cerr<<"can't get a tree with clusters !\n";
1495 TObjArray *ClusterArray = new TObjArray(400);
1497 ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray);
1499 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1500 cout<<"found "<<nEntries<<" in ClusterTree"<<endl;
1502 // Loop through all entries in the tree
1504 AliTRDcluster *c = 0;
1508 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1511 nbytes += ClusterTree->GetEvent(iEntry);
1513 // Get the number of points in the detector
1514 Int_t nCluster = ClusterArray->GetEntriesFast();
1515 printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1517 // Loop through all TRD digits
1518 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1519 c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1520 AliTRDcluster *co = new AliTRDcluster(*c);
1521 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1522 Int_t ltb = co->GetLocalTimeBin();
1523 if(ltb != 0) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1525 delete ClusterArray->RemoveAt(iCluster);
1530 delete ClusterArray;
1536 //__________________________________________________________________
1537 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const {
1539 Int_t label=123456789, index, i, j;
1540 Int_t ncl=pt->GetNumberOfClusters();
1541 const Int_t range = fTrSec[0]->GetOuterTimeBin()+1;
1545 // Int_t s[range][2];
1546 Int_t **s = new Int_t* [range];
1547 for (i=0; i<range; i++) {
1548 s[i] = new Int_t[2];
1550 for (i=0; i<range; i++) {
1556 for (i=0; i<ncl; i++) {
1557 index=pt->GetClusterIndex(i);
1558 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1564 for (i=0; i<ncl; i++) {
1565 index=pt->GetClusterIndex(i);
1566 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1567 for (Int_t k=0; k<3; k++) {
1568 label=c->GetLabel(k);
1569 label_added=kFALSE; j=0;
1571 while ( (!label_added) && ( j < range ) ) {
1572 if (s[j][0]==label || s[j][1]==0) {
1586 for (i=0; i<range; i++) {
1588 max=s[i][1]; label=s[i][0];
1592 for (i=0; i<range; i++) {
1598 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1600 pt->SetLabel(label);
1605 //__________________________________________________________________
1606 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const {
1607 Int_t ncl=t->GetNumberOfClusters();
1608 for (Int_t i=from; i<ncl; i++) {
1609 Int_t index = t->GetClusterIndex(i);
1610 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1616 //_____________________________________________________________________
1617 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t r, Double_t tgl, Double_t pt)
1619 // Parametrised "expected" error of the cluster reconstruction in Y
1621 Double_t s = 0.08 * 0.08;
1625 //_____________________________________________________________________
1626 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t r, Double_t tgl)
1628 // Parametrised "expected" error of the cluster reconstruction in Z
1630 Double_t s = 6 * 6 /12.;
1635 //_____________________________________________________________________
1636 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t local_tb) const
1639 // Returns radial position which corresponds to time bin <local_tb>
1640 // in tracking sector <sector> and plane <plane>
1643 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, local_tb);
1644 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1645 return fTrSec[sector]->GetLayer(pl)->GetX();
1650 //_______________________________________________________
1651 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1652 Double_t dx, Double_t rho, Double_t x0, Int_t tb_index)
1655 // AliTRDpropagationLayer constructor
1658 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = x0;
1659 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tb_index;
1662 for(Int_t i=0; i < (Int_t) kZONES; i++) {
1663 fZc[i]=0; fZmax[i] = 0;
1668 if(fTimeBinIndex >= 0) {
1669 fClusters = new AliTRDcluster*[kMAX_CLUSTER_PER_TIME_BIN];
1670 fIndex = new UInt_t[kMAX_CLUSTER_PER_TIME_BIN];
1683 //_______________________________________________________
1684 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1685 Double_t Zmax, Double_t Ymax, Double_t rho,
1686 Double_t x0, Double_t Yc, Double_t Zc)
1689 // Sets hole in the layer
1702 //_______________________________________________________
1703 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1706 // AliTRDtrackingSector Constructor
1715 for(UInt_t i=0; i < kMAX_TIME_BIN_INDEX; i++) fTimeBinIndex[i] = -1;
1718 AliTRDpropagationLayer* ppl;
1720 Double_t x, xin, xout, dx, rho, x0;
1723 // set time bins in the gas of the TPC
1725 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
1726 rho = 0.9e-3; x0 = 28.94;
1728 for(Int_t i=0; i<steps; i++) {
1729 x = xin + i*dx + dx/2;
1730 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1734 // set time bins in the outer field cage vessel
1736 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1737 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1740 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1741 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1744 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; x0 = 41.28; // Nomex
1745 steps = 5; dx = (xout - xin)/steps;
1746 for(Int_t i=0; i<steps; i++) {
1747 x = xin + i*dx + dx/2;
1748 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1752 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1753 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1756 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1757 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1761 // set time bins in CO2
1763 xin = xout; xout = 275.0;
1764 steps = 50; dx = (xout - xin)/steps;
1765 rho = 1.977e-3; x0 = 36.2;
1767 for(Int_t i=0; i<steps; i++) {
1768 x = xin + i*dx + dx/2;
1769 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1773 // set time bins in the outer containment vessel
1775 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; x0 = 24.01; // Al
1776 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1779 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1780 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1783 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1784 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1787 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; x0 = 41.28; // Nomex
1788 steps = 10; dx = (xout - xin)/steps;
1789 for(Int_t i=0; i<steps; i++) {
1790 x = xin + i*dx + dx/2;
1791 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1795 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; x0 = 44.86; // prepreg
1796 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1799 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; x0 = 44.77; // Tedlar
1800 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1803 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; x0 = 24.01; // Al
1804 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1807 Double_t xtrd = (Double_t) fGeom->Rmin();
1809 // add layers between TPC and TRD (Air temporarily)
1810 xin = xout; xout = xtrd;
1811 steps = 50; dx = (xout - xin)/steps;
1812 rho = 1.2e-3; x0 = 36.66;
1814 for(Int_t i=0; i<steps; i++) {
1815 x = xin + i*dx + dx/2;
1816 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1821 Double_t alpha=AliTRDgeometry::GetAlpha();
1823 // add layers for each of the planes
1825 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
1826 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
1827 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
1828 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
1829 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
1830 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
1831 Double_t dxPlane = dxTEC + dxSpace;
1833 Int_t tbBefore = (Int_t) (dxAmp/fPar->GetTimeBinSize());
1834 if(tbBefore > fPar->GetTimeBefore()) tbBefore = fPar->GetTimeBefore();
1837 const Int_t nChambers = AliTRDgeometry::Ncham();
1838 Double_t Ymax = 0, holeYmax = 0;
1839 Double_t * Zc = new Double_t[nChambers];
1840 Double_t * Zmax = new Double_t[nChambers];
1841 Double_t holeZmax = 1000.; // the whole sector is missing
1844 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
1847 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
1848 steps = 12; dx = (xout - xin)/steps; rho = 0.074; x0 = 40.6;
1849 for(Int_t i=0; i<steps; i++) {
1850 x = xin + i*dx + dx/2;
1851 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1852 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1853 holeYmax = x*TMath::Tan(0.5*alpha);
1854 ppl->SetHole(holeYmax, holeZmax);
1856 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1857 holeYmax = x*TMath::Tan(0.5*alpha);
1858 ppl->SetHole(holeYmax, holeZmax);
1863 Ymax = fGeom->GetChamberWidth(plane)/2;
1864 for(Int_t ch = 0; ch < nChambers; ch++) {
1865 Zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
1866 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
1867 Float_t row0 = fPar->GetRow0(plane,ch,0);
1868 Int_t nPads = fPar->GetRowMax(plane,ch,0);
1869 Zc[ch] = (pad * nPads)/2 + row0 - pad/2;
1872 dx = fPar->GetTimeBinSize();
1873 rho = 0.00295 * 0.85; x0 = 11.0;
1875 Double_t x0 = (Double_t) fPar->GetTime0(plane);
1876 Double_t xbottom = x0 - dxDrift;
1877 Double_t xtop = x0 + dxAmp;
1879 // Amplification region
1881 steps = (Int_t) (dxAmp/dx);
1883 for(tb = 0; tb < steps; tb++) {
1884 x = x0 + tb * dx + dx/2;
1885 tb_index = CookTimeBinIndex(plane, -tb-1);
1886 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1888 for(Int_t ch = 0; ch < nChambers; ch++) {
1889 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1891 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1892 holeYmax = x*TMath::Tan(0.5*alpha);
1893 ppl->SetHole(holeYmax, holeZmax);
1895 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1896 holeYmax = x*TMath::Tan(0.5*alpha);
1897 ppl->SetHole(holeYmax, holeZmax);
1901 tb_index = CookTimeBinIndex(plane, -steps);
1902 x = (x + dx/2 + xtop)/2;
1904 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1906 for(Int_t ch = 0; ch < nChambers; ch++) {
1907 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1909 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1910 holeYmax = x*TMath::Tan(0.5*alpha);
1911 ppl->SetHole(holeYmax, holeZmax);
1913 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1914 holeYmax = x*TMath::Tan(0.5*alpha);
1915 ppl->SetHole(holeYmax, holeZmax);
1920 dx = fPar->GetTimeBinSize();
1921 steps = (Int_t) (dxDrift/dx);
1923 for(tb = 0; tb < steps; tb++) {
1924 x = x0 - tb * dx - dx/2;
1925 tb_index = CookTimeBinIndex(plane, tb);
1927 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1929 for(Int_t ch = 0; ch < nChambers; ch++) {
1930 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1932 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1933 holeYmax = x*TMath::Tan(0.5*alpha);
1934 ppl->SetHole(holeYmax, holeZmax);
1936 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1937 holeYmax = x*TMath::Tan(0.5*alpha);
1938 ppl->SetHole(holeYmax, holeZmax);
1942 tb_index = CookTimeBinIndex(plane, steps);
1943 x = (x - dx/2 + xbottom)/2;
1945 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,tb_index);
1947 for(Int_t ch = 0; ch < nChambers; ch++) {
1948 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1950 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1951 holeYmax = x*TMath::Tan(0.5*alpha);
1952 ppl->SetHole(holeYmax, holeZmax);
1954 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1955 holeYmax = x*TMath::Tan(0.5*alpha);
1956 ppl->SetHole(holeYmax, holeZmax);
1961 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; x0 = 33.0;
1962 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,x0,-1);
1963 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1964 holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
1965 ppl->SetHole(holeYmax, holeZmax);
1967 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1968 holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
1969 ppl->SetHole(holeYmax, holeZmax);
1974 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
1975 steps = 5; dx = (xout - xin)/steps; rho = 0.074; x0 = 40.6;
1976 for(Int_t i=0; i<steps; i++) {
1977 x = xin + i*dx + dx/2;
1978 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1979 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1980 holeYmax = x*TMath::Tan(0.5*alpha);
1981 ppl->SetHole(holeYmax, holeZmax);
1983 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1984 holeYmax = x*TMath::Tan(0.5*alpha);
1985 ppl->SetHole(holeYmax, holeZmax);
1990 // Space between the chambers, air
1991 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
1992 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; x0 = 36.66;
1993 for(Int_t i=0; i<steps; i++) {
1994 x = xin + i*dx + dx/2;
1995 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
1996 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1997 holeYmax = x*TMath::Tan(0.5*alpha);
1998 ppl->SetHole(holeYmax, holeZmax);
2000 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2001 holeYmax = x*TMath::Tan(0.5*alpha);
2002 ppl->SetHole(holeYmax, holeZmax);
2008 // Space between the TRD and RICH
2009 Double_t xRICH = 500.;
2010 xin = xout; xout = xRICH;
2011 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; x0 = 36.66;
2012 for(Int_t i=0; i<steps; i++) {
2013 x = xin + i*dx + dx/2;
2014 ppl = new AliTRDpropagationLayer(x,dx,rho,x0,-1);
2024 //______________________________________________________
2026 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t local_tb) const
2029 // depending on the digitization parameters calculates "global"
2030 // time bin index for timebin <local_tb> in plane <plane>
2033 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2034 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2035 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2037 Int_t tbAmp = fPar->GetTimeBefore();
2038 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2039 Int_t tbDrift = fPar->GetTimeMax();
2040 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2042 Int_t tb_per_plane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2044 Int_t gtb = (plane+1) * tb_per_plane - local_tb - 1;
2046 if((local_tb < 0) &&
2047 (TMath::Abs(local_tb) > TMath::Min(tbAmp,maxAmp))) return -1;
2048 if(local_tb >= TMath::Min(tbDrift,maxDrift)) return -1;
2055 //______________________________________________________
2057 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2060 // For all sensitive time bins sets corresponding layer index
2061 // in the array fTimeBins
2066 for(Int_t i = 0; i < fN; i++) {
2067 index = fLayers[i]->GetTimeBinIndex();
2069 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2071 if(index < 0) continue;
2072 if(index >= (Int_t) kMAX_TIME_BIN_INDEX) {
2073 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2074 printf(" index %d exceeds allowed maximum of %d!\n",
2075 index, kMAX_TIME_BIN_INDEX-1);
2078 fTimeBinIndex[index] = i;
2081 Double_t x1, dx1, x2, dx2, gap;
2083 for(Int_t i = 0; i < fN-1; i++) {
2084 x1 = fLayers[i]->GetX();
2085 dx1 = fLayers[i]->GetdX();
2086 x2 = fLayers[i+1]->GetX();
2087 dx2 = fLayers[i+1]->GetdX();
2088 gap = (x2 - dx2/2) - (x1 + dx1/2);
2090 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2091 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2094 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2095 printf(" (%f - %f) - (%f + %f) = %f\n",
2096 x2, dx2/2, x1, dx1, gap);
2102 //______________________________________________________
2105 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2108 // Returns the number of time bin which in radial position is closest to <x>
2111 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2112 if(x <= fLayers[0]->GetX()) return 0;
2114 Int_t b=0, e=fN-1, m=(b+e)/2;
2115 for (; b<e; m=(b+e)/2) {
2116 if (x > fLayers[m]->GetX()) b=m+1;
2119 if(TMath::Abs(x - fLayers[m]->GetX()) >
2120 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2125 //______________________________________________________
2127 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2130 // Returns number of the innermost SENSITIVE propagation layer
2133 return GetLayerNumber(0);
2136 //______________________________________________________
2138 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2141 // Returns number of the outermost SENSITIVE time bin
2144 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2147 //______________________________________________________
2149 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2152 // Returns number of SENSITIVE time bins
2156 for(tb = kMAX_TIME_BIN_INDEX-1; tb >=0; tb--) {
2157 layer = GetLayerNumber(tb);
2163 //______________________________________________________
2165 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2168 // Insert layer <pl> in fLayers array.
2169 // Layers are sorted according to X coordinate.
2171 if ( fN == ((Int_t) kMAX_LAYERS_PER_SECTOR)) {
2172 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2175 if (fN==0) {fLayers[fN++] = pl; return;}
2176 Int_t i=Find(pl->GetX());
2178 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2179 fLayers[i]=pl; fN++;
2183 //______________________________________________________
2185 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2188 // Returns index of the propagation layer nearest to X
2191 if (x <= fLayers[0]->GetX()) return 0;
2192 if (x > fLayers[fN-1]->GetX()) return fN;
2193 Int_t b=0, e=fN-1, m=(b+e)/2;
2194 for (; b<e; m=(b+e)/2) {
2195 if (x > fLayers[m]->GetX()) b=m+1;
2201 //______________________________________________________
2203 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2204 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &x0,
2205 Bool_t &lookForCluster) const
2208 // Returns radial step <dx>, density <rho>, rad. length <x0>,
2209 // and sensitivity <lookForCluster> in point <y,z>
2215 lookForCluster = kFALSE;
2217 // check dead regions
2218 if(fTimeBinIndex >= 0) {
2219 for(Int_t ch = 0; ch < (Int_t) kZONES; ch++) {
2220 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2221 lookForCluster = kTRUE;
2223 if(TMath::Abs(y) > fYmax) lookForCluster = kFALSE;
2224 if(!lookForCluster) {
2225 // rho = 1.7; x0 = 33.0; // G10
2230 if(fHole && (TMath::Abs(y - fHoleYc) < fHoleYmax) &&
2231 (TMath::Abs(z - fHoleZc) < fHoleZmax)) {
2232 lookForCluster = kFALSE;
2240 //______________________________________________________
2242 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2245 // Insert cluster in cluster array.
2246 // Clusters are sorted according to Y coordinate.
2248 if(fTimeBinIndex < 0) {
2249 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2253 if (fN== (Int_t) kMAX_CLUSTER_PER_TIME_BIN) {
2254 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2257 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2258 Int_t i=Find(c->GetY());
2259 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2260 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2261 fIndex[i]=index; fClusters[i]=c; fN++;
2264 //______________________________________________________
2266 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2268 // Returns index of the cluster nearest in Y
2270 if (y <= fClusters[0]->GetY()) return 0;
2271 if (y > fClusters[fN-1]->GetY()) return fN;
2272 Int_t b=0, e=fN-1, m=(b+e)/2;
2273 for (; b<e; m=(b+e)/2) {
2274 if (y > fClusters[m]->GetY()) b=m+1;