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.23 2003/02/10 14:06:10 cblume
19 Add tracking without tilted pads as option
21 Revision 1.22 2003/01/30 15:19:58 cblume
24 Revision 1.21 2003/01/27 16:34:49 cblume
25 Update of tracking by Sergei and Chuncheng
27 Revision 1.20 2002/11/07 15:52:09 cblume
28 Update of tracking code for tilted pads
30 Revision 1.19 2002/10/22 15:53:08 alibrary
31 Introducing Riostream.h
33 Revision 1.18 2002/10/14 14:57:44 hristov
34 Merging the VirtualMC branch to the main development branch (HEAD)
36 Revision 1.14.6.2 2002/07/24 10:09:31 alibrary
39 Revision 1.17 2002/06/13 12:09:58 hristov
42 Revision 1.16 2002/06/12 09:54:36 cblume
43 Update of tracking code provided by Sergei
45 Revision 1.14 2001/11/14 10:50:46 cblume
46 Changes in digits IO. Add merging of summable digits
48 Revision 1.13 2001/05/30 12:17:47 hristov
49 Loop variables declared once
51 Revision 1.12 2001/05/28 17:07:58 hristov
52 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)
54 Revision 1.8 2000/12/20 13:00:44 cblume
55 Modifications for the HP-compiler
57 Revision 1.7 2000/12/08 16:07:02 cblume
58 Update of the tracking by Sergei
60 Revision 1.6 2000/11/30 17:38:08 cblume
61 Changes to get in line with new STEER and EVGEN
63 Revision 1.5 2000/11/14 14:40:27 cblume
64 Correction for the Sun compiler (kTRUE and kFALSE)
66 Revision 1.4 2000/11/10 14:57:52 cblume
67 Changes in the geometry constants for the DEC compiler
69 Revision 1.3 2000/10/15 23:40:01 cblume
72 Revision 1.2 2000/10/06 16:49:46 cblume
75 Revision 1.1.2.2 2000/10/04 16:34:58 cblume
76 Replace include files by forward declarations
78 Revision 1.1.2.1 2000/09/22 14:47:52 cblume
83 #include <Riostream.h>
88 #include <TObjArray.h>
90 #include "AliTRDgeometry.h"
91 #include "AliTRDparameter.h"
92 #include "AliTRDgeometryDetail.h"
93 #include "AliTRDcluster.h"
94 #include "AliTRDtrack.h"
95 #include "../TPC/AliTPCtrack.h"
97 #include "AliTRDtracker.h"
99 ClassImp(AliTRDtracker)
101 const Float_t AliTRDtracker::fSeedDepth = 0.5;
102 const Float_t AliTRDtracker::fSeedStep = 0.10;
103 const Float_t AliTRDtracker::fSeedGap = 0.25;
105 const Float_t AliTRDtracker::fMaxSeedDeltaZ12 = 40.;
106 const Float_t AliTRDtracker::fMaxSeedDeltaZ = 25.;
107 const Float_t AliTRDtracker::fMaxSeedC = 0.0052;
108 const Float_t AliTRDtracker::fMaxSeedTan = 1.2;
109 const Float_t AliTRDtracker::fMaxSeedVertexZ = 150.;
111 const Double_t AliTRDtracker::fSeedErrorSY = 0.2;
112 const Double_t AliTRDtracker::fSeedErrorSY3 = 2.5;
113 const Double_t AliTRDtracker::fSeedErrorSZ = 0.1;
115 const Float_t AliTRDtracker::fMinClustersInSeed = 0.7;
117 const Float_t AliTRDtracker::fMinClustersInTrack = 0.5;
118 const Float_t AliTRDtracker::fMinFractionOfFoundClusters = 0.8;
120 const Float_t AliTRDtracker::fSkipDepth = 0.3;
121 const Float_t AliTRDtracker::fLabelFraction = 0.8;
122 const Float_t AliTRDtracker::fWideRoad = 20.;
124 const Double_t AliTRDtracker::fMaxChi2 = 12.;
127 //____________________________________________________________________
128 AliTRDtracker::AliTRDtracker(const TFile *geomfile)
136 fAddTRDseeds = kFALSE;
140 TDirectory *savedir=gDirectory;
141 TFile *in=(TFile*)geomfile;
143 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
144 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
149 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
150 fPar = (AliTRDparameter*) in->Get("TRDparameter");
155 // fTzero = geo->GetT0();
156 printf("Found geometry version %d on file \n", fGeom->IsVersion());
159 printf("AliTRDtracker::AliTRDtracker(): cann't find TRD geometry!\n");
160 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
161 fGeom = new AliTRDgeometryDetail();
162 fPar = new AliTRDparameter();
168 // fGeom->SetT0(fTzero);
173 fClusters = new TObjArray(2000);
175 fSeeds = new TObjArray(2000);
177 fTracks = new TObjArray(1000);
179 for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
180 Int_t tr_s = CookSectorIndex(geom_s);
181 fTrSec[tr_s] = new AliTRDtrackingSector(fGeom, geom_s, fPar);
184 Float_t tilt_angle = TMath::Abs(fPar->GetTiltingAngle());
185 if(tilt_angle < 0.1) {
192 if(fNoTilt && (tilt_angle > 0.1)) fSY2corr = fSY2corr + tilt_angle * 0.05;
195 // calculate max gap on track
197 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
198 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
200 Double_t dx = (Double_t) fPar->GetTimeBinSize();
201 Int_t tbAmp = fPar->GetTimeBefore();
202 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
203 if(kTRUE) maxAmp = 0; // intentional until we change the parameter class
204 Int_t tbDrift = fPar->GetTimeMax();
205 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
207 tbDrift = TMath::Min(tbDrift,maxDrift);
208 tbAmp = TMath::Min(tbAmp,maxAmp);
210 fTimeBinsPerPlane = tbAmp + tbDrift;
211 fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fSkipDepth);
217 //___________________________________________________________________
218 AliTRDtracker::~AliTRDtracker()
226 for(Int_t geom_s = 0; geom_s < kTRACKING_SECTORS; geom_s++) {
227 delete fTrSec[geom_s];
231 //_____________________________________________________________________
232 inline Double_t f1trd(Double_t x1,Double_t y1,
233 Double_t x2,Double_t y2,
234 Double_t x3,Double_t y3)
237 // Initial approximation of the track curvature
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));
245 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
247 return -xr*yr/sqrt(xr*xr+yr*yr);
250 //_____________________________________________________________________
251 inline Double_t f2trd(Double_t x1,Double_t y1,
252 Double_t x2,Double_t y2,
253 Double_t x3,Double_t y3)
256 // Initial approximation of the track curvature times X coordinate
257 // of the center of curvature
260 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
261 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
262 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
263 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
264 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
266 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
268 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
271 //_____________________________________________________________________
272 inline Double_t f3trd(Double_t x1,Double_t y1,
273 Double_t x2,Double_t y2,
274 Double_t z1,Double_t z2)
277 // Initial approximation of the tangent of the track dip angle
280 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
283 //___________________________________________________________________
284 Int_t AliTRDtracker::Clusters2Tracks(const TFile *inp, TFile *out)
287 // Finds tracks within the TRD. File <inp> is expected to contain seeds
288 // at the outer part of the TRD. If <inp> is NULL, the seeds
289 // are found within the TRD if fAddTRDseeds is TRUE.
290 // The tracks are propagated to the innermost time bin
291 // of the TRD and stored in file <out>.
296 TDirectory *savedir=gDirectory;
300 if (!out->IsOpen()) {
301 cerr<<"AliTRDtracker::Clusters2Tracks(): output file is not open !\n";
305 sprintf(tname,"seedTRDtoTPC_%d",fEvent);
306 TTree tpc_tree(tname,"Tree with seeds from TRD at outer TPC pad row");
307 AliTPCtrack *iotrack=0;
308 tpc_tree.Branch("tracks","AliTPCtrack",&iotrack,32000,0);
310 sprintf(tname,"TreeT%d_TRD",fEvent);
311 TTree trd_tree(tname,"TRD tracks at inner TRD time bin");
312 AliTRDtrack *iotrack_trd=0;
313 trd_tree.Branch("tracks","AliTRDtrack",&iotrack_trd,32000,0);
315 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
316 Float_t foundMin = fMinClustersInTrack * timeBins;
319 TFile *in=(TFile*)inp;
321 cerr<<"AliTRDtracker::Clusters2Tracks(): file with seeds is not open !\n";
322 cerr<<" ... going for seeds finding inside the TRD\n";
326 sprintf(tname,"TRDb_%d",fEvent);
327 TTree *seedTree=(TTree*)in->Get(tname);
329 cerr<<"AliTRDtracker::Clusters2Tracks(): ";
330 cerr<<"can't get a tree with track seeds !\n";
333 AliTRDtrack *seed=new AliTRDtrack;
334 seedTree->SetBranchAddress("tracks",&seed);
336 Int_t n=(Int_t)seedTree->GetEntries();
337 for (Int_t i=0; i<n; i++) {
338 seedTree->GetEvent(i);
339 seed->ResetCovariance();
340 AliTRDtrack *tr = new AliTRDtrack(*seed,seed->GetAlpha());
352 // find tracks from loaded seeds
354 Int_t nseed=fSeeds->GetEntriesFast();
356 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
358 for (i=0; i<nseed; i++) {
359 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
360 FollowProlongation(t, innerTB);
361 if (t.GetNumberOfClusters() >= foundMin) {
363 CookLabel(pt, 1-fLabelFraction);
370 if(PropagateToTPC(t)) {
371 AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
376 delete fSeeds->RemoveAt(i);
380 cout<<"Number of loaded seeds: "<<nseed<<endl;
381 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
383 // after tracks from loaded seeds are found and the corresponding
384 // clusters are used, look for additional seeds from TRD
387 // Find tracks for the seeds in the TRD
388 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
390 Int_t nSteps = (Int_t) (fSeedDepth / fSeedStep);
391 Int_t gap = (Int_t) (timeBins * fSeedGap);
392 Int_t step = (Int_t) (timeBins * fSeedStep);
394 // make a first turn with tight cut on initial curvature
395 for(Int_t turn = 1; turn <= 2; turn++) {
397 nSteps = (Int_t) (fSeedDepth / (3*fSeedStep));
398 step = (Int_t) (timeBins * (3*fSeedStep));
400 for(Int_t i=0; i<nSteps; i++) {
401 Int_t outer=timeBins-1-i*step;
402 Int_t inner=outer-gap;
404 nseed=fSeeds->GetEntriesFast();
406 MakeSeeds(inner, outer, turn);
408 nseed=fSeeds->GetEntriesFast();
409 printf("\n turn %d, step %d: number of seeds for TRD inward %d\n",
412 for (Int_t i=0; i<nseed; i++) {
413 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
414 FollowProlongation(t,innerTB);
415 if (t.GetNumberOfClusters() >= foundMin) {
417 CookLabel(pt, 1-fLabelFraction);
422 if(PropagateToTPC(t)) {
423 AliTPCtrack *tpc = new AliTPCtrack(*pt,pt->GetAlpha());
429 delete fSeeds->RemoveAt(i);
438 cout<<"Total number of found tracks: "<<found<<endl;
449 //_____________________________________________________________________________
450 Int_t AliTRDtracker::PropagateBack(const TFile *inp, TFile *out) {
452 // Reads seeds from file <inp>. The seeds are AliTPCtrack's found and
453 // backpropagated by the TPC tracker. Each seed is first propagated
454 // to the TRD, and then its prolongation is searched in the TRD.
455 // If sufficiently long continuation of the track is found in the TRD
456 // the track is updated, otherwise it's stored as originaly defined
457 // by the TPC tracker.
462 TDirectory *savedir=gDirectory;
464 TFile *in=(TFile*)inp;
467 cerr<<"AliTRDtracker::PropagateBack(): ";
468 cerr<<"file with back propagated TPC tracks is not open !\n";
472 if (!out->IsOpen()) {
473 cerr<<"AliTRDtracker::PropagateBack(): ";
474 cerr<<"file for back propagated TRD tracks is not open !\n";
480 sprintf(tname,"seedsTPCtoTRD_%d",fEvent);
481 TTree *seedTree=(TTree*)in->Get(tname);
483 cerr<<"AliTRDtracker::PropagateBack(): ";
484 cerr<<"can't get a tree with seeds from TPC !\n";
485 cerr<<"check if your version of TPC tracker creates tree "<<tname<<"\n";
489 AliTPCtrack *seed=new AliTPCtrack;
490 seedTree->SetBranchAddress("tracks",&seed);
492 Int_t n=(Int_t)seedTree->GetEntries();
493 for (Int_t i=0; i<n; i++) {
494 seedTree->GetEvent(i);
495 Int_t lbl = seed->GetLabel();
496 AliTRDtrack *tr = new AliTRDtrack(*seed,seed->GetAlpha());
497 tr->SetSeedLabel(lbl);
507 AliTPCtrack *otrack=0;
509 sprintf(tname,"seedsTRDtoTOF1_%d",fEvent);
510 TTree tofTree1(tname,"Tracks back propagated through TPC and TRD");
511 tofTree1.Branch("tracks","AliTPCtrack",&otrack,32000,0);
513 sprintf(tname,"seedsTRDtoTOF2_%d",fEvent);
514 TTree tofTree2(tname,"Tracks back propagated through TPC and TRD");
515 tofTree2.Branch("tracks","AliTPCtrack",&otrack,32000,0);
517 sprintf(tname,"seedsTRDtoPHOS_%d",fEvent);
518 TTree phosTree(tname,"Tracks back propagated through TPC and TRD");
519 phosTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
521 sprintf(tname,"seedsTRDtoRICH_%d",fEvent);
522 TTree richTree(tname,"Tracks back propagated through TPC and TRD");
523 richTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
525 sprintf(tname,"TRDb_%d",fEvent);
526 TTree trdTree(tname,"Back propagated TRD tracks at outer TRD time bin");
527 AliTRDtrack *otrack_trd=0;
528 trdTree.Branch("tracks","AliTRDtrack",&otrack_trd,32000,0);
532 Int_t nseed=fSeeds->GetEntriesFast();
534 // Float_t foundMin = fMinClustersInTrack * fTimeBinsPerPlane * fGeom->Nplan();
535 Float_t foundMin = 0;
537 Int_t outermost_tb = fTrSec[0]->GetOuterTimeBin();
539 for (Int_t i=0; i<nseed; i++) {
541 AliTRDtrack *ps=(AliTRDtrack*)fSeeds->UncheckedAt(i), &s=*ps;
542 Int_t expectedClr = FollowBackProlongation(s);
543 Int_t foundClr = s.GetNumberOfClusters();
544 Int_t last_tb = fTrSec[0]->GetLayerNumber(s.GetX());
546 // printf("seed %d: found %d out of %d expected clusters, Min is %f\n",
547 // i, foundClr, expectedClr, foundMin);
549 if (foundClr >= foundMin) {
552 CookLabel(ps, 1-fLabelFraction);
556 // Propagate to outer reference plane [SR, GSI, 18.02.2003]
557 ps->PropagateTo(364.8);
563 if(((expectedClr < 10) && (last_tb == outermost_tb)) ||
564 ((expectedClr >= 10) &&
565 (((Float_t) foundClr) / ((Float_t) expectedClr) >=
566 fMinFractionOfFoundClusters) && (last_tb == outermost_tb))) {
568 Double_t x_tof = 375.5;
570 if(PropagateToOuterPlane(s,x_tof)) {
571 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
578 if(PropagateToOuterPlane(s,x_tof)) {
579 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
584 Double_t x_phos = 460.;
586 if(PropagateToOuterPlane(s,x_phos)) {
587 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
592 Double_t x_rich = 490+1.267;
594 if(PropagateToOuterPlane(s,x_rich)) {
595 AliTPCtrack *pt = new AliTPCtrack(*ps,ps->GetAlpha());
613 cerr<<"Number of seeds: "<<nseed<<endl;
614 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
623 //---------------------------------------------------------------------------
624 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
626 // Starting from current position on track=t this function tries
627 // to extrapolate the track up to timeBin=0 and to confirm prolongation
628 // if a close cluster is found. Returns the number of clusters
629 // expected to be found in sensitive layers
631 Float_t wIndex, wTB, wChi2;
632 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
633 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
634 Float_t wPx, wPy, wPz, wC;
636 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
638 Int_t trackIndex = t.GetLabel();
640 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
642 Int_t try_again=fMaxGap;
644 Double_t alpha=t.GetAlpha();
646 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
647 if (alpha < 0. ) alpha += 2.*TMath::Pi();
649 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
651 Double_t rad_length, rho, x, dx, y, ymax, z;
653 Int_t expectedNumberOfClusters = 0;
654 Bool_t lookForCluster;
656 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
659 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
661 y = t.GetY(); z = t.GetZ();
663 // first propagate to the inner surface of the current time bin
664 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
665 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
666 if(!t.PropagateTo(x,rad_length,rho)) break;
668 ymax = x*TMath::Tan(0.5*alpha);
671 if (!t.Rotate(alpha)) break;
672 if(!t.PropagateTo(x,rad_length,rho)) break;
673 } else if (y <-ymax) {
675 if (!t.Rotate(-alpha)) break;
676 if(!t.PropagateTo(x,rad_length,rho)) break;
679 y = t.GetY(); z = t.GetZ();
681 // now propagate to the middle plane of the next time bin
682 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
683 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
684 if(!t.PropagateTo(x,rad_length,rho)) break;
686 ymax = x*TMath::Tan(0.5*alpha);
689 if (!t.Rotate(alpha)) break;
690 if(!t.PropagateTo(x,rad_length,rho)) break;
691 } else if (y <-ymax) {
693 if (!t.Rotate(-alpha)) break;
694 if(!t.PropagateTo(x,rad_length,rho)) break;
700 expectedNumberOfClusters++;
701 wIndex = (Float_t) t.GetLabel();
704 AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr-1));
706 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
707 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
710 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
711 else return expectedNumberOfClusters;
715 wYwindow = (Float_t) road;
716 t.GetPxPyPz(Px,Py,Pz);
720 wC = (Float_t) t.GetC();
721 wSigmaC2 = (Float_t) t.GetSigmaC2();
722 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
723 wSigmaY2 = (Float_t) t.GetSigmaY2();
724 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
727 if (road>fWideRoad) {
728 if (t.GetNumberOfClusters()>4)
729 cerr<<t.GetNumberOfClusters()
730 <<"FindProlongation warning: Too broad road !\n";
737 Double_t max_chi2=fMaxChi2;
739 wYclosest = 12345678;
740 wYcorrect = 12345678;
741 wZclosest = 12345678;
742 wZcorrect = 12345678;
743 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
745 // Find the closest correct cluster for debugging purposes
747 Float_t minDY = 1000000;
748 for (Int_t i=0; i<time_bin; i++) {
749 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
750 if((c->GetLabel(0) != trackIndex) &&
751 (c->GetLabel(1) != trackIndex) &&
752 (c->GetLabel(2) != trackIndex)) continue;
753 if(TMath::Abs(c->GetY() - y) > minDY) continue;
754 minDY = TMath::Abs(c->GetY() - y);
755 wYcorrect = c->GetY();
756 wZcorrect = c->GetZ();
758 Double_t h01 = GetTiltFactor(c);
759 wChi2 = t.GetPredictedChi2(c, h01);
763 // Now go for the real cluster search
767 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
768 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
769 if (c->GetY() > y+road) break;
770 if (c->IsUsed() > 0) continue;
771 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
773 Double_t h01 = GetTiltFactor(c);
774 Double_t chi2=t.GetPredictedChi2(c,h01);
776 if (chi2 > max_chi2) continue;
779 index=time_bin.GetIndex(i);
784 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
785 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
787 if (c->GetY() > y+road) break;
788 if (c->IsUsed() > 0) continue;
789 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
791 Double_t h01 = GetTiltFactor(c);
792 Double_t chi2=t.GetPredictedChi2(c, h01);
794 if (chi2 > max_chi2) continue;
797 index=time_bin.GetIndex(i);
803 wYclosest = cl->GetY();
804 wZclosest = cl->GetZ();
805 Double_t h01 = GetTiltFactor(cl);
807 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
808 if(!t.Update(cl,max_chi2,index,h01)) {
809 if(!try_again--) return 0;
811 else try_again=fMaxGap;
814 if (try_again==0) break;
819 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
821 printf(" %f", wIndex); //1
822 printf(" %f", wTB); //2
823 printf(" %f", wYrt); //3
824 printf(" %f", wYclosest); //4
825 printf(" %f", wYcorrect); //5
826 printf(" %f", wYwindow); //6
827 printf(" %f", wZrt); //7
828 printf(" %f", wZclosest); //8
829 printf(" %f", wZcorrect); //9
830 printf(" %f", wZwindow); //10
831 printf(" %f", wPx); //11
832 printf(" %f", wPy); //12
833 printf(" %f", wPz); //13
834 printf(" %f", wSigmaC2*1000000); //14
835 printf(" %f", wSigmaTgl2*1000); //15
836 printf(" %f", wSigmaY2); //16
837 // printf(" %f", wSigmaZ2); //17
838 printf(" %f", wChi2); //17
839 printf(" %f", wC); //18
846 return expectedNumberOfClusters;
851 //___________________________________________________________________
853 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
855 // Starting from current radial position of track <t> this function
856 // extrapolates the track up to outer timebin and in the sensitive
857 // layers confirms prolongation if a close cluster is found.
858 // Returns the number of clusters expected to be found in sensitive layers
860 Float_t wIndex, wTB, wChi2;
861 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
862 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
863 Float_t wPx, wPy, wPz, wC;
865 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
867 Int_t trackIndex = t.GetLabel();
869 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
871 Int_t try_again=fMaxGap;
873 Double_t alpha=t.GetAlpha();
875 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
876 if (alpha < 0. ) alpha += 2.*TMath::Pi();
878 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
880 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
882 Double_t rad_length, rho, x, dx, y, ymax, z;
884 Bool_t lookForCluster;
886 Int_t expectedNumberOfClusters = 0;
889 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
892 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB; nr++) {
894 y = t.GetY(); z = t.GetZ();
896 // first propagate to the outer surface of the current time bin
898 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
899 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
901 if(!t.PropagateTo(x,rad_length,rho)) break;
904 ymax = x*TMath::Tan(0.5*alpha);
908 if (!t.Rotate(alpha)) break;
909 if(!t.PropagateTo(x,rad_length,rho)) break;
910 } else if (y <-ymax) {
912 if (!t.Rotate(-alpha)) break;
913 if(!t.PropagateTo(x,rad_length,rho)) break;
915 y = t.GetY(); z = t.GetZ();
917 // now propagate to the middle plane of the next time bin
918 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
920 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
922 if(!t.PropagateTo(x,rad_length,rho)) break;
926 ymax = x*TMath::Tan(0.5*alpha);
928 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
932 if (!t.Rotate(alpha)) break;
933 if(!t.PropagateTo(x,rad_length,rho)) break;
934 } else if (y <-ymax) {
936 if (!t.Rotate(-alpha)) break;
937 if(!t.PropagateTo(x,rad_length,rho)) break;
940 // printf("label %d, pl %d, lookForCluster %d \n",
941 // trackIndex, nr+1, lookForCluster);
944 expectedNumberOfClusters++;
946 wIndex = (Float_t) t.GetLabel();
947 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
949 AliTRDpropagationLayer& time_bin=*(fTrSec[s]->GetLayer(nr+1));
950 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
951 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
952 if((t.GetSigmaY2() + sy2) < 0) break;
953 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
954 Double_t y=t.GetY(), z=t.GetZ();
958 wYwindow = (Float_t) road;
959 t.GetPxPyPz(Px,Py,Pz);
963 wC = (Float_t) t.GetC();
964 wSigmaC2 = (Float_t) t.GetSigmaC2();
965 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
966 wSigmaY2 = (Float_t) t.GetSigmaY2();
967 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
970 if (road>fWideRoad) {
971 if (t.GetNumberOfClusters()>4)
972 cerr<<t.GetNumberOfClusters()
973 <<"FindProlongation warning: Too broad road !\n";
980 Double_t max_chi2=fMaxChi2;
982 wYclosest = 12345678;
983 wYcorrect = 12345678;
984 wZclosest = 12345678;
985 wZcorrect = 12345678;
986 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
988 // Find the closest correct cluster for debugging purposes
990 Float_t minDY = 1000000;
991 for (Int_t i=0; i<time_bin; i++) {
992 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
993 if((c->GetLabel(0) != trackIndex) &&
994 (c->GetLabel(1) != trackIndex) &&
995 (c->GetLabel(2) != trackIndex)) continue;
996 if(TMath::Abs(c->GetY() - y) > minDY) continue;
997 minDY = TMath::Abs(c->GetY() - y);
998 wYcorrect = c->GetY();
999 wZcorrect = c->GetZ();
1001 Double_t h01 = GetTiltFactor(c);
1002 wChi2 = t.GetPredictedChi2(c, h01);
1006 // Now go for the real cluster search
1010 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
1011 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
1012 if (c->GetY() > y+road) break;
1013 // if (c->IsUsed() > 0) continue;
1014 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1016 Double_t h01 = GetTiltFactor(c);
1017 Double_t chi2=t.GetPredictedChi2(c,h01);
1019 if (chi2 > max_chi2) continue;
1022 index=time_bin.GetIndex(i);
1027 for (Int_t i=time_bin.Find(y-road); i<time_bin; i++) {
1028 AliTRDcluster* c=(AliTRDcluster*)(time_bin[i]);
1030 if (c->GetY() > y+road) break;
1031 // if (c->IsUsed() > 0) continue;
1032 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1034 Double_t h01 = GetTiltFactor(c);
1035 Double_t chi2=t.GetPredictedChi2(c,h01);
1037 if (chi2 > max_chi2) continue;
1040 index=time_bin.GetIndex(i);
1045 wYclosest = cl->GetY();
1046 wZclosest = cl->GetZ();
1048 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1049 Double_t h01 = GetTiltFactor(cl);
1050 if(!t.Update(cl,max_chi2,index,h01)) {
1051 if(!try_again--) return 0;
1053 else try_again=fMaxGap;
1056 if (try_again==0) break;
1061 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1063 printf(" %f", wIndex); //1
1064 printf(" %f", wTB); //2
1065 printf(" %f", wYrt); //3
1066 printf(" %f", wYclosest); //4
1067 printf(" %f", wYcorrect); //5
1068 printf(" %f", wYwindow); //6
1069 printf(" %f", wZrt); //7
1070 printf(" %f", wZclosest); //8
1071 printf(" %f", wZcorrect); //9
1072 printf(" %f", wZwindow); //10
1073 printf(" %f", wPx); //11
1074 printf(" %f", wPy); //12
1075 printf(" %f", wPz); //13
1076 printf(" %f", wSigmaC2*1000000); //14
1077 printf(" %f", wSigmaTgl2*1000); //15
1078 printf(" %f", wSigmaY2); //16
1079 // printf(" %f", wSigmaZ2); //17
1080 printf(" %f", wChi2); //17
1081 printf(" %f", wC); //18
1088 return expectedNumberOfClusters;
1091 //___________________________________________________________________
1093 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1095 // Starting from current radial position of track <t> this function
1096 // extrapolates the track up to radial position <xToGo>.
1097 // Returns 1 if track reaches the plane, and 0 otherwise
1099 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1101 Double_t alpha=t.GetAlpha();
1103 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1104 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1106 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1108 Bool_t lookForCluster;
1109 Double_t rad_length, rho, x, dx, y, ymax, z;
1113 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1115 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1117 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1119 y = t.GetY(); z = t.GetZ();
1121 // first propagate to the outer surface of the current time bin
1122 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
1123 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1124 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1126 ymax = x*TMath::Tan(0.5*alpha);
1129 if (!t.Rotate(alpha)) return 0;
1130 } else if (y <-ymax) {
1132 if (!t.Rotate(-alpha)) return 0;
1134 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1136 y = t.GetY(); z = t.GetZ();
1138 // now propagate to the middle plane of the next time bin
1139 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
1140 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1141 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1143 ymax = x*TMath::Tan(0.5*alpha);
1146 if (!t.Rotate(alpha)) return 0;
1147 } else if (y <-ymax) {
1149 if (!t.Rotate(-alpha)) return 0;
1151 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1156 //___________________________________________________________________
1158 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1160 // Starting from current radial position of track <t> this function
1161 // extrapolates the track up to radial position of the outermost
1162 // padrow of the TPC.
1163 // Returns 1 if track reaches the TPC, and 0 otherwise
1165 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1167 Double_t alpha=t.GetAlpha();
1169 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1170 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1172 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1174 Bool_t lookForCluster;
1175 Double_t rad_length, rho, x, dx, y, ymax, z;
1179 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1181 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1183 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1185 y = t.GetY(); z = t.GetZ();
1187 // first propagate to the outer surface of the current time bin
1188 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
1189 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1190 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1192 ymax = x*TMath::Tan(0.5*alpha);
1195 if (!t.Rotate(alpha)) return 0;
1196 } else if (y <-ymax) {
1198 if (!t.Rotate(-alpha)) return 0;
1200 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1202 y = t.GetY(); z = t.GetZ();
1204 // now propagate to the middle plane of the next time bin
1205 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,rad_length,lookForCluster);
1206 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1207 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1209 ymax = x*TMath::Tan(0.5*alpha);
1212 if (!t.Rotate(alpha)) return 0;
1213 } else if (y <-ymax) {
1215 if (!t.Rotate(-alpha)) return 0;
1217 if(!t.PropagateTo(x,rad_length,rho)) return 0;
1223 //_____________________________________________________________________________
1224 void AliTRDtracker::LoadEvent()
1226 // Fills clusters into TRD tracking_sectors
1227 // Note that the numbering scheme for the TRD tracking_sectors
1228 // differs from that of TRD sectors
1230 ReadClusters(fClusters);
1231 Int_t ncl=fClusters->GetEntriesFast();
1232 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1236 printf("\r %d left ",ncl);
1237 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1238 Int_t detector=c->GetDetector(), local_time_bin=c->GetLocalTimeBin();
1239 Int_t sector=fGeom->GetSector(detector);
1240 Int_t plane=fGeom->GetPlane(detector);
1242 Int_t tracking_sector = CookSectorIndex(sector);
1244 Int_t gtb = fTrSec[tracking_sector]->CookTimeBinIndex(plane,local_time_bin);
1245 if(gtb < 0) continue;
1246 Int_t layer = fTrSec[tracking_sector]->GetLayerNumber(gtb);
1249 fTrSec[tracking_sector]->GetLayer(layer)->InsertCluster(c,index);
1255 //_____________________________________________________________________________
1256 void AliTRDtracker::UnloadEvent()
1259 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1264 nentr = fClusters->GetEntriesFast();
1265 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1267 nentr = fSeeds->GetEntriesFast();
1268 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1270 nentr = fTracks->GetEntriesFast();
1271 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1273 Int_t nsec = AliTRDgeometry::kNsect;
1275 for (i = 0; i < nsec; i++) {
1276 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1277 fTrSec[i]->GetLayer(pl)->Clear();
1283 //__________________________________________________________________________
1284 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1286 // Creates track seeds using clusters in timeBins=i1,i2
1289 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1293 Double_t x[5], c[15];
1294 Int_t max_sec=AliTRDgeometry::kNsect;
1296 Double_t alpha=AliTRDgeometry::GetAlpha();
1297 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1298 Double_t cs=cos(alpha), sn=sin(alpha);
1299 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1302 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1303 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1305 Double_t x1 =fTrSec[0]->GetX(i1);
1306 Double_t xx2=fTrSec[0]->GetX(i2);
1308 for (Int_t ns=0; ns<max_sec; ns++) {
1310 Int_t nl2 = *(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1311 Int_t nl=(*fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1312 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1313 Int_t nu=(*fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1314 Int_t nu2=(*fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1316 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1318 for (Int_t is=0; is < r1; is++) {
1319 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1321 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1323 const AliTRDcluster *cl;
1324 Double_t x2, y2, z2;
1325 Double_t x3=0., y3=0.;
1328 if(turn != 2) continue;
1329 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+max_sec)%max_sec]->GetLayer(i2));
1331 y2=cl->GetY(); z2=cl->GetZ();
1336 else if (js<nl2+nl) {
1337 if(turn != 1) continue;
1338 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+max_sec)%max_sec]->GetLayer(i2));
1340 y2=cl->GetY(); z2=cl->GetZ();
1345 else if (js<nl2+nl+nm) {
1346 if(turn != 1) continue;
1347 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1349 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1351 else if (js<nl2+nl+nm+nu) {
1352 if(turn != 1) continue;
1353 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%max_sec]->GetLayer(i2));
1354 cl=r2[js-nl2-nl-nm];
1355 y2=cl->GetY(); z2=cl->GetZ();
1361 if(turn != 2) continue;
1362 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%max_sec]->GetLayer(i2));
1363 cl=r2[js-nl2-nl-nm-nu];
1364 y2=cl->GetY(); z2=cl->GetZ();
1370 if(TMath::Abs(z1-z2) > fMaxSeedDeltaZ12) continue;
1372 Double_t zz=z1 - z1/x1*(x1-x2);
1374 if (TMath::Abs(zz-z2)>fMaxSeedDeltaZ) continue;
1376 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1377 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1381 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1383 if (TMath::Abs(x[4]) > fMaxSeedC) continue;
1385 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1387 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1389 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1391 if (TMath::Abs(x[3]) > fMaxSeedTan) continue;
1393 Double_t a=asin(x[2]);
1394 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1396 if (TMath::Abs(zv)>fMaxSeedVertexZ) continue;
1398 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1399 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1400 Double_t sy3=fSeedErrorSY3, sy=fSeedErrorSY, sz=fSeedErrorSZ;
1403 Double_t h01 = GetTiltFactor(r1[is]);
1404 Double_t xu_factor = 100.;
1410 sy1=sy1+sz1*h01*h01;
1411 Double_t syz=sz1*(-h01);
1412 // end of tilt changes
1414 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1415 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1416 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1417 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1418 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1419 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1420 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1421 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1422 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1423 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1427 // c[1]=0.; c[2]=sz1;
1428 c[1]=syz; c[2]=sz1*xu_factor;
1429 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1430 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1431 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1432 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1433 c[13]=f30*sy1*f40+f32*sy2*f42;
1434 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1436 UInt_t index=r1.GetIndex(is);
1438 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1440 Int_t rc=FollowProlongation(*track, i2);
1443 (track->GetNumberOfClusters() <
1444 (outer-inner)*fMinClustersInSeed)) delete track;
1446 fSeeds->AddLast(track); fNseeds++;
1447 cerr<<"\r found seed "<<fNseeds;
1454 //_____________________________________________________________________________
1455 void AliTRDtracker::ReadClusters(TObjArray *array, const TFile *inp)
1458 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1459 // from the file. The names of the cluster tree and branches
1460 // should match the ones used in AliTRDclusterizer::WriteClusters()
1463 TDirectory *savedir=gDirectory;
1466 TFile *in=(TFile*)inp;
1467 if (!in->IsOpen()) {
1468 cerr<<"AliTRDtracker::ReadClusters(): input file is not open !\n";
1476 Char_t treeName[12];
1477 sprintf(treeName,"TreeR%d_TRD",fEvent);
1478 TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1480 TObjArray *ClusterArray = new TObjArray(400);
1482 ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray);
1484 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1485 printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1487 // Loop through all entries in the tree
1489 AliTRDcluster *c = 0;
1492 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1495 nbytes += ClusterTree->GetEvent(iEntry);
1497 // Get the number of points in the detector
1498 Int_t nCluster = ClusterArray->GetEntriesFast();
1499 printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1501 // Loop through all TRD digits
1502 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1503 c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1504 AliTRDcluster *co = new AliTRDcluster(*c);
1505 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1506 Int_t ltb = co->GetLocalTimeBin();
1507 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1508 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1510 delete ClusterArray->RemoveAt(iCluster);
1514 delete ClusterArray;
1519 //______________________________________________________________________
1520 void AliTRDtracker::ReadClusters(TObjArray *array, const Char_t *filename)
1523 // Reads AliTRDclusters from file <filename>. The names of the cluster
1524 // tree and branches should match the ones used in
1525 // AliTRDclusterizer::WriteClusters()
1526 // if <array> == 0, clusters are added into AliTRDtracker fCluster array
1529 TDirectory *savedir=gDirectory;
1531 TFile *file = TFile::Open(filename);
1532 if (!file->IsOpen()) {
1533 cerr<<"Can't open file with TRD clusters"<<endl;
1537 Char_t treeName[12];
1538 sprintf(treeName,"TreeR%d_TRD",fEvent);
1539 TTree *ClusterTree = (TTree*) gDirectory->Get(treeName);
1542 cerr<<"AliTRDtracker::ReadClusters(): ";
1543 cerr<<"can't get a tree with clusters !\n";
1547 TObjArray *ClusterArray = new TObjArray(400);
1549 ClusterTree->GetBranch("TRDcluster")->SetAddress(&ClusterArray);
1551 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1552 cout<<"found "<<nEntries<<" in ClusterTree"<<endl;
1554 // Loop through all entries in the tree
1556 AliTRDcluster *c = 0;
1560 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1563 nbytes += ClusterTree->GetEvent(iEntry);
1565 // Get the number of points in the detector
1566 Int_t nCluster = ClusterArray->GetEntriesFast();
1567 printf("\n Read %d clusters from entry %d", nCluster, iEntry);
1569 // Loop through all TRD digits
1570 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1571 c = (AliTRDcluster*)ClusterArray->UncheckedAt(iCluster);
1572 AliTRDcluster *co = new AliTRDcluster(*c);
1573 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1574 Int_t ltb = co->GetLocalTimeBin();
1575 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1576 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1578 delete ClusterArray->RemoveAt(iCluster);
1583 delete ClusterArray;
1589 //__________________________________________________________________
1590 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const {
1592 Int_t label=123456789, index, i, j;
1593 Int_t ncl=pt->GetNumberOfClusters();
1594 const Int_t range = fTrSec[0]->GetOuterTimeBin()+1;
1598 // Int_t s[range][2];
1599 Int_t **s = new Int_t* [range];
1600 for (i=0; i<range; i++) {
1601 s[i] = new Int_t[2];
1603 for (i=0; i<range; i++) {
1609 for (i=0; i<ncl; i++) {
1610 index=pt->GetClusterIndex(i);
1611 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1617 for (i=0; i<ncl; i++) {
1618 index=pt->GetClusterIndex(i);
1619 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1620 for (Int_t k=0; k<3; k++) {
1621 label=c->GetLabel(k);
1622 label_added=kFALSE; j=0;
1624 while ( (!label_added) && ( j < range ) ) {
1625 if (s[j][0]==label || s[j][1]==0) {
1639 for (i=0; i<range; i++) {
1641 max=s[i][1]; label=s[i][0];
1645 for (i=0; i<range; i++) {
1651 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1653 pt->SetLabel(label);
1658 //__________________________________________________________________
1659 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const {
1660 Int_t ncl=t->GetNumberOfClusters();
1661 for (Int_t i=from; i<ncl; i++) {
1662 Int_t index = t->GetClusterIndex(i);
1663 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1669 //_____________________________________________________________________
1670 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t r, Double_t tgl, Double_t pt)
1672 // Parametrised "expected" error of the cluster reconstruction in Y
1674 Double_t s = 0.08 * 0.08;
1678 //_____________________________________________________________________
1679 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t r, Double_t tgl)
1681 // Parametrised "expected" error of the cluster reconstruction in Z
1683 Double_t s = 9 * 9 /12.;
1688 //_____________________________________________________________________
1689 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t local_tb) const
1692 // Returns radial position which corresponds to time bin <local_tb>
1693 // in tracking sector <sector> and plane <plane>
1696 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, local_tb);
1697 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1698 return fTrSec[sector]->GetLayer(pl)->GetX();
1703 //_______________________________________________________
1704 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1705 Double_t dx, Double_t rho, Double_t rad_length, Int_t tb_index)
1708 // AliTRDpropagationLayer constructor
1711 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = rad_length;
1712 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tb_index;
1715 for(Int_t i=0; i < (Int_t) kZONES; i++) {
1716 fZc[i]=0; fZmax[i] = 0;
1721 if(fTimeBinIndex >= 0) {
1722 fClusters = new AliTRDcluster*[kMAX_CLUSTER_PER_TIME_BIN];
1723 fIndex = new UInt_t[kMAX_CLUSTER_PER_TIME_BIN];
1736 //_______________________________________________________
1737 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1738 Double_t Zmax, Double_t Ymax, Double_t rho,
1739 Double_t rad_length, Double_t Yc, Double_t Zc)
1742 // Sets hole in the layer
1751 fHoleX0 = rad_length;
1755 //_______________________________________________________
1756 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1759 // AliTRDtrackingSector Constructor
1768 for(UInt_t i=0; i < kMAX_TIME_BIN_INDEX; i++) fTimeBinIndex[i] = -1;
1771 AliTRDpropagationLayer* ppl;
1773 Double_t x, xin, xout, dx, rho, rad_length;
1776 // set time bins in the gas of the TPC
1778 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
1779 rho = 0.9e-3; rad_length = 28.94;
1781 for(Int_t i=0; i<steps; i++) {
1782 x = xin + i*dx + dx/2;
1783 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1787 // set time bins in the outer field cage vessel
1789 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; rad_length = 44.77; // Tedlar
1790 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1793 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; rad_length = 44.86; // prepreg
1794 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1797 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; rad_length = 41.28; // Nomex
1798 steps = 5; dx = (xout - xin)/steps;
1799 for(Int_t i=0; i<steps; i++) {
1800 x = xin + i*dx + dx/2;
1801 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1805 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; rad_length = 44.86; // prepreg
1806 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1809 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; rad_length = 44.77; // Tedlar
1810 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1814 // set time bins in CO2
1816 xin = xout; xout = 275.0;
1817 steps = 50; dx = (xout - xin)/steps;
1818 rho = 1.977e-3; rad_length = 36.2;
1820 for(Int_t i=0; i<steps; i++) {
1821 x = xin + i*dx + dx/2;
1822 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1826 // set time bins in the outer containment vessel
1828 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; rad_length = 24.01; // Al
1829 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1832 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; rad_length = 44.77; // Tedlar
1833 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1836 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; rad_length = 44.86; // prepreg
1837 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1840 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; rad_length = 41.28; // Nomex
1841 steps = 10; dx = (xout - xin)/steps;
1842 for(Int_t i=0; i<steps; i++) {
1843 x = xin + i*dx + dx/2;
1844 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1848 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; rad_length = 44.86; // prepreg
1849 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1852 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; rad_length = 44.77; // Tedlar
1853 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1856 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; rad_length = 24.01; // Al
1857 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
1860 Double_t xtrd = (Double_t) fGeom->Rmin();
1862 // add layers between TPC and TRD (Air temporarily)
1863 xin = xout; xout = xtrd;
1864 steps = 50; dx = (xout - xin)/steps;
1865 rho = 1.2e-3; rad_length = 36.66;
1867 for(Int_t i=0; i<steps; i++) {
1868 x = xin + i*dx + dx/2;
1869 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1874 Double_t alpha=AliTRDgeometry::GetAlpha();
1876 // add layers for each of the planes
1878 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
1879 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
1880 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
1881 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
1882 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
1883 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
1884 Double_t dxPlane = dxTEC + dxSpace;
1887 const Int_t nChambers = AliTRDgeometry::Ncham();
1888 Double_t Ymax = 0, holeYmax = 0;
1889 Double_t * Zc = new Double_t[nChambers];
1890 Double_t * Zmax = new Double_t[nChambers];
1891 Double_t holeZmax = 1000.; // the whole sector is missing
1893 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
1896 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
1897 steps = 12; dx = (xout - xin)/steps; rho = 0.074; rad_length = 40.6;
1898 for(Int_t i=0; i<steps; i++) {
1899 x = xin + i*dx + dx/2;
1900 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
1901 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1902 holeYmax = x*TMath::Tan(0.5*alpha);
1903 ppl->SetHole(holeYmax, holeZmax);
1905 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1906 holeYmax = x*TMath::Tan(0.5*alpha);
1907 ppl->SetHole(holeYmax, holeZmax);
1912 Ymax = fGeom->GetChamberWidth(plane)/2;
1913 for(Int_t ch = 0; ch < nChambers; ch++) {
1914 Zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
1915 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
1916 Float_t row0 = fPar->GetRow0(plane,ch,0);
1917 Int_t nPads = fPar->GetRowMax(plane,ch,0);
1918 Zc[ch] = (pad * nPads)/2 + row0 - pad/2;
1921 dx = fPar->GetTimeBinSize();
1922 rho = 0.00295 * 0.85; rad_length = 11.0;
1924 Double_t x0 = (Double_t) fPar->GetTime0(plane);
1925 Double_t xbottom = x0 - dxDrift;
1926 Double_t xtop = x0 + dxAmp;
1928 // Amplification region
1930 steps = (Int_t) (dxAmp/dx);
1932 for(tb = 0; tb < steps; tb++) {
1933 x = x0 + tb * dx + dx/2;
1934 tb_index = CookTimeBinIndex(plane, -tb-1);
1935 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,tb_index);
1937 for(Int_t ch = 0; ch < nChambers; ch++) {
1938 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1940 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1941 holeYmax = x*TMath::Tan(0.5*alpha);
1942 ppl->SetHole(holeYmax, holeZmax);
1944 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1945 holeYmax = x*TMath::Tan(0.5*alpha);
1946 ppl->SetHole(holeYmax, holeZmax);
1950 tb_index = CookTimeBinIndex(plane, -steps);
1951 x = (x + dx/2 + xtop)/2;
1953 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,tb_index);
1955 for(Int_t ch = 0; ch < nChambers; ch++) {
1956 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1958 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1959 holeYmax = x*TMath::Tan(0.5*alpha);
1960 ppl->SetHole(holeYmax, holeZmax);
1962 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1963 holeYmax = x*TMath::Tan(0.5*alpha);
1964 ppl->SetHole(holeYmax, holeZmax);
1969 dx = fPar->GetTimeBinSize();
1970 steps = (Int_t) (dxDrift/dx);
1972 for(tb = 0; tb < steps; tb++) {
1973 x = x0 - tb * dx - dx/2;
1974 tb_index = CookTimeBinIndex(plane, tb);
1976 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,tb_index);
1978 for(Int_t ch = 0; ch < nChambers; ch++) {
1979 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1981 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
1982 holeYmax = x*TMath::Tan(0.5*alpha);
1983 ppl->SetHole(holeYmax, holeZmax);
1985 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
1986 holeYmax = x*TMath::Tan(0.5*alpha);
1987 ppl->SetHole(holeYmax, holeZmax);
1991 tb_index = CookTimeBinIndex(plane, steps);
1992 x = (x - dx/2 + xbottom)/2;
1994 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,tb_index);
1996 for(Int_t ch = 0; ch < nChambers; ch++) {
1997 ppl->SetZmax(ch, Zc[ch], Zmax[ch]);
1999 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
2000 holeYmax = x*TMath::Tan(0.5*alpha);
2001 ppl->SetHole(holeYmax, holeZmax);
2003 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2004 holeYmax = x*TMath::Tan(0.5*alpha);
2005 ppl->SetHole(holeYmax, holeZmax);
2010 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; rad_length = 33.0;
2011 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,rad_length,-1);
2012 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
2013 holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
2014 ppl->SetHole(holeYmax, holeZmax);
2016 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2017 holeYmax = (xin+dx/2)*TMath::Tan(0.5*alpha);
2018 ppl->SetHole(holeYmax, holeZmax);
2023 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2024 steps = 5; dx = (xout - xin)/steps; rho = 0.074; rad_length = 40.6;
2025 for(Int_t i=0; i<steps; i++) {
2026 x = xin + i*dx + dx/2;
2027 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
2028 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
2029 holeYmax = x*TMath::Tan(0.5*alpha);
2030 ppl->SetHole(holeYmax, holeZmax);
2032 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2033 holeYmax = x*TMath::Tan(0.5*alpha);
2034 ppl->SetHole(holeYmax, holeZmax);
2039 // Space between the chambers, air
2040 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2041 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; rad_length = 36.66;
2042 for(Int_t i=0; i<steps; i++) {
2043 x = xin + i*dx + dx/2;
2044 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
2045 if((fGeom->GetPHOShole()) && (fGeomSector >= 2) && (fGeomSector <= 6)) {
2046 holeYmax = x*TMath::Tan(0.5*alpha);
2047 ppl->SetHole(holeYmax, holeZmax);
2049 if((fGeom->GetRICHhole()) && (fGeomSector >= 12) && (fGeomSector <= 14)) {
2050 holeYmax = x*TMath::Tan(0.5*alpha);
2051 ppl->SetHole(holeYmax, holeZmax);
2057 // Space between the TRD and RICH
2058 Double_t xRICH = 500.;
2059 xin = xout; xout = xRICH;
2060 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; rad_length = 36.66;
2061 for(Int_t i=0; i<steps; i++) {
2062 x = xin + i*dx + dx/2;
2063 ppl = new AliTRDpropagationLayer(x,dx,rho,rad_length,-1);
2073 //______________________________________________________
2075 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t local_tb) const
2078 // depending on the digitization parameters calculates "global"
2079 // time bin index for timebin <local_tb> in plane <plane>
2082 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2083 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2084 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2086 Int_t tbAmp = fPar->GetTimeBefore();
2087 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2088 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2089 Int_t tbDrift = fPar->GetTimeMax();
2090 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2092 Int_t tb_per_plane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2094 Int_t gtb = (plane+1) * tb_per_plane - local_tb - 1 - TMath::Min(tbAmp,maxAmp);
2096 if((local_tb < 0) &&
2097 (TMath::Abs(local_tb) > TMath::Min(tbAmp,maxAmp))) return -1;
2098 if(local_tb >= TMath::Min(tbDrift,maxDrift)) return -1;
2105 //______________________________________________________
2107 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2110 // For all sensitive time bins sets corresponding layer index
2111 // in the array fTimeBins
2116 for(Int_t i = 0; i < fN; i++) {
2117 index = fLayers[i]->GetTimeBinIndex();
2119 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2121 if(index < 0) continue;
2122 if(index >= (Int_t) kMAX_TIME_BIN_INDEX) {
2123 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2124 printf(" index %d exceeds allowed maximum of %d!\n",
2125 index, kMAX_TIME_BIN_INDEX-1);
2128 fTimeBinIndex[index] = i;
2131 Double_t x1, dx1, x2, dx2, gap;
2133 for(Int_t i = 0; i < fN-1; i++) {
2134 x1 = fLayers[i]->GetX();
2135 dx1 = fLayers[i]->GetdX();
2136 x2 = fLayers[i+1]->GetX();
2137 dx2 = fLayers[i+1]->GetdX();
2138 gap = (x2 - dx2/2) - (x1 + dx1/2);
2140 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2141 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2144 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2145 printf(" (%f - %f) - (%f + %f) = %f\n",
2146 x2, dx2/2, x1, dx1, gap);
2152 //______________________________________________________
2155 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2158 // Returns the number of time bin which in radial position is closest to <x>
2161 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2162 if(x <= fLayers[0]->GetX()) return 0;
2164 Int_t b=0, e=fN-1, m=(b+e)/2;
2165 for (; b<e; m=(b+e)/2) {
2166 if (x > fLayers[m]->GetX()) b=m+1;
2169 if(TMath::Abs(x - fLayers[m]->GetX()) >
2170 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2175 //______________________________________________________
2177 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2180 // Returns number of the innermost SENSITIVE propagation layer
2183 return GetLayerNumber(0);
2186 //______________________________________________________
2188 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2191 // Returns number of the outermost SENSITIVE time bin
2194 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2197 //______________________________________________________
2199 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2202 // Returns number of SENSITIVE time bins
2206 for(tb = kMAX_TIME_BIN_INDEX-1; tb >=0; tb--) {
2207 layer = GetLayerNumber(tb);
2213 //______________________________________________________
2215 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2218 // Insert layer <pl> in fLayers array.
2219 // Layers are sorted according to X coordinate.
2221 if ( fN == ((Int_t) kMAX_LAYERS_PER_SECTOR)) {
2222 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2225 if (fN==0) {fLayers[fN++] = pl; return;}
2226 Int_t i=Find(pl->GetX());
2228 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2229 fLayers[i]=pl; fN++;
2233 //______________________________________________________
2235 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2238 // Returns index of the propagation layer nearest to X
2241 if (x <= fLayers[0]->GetX()) return 0;
2242 if (x > fLayers[fN-1]->GetX()) return fN;
2243 Int_t b=0, e=fN-1, m=(b+e)/2;
2244 for (; b<e; m=(b+e)/2) {
2245 if (x > fLayers[m]->GetX()) b=m+1;
2251 //______________________________________________________
2253 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2254 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &rad_length,
2255 Bool_t &lookForCluster) const
2258 // Returns radial step <dx>, density <rho>, rad. length <rad_length>,
2259 // and sensitivity <lookForCluster> in point <y,z>
2265 lookForCluster = kFALSE;
2267 // check dead regions
2268 if(fTimeBinIndex >= 0) {
2269 for(Int_t ch = 0; ch < (Int_t) kZONES; ch++) {
2270 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2271 lookForCluster = kTRUE;
2272 // else { rho = 1.7; rad_length = 33.0; } // G10
2274 if(TMath::Abs(y) > fYmax) lookForCluster = kFALSE;
2275 if(!lookForCluster) {
2276 // rho = 1.7; rad_length = 33.0; // G10
2281 if(fHole && (TMath::Abs(y - fHoleYc) < fHoleYmax) &&
2282 (TMath::Abs(z - fHoleZc) < fHoleZmax)) {
2283 lookForCluster = kFALSE;
2285 rad_length = fHoleX0;
2291 //______________________________________________________
2293 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2296 // Insert cluster in cluster array.
2297 // Clusters are sorted according to Y coordinate.
2299 if(fTimeBinIndex < 0) {
2300 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2304 if (fN== (Int_t) kMAX_CLUSTER_PER_TIME_BIN) {
2305 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2308 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2309 Int_t i=Find(c->GetY());
2310 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2311 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2312 fIndex[i]=index; fClusters[i]=c; fN++;
2315 //______________________________________________________
2317 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2319 // Returns index of the cluster nearest in Y
2321 if (y <= fClusters[0]->GetY()) return 0;
2322 if (y > fClusters[fN-1]->GetY()) return fN;
2323 Int_t b=0, e=fN-1, m=(b+e)/2;
2324 for (; b<e; m=(b+e)/2) {
2325 if (y > fClusters[m]->GetY()) b=m+1;
2331 //---------------------------------------------------------
2333 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2335 // Returns correction factor for tilted pads geometry
2338 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2339 Int_t det = c->GetDetector();
2340 Int_t plane = fGeom->GetPlane(det);
2342 if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2344 if(fNoTilt) h01 = 0;