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 ///////////////////////////////////////////////////////////////////////////////
20 // The standard TRD tracker //
22 ///////////////////////////////////////////////////////////////////////////////
24 #include <Riostream.h>
28 #include <TObjArray.h>
30 #include "AliTRDgeometry.h"
31 #include "AliTRDparameter.h"
32 #include "AliTRDgeometryHole.h"
33 #include "AliTRDcluster.h"
34 #include "AliTRDtrack.h"
35 #include "AliBarrelTrack.h"
38 #include "AliTRDtracker.h"
40 ClassImp(AliTRDtracker)
42 const Float_t AliTRDtracker::fgkSeedDepth = 0.5;
43 const Float_t AliTRDtracker::fgkSeedStep = 0.10;
44 const Float_t AliTRDtracker::fgkSeedGap = 0.25;
46 const Float_t AliTRDtracker::fgkMaxSeedDeltaZ12 = 40.;
47 const Float_t AliTRDtracker::fgkMaxSeedDeltaZ = 25.;
48 const Float_t AliTRDtracker::fgkMaxSeedC = 0.0052;
49 const Float_t AliTRDtracker::fgkMaxSeedTan = 1.2;
50 const Float_t AliTRDtracker::fgkMaxSeedVertexZ = 150.;
52 const Double_t AliTRDtracker::fgkSeedErrorSY = 0.2;
53 const Double_t AliTRDtracker::fgkSeedErrorSY3 = 2.5;
54 const Double_t AliTRDtracker::fgkSeedErrorSZ = 0.1;
56 const Float_t AliTRDtracker::fgkMinClustersInSeed = 0.7;
58 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5;
59 const Float_t AliTRDtracker::fgkMinFractionOfFoundClusters = 0.8;
61 const Float_t AliTRDtracker::fgkSkipDepth = 0.3;
62 const Float_t AliTRDtracker::fgkLabelFraction = 0.8;
63 const Float_t AliTRDtracker::fgkWideRoad = 20.;
65 const Double_t AliTRDtracker::fgkMaxChi2 = 12.;
67 const Int_t AliTRDtracker::fgkFirstPlane = 5;
68 const Int_t AliTRDtracker::fgkLastPlane = 17;
71 //____________________________________________________________________
72 AliTRDtracker::AliTRDtracker():AliTracker(),
89 // Default constructor
91 for(Int_t i=0;i<kTrackingSectors;i++) fTrSec[i]=0;
92 for(Int_t j=0;j<5;j++)
93 for(Int_t k=0;k<18;k++) fHoles[j][k]=kFALSE;
95 //____________________________________________________________________
96 AliTRDtracker::AliTRDtracker(const TFile *geomfile):AliTracker()
102 //Float_t fTzero = 0;
104 fAddTRDseeds = kFALSE;
108 TDirectory *savedir=gDirectory;
109 TFile *in=(TFile*)geomfile;
111 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
112 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
117 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
118 fPar = (AliTRDparameter*) in->Get("TRDparameter");
123 // fTzero = geo->GetT0();
124 printf("Found geometry version %d on file \n", fGeom->IsVersion());
127 printf("AliTRDtracker::AliTRDtracker(): can't find TRD geometry!\n");
128 //printf("The DETAIL TRD geometry will be used\n");
129 //fGeom = new AliTRDgeometryDetail();
130 fGeom = new AliTRDgeometryHole();
131 fGeom->SetPHOShole();
132 fGeom->SetRICHhole();
136 printf("AliTRDtracker::AliTRDtracker(): can't find TRD parameter!\n");
137 printf("The DEFAULT TRD parameter will be used\n");
138 fPar = new AliTRDparameter();
145 // fGeom->SetT0(fTzero);
148 fClusters = new TObjArray(2000);
150 fSeeds = new TObjArray(2000);
152 fTracks = new TObjArray(1000);
154 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
155 Int_t trS = CookSectorIndex(geomS);
156 fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS, fPar);
157 for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
158 fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
162 Float_t tiltAngle = TMath::Abs(fPar->GetTiltingAngle());
163 if(tiltAngle < 0.1) {
170 if(fNoTilt && (tiltAngle > 0.1)) fSY2corr = fSY2corr + tiltAngle * 0.05;
173 // calculate max gap on track
175 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
176 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
178 Double_t dx = (Double_t) fPar->GetTimeBinSize();
179 Int_t tbAmp = fPar->GetTimeBefore();
180 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
181 if(kTRUE) maxAmp = 0; // intentional until we change the parameter class
182 Int_t tbDrift = fPar->GetTimeMax();
183 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
185 tbDrift = TMath::Min(tbDrift,maxDrift);
186 tbAmp = TMath::Min(tbAmp,maxAmp);
188 fTimeBinsPerPlane = tbAmp + tbDrift;
189 fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fgkSkipDepth);
194 // Barrel Tracks [SR, 03.04.2003]
204 //___________________________________________________________________
205 AliTRDtracker::~AliTRDtracker()
208 // Destructor of AliTRDtracker
226 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
227 delete fTrSec[geomS];
231 //_____________________________________________________________________
233 void AliTRDtracker::SetBarrelTree(const char *mode) {
238 if (!IsStoringBarrel()) return;
240 TDirectory *sav = gDirectory;
241 if (!fBarrelFile) fBarrelFile = new TFile("AliBarrelTracks.root", "UPDATE");
244 sprintf(buff, "BarrelTRD_%d_%s", GetEventNumber(), mode);
247 fBarrelTree = new TTree(buff, "Barrel TPC tracks");
249 Int_t nRefs = fgkLastPlane - fgkFirstPlane + 1;
251 if (!fBarrelArray) fBarrelArray = new TClonesArray("AliBarrelTrack", nRefs);
252 for(Int_t i=0; i<nRefs; i++) new((*fBarrelArray)[i]) AliBarrelTrack();
254 fBarrelTree->Branch("tracks", &fBarrelArray);
258 //_____________________________________________________________________
260 void AliTRDtracker::StoreBarrelTrack(AliTRDtrack *ps, Int_t refPlane, Int_t isIn) {
265 if (!IsStoringBarrel()) return;
267 static Int_t nClusters;
269 static Double_t chi2;
271 static Bool_t wasLast = kTRUE;
273 Int_t newClusters, newWrong;
278 fBarrelArray->Clear();
279 nClusters = nWrong = 0;
285 fBarrelTrack = (AliBarrelTrack*)(*fBarrelArray)[index++];
286 ps->GetBarrelTrack(fBarrelTrack);
288 newClusters = ps->GetNumberOfClusters() - nClusters;
289 newWrong = ps->GetNWrong() - nWrong;
290 newChi2 = ps->GetChi2() - chi2;
292 nClusters = ps->GetNumberOfClusters();
293 nWrong = ps->GetNWrong();
294 chi2 = ps->GetChi2();
296 if (refPlane != fgkLastPlane) {
297 fBarrelTrack->SetNClusters(newClusters, newChi2);
298 fBarrelTrack->SetNWrongClusters(newWrong);
303 fBarrelTrack->SetRefPlane(refPlane, isIn);
306 //_____________________________________________________________________
308 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track) {
310 // Rotates the track when necessary
313 Double_t alpha = AliTRDgeometry::GetAlpha();
314 Double_t y = track->GetY();
315 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
317 //Int_t ns = AliTRDgeometry::kNsect;
318 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
322 if (!track->Rotate(alpha)) return kFALSE;
323 } else if (y <-ymax) {
325 if (!track->Rotate(-alpha)) return kFALSE;
331 //_____________________________________________________________________
332 inline Double_t f1trd(Double_t x1,Double_t y1,
333 Double_t x2,Double_t y2,
334 Double_t x3,Double_t y3)
337 // Initial approximation of the track curvature
339 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
340 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
341 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
342 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
343 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
345 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
347 return -xr*yr/sqrt(xr*xr+yr*yr);
350 //_____________________________________________________________________
351 inline Double_t f2trd(Double_t x1,Double_t y1,
352 Double_t x2,Double_t y2,
353 Double_t x3,Double_t y3)
356 // Initial approximation of the track curvature times X coordinate
357 // of the center of curvature
360 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
361 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
362 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
363 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
364 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
366 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
368 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
371 //_____________________________________________________________________
372 inline Double_t f3trd(Double_t x1,Double_t y1,
373 Double_t x2,Double_t y2,
374 Double_t z1,Double_t z2)
377 // Initial approximation of the tangent of the track dip angle
380 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
384 AliTRDcluster * AliTRDtracker::GetCluster(AliTRDtrack * track, Int_t plane, Int_t timebin){
386 //try to find cluster in the backup list
388 AliTRDcluster * cl =0;
389 UInt_t *indexes = track->GetBackupIndexes();
390 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
391 if (indexes[i]==0) break;
392 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
394 if (cli->GetLocalTimeBin()!=timebin) continue;
395 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
405 Int_t AliTRDtracker::GetLastPlane(AliTRDtrack * track){
407 //return last updated plane
409 UInt_t *indexes = track->GetBackupIndexes();
410 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
411 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
413 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
414 if (iplane>lastplane) {
420 //___________________________________________________________________
421 Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
424 // Finds tracks within the TRD. The ESD event is expected to contain seeds
425 // at the outer part of the TRD. The seeds
426 // are found within the TRD if fAddTRDseeds is TRUE.
427 // The tracks are propagated to the innermost time bin
428 // of the TRD and the ESD event is updated
431 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
432 Float_t foundMin = fgkMinClustersInTrack * timeBins;
435 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
437 Int_t n = event->GetNumberOfTracks();
438 for (Int_t i=0; i<n; i++) {
439 AliESDtrack* seed=event->GetTrack(i);
440 ULong_t status=seed->GetStatus();
441 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
442 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
445 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
446 //seed2->ResetCovariance();
447 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
449 FollowProlongation(t, innerTB);
450 if (t.GetNumberOfClusters() >= foundMin) {
452 CookLabel(pt, 1-fgkLabelFraction);
456 // cout<<found<<'\r';
458 if(PropagateToTPC(t)) {
459 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
465 cout<<"Number of loaded seeds: "<<nseed<<endl;
466 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
468 // after tracks from loaded seeds are found and the corresponding
469 // clusters are used, look for additional seeds from TRD
472 // Find tracks for the seeds in the TRD
473 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
475 Int_t nSteps = (Int_t) (fgkSeedDepth / fgkSeedStep);
476 Int_t gap = (Int_t) (timeBins * fgkSeedGap);
477 Int_t step = (Int_t) (timeBins * fgkSeedStep);
479 // make a first turn with tight cut on initial curvature
480 for(Int_t turn = 1; turn <= 2; turn++) {
482 nSteps = (Int_t) (fgkSeedDepth / (3*fgkSeedStep));
483 step = (Int_t) (timeBins * (3*fgkSeedStep));
485 for(Int_t i=0; i<nSteps; i++) {
486 Int_t outer=timeBins-1-i*step;
487 Int_t inner=outer-gap;
489 nseed=fSeeds->GetEntriesFast();
491 MakeSeeds(inner, outer, turn);
493 nseed=fSeeds->GetEntriesFast();
494 // printf("\n turn %d, step %d: number of seeds for TRD inward %d\n",
497 for (Int_t i=0; i<nseed; i++) {
498 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
499 FollowProlongation(t,innerTB);
500 if (t.GetNumberOfClusters() >= foundMin) {
502 CookLabel(pt, 1-fgkLabelFraction);
505 // cout<<found<<'\r';
506 if(PropagateToTPC(t)) {
508 track.UpdateTrackParams(pt,AliESDtrack::kTRDin);
509 event->AddTrack(&track);
510 // track.SetTRDtrack(new AliTRDtrack(*pt));
513 delete fSeeds->RemoveAt(i);
520 cout<<"Total number of found tracks: "<<found<<endl;
527 //_____________________________________________________________________________
528 Int_t AliTRDtracker::PropagateBack(AliESD* event) {
530 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
531 // backpropagated by the TPC tracker. Each seed is first propagated
532 // to the TRD, and then its prolongation is searched in the TRD.
533 // If sufficiently long continuation of the track is found in the TRD
534 // the track is updated, otherwise it's stored as originaly defined
535 // by the TPC tracker.
539 Float_t foundMin = 20;
541 Int_t n = event->GetNumberOfTracks();
542 for (Int_t i=0; i<n; i++) {
543 AliESDtrack* seed=event->GetTrack(i);
544 ULong_t status=seed->GetStatus();
545 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
546 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
548 Int_t lbl = seed->GetLabel();
549 AliTRDtrack *track = new AliTRDtrack(*seed);
550 track->SetSeedLabel(lbl);
552 Float_t p4 = track->GetC();
554 Int_t expectedClr = FollowBackProlongation(*track);
555 if (track->GetNumberOfClusters()<expectedClr/3){
556 AliTRDtrack *track1 = new AliTRDtrack(*seed);
557 track1->SetSeedLabel(lbl);
558 FollowBackProlongation(*track1);
559 AliTRDtrack *track2= new AliTRDtrack(*seed);
560 track->SetSeedLabel(lbl);
561 FollowBackProlongation(*track2);
565 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)>0.2) {
567 continue; //too big change of curvature - to be checked
569 Int_t foundClr = track->GetNumberOfClusters();
570 if (foundClr >= foundMin) {
573 CookLabel(track, 1-fgkLabelFraction);
577 // Propagate to outer reference plane [SR, GSI, 18.02.2003]
578 // track->PropagateTo(364.8); why?
580 //seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
582 if (track->GetNCross()==0) seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
584 if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
588 //Propagation to the TOF (I.Belikov)
590 if (track->GetStop()==kFALSE){
593 Double_t c2=track->GetC()*xtof - track->GetEta();
594 if (TMath::Abs(c2)>=0.85) {
598 Double_t xTOF0 = 371. ;
599 PropagateToOuterPlane(*track,xTOF0);
601 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
602 Double_t y=track->GetYat(xtof);
604 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
608 } else if (y <-ymax) {
609 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
615 if (track->PropagateTo(xtof)) {
616 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
617 seed->SetTRDtrack(new AliTRDtrack(*track));
618 if (track->GetNumberOfClusters()>foundMin) found++;
621 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
622 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
623 //seed->SetStatus(AliESDtrack::kTRDStop);
624 seed->SetTRDtrack(new AliTRDtrack(*track));
631 //End of propagation to the TOF
632 //if (foundClr>foundMin)
633 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
638 cerr<<"Number of seeds: "<<fNseeds<<endl;
639 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
641 fSeeds->Clear(); fNseeds=0;
647 //_____________________________________________________________________________
648 Int_t AliTRDtracker::RefitInward(AliESD* event)
651 // Refits tracks within the TRD. The ESD event is expected to contain seeds
652 // at the outer part of the TRD.
653 // The tracks are propagated to the innermost time bin
654 // of the TRD and the ESD event is updated
655 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
658 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
659 Float_t foundMin = fgkMinClustersInTrack * timeBins;
662 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
664 Int_t n = event->GetNumberOfTracks();
665 for (Int_t i=0; i<n; i++) {
666 AliESDtrack* seed=event->GetTrack(i);
667 ULong_t status=seed->GetStatus();
668 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
669 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
672 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
673 seed2->ResetCovariance(5.);
674 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
675 UInt_t * indexes2 = seed2->GetIndexes();
676 UInt_t * indexes3 = pt->GetBackupIndexes();
677 for (Int_t i=0;i<200;i++) {
678 if (indexes2[i]==0) break;
679 indexes3[i] = indexes2[i];
681 //AliTRDtrack *pt = seed2;
683 FollowProlongation(t, innerTB);
685 if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
686 // debug - why we dont go back?
687 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
688 UInt_t * indexes2 = seed2->GetIndexes();
689 UInt_t * indexes3 = pt2->GetBackupIndexes();
690 for (Int_t i=0;i<200;i++) {
691 if (indexes2[i]==0) break;
692 indexes3[i] = indexes2[i];
694 FollowProlongation(*pt2, innerTB);
698 if (t.GetNumberOfClusters() >= foundMin) {
700 //CookLabel(pt, 1-fgkLabelFraction);
704 // cout<<found<<'\r';
706 if(PropagateToTPC(t)) {
707 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
709 //if not prolongation to TPC - propagate without update
710 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
711 seed2->ResetCovariance(5.);
712 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
714 if (PropagateToTPC(*pt2))
715 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
723 cout<<"Number of loaded seeds: "<<nseed<<endl;
724 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
731 //---------------------------------------------------------------------------
732 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
734 // Starting from current position on track=t this function tries
735 // to extrapolate the track up to timeBin=0 and to confirm prolongation
736 // if a close cluster is found. Returns the number of clusters
737 // expected to be found in sensitive layers
739 Float_t wIndex, wTB, wChi2;
740 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
741 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
742 Float_t wPx, wPy, wPz, wC;
744 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
745 Int_t lastplane = GetLastPlane(&t);
747 Int_t trackIndex = t.GetLabel();
749 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
751 Int_t tryAgain=fMaxGap;
753 Double_t alpha=t.GetAlpha();
754 alpha = TVector2::Phi_0_2pi(alpha);
756 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
757 Double_t radLength, rho, x, dx, y, ymax, z;
759 Int_t expectedNumberOfClusters = 0;
760 Bool_t lookForCluster;
762 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
765 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
767 y = t.GetY(); z = t.GetZ();
769 // first propagate to the inner surface of the current time bin
770 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
771 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
772 if(!t.PropagateTo(x,radLength,rho)) break;
774 ymax = x*TMath::Tan(0.5*alpha);
777 if (!t.Rotate(alpha)) break;
778 if(!t.PropagateTo(x,radLength,rho)) break;
779 } else if (y <-ymax) {
781 if (!t.Rotate(-alpha)) break;
782 if(!t.PropagateTo(x,radLength,rho)) break;
785 y = t.GetY(); z = t.GetZ();
787 // now propagate to the middle plane of the next time bin
788 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
789 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
790 if(!t.PropagateTo(x,radLength,rho)) break;
792 ymax = x*TMath::Tan(0.5*alpha);
795 if (!t.Rotate(alpha)) break;
796 if(!t.PropagateTo(x,radLength,rho)) break;
797 } else if (y <-ymax) {
799 if (!t.Rotate(-alpha)) break;
800 if(!t.PropagateTo(x,radLength,rho)) break;
806 expectedNumberOfClusters++;
807 wIndex = (Float_t) t.GetLabel();
810 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
812 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
813 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
816 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
817 else return expectedNumberOfClusters;
821 wYwindow = (Float_t) road;
822 t.GetPxPyPz(px,py,pz);
826 wC = (Float_t) t.GetC();
827 wSigmaC2 = (Float_t) t.GetSigmaC2();
828 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
829 wSigmaY2 = (Float_t) t.GetSigmaY2();
830 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
837 Double_t maxChi2=fgkMaxChi2;
839 wYclosest = 12345678;
840 wYcorrect = 12345678;
841 wZclosest = 12345678;
842 wZcorrect = 12345678;
843 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
845 // Find the closest correct cluster for debugging purposes
847 Float_t minDY = 1000000;
848 for (Int_t i=0; i<timeBin; i++) {
849 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
850 if((c->GetLabel(0) != trackIndex) &&
851 (c->GetLabel(1) != trackIndex) &&
852 (c->GetLabel(2) != trackIndex)) continue;
853 if(TMath::Abs(c->GetY() - y) > minDY) continue;
854 minDY = TMath::Abs(c->GetY() - y);
855 wYcorrect = c->GetY();
856 wZcorrect = c->GetZ();
858 Double_t h01 = GetTiltFactor(c);
859 wChi2 = t.GetPredictedChi2(c, h01);
863 // Now go for the real cluster search
867 //find cluster in history
870 AliTRDcluster * cl0 = timeBin[0];
874 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
875 if (plane>lastplane) continue;
876 Int_t timebin = cl0->GetLocalTimeBin();
877 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
880 Double_t h01 = GetTiltFactor(cl);
881 maxChi2=t.GetPredictedChi2(cl,h01);
883 if ((!cl) && road>fgkWideRoad) {
884 //if (t.GetNumberOfClusters()>4)
885 // cerr<<t.GetNumberOfClusters()
886 // <<"FindProlongation warning: Too broad road !\n";
893 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
894 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
895 if (c->GetY() > y+road) break;
896 if (c->IsUsed() > 0) continue;
897 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
899 Double_t h01 = GetTiltFactor(c);
900 Double_t chi2=t.GetPredictedChi2(c,h01);
902 if (chi2 > maxChi2) continue;
905 index=timeBin.GetIndex(i);
911 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
912 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
914 if (c->GetY() > y+road) break;
915 if (c->IsUsed() > 0) continue;
916 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
918 Double_t h01 = GetTiltFactor(c);
919 Double_t chi2=t.GetPredictedChi2(c, h01);
921 if (chi2 > maxChi2) continue;
924 index=timeBin.GetIndex(i);
928 wYclosest = cl->GetY();
929 wZclosest = cl->GetZ();
930 Double_t h01 = GetTiltFactor(cl);
932 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
933 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
934 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
935 Int_t det = cl->GetDetector();
936 Int_t plane = fGeom->GetPlane(det);
938 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
939 //if(!t.Update(cl,maxChi2,index,h01)) {
940 //if(!tryAgain--) return 0;
942 else tryAgain=fMaxGap;
945 //if (tryAgain==0) break;
950 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
952 printf(" %f", wIndex); //1
953 printf(" %f", wTB); //2
954 printf(" %f", wYrt); //3
955 printf(" %f", wYclosest); //4
956 printf(" %f", wYcorrect); //5
957 printf(" %f", wYwindow); //6
958 printf(" %f", wZrt); //7
959 printf(" %f", wZclosest); //8
960 printf(" %f", wZcorrect); //9
961 printf(" %f", wZwindow); //10
962 printf(" %f", wPx); //11
963 printf(" %f", wPy); //12
964 printf(" %f", wPz); //13
965 printf(" %f", wSigmaC2*1000000); //14
966 printf(" %f", wSigmaTgl2*1000); //15
967 printf(" %f", wSigmaY2); //16
968 // printf(" %f", wSigmaZ2); //17
969 printf(" %f", wChi2); //17
970 printf(" %f", wC); //18
977 return expectedNumberOfClusters;
982 //___________________________________________________________________
984 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
986 // Starting from current radial position of track <t> this function
987 // extrapolates the track up to outer timebin and in the sensitive
988 // layers confirms prolongation if a close cluster is found.
989 // Returns the number of clusters expected to be found in sensitive layers
992 Float_t wIndex, wTB, wChi2;
993 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
994 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
995 Float_t wPx, wPy, wPz, wC;
997 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
999 Int_t trackIndex = t.GetLabel();
1000 Int_t tryAgain=fMaxGap;
1002 Double_t alpha=t.GetAlpha();
1003 TVector2::Phi_0_2pi(alpha);
1007 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
1008 Double_t radLength, rho, x, dx, y, ymax = 0, z;
1009 Bool_t lookForCluster;
1011 Int_t expectedNumberOfClusters = 0;
1014 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1016 Int_t nRefPlane = fgkFirstPlane;
1017 Bool_t isNewLayer = kFALSE;
1023 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
1028 // first propagate to the outer surface of the current time bin
1031 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1032 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
1036 if(!t.PropagateTo(x,radLength,rho)) break;
1037 // if (!AdjustSector(&t)) break;
1039 // MI -fix untill correct material desription will be implemented
1041 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material
1042 if (!AdjustSector(&t)) break;
1043 Int_t cross = kFALSE;
1045 if (TMath::Abs(angle - t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
1046 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1047 if (currentzone==-10) cross = kTRUE; // we are in the frame
1048 if (currentzone>-10){ // layer knows where we are
1049 if (zone==-10) zone = currentzone;
1050 if (zone!=currentzone) cross=kTRUE;
1054 if (t.GetNCross()==1) t.MakeBackupTrack();
1055 if (t.GetNCross()>2) break;
1061 if (!t.PropagateTo(x,radLength,rho)) break;
1066 // Barrel Tracks [SR, 04.04.2003]
1069 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1070 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1072 // if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1075 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1076 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1078 } else {isNewLayer = kFALSE;}
1083 // now propagate to the middle plane of the next time bin
1084 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1086 x = fTrSec[s]->GetLayer(nr+1)->GetX();
1087 if(!t.PropagateTo(x,radLength,rho)) break;
1088 if (!AdjustSector(&t)) break;
1090 if(!t.PropagateTo(x,radLength,rho)) break;
1095 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1096 // printf("label %d, pl %d, lookForCluster %d \n",
1097 // trackIndex, nr+1, lookForCluster);
1099 if(lookForCluster) {
1100 expectedNumberOfClusters++;
1102 wIndex = (Float_t) t.GetLabel();
1103 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1105 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1106 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1107 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1108 if((t.GetSigmaY2() + sy2) < 0) break;
1109 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1110 Double_t y=t.GetY(), z=t.GetZ();
1114 wYwindow = (Float_t) road;
1115 t.GetPxPyPz(px,py,pz);
1119 wC = (Float_t) t.GetC();
1120 wSigmaC2 = (Float_t) t.GetSigmaC2();
1121 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1122 wSigmaY2 = (Float_t) t.GetSigmaY2();
1123 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1126 if (road>fgkWideRoad) {
1127 if (t.GetNumberOfClusters()>4)
1128 cerr<<t.GetNumberOfClusters()
1129 <<"FindProlongation warning: Too broad road !\n";
1133 AliTRDcluster *cl=0;
1136 Double_t maxChi2=fgkMaxChi2;
1141 maxChi2 = 10 * fgkMaxChi2;
1144 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1145 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1146 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1149 wYclosest = 12345678;
1150 wYcorrect = 12345678;
1151 wZclosest = 12345678;
1152 wZcorrect = 12345678;
1153 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1155 // Find the closest correct cluster for debugging purposes
1158 for (Int_t i=0; i<timeBin; i++) {
1159 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1160 if((c->GetLabel(0) != trackIndex) &&
1161 (c->GetLabel(1) != trackIndex) &&
1162 (c->GetLabel(2) != trackIndex)) continue;
1163 if(TMath::Abs(c->GetY() - y) > minDY) continue;
1164 //minDY = TMath::Abs(c->GetY() - y);
1165 minDY = c->GetY() - y;
1166 wYcorrect = c->GetY();
1167 wZcorrect = c->GetZ();
1169 Double_t h01 = GetTiltFactor(c);
1170 wChi2 = t.GetPredictedChi2(c, h01);
1174 // Now go for the real cluster search
1178 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1179 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1180 if (c->GetY() > y+road) break;
1181 if (c->IsUsed() > 0) continue;
1182 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1184 Double_t h01 = GetTiltFactor(c);
1185 chi2=t.GetPredictedChi2(c,h01);
1187 if (chi2 > maxChi2) continue;
1190 index=timeBin.GetIndex(i);
1193 if((c->GetLabel(0) != trackIndex) &&
1194 (c->GetLabel(1) != trackIndex) &&
1195 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1200 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1201 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1203 if (c->GetY() > y+road) break;
1204 if (c->IsUsed() > 0) continue;
1205 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1207 Double_t h01 = GetTiltFactor(c);
1208 chi2=t.GetPredictedChi2(c,h01);
1210 if (chi2 > maxChi2) continue;
1213 index=timeBin.GetIndex(i);
1218 wYclosest = cl->GetY();
1219 wZclosest = cl->GetZ();
1221 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1222 Double_t h01 = GetTiltFactor(cl);
1223 Int_t det = cl->GetDetector();
1224 Int_t plane = fGeom->GetPlane(det);
1226 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1227 //if(!t.Update(cl,maxChi2,index,h01)) {
1228 if(!tryAgain--) return 0;
1230 else tryAgain=fMaxGap;
1233 if (tryAgain==0) break;
1236 //if (minDY < 1000000 && isNewLayer)
1237 //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() << "\t" <<
1238 // road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
1242 isNewLayer = kFALSE;
1245 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1247 printf(" %f", wIndex); //1
1248 printf(" %f", wTB); //2
1249 printf(" %f", wYrt); //3
1250 printf(" %f", wYclosest); //4
1251 printf(" %f", wYcorrect); //5
1252 printf(" %f", wYwindow); //6
1253 printf(" %f", wZrt); //7
1254 printf(" %f", wZclosest); //8
1255 printf(" %f", wZcorrect); //9
1256 printf(" %f", wZwindow); //10
1257 printf(" %f", wPx); //11
1258 printf(" %f", wPy); //12
1259 printf(" %f", wPz); //13
1260 printf(" %f", wSigmaC2*1000000); //14
1261 printf(" %f", wSigmaTgl2*1000); //15
1262 printf(" %f", wSigmaY2); //16
1263 // printf(" %f", wSigmaZ2); //17
1264 printf(" %f", wChi2); //17
1265 printf(" %f", wC); //18
1276 return expectedNumberOfClusters;
1281 //---------------------------------------------------------------------------
1282 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1284 // Starting from current position on track=t this function tries
1285 // to extrapolate the track up to timeBin=0 and to reuse already
1286 // assigned clusters. Returns the number of clusters
1287 // expected to be found in sensitive layers
1288 // get indices of assigned clusters for each layer
1289 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1292 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1293 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1294 Int_t index = t.GetClusterIndex(i);
1295 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1297 Int_t detector=cl->GetDetector();
1298 Int_t localTimeBin=cl->GetLocalTimeBin();
1299 Int_t sector=fGeom->GetSector(detector);
1300 Int_t plane=fGeom->GetPlane(detector);
1302 Int_t trackingSector = CookSectorIndex(sector);
1304 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1305 if(gtb < 0) continue;
1306 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1307 iCluster[layer] = index;
1311 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1313 Double_t alpha=t.GetAlpha();
1314 alpha = TVector2::Phi_0_2pi(alpha);
1316 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1317 Double_t radLength, rho, x, dx, y, ymax, z;
1319 Int_t expectedNumberOfClusters = 0;
1320 Bool_t lookForCluster;
1322 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1325 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1327 y = t.GetY(); z = t.GetZ();
1329 // first propagate to the inner surface of the current time bin
1330 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1331 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1332 if(!t.PropagateTo(x,radLength,rho)) break;
1334 ymax = x*TMath::Tan(0.5*alpha);
1337 if (!t.Rotate(alpha)) break;
1338 if(!t.PropagateTo(x,radLength,rho)) break;
1339 } else if (y <-ymax) {
1341 if (!t.Rotate(-alpha)) break;
1342 if(!t.PropagateTo(x,radLength,rho)) break;
1345 y = t.GetY(); z = t.GetZ();
1347 // now propagate to the middle plane of the next time bin
1348 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1349 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1350 if(!t.PropagateTo(x,radLength,rho)) break;
1352 ymax = x*TMath::Tan(0.5*alpha);
1355 if (!t.Rotate(alpha)) break;
1356 if(!t.PropagateTo(x,radLength,rho)) break;
1357 } else if (y <-ymax) {
1359 if (!t.Rotate(-alpha)) break;
1360 if(!t.PropagateTo(x,radLength,rho)) break;
1363 if(lookForCluster) expectedNumberOfClusters++;
1365 // use assigned cluster
1366 if (!iCluster[nr-1]) continue;
1367 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1368 Double_t h01 = GetTiltFactor(cl);
1369 Double_t chi2=t.GetPredictedChi2(cl, h01);
1370 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1371 t.Update(cl,chi2,iCluster[nr-1],h01);
1374 return expectedNumberOfClusters;
1377 //___________________________________________________________________
1379 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1381 // Starting from current radial position of track <t> this function
1382 // extrapolates the track up to radial position <xToGo>.
1383 // Returns 1 if track reaches the plane, and 0 otherwise
1385 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1387 Double_t alpha=t.GetAlpha();
1389 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1390 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1392 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1394 Bool_t lookForCluster;
1395 Double_t radLength, rho, x, dx, y, ymax, z;
1399 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1401 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1403 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1405 y = t.GetY(); z = t.GetZ();
1407 // first propagate to the outer surface of the current time bin
1408 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1409 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1410 if(!t.PropagateTo(x,radLength,rho)) return 0;
1412 ymax = x*TMath::Tan(0.5*alpha);
1415 if (!t.Rotate(alpha)) return 0;
1416 } else if (y <-ymax) {
1418 if (!t.Rotate(-alpha)) return 0;
1420 if(!t.PropagateTo(x,radLength,rho)) return 0;
1422 y = t.GetY(); z = t.GetZ();
1424 // now propagate to the middle plane of the next time bin
1425 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1426 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1427 if(!t.PropagateTo(x,radLength,rho)) return 0;
1429 ymax = x*TMath::Tan(0.5*alpha);
1432 if (!t.Rotate(alpha)) return 0;
1433 } else if (y <-ymax) {
1435 if (!t.Rotate(-alpha)) return 0;
1437 if(!t.PropagateTo(x,radLength,rho)) return 0;
1442 //___________________________________________________________________
1444 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1446 // Starting from current radial position of track <t> this function
1447 // extrapolates the track up to radial position of the outermost
1448 // padrow of the TPC.
1449 // Returns 1 if track reaches the TPC, and 0 otherwise
1451 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1453 Double_t alpha=t.GetAlpha();
1454 alpha = TVector2::Phi_0_2pi(alpha);
1456 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1458 Bool_t lookForCluster;
1459 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1463 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1464 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1466 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1471 // first propagate to the outer surface of the current time bin
1472 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1473 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1475 if(!t.PropagateTo(x,radLength,rho)) return 0;
1477 if(!t.PropagateTo(x,radLength,rho)) return 0;
1482 // now propagate to the middle plane of the next time bin
1483 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1484 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1486 if(!t.PropagateTo(x,radLength,rho)) return 0;
1488 if(!t.PropagateTo(x,radLength,rho)) return 0;
1493 //_____________________________________________________________________________
1494 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1496 // Fills clusters into TRD tracking_sectors
1497 // Note that the numbering scheme for the TRD tracking_sectors
1498 // differs from that of TRD sectors
1500 if (ReadClusters(fClusters,cTree)) {
1501 Error("LoadClusters","Problem with reading the clusters !");
1504 Int_t ncl=fClusters->GetEntriesFast();
1506 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1509 for (Int_t ichamber=0;ichamber<5;ichamber++)
1510 for (Int_t isector=0;isector<18;isector++){
1511 fHoles[ichamber][isector]=kTRUE;
1516 // printf("\r %d left ",ncl);
1517 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1518 Int_t detector=c->GetDetector();
1519 Int_t localTimeBin=c->GetLocalTimeBin();
1520 Int_t sector=fGeom->GetSector(detector);
1521 Int_t plane=fGeom->GetPlane(detector);
1523 Int_t trackingSector = CookSectorIndex(sector);
1524 if (c->GetLabel(0)>0){
1525 Int_t chamber = fGeom->GetChamber(detector);
1526 fHoles[chamber][trackingSector]=kFALSE;
1529 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1530 if(gtb < 0) continue;
1531 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1534 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1540 for (Int_t isector=0;isector<18;isector++){
1541 for (Int_t ichamber=0;ichamber<5;ichamber++)
1542 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1543 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1544 fGeom->IsHole(0,ichamber,17-isector));
1550 //_____________________________________________________________________________
1551 void AliTRDtracker::UnloadClusters()
1554 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1559 nentr = fClusters->GetEntriesFast();
1560 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1563 nentr = fSeeds->GetEntriesFast();
1564 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1566 nentr = fTracks->GetEntriesFast();
1567 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1569 Int_t nsec = AliTRDgeometry::kNsect;
1571 for (i = 0; i < nsec; i++) {
1572 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1573 fTrSec[i]->GetLayer(pl)->Clear();
1579 //__________________________________________________________________________
1580 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1582 // Creates track seeds using clusters in timeBins=i1,i2
1585 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1589 Double_t x[5], c[15];
1590 Int_t maxSec=AliTRDgeometry::kNsect;
1592 Double_t alpha=AliTRDgeometry::GetAlpha();
1593 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1594 Double_t cs=cos(alpha), sn=sin(alpha);
1595 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1598 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1599 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1601 Double_t x1 =fTrSec[0]->GetX(i1);
1602 Double_t xx2=fTrSec[0]->GetX(i2);
1604 for (Int_t ns=0; ns<maxSec; ns++) {
1606 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1607 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1608 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1609 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1610 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1612 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1614 for (Int_t is=0; is < r1; is++) {
1615 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1617 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1619 const AliTRDcluster *cl;
1620 Double_t x2, y2, z2;
1621 Double_t x3=0., y3=0.;
1624 if(turn != 2) continue;
1625 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1627 y2=cl->GetY(); z2=cl->GetZ();
1632 else if (js<nl2+nl) {
1633 if(turn != 1) continue;
1634 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1636 y2=cl->GetY(); z2=cl->GetZ();
1641 else if (js<nl2+nl+nm) {
1642 if(turn != 1) continue;
1643 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1645 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1647 else if (js<nl2+nl+nm+nu) {
1648 if(turn != 1) continue;
1649 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1650 cl=r2[js-nl2-nl-nm];
1651 y2=cl->GetY(); z2=cl->GetZ();
1657 if(turn != 2) continue;
1658 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1659 cl=r2[js-nl2-nl-nm-nu];
1660 y2=cl->GetY(); z2=cl->GetZ();
1666 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1668 Double_t zz=z1 - z1/x1*(x1-x2);
1670 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1672 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1673 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1677 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1679 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
1681 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1683 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1685 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1687 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1689 Double_t a=asin(x[2]);
1690 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1692 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1694 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1695 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1696 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
1699 Double_t h01 = GetTiltFactor(r1[is]);
1700 Double_t xuFactor = 100.;
1706 sy1=sy1+sz1*h01*h01;
1707 Double_t syz=sz1*(-h01);
1708 // end of tilt changes
1710 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1711 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1712 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1713 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1714 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1715 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1716 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1717 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1718 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1719 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1723 // c[1]=0.; c[2]=sz1;
1724 c[1]=syz; c[2]=sz1*xuFactor;
1725 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1726 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1727 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1728 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1729 c[13]=f30*sy1*f40+f32*sy2*f42;
1730 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1732 UInt_t index=r1.GetIndex(is);
1734 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1736 Int_t rc=FollowProlongation(*track, i2);
1739 (track->GetNumberOfClusters() <
1740 (outer-inner)*fgkMinClustersInSeed)) delete track;
1742 fSeeds->AddLast(track); fNseeds++;
1743 // cerr<<"\r found seed "<<fNseeds;
1750 //_____________________________________________________________________________
1751 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1754 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1755 // from the file. The names of the cluster tree and branches
1756 // should match the ones used in AliTRDclusterizer::WriteClusters()
1758 TObjArray *clusterArray = new TObjArray(400);
1760 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1762 Error("ReadClusters","Can't get the branch !");
1765 branch->SetAddress(&clusterArray);
1767 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1768 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1770 // Loop through all entries in the tree
1772 AliTRDcluster *c = 0;
1775 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1778 nbytes += ClusterTree->GetEvent(iEntry);
1780 // Get the number of points in the detector
1781 Int_t nCluster = clusterArray->GetEntriesFast();
1782 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1784 // Loop through all TRD digits
1785 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1786 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1787 AliTRDcluster *co = new AliTRDcluster(*c);
1788 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1789 Int_t ltb = co->GetLocalTimeBin();
1790 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1791 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1793 delete clusterArray->RemoveAt(iCluster);
1797 delete clusterArray;
1802 //__________________________________________________________________
1803 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1806 // This cooks a label. Mmmmh, smells good...
1809 Int_t label=123456789, index, i, j;
1810 Int_t ncl=pt->GetNumberOfClusters();
1811 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1815 // Int_t s[kRange][2];
1816 Int_t **s = new Int_t* [kRange];
1817 for (i=0; i<kRange; i++) {
1818 s[i] = new Int_t[2];
1820 for (i=0; i<kRange; i++) {
1826 for (i=0; i<ncl; i++) {
1827 index=pt->GetClusterIndex(i);
1828 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1834 for (i=0; i<ncl; i++) {
1835 index=pt->GetClusterIndex(i);
1836 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1837 for (Int_t k=0; k<3; k++) {
1838 label=c->GetLabel(k);
1839 labelAdded=kFALSE; j=0;
1841 while ( (!labelAdded) && ( j < kRange ) ) {
1842 if (s[j][0]==label || s[j][1]==0) {
1856 for (i=0; i<kRange; i++) {
1858 max=s[i][1]; label=s[i][0];
1862 for (i=0; i<kRange; i++) {
1868 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1870 pt->SetLabel(label);
1875 //__________________________________________________________________
1876 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1879 // Use clusters, but don't abuse them!
1882 Int_t ncl=t->GetNumberOfClusters();
1883 for (Int_t i=from; i<ncl; i++) {
1884 Int_t index = t->GetClusterIndex(i);
1885 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1891 //_____________________________________________________________________
1892 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1894 // Parametrised "expected" error of the cluster reconstruction in Y
1896 Double_t s = 0.08 * 0.08;
1900 //_____________________________________________________________________
1901 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1903 // Parametrised "expected" error of the cluster reconstruction in Z
1905 Double_t s = 9 * 9 /12.;
1909 //_____________________________________________________________________
1910 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
1913 // Returns radial position which corresponds to time bin <localTB>
1914 // in tracking sector <sector> and plane <plane>
1917 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
1918 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1919 return fTrSec[sector]->GetLayer(pl)->GetX();
1924 //_______________________________________________________
1925 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1926 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1929 // AliTRDpropagationLayer constructor
1932 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1933 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1936 for(Int_t i=0; i < (Int_t) kZones; i++) {
1937 fZc[i]=0; fZmax[i] = 0;
1942 if(fTimeBinIndex >= 0) {
1943 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1944 fIndex = new UInt_t[kMaxClusterPerTimeBin];
1947 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
1958 //_______________________________________________________
1959 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1960 Double_t Zmax, Double_t Ymax, Double_t rho,
1961 Double_t radLength, Double_t Yc, Double_t Zc)
1964 // Sets hole in the layer
1972 fHoleX0 = radLength;
1976 //_______________________________________________________
1977 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1980 // AliTRDtrackingSector Constructor
1989 // get holes description from geometry
1990 Bool_t holes[AliTRDgeometry::kNcham];
1991 //printf("sector\t%d\t",gs);
1992 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
1993 holes[icham] = fGeom->IsHole(0,icham,gs);
1994 //printf("%d",holes[icham]);
1998 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
2001 AliTRDpropagationLayer* ppl;
2003 Double_t x, xin, xout, dx, rho, radLength;
2006 // set time bins in the gas of the TPC
2008 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
2009 rho = 0.9e-3; radLength = 28.94;
2011 for(Int_t i=0; i<steps; i++) {
2012 x = xin + i*dx + dx/2;
2013 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2017 // set time bins in the outer field cage vessel
2019 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2020 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2023 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2024 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2027 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2028 steps = 5; dx = (xout - xin)/steps;
2029 for(Int_t i=0; i<steps; i++) {
2030 x = xin + i*dx + dx/2;
2031 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2035 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2036 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2039 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2040 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2044 // set time bins in CO2
2046 xin = xout; xout = 275.0;
2047 steps = 50; dx = (xout - xin)/steps;
2048 rho = 1.977e-3; radLength = 36.2;
2050 for(Int_t i=0; i<steps; i++) {
2051 x = xin + i*dx + dx/2;
2052 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2056 // set time bins in the outer containment vessel
2058 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2059 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2062 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2063 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2066 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2067 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2070 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2071 steps = 10; dx = (xout - xin)/steps;
2072 for(Int_t i=0; i<steps; i++) {
2073 x = xin + i*dx + dx/2;
2074 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2078 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2079 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2082 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2083 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2086 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2087 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2090 Double_t xtrd = (Double_t) fGeom->Rmin();
2092 // add layers between TPC and TRD (Air temporarily)
2093 xin = xout; xout = xtrd;
2094 steps = 50; dx = (xout - xin)/steps;
2095 rho = 1.2e-3; radLength = 36.66;
2097 for(Int_t i=0; i<steps; i++) {
2098 x = xin + i*dx + dx/2;
2099 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2104 // Double_t alpha=AliTRDgeometry::GetAlpha();
2106 // add layers for each of the planes
2108 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2109 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2110 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2111 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2112 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2113 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2114 Double_t dxPlane = dxTEC + dxSpace;
2117 const Int_t kNchambers = AliTRDgeometry::Ncham();
2120 Double_t ymaxsensitive=0;
2121 Double_t *zc = new Double_t[kNchambers];
2122 Double_t *zmax = new Double_t[kNchambers];
2123 Double_t *zmaxsensitive = new Double_t[kNchambers];
2124 // Double_t holeZmax = 1000.; // the whole sector is missing
2126 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2129 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2130 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2131 for(Int_t i=0; i<steps; i++) {
2132 x = xin + i*dx + dx/2;
2133 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2137 ymax = fGeom->GetChamberWidth(plane)/2.;
2138 ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2140 for(Int_t ch = 0; ch < kNchambers; ch++) {
2141 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2142 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2143 Float_t row0 = fPar->GetRow0(plane,ch,0);
2144 Int_t nPads = fPar->GetRowMax(plane,ch,0);
2145 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2146 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2147 zc[ch] = (pad * nPads)/2 + row0;
2148 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2152 dx = fPar->GetTimeBinSize();
2153 rho = 0.00295 * 0.85; radLength = 11.0;
2155 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2156 Double_t xbottom = x0 - dxDrift;
2157 Double_t xtop = x0 + dxAmp;
2159 // Amplification region
2160 steps = (Int_t) (dxAmp/dx);
2162 for(tb = 0; tb < steps; tb++) {
2163 x = x0 + tb * dx + dx/2;
2164 tbIndex = CookTimeBinIndex(plane, -tb-1);
2165 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2166 ppl->SetYmax(ymax,ymaxsensitive);
2167 ppl->SetZ(zc, zmax, zmaxsensitive);
2168 ppl->SetHoles(holes);
2171 tbIndex = CookTimeBinIndex(plane, -steps);
2172 x = (x + dx/2 + xtop)/2;
2174 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2175 ppl->SetYmax(ymax,ymaxsensitive);
2176 ppl->SetZ(zc, zmax,zmaxsensitive);
2177 ppl->SetHoles(holes);
2181 dx = fPar->GetTimeBinSize();
2182 steps = (Int_t) (dxDrift/dx);
2184 for(tb = 0; tb < steps; tb++) {
2185 x = x0 - tb * dx - dx/2;
2186 tbIndex = CookTimeBinIndex(plane, tb);
2188 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2189 ppl->SetYmax(ymax,ymaxsensitive);
2190 ppl->SetZ(zc, zmax, zmaxsensitive);
2191 ppl->SetHoles(holes);
2194 tbIndex = CookTimeBinIndex(plane, steps);
2195 x = (x - dx/2 + xbottom)/2;
2197 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2198 ppl->SetYmax(ymax,ymaxsensitive);
2199 ppl->SetZ(zc, zmax, zmaxsensitive);
2200 ppl->SetHoles(holes);
2204 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2205 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2206 ppl->SetYmax(ymax,ymaxsensitive);
2207 ppl->SetZ(zc, zmax,zmax);
2208 ppl->SetHoles(holes);
2212 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2213 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2214 for(Int_t i=0; i<steps; i++) {
2215 x = xin + i*dx + dx/2;
2216 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2217 ppl->SetYmax(ymax,ymaxsensitive);
2218 ppl->SetZ(zc, zmax,zmax);
2219 ppl->SetHoles(holes);
2223 // Space between the chambers, air
2224 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2225 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2226 for(Int_t i=0; i<steps; i++) {
2227 x = xin + i*dx + dx/2;
2228 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2233 // Space between the TRD and RICH
2234 Double_t xRICH = 500.;
2235 xin = xout; xout = xRICH;
2236 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2237 for(Int_t i=0; i<steps; i++) {
2238 x = xin + i*dx + dx/2;
2239 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2249 //______________________________________________________
2251 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2254 // depending on the digitization parameters calculates "global"
2255 // time bin index for timebin <localTB> in plane <plane>
2258 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2259 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2260 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2262 Int_t tbAmp = fPar->GetTimeBefore();
2263 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2264 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2265 Int_t tbDrift = fPar->GetTimeMax();
2266 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2268 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2270 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2273 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2274 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2281 //______________________________________________________
2283 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2286 // For all sensitive time bins sets corresponding layer index
2287 // in the array fTimeBins
2292 for(Int_t i = 0; i < fN; i++) {
2293 index = fLayers[i]->GetTimeBinIndex();
2295 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2297 if(index < 0) continue;
2298 if(index >= (Int_t) kMaxTimeBinIndex) {
2299 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2300 printf(" index %d exceeds allowed maximum of %d!\n",
2301 index, kMaxTimeBinIndex-1);
2304 fTimeBinIndex[index] = i;
2307 Double_t x1, dx1, x2, dx2, gap;
2309 for(Int_t i = 0; i < fN-1; i++) {
2310 x1 = fLayers[i]->GetX();
2311 dx1 = fLayers[i]->GetdX();
2312 x2 = fLayers[i+1]->GetX();
2313 dx2 = fLayers[i+1]->GetdX();
2314 gap = (x2 - dx2/2) - (x1 + dx1/2);
2316 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2317 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2320 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2321 printf(" (%f - %f) - (%f + %f) = %f\n",
2322 x2, dx2/2, x1, dx1, gap);
2328 //______________________________________________________
2331 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2334 // Returns the number of time bin which in radial position is closest to <x>
2337 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2338 if(x <= fLayers[0]->GetX()) return 0;
2340 Int_t b=0, e=fN-1, m=(b+e)/2;
2341 for (; b<e; m=(b+e)/2) {
2342 if (x > fLayers[m]->GetX()) b=m+1;
2345 if(TMath::Abs(x - fLayers[m]->GetX()) >
2346 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2351 //______________________________________________________
2353 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2356 // Returns number of the innermost SENSITIVE propagation layer
2359 return GetLayerNumber(0);
2362 //______________________________________________________
2364 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2367 // Returns number of the outermost SENSITIVE time bin
2370 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2373 //______________________________________________________
2375 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2378 // Returns number of SENSITIVE time bins
2382 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2383 layer = GetLayerNumber(tb);
2389 //______________________________________________________
2391 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2394 // Insert layer <pl> in fLayers array.
2395 // Layers are sorted according to X coordinate.
2397 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2398 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2401 if (fN==0) {fLayers[fN++] = pl; return;}
2402 Int_t i=Find(pl->GetX());
2404 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2405 fLayers[i]=pl; fN++;
2409 //______________________________________________________
2411 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2414 // Returns index of the propagation layer nearest to X
2417 if (x <= fLayers[0]->GetX()) return 0;
2418 if (x > fLayers[fN-1]->GetX()) return fN;
2419 Int_t b=0, e=fN-1, m=(b+e)/2;
2420 for (; b<e; m=(b+e)/2) {
2421 if (x > fLayers[m]->GetX()) b=m+1;
2427 //______________________________________________________
2428 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2431 // set centers and the width of sectors
2432 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2433 fZc[icham] = center[icham];
2434 fZmax[icham] = w[icham];
2435 fZmaxSensitive[icham] = wsensitive[icham];
2436 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2439 //______________________________________________________
2441 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2444 // set centers and the width of sectors
2446 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2447 fIsHole[icham] = holes[icham];
2448 if (holes[icham]) fHole = kTRUE;
2454 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2455 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
2456 Bool_t &lookForCluster) const
2459 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2460 // and sensitivity <lookForCluster> in point <y,z>
2466 lookForCluster = kFALSE;
2468 // check dead regions in sensitive volume
2469 if(fTimeBinIndex >= 0) {
2472 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2473 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2475 lookForCluster = !(fIsHole[zone]);
2476 if(TMath::Abs(y) > fYmaxSensitive){
2477 lookForCluster = kFALSE;
2479 if (fIsHole[zone]) {
2491 if (fHole==kFALSE) return;
2493 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2494 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2505 Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2509 if (fTimeBinIndex < 0) return -20; //unknown
2510 Int_t zone=-10; // dead zone
2511 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2512 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2519 //______________________________________________________
2521 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2524 // Insert cluster in cluster array.
2525 // Clusters are sorted according to Y coordinate.
2527 if(fTimeBinIndex < 0) {
2528 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2532 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2533 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2536 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2537 Int_t i=Find(c->GetY());
2538 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2539 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2540 fIndex[i]=index; fClusters[i]=c; fN++;
2543 //______________________________________________________
2545 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2547 // Returns index of the cluster nearest in Y
2549 if (y <= fClusters[0]->GetY()) return 0;
2550 if (y > fClusters[fN-1]->GetY()) return fN;
2551 Int_t b=0, e=fN-1, m=(b+e)/2;
2552 for (; b<e; m=(b+e)/2) {
2553 if (y > fClusters[m]->GetY()) b=m+1;
2559 //---------------------------------------------------------
2561 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2563 // Returns correction factor for tilted pads geometry
2566 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2567 Int_t det = c->GetDetector();
2568 Int_t plane = fGeom->GetPlane(det);
2570 //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2571 if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
2573 if(fNoTilt) h01 = 0;