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);
512 delete fSeeds->RemoveAt(i);
519 cout<<"Total number of found tracks: "<<found<<endl;
526 //_____________________________________________________________________________
527 Int_t AliTRDtracker::PropagateBack(AliESD* event) {
529 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
530 // backpropagated by the TPC tracker. Each seed is first propagated
531 // to the TRD, and then its prolongation is searched in the TRD.
532 // If sufficiently long continuation of the track is found in the TRD
533 // the track is updated, otherwise it's stored as originaly defined
534 // by the TPC tracker.
538 Float_t foundMin = 40;
540 Int_t n = event->GetNumberOfTracks();
541 for (Int_t i=0; i<n; i++) {
542 AliESDtrack* seed=event->GetTrack(i);
543 ULong_t status=seed->GetStatus();
544 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
545 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
547 Int_t lbl = seed->GetLabel();
548 AliTRDtrack *track = new AliTRDtrack(*seed);
549 track->SetSeedLabel(lbl);
551 Float_t p4 = track->GetC();
553 Int_t expectedClr = FollowBackProlongation(*track);
554 if (track->GetNumberOfClusters()<expectedClr/3){
555 AliTRDtrack *track1 = new AliTRDtrack(*seed);
556 track1->SetSeedLabel(lbl);
557 FollowBackProlongation(*track1);
558 AliTRDtrack *track2= new AliTRDtrack(*seed);
559 track->SetSeedLabel(lbl);
560 FollowBackProlongation(*track2);
564 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)>0.2) {
566 continue; //too big change of curvature - to be checked
568 Int_t foundClr = track->GetNumberOfClusters();
569 if (foundClr >= foundMin) {
572 // CookLabel(track, 1-fgkLabelFraction);
576 // Propagate to outer reference plane [SR, GSI, 18.02.2003]
577 // track->PropagateTo(364.8); why?
579 //seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
583 //Propagation to the TOF (I.Belikov)
585 if (track->GetStop()==kFALSE){
588 Double_t c2=track->GetC()*xtof - track->GetEta();
589 if (TMath::Abs(c2)>=0.85) continue;
590 Double_t xTOF0 = 371. ;
591 PropagateToOuterPlane(*track,xTOF0);
593 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
594 Double_t y=track->GetYat(xtof);
596 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
600 } else if (y <-ymax) {
601 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
607 if (track->PropagateTo(xtof)) {
608 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
609 if (track->GetNumberOfClusters()>foundMin) found++;
612 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
613 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
614 seed->SetStatus(AliESDtrack::kTRDStop);
621 //End of propagation to the TOF
622 //if (foundClr>foundMin)
623 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
628 cerr<<"Number of seeds: "<<fNseeds<<endl;
629 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
631 fSeeds->Clear(); fNseeds=0;
637 //_____________________________________________________________________________
638 Int_t AliTRDtracker::RefitInward(AliESD* event)
641 // Refits tracks within the TRD. The ESD event is expected to contain seeds
642 // at the outer part of the TRD.
643 // The tracks are propagated to the innermost time bin
644 // of the TRD and the ESD event is updated
645 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
648 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
649 Float_t foundMin = fgkMinClustersInTrack * timeBins;
652 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
654 Int_t n = event->GetNumberOfTracks();
655 for (Int_t i=0; i<n; i++) {
656 AliESDtrack* seed=event->GetTrack(i);
657 ULong_t status=seed->GetStatus();
658 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
659 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
662 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
663 seed2->ResetCovariance(5.);
664 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
665 UInt_t * indexes2 = seed2->GetIndexes();
666 UInt_t * indexes3 = pt->GetBackupIndexes();
667 for (Int_t i=0;i<200;i++) {
668 if (indexes2[i]==0) break;
669 indexes3[i] = indexes2[i];
671 //AliTRDtrack *pt = seed2;
673 FollowProlongation(t, innerTB);
675 if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
676 // debug - why we dont go back?
677 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
678 UInt_t * indexes2 = seed2->GetIndexes();
679 UInt_t * indexes3 = pt2->GetBackupIndexes();
680 for (Int_t i=0;i<200;i++) {
681 if (indexes2[i]==0) break;
682 indexes3[i] = indexes2[i];
684 FollowProlongation(*pt2, innerTB);
688 if (t.GetNumberOfClusters() >= foundMin) {
690 //CookLabel(pt, 1-fgkLabelFraction);
694 // cout<<found<<'\r';
696 if(PropagateToTPC(t)) {
697 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
699 //if not prolongation to TPC - propagate without update
700 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
701 seed2->ResetCovariance(5.);
702 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
704 if (PropagateToTPC(*pt2))
705 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
713 cout<<"Number of loaded seeds: "<<nseed<<endl;
714 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
721 //---------------------------------------------------------------------------
722 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
724 // Starting from current position on track=t this function tries
725 // to extrapolate the track up to timeBin=0 and to confirm prolongation
726 // if a close cluster is found. Returns the number of clusters
727 // expected to be found in sensitive layers
729 Float_t wIndex, wTB, wChi2;
730 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
731 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
732 Float_t wPx, wPy, wPz, wC;
734 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
735 Int_t lastplane = GetLastPlane(&t);
737 Int_t trackIndex = t.GetLabel();
739 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
741 Int_t tryAgain=fMaxGap;
743 Double_t alpha=t.GetAlpha();
744 alpha = TVector2::Phi_0_2pi(alpha);
746 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
747 Double_t radLength, rho, x, dx, y, ymax, z;
749 Int_t expectedNumberOfClusters = 0;
750 Bool_t lookForCluster;
752 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
755 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
757 y = t.GetY(); z = t.GetZ();
759 // first propagate to the inner surface of the current time bin
760 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
761 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
762 if(!t.PropagateTo(x,radLength,rho)) break;
764 ymax = x*TMath::Tan(0.5*alpha);
767 if (!t.Rotate(alpha)) break;
768 if(!t.PropagateTo(x,radLength,rho)) break;
769 } else if (y <-ymax) {
771 if (!t.Rotate(-alpha)) break;
772 if(!t.PropagateTo(x,radLength,rho)) break;
775 y = t.GetY(); z = t.GetZ();
777 // now propagate to the middle plane of the next time bin
778 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
779 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
780 if(!t.PropagateTo(x,radLength,rho)) break;
782 ymax = x*TMath::Tan(0.5*alpha);
785 if (!t.Rotate(alpha)) break;
786 if(!t.PropagateTo(x,radLength,rho)) break;
787 } else if (y <-ymax) {
789 if (!t.Rotate(-alpha)) break;
790 if(!t.PropagateTo(x,radLength,rho)) break;
796 expectedNumberOfClusters++;
797 wIndex = (Float_t) t.GetLabel();
800 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
802 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
803 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
806 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
807 else return expectedNumberOfClusters;
811 wYwindow = (Float_t) road;
812 t.GetPxPyPz(px,py,pz);
816 wC = (Float_t) t.GetC();
817 wSigmaC2 = (Float_t) t.GetSigmaC2();
818 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
819 wSigmaY2 = (Float_t) t.GetSigmaY2();
820 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
827 Double_t maxChi2=fgkMaxChi2;
829 wYclosest = 12345678;
830 wYcorrect = 12345678;
831 wZclosest = 12345678;
832 wZcorrect = 12345678;
833 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
835 // Find the closest correct cluster for debugging purposes
837 Float_t minDY = 1000000;
838 for (Int_t i=0; i<timeBin; i++) {
839 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
840 if((c->GetLabel(0) != trackIndex) &&
841 (c->GetLabel(1) != trackIndex) &&
842 (c->GetLabel(2) != trackIndex)) continue;
843 if(TMath::Abs(c->GetY() - y) > minDY) continue;
844 minDY = TMath::Abs(c->GetY() - y);
845 wYcorrect = c->GetY();
846 wZcorrect = c->GetZ();
848 Double_t h01 = GetTiltFactor(c);
849 wChi2 = t.GetPredictedChi2(c, h01);
853 // Now go for the real cluster search
857 //find cluster in history
860 AliTRDcluster * cl0 = timeBin[0];
864 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
865 if (plane>lastplane) continue;
866 Int_t timebin = cl0->GetLocalTimeBin();
867 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
870 Double_t h01 = GetTiltFactor(cl);
871 maxChi2=t.GetPredictedChi2(cl,h01);
873 if ((!cl) && road>fgkWideRoad) {
874 //if (t.GetNumberOfClusters()>4)
875 // cerr<<t.GetNumberOfClusters()
876 // <<"FindProlongation warning: Too broad road !\n";
883 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
884 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
885 if (c->GetY() > y+road) break;
886 if (c->IsUsed() > 0) continue;
887 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
889 Double_t h01 = GetTiltFactor(c);
890 Double_t chi2=t.GetPredictedChi2(c,h01);
892 if (chi2 > maxChi2) continue;
895 index=timeBin.GetIndex(i);
901 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
902 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
904 if (c->GetY() > y+road) break;
905 if (c->IsUsed() > 0) continue;
906 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
908 Double_t h01 = GetTiltFactor(c);
909 Double_t chi2=t.GetPredictedChi2(c, h01);
911 if (chi2 > maxChi2) continue;
914 index=timeBin.GetIndex(i);
918 wYclosest = cl->GetY();
919 wZclosest = cl->GetZ();
920 Double_t h01 = GetTiltFactor(cl);
922 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
923 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
924 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
925 Int_t det = cl->GetDetector();
926 Int_t plane = fGeom->GetPlane(det);
928 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
929 //if(!t.Update(cl,maxChi2,index,h01)) {
930 //if(!tryAgain--) return 0;
932 else tryAgain=fMaxGap;
935 //if (tryAgain==0) break;
940 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
942 printf(" %f", wIndex); //1
943 printf(" %f", wTB); //2
944 printf(" %f", wYrt); //3
945 printf(" %f", wYclosest); //4
946 printf(" %f", wYcorrect); //5
947 printf(" %f", wYwindow); //6
948 printf(" %f", wZrt); //7
949 printf(" %f", wZclosest); //8
950 printf(" %f", wZcorrect); //9
951 printf(" %f", wZwindow); //10
952 printf(" %f", wPx); //11
953 printf(" %f", wPy); //12
954 printf(" %f", wPz); //13
955 printf(" %f", wSigmaC2*1000000); //14
956 printf(" %f", wSigmaTgl2*1000); //15
957 printf(" %f", wSigmaY2); //16
958 // printf(" %f", wSigmaZ2); //17
959 printf(" %f", wChi2); //17
960 printf(" %f", wC); //18
967 return expectedNumberOfClusters;
972 //___________________________________________________________________
974 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
976 // Starting from current radial position of track <t> this function
977 // extrapolates the track up to outer timebin and in the sensitive
978 // layers confirms prolongation if a close cluster is found.
979 // Returns the number of clusters expected to be found in sensitive layers
981 Float_t wIndex, wTB, wChi2;
982 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
983 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
984 Float_t wPx, wPy, wPz, wC;
986 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
988 Int_t trackIndex = t.GetLabel();
989 Int_t tryAgain=fMaxGap;
991 Double_t alpha=t.GetAlpha();
992 TVector2::Phi_0_2pi(alpha);
996 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
997 Double_t radLength, rho, x, dx, y, ymax = 0, z;
998 Bool_t lookForCluster;
1000 Int_t expectedNumberOfClusters = 0;
1003 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1005 Int_t nRefPlane = fgkFirstPlane;
1006 Bool_t isNewLayer = kFALSE;
1012 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
1017 // first propagate to the outer surface of the current time bin
1020 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1021 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
1025 if(!t.PropagateTo(x,radLength,rho)) break;
1026 // if (!AdjustSector(&t)) break;
1028 // MI -fix untill correct material desription will be implemented
1030 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material -
1031 if (!AdjustSector(&t)) break;
1032 if (TMath::Abs(angle - t.GetAlpha())>0.000001) break; //better to stop track
1033 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1034 if (currentzone==-10) break; // we are in the frame
1035 if (currentzone>-10){ // layer knows where we are
1036 if (zone==-10) zone = currentzone;
1037 if (zone!=currentzone) break;
1042 if (!t.PropagateTo(x,radLength,rho)) break;
1047 // Barrel Tracks [SR, 04.04.2003]
1050 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1051 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1053 // if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1056 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1057 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1059 } else {isNewLayer = kFALSE;}
1064 // now propagate to the middle plane of the next time bin
1065 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1067 x = fTrSec[s]->GetLayer(nr+1)->GetX();
1068 if(!t.PropagateTo(x,radLength,rho)) break;
1069 if (!AdjustSector(&t)) break;
1071 if(!t.PropagateTo(x,radLength,rho)) break;
1076 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1077 // printf("label %d, pl %d, lookForCluster %d \n",
1078 // trackIndex, nr+1, lookForCluster);
1080 if(lookForCluster) {
1081 expectedNumberOfClusters++;
1083 wIndex = (Float_t) t.GetLabel();
1084 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1086 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1087 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1088 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1089 if((t.GetSigmaY2() + sy2) < 0) break;
1090 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1091 Double_t y=t.GetY(), z=t.GetZ();
1095 wYwindow = (Float_t) road;
1096 t.GetPxPyPz(px,py,pz);
1100 wC = (Float_t) t.GetC();
1101 wSigmaC2 = (Float_t) t.GetSigmaC2();
1102 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1103 wSigmaY2 = (Float_t) t.GetSigmaY2();
1104 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1107 if (road>fgkWideRoad) {
1108 if (t.GetNumberOfClusters()>4)
1109 cerr<<t.GetNumberOfClusters()
1110 <<"FindProlongation warning: Too broad road !\n";
1114 AliTRDcluster *cl=0;
1117 Double_t maxChi2=fgkMaxChi2;
1122 maxChi2 = 10 * fgkMaxChi2;
1125 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1126 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1127 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1130 wYclosest = 12345678;
1131 wYcorrect = 12345678;
1132 wZclosest = 12345678;
1133 wZcorrect = 12345678;
1134 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1136 // Find the closest correct cluster for debugging purposes
1139 for (Int_t i=0; i<timeBin; i++) {
1140 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1141 if((c->GetLabel(0) != trackIndex) &&
1142 (c->GetLabel(1) != trackIndex) &&
1143 (c->GetLabel(2) != trackIndex)) continue;
1144 if(TMath::Abs(c->GetY() - y) > minDY) continue;
1145 //minDY = TMath::Abs(c->GetY() - y);
1146 minDY = c->GetY() - y;
1147 wYcorrect = c->GetY();
1148 wZcorrect = c->GetZ();
1150 Double_t h01 = GetTiltFactor(c);
1151 wChi2 = t.GetPredictedChi2(c, h01);
1155 // Now go for the real cluster search
1159 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1160 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1161 if (c->GetY() > y+road) break;
1162 if (c->IsUsed() > 0) continue;
1163 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1165 Double_t h01 = GetTiltFactor(c);
1166 chi2=t.GetPredictedChi2(c,h01);
1168 if (chi2 > maxChi2) continue;
1171 index=timeBin.GetIndex(i);
1174 if((c->GetLabel(0) != trackIndex) &&
1175 (c->GetLabel(1) != trackIndex) &&
1176 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1181 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1182 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1184 if (c->GetY() > y+road) break;
1185 if (c->IsUsed() > 0) continue;
1186 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1188 Double_t h01 = GetTiltFactor(c);
1189 chi2=t.GetPredictedChi2(c,h01);
1191 if (chi2 > maxChi2) continue;
1194 index=timeBin.GetIndex(i);
1199 wYclosest = cl->GetY();
1200 wZclosest = cl->GetZ();
1202 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1203 Double_t h01 = GetTiltFactor(cl);
1204 Int_t det = cl->GetDetector();
1205 Int_t plane = fGeom->GetPlane(det);
1207 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1208 //if(!t.Update(cl,maxChi2,index,h01)) {
1209 if(!tryAgain--) return 0;
1211 else tryAgain=fMaxGap;
1214 if (tryAgain==0) break;
1217 //if (minDY < 1000000 && isNewLayer)
1218 //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() << "\t" <<
1219 // road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
1223 isNewLayer = kFALSE;
1226 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1228 printf(" %f", wIndex); //1
1229 printf(" %f", wTB); //2
1230 printf(" %f", wYrt); //3
1231 printf(" %f", wYclosest); //4
1232 printf(" %f", wYcorrect); //5
1233 printf(" %f", wYwindow); //6
1234 printf(" %f", wZrt); //7
1235 printf(" %f", wZclosest); //8
1236 printf(" %f", wZcorrect); //9
1237 printf(" %f", wZwindow); //10
1238 printf(" %f", wPx); //11
1239 printf(" %f", wPy); //12
1240 printf(" %f", wPz); //13
1241 printf(" %f", wSigmaC2*1000000); //14
1242 printf(" %f", wSigmaTgl2*1000); //15
1243 printf(" %f", wSigmaY2); //16
1244 // printf(" %f", wSigmaZ2); //17
1245 printf(" %f", wChi2); //17
1246 printf(" %f", wC); //18
1257 return expectedNumberOfClusters;
1262 //---------------------------------------------------------------------------
1263 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1265 // Starting from current position on track=t this function tries
1266 // to extrapolate the track up to timeBin=0 and to reuse already
1267 // assigned clusters. Returns the number of clusters
1268 // expected to be found in sensitive layers
1269 // get indices of assigned clusters for each layer
1270 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1273 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1274 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1275 Int_t index = t.GetClusterIndex(i);
1276 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1278 Int_t detector=cl->GetDetector();
1279 Int_t localTimeBin=cl->GetLocalTimeBin();
1280 Int_t sector=fGeom->GetSector(detector);
1281 Int_t plane=fGeom->GetPlane(detector);
1283 Int_t trackingSector = CookSectorIndex(sector);
1285 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1286 if(gtb < 0) continue;
1287 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1288 iCluster[layer] = index;
1292 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1294 Double_t alpha=t.GetAlpha();
1295 alpha = TVector2::Phi_0_2pi(alpha);
1297 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1298 Double_t radLength, rho, x, dx, y, ymax, z;
1300 Int_t expectedNumberOfClusters = 0;
1301 Bool_t lookForCluster;
1303 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1306 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1308 y = t.GetY(); z = t.GetZ();
1310 // first propagate to the inner surface of the current time bin
1311 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1312 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1313 if(!t.PropagateTo(x,radLength,rho)) break;
1315 ymax = x*TMath::Tan(0.5*alpha);
1318 if (!t.Rotate(alpha)) break;
1319 if(!t.PropagateTo(x,radLength,rho)) break;
1320 } else if (y <-ymax) {
1322 if (!t.Rotate(-alpha)) break;
1323 if(!t.PropagateTo(x,radLength,rho)) break;
1326 y = t.GetY(); z = t.GetZ();
1328 // now propagate to the middle plane of the next time bin
1329 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1330 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1331 if(!t.PropagateTo(x,radLength,rho)) break;
1333 ymax = x*TMath::Tan(0.5*alpha);
1336 if (!t.Rotate(alpha)) break;
1337 if(!t.PropagateTo(x,radLength,rho)) break;
1338 } else if (y <-ymax) {
1340 if (!t.Rotate(-alpha)) break;
1341 if(!t.PropagateTo(x,radLength,rho)) break;
1344 if(lookForCluster) expectedNumberOfClusters++;
1346 // use assigned cluster
1347 if (!iCluster[nr-1]) continue;
1348 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1349 Double_t h01 = GetTiltFactor(cl);
1350 Double_t chi2=t.GetPredictedChi2(cl, h01);
1351 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1352 t.Update(cl,chi2,iCluster[nr-1],h01);
1355 return expectedNumberOfClusters;
1358 //___________________________________________________________________
1360 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1362 // Starting from current radial position of track <t> this function
1363 // extrapolates the track up to radial position <xToGo>.
1364 // Returns 1 if track reaches the plane, and 0 otherwise
1366 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1368 Double_t alpha=t.GetAlpha();
1370 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1371 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1373 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1375 Bool_t lookForCluster;
1376 Double_t radLength, rho, x, dx, y, ymax, z;
1380 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1382 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1384 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1386 y = t.GetY(); z = t.GetZ();
1388 // first propagate to the outer surface of the current time bin
1389 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1390 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1391 if(!t.PropagateTo(x,radLength,rho)) return 0;
1393 ymax = x*TMath::Tan(0.5*alpha);
1396 if (!t.Rotate(alpha)) return 0;
1397 } else if (y <-ymax) {
1399 if (!t.Rotate(-alpha)) return 0;
1401 if(!t.PropagateTo(x,radLength,rho)) return 0;
1403 y = t.GetY(); z = t.GetZ();
1405 // now propagate to the middle plane of the next time bin
1406 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1407 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1408 if(!t.PropagateTo(x,radLength,rho)) return 0;
1410 ymax = x*TMath::Tan(0.5*alpha);
1413 if (!t.Rotate(alpha)) return 0;
1414 } else if (y <-ymax) {
1416 if (!t.Rotate(-alpha)) return 0;
1418 if(!t.PropagateTo(x,radLength,rho)) return 0;
1423 //___________________________________________________________________
1425 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1427 // Starting from current radial position of track <t> this function
1428 // extrapolates the track up to radial position of the outermost
1429 // padrow of the TPC.
1430 // Returns 1 if track reaches the TPC, and 0 otherwise
1432 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1434 Double_t alpha=t.GetAlpha();
1435 alpha = TVector2::Phi_0_2pi(alpha);
1437 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1439 Bool_t lookForCluster;
1440 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1444 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1445 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1447 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1452 // first propagate to the outer surface of the current time bin
1453 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1454 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1456 if(!t.PropagateTo(x,radLength,rho)) return 0;
1458 if(!t.PropagateTo(x,radLength,rho)) return 0;
1463 // now propagate to the middle plane of the next time bin
1464 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1465 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1467 if(!t.PropagateTo(x,radLength,rho)) return 0;
1469 if(!t.PropagateTo(x,radLength,rho)) return 0;
1474 //_____________________________________________________________________________
1475 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1477 // Fills clusters into TRD tracking_sectors
1478 // Note that the numbering scheme for the TRD tracking_sectors
1479 // differs from that of TRD sectors
1481 if (ReadClusters(fClusters,cTree)) {
1482 Error("LoadClusters","Problem with reading the clusters !");
1485 Int_t ncl=fClusters->GetEntriesFast();
1487 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1490 for (Int_t ichamber=0;ichamber<5;ichamber++)
1491 for (Int_t isector=0;isector<18;isector++){
1492 fHoles[ichamber][isector]=kTRUE;
1497 // printf("\r %d left ",ncl);
1498 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1499 Int_t detector=c->GetDetector();
1500 Int_t localTimeBin=c->GetLocalTimeBin();
1501 Int_t sector=fGeom->GetSector(detector);
1502 Int_t plane=fGeom->GetPlane(detector);
1504 Int_t trackingSector = CookSectorIndex(sector);
1505 if (c->GetLabel(0)>0){
1506 Int_t chamber = fGeom->GetChamber(detector);
1507 fHoles[chamber][trackingSector]=kFALSE;
1510 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1511 if(gtb < 0) continue;
1512 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1515 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1521 for (Int_t isector=0;isector<18;isector++){
1522 for (Int_t ichamber=0;ichamber<5;ichamber++)
1523 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1524 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1525 fGeom->IsHole(0,ichamber,17-isector));
1531 //_____________________________________________________________________________
1532 void AliTRDtracker::UnloadClusters()
1535 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1540 nentr = fClusters->GetEntriesFast();
1541 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1544 nentr = fSeeds->GetEntriesFast();
1545 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1547 nentr = fTracks->GetEntriesFast();
1548 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1550 Int_t nsec = AliTRDgeometry::kNsect;
1552 for (i = 0; i < nsec; i++) {
1553 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1554 fTrSec[i]->GetLayer(pl)->Clear();
1560 //__________________________________________________________________________
1561 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1563 // Creates track seeds using clusters in timeBins=i1,i2
1566 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1570 Double_t x[5], c[15];
1571 Int_t maxSec=AliTRDgeometry::kNsect;
1573 Double_t alpha=AliTRDgeometry::GetAlpha();
1574 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1575 Double_t cs=cos(alpha), sn=sin(alpha);
1576 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1579 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1580 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1582 Double_t x1 =fTrSec[0]->GetX(i1);
1583 Double_t xx2=fTrSec[0]->GetX(i2);
1585 for (Int_t ns=0; ns<maxSec; ns++) {
1587 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1588 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1589 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1590 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1591 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1593 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1595 for (Int_t is=0; is < r1; is++) {
1596 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1598 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1600 const AliTRDcluster *cl;
1601 Double_t x2, y2, z2;
1602 Double_t x3=0., y3=0.;
1605 if(turn != 2) continue;
1606 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1608 y2=cl->GetY(); z2=cl->GetZ();
1613 else if (js<nl2+nl) {
1614 if(turn != 1) continue;
1615 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1617 y2=cl->GetY(); z2=cl->GetZ();
1622 else if (js<nl2+nl+nm) {
1623 if(turn != 1) continue;
1624 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1626 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1628 else if (js<nl2+nl+nm+nu) {
1629 if(turn != 1) continue;
1630 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1631 cl=r2[js-nl2-nl-nm];
1632 y2=cl->GetY(); z2=cl->GetZ();
1638 if(turn != 2) continue;
1639 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1640 cl=r2[js-nl2-nl-nm-nu];
1641 y2=cl->GetY(); z2=cl->GetZ();
1647 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1649 Double_t zz=z1 - z1/x1*(x1-x2);
1651 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1653 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1654 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1658 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1660 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
1662 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1664 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1666 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1668 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1670 Double_t a=asin(x[2]);
1671 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1673 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1675 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1676 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1677 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
1680 Double_t h01 = GetTiltFactor(r1[is]);
1681 Double_t xuFactor = 100.;
1687 sy1=sy1+sz1*h01*h01;
1688 Double_t syz=sz1*(-h01);
1689 // end of tilt changes
1691 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1692 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1693 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1694 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1695 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1696 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1697 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1698 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1699 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1700 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1704 // c[1]=0.; c[2]=sz1;
1705 c[1]=syz; c[2]=sz1*xuFactor;
1706 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1707 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1708 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1709 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1710 c[13]=f30*sy1*f40+f32*sy2*f42;
1711 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1713 UInt_t index=r1.GetIndex(is);
1715 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1717 Int_t rc=FollowProlongation(*track, i2);
1720 (track->GetNumberOfClusters() <
1721 (outer-inner)*fgkMinClustersInSeed)) delete track;
1723 fSeeds->AddLast(track); fNseeds++;
1724 // cerr<<"\r found seed "<<fNseeds;
1731 //_____________________________________________________________________________
1732 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1735 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1736 // from the file. The names of the cluster tree and branches
1737 // should match the ones used in AliTRDclusterizer::WriteClusters()
1739 TObjArray *clusterArray = new TObjArray(400);
1741 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1743 Error("ReadClusters","Can't get the branch !");
1746 branch->SetAddress(&clusterArray);
1748 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1749 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1751 // Loop through all entries in the tree
1753 AliTRDcluster *c = 0;
1756 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1759 nbytes += ClusterTree->GetEvent(iEntry);
1761 // Get the number of points in the detector
1762 Int_t nCluster = clusterArray->GetEntriesFast();
1763 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1765 // Loop through all TRD digits
1766 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1767 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1768 AliTRDcluster *co = new AliTRDcluster(*c);
1769 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1770 Int_t ltb = co->GetLocalTimeBin();
1771 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1772 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1774 delete clusterArray->RemoveAt(iCluster);
1778 delete clusterArray;
1783 //__________________________________________________________________
1784 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1787 // This cooks a label. Mmmmh, smells good...
1790 Int_t label=123456789, index, i, j;
1791 Int_t ncl=pt->GetNumberOfClusters();
1792 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1796 // Int_t s[kRange][2];
1797 Int_t **s = new Int_t* [kRange];
1798 for (i=0; i<kRange; i++) {
1799 s[i] = new Int_t[2];
1801 for (i=0; i<kRange; i++) {
1807 for (i=0; i<ncl; i++) {
1808 index=pt->GetClusterIndex(i);
1809 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1815 for (i=0; i<ncl; i++) {
1816 index=pt->GetClusterIndex(i);
1817 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1818 for (Int_t k=0; k<3; k++) {
1819 label=c->GetLabel(k);
1820 labelAdded=kFALSE; j=0;
1822 while ( (!labelAdded) && ( j < kRange ) ) {
1823 if (s[j][0]==label || s[j][1]==0) {
1837 for (i=0; i<kRange; i++) {
1839 max=s[i][1]; label=s[i][0];
1843 for (i=0; i<kRange; i++) {
1849 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1851 pt->SetLabel(label);
1856 //__________________________________________________________________
1857 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1860 // Use clusters, but don't abuse them!
1863 Int_t ncl=t->GetNumberOfClusters();
1864 for (Int_t i=from; i<ncl; i++) {
1865 Int_t index = t->GetClusterIndex(i);
1866 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1872 //_____________________________________________________________________
1873 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1875 // Parametrised "expected" error of the cluster reconstruction in Y
1877 Double_t s = 0.08 * 0.08;
1881 //_____________________________________________________________________
1882 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1884 // Parametrised "expected" error of the cluster reconstruction in Z
1886 Double_t s = 9 * 9 /12.;
1890 //_____________________________________________________________________
1891 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
1894 // Returns radial position which corresponds to time bin <localTB>
1895 // in tracking sector <sector> and plane <plane>
1898 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
1899 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1900 return fTrSec[sector]->GetLayer(pl)->GetX();
1905 //_______________________________________________________
1906 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1907 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1910 // AliTRDpropagationLayer constructor
1913 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1914 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1917 for(Int_t i=0; i < (Int_t) kZones; i++) {
1918 fZc[i]=0; fZmax[i] = 0;
1923 if(fTimeBinIndex >= 0) {
1924 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1925 fIndex = new UInt_t[kMaxClusterPerTimeBin];
1928 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
1939 //_______________________________________________________
1940 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1941 Double_t Zmax, Double_t Ymax, Double_t rho,
1942 Double_t radLength, Double_t Yc, Double_t Zc)
1945 // Sets hole in the layer
1953 fHoleX0 = radLength;
1957 //_______________________________________________________
1958 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1961 // AliTRDtrackingSector Constructor
1970 // get holes description from geometry
1971 Bool_t holes[AliTRDgeometry::kNcham];
1972 //printf("sector\t%d\t",gs);
1973 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
1974 holes[icham] = fGeom->IsHole(0,icham,gs);
1975 //printf("%d",holes[icham]);
1979 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
1982 AliTRDpropagationLayer* ppl;
1984 Double_t x, xin, xout, dx, rho, radLength;
1987 // set time bins in the gas of the TPC
1989 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
1990 rho = 0.9e-3; radLength = 28.94;
1992 for(Int_t i=0; i<steps; i++) {
1993 x = xin + i*dx + dx/2;
1994 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
1998 // set time bins in the outer field cage vessel
2000 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2001 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2004 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2005 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2008 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2009 steps = 5; dx = (xout - xin)/steps;
2010 for(Int_t i=0; i<steps; i++) {
2011 x = xin + i*dx + dx/2;
2012 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2016 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2017 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2020 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2021 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2025 // set time bins in CO2
2027 xin = xout; xout = 275.0;
2028 steps = 50; dx = (xout - xin)/steps;
2029 rho = 1.977e-3; radLength = 36.2;
2031 for(Int_t i=0; i<steps; i++) {
2032 x = xin + i*dx + dx/2;
2033 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2037 // set time bins in the outer containment vessel
2039 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2040 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2043 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2044 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2047 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2048 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2051 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2052 steps = 10; dx = (xout - xin)/steps;
2053 for(Int_t i=0; i<steps; i++) {
2054 x = xin + i*dx + dx/2;
2055 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2059 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2060 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2063 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2064 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2067 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2068 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2071 Double_t xtrd = (Double_t) fGeom->Rmin();
2073 // add layers between TPC and TRD (Air temporarily)
2074 xin = xout; xout = xtrd;
2075 steps = 50; dx = (xout - xin)/steps;
2076 rho = 1.2e-3; radLength = 36.66;
2078 for(Int_t i=0; i<steps; i++) {
2079 x = xin + i*dx + dx/2;
2080 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2085 // Double_t alpha=AliTRDgeometry::GetAlpha();
2087 // add layers for each of the planes
2089 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2090 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2091 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2092 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2093 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2094 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2095 Double_t dxPlane = dxTEC + dxSpace;
2098 const Int_t kNchambers = AliTRDgeometry::Ncham();
2101 Double_t ymaxsensitive=0;
2102 Double_t *zc = new Double_t[kNchambers];
2103 Double_t *zmax = new Double_t[kNchambers];
2104 Double_t *zmaxsensitive = new Double_t[kNchambers];
2105 // Double_t holeZmax = 1000.; // the whole sector is missing
2107 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2110 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2111 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2112 for(Int_t i=0; i<steps; i++) {
2113 x = xin + i*dx + dx/2;
2114 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2118 ymax = fGeom->GetChamberWidth(plane)/2.;
2119 ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2121 for(Int_t ch = 0; ch < kNchambers; ch++) {
2122 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2123 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2124 Float_t row0 = fPar->GetRow0(plane,ch,0);
2125 Int_t nPads = fPar->GetRowMax(plane,ch,0);
2126 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2127 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2128 zc[ch] = (pad * nPads)/2 + row0;
2129 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2133 dx = fPar->GetTimeBinSize();
2134 rho = 0.00295 * 0.85; radLength = 11.0;
2136 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2137 Double_t xbottom = x0 - dxDrift;
2138 Double_t xtop = x0 + dxAmp;
2140 // Amplification region
2141 steps = (Int_t) (dxAmp/dx);
2143 for(tb = 0; tb < steps; tb++) {
2144 x = x0 + tb * dx + dx/2;
2145 tbIndex = CookTimeBinIndex(plane, -tb-1);
2146 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2147 ppl->SetYmax(ymax,ymaxsensitive);
2148 ppl->SetZ(zc, zmax, zmaxsensitive);
2149 ppl->SetHoles(holes);
2152 tbIndex = CookTimeBinIndex(plane, -steps);
2153 x = (x + dx/2 + xtop)/2;
2155 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2156 ppl->SetYmax(ymax,ymaxsensitive);
2157 ppl->SetZ(zc, zmax,zmaxsensitive);
2158 ppl->SetHoles(holes);
2162 dx = fPar->GetTimeBinSize();
2163 steps = (Int_t) (dxDrift/dx);
2165 for(tb = 0; tb < steps; tb++) {
2166 x = x0 - tb * dx - dx/2;
2167 tbIndex = CookTimeBinIndex(plane, tb);
2169 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2170 ppl->SetYmax(ymax,ymaxsensitive);
2171 ppl->SetZ(zc, zmax, zmaxsensitive);
2172 ppl->SetHoles(holes);
2175 tbIndex = CookTimeBinIndex(plane, steps);
2176 x = (x - dx/2 + xbottom)/2;
2178 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2179 ppl->SetYmax(ymax,ymaxsensitive);
2180 ppl->SetZ(zc, zmax, zmaxsensitive);
2181 ppl->SetHoles(holes);
2185 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2186 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2187 ppl->SetYmax(ymax,ymaxsensitive);
2188 ppl->SetZ(zc, zmax,zmax);
2189 ppl->SetHoles(holes);
2193 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2194 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2195 for(Int_t i=0; i<steps; i++) {
2196 x = xin + i*dx + dx/2;
2197 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2198 ppl->SetYmax(ymax,ymaxsensitive);
2199 ppl->SetZ(zc, zmax,zmax);
2200 ppl->SetHoles(holes);
2204 // Space between the chambers, air
2205 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2206 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2207 for(Int_t i=0; i<steps; i++) {
2208 x = xin + i*dx + dx/2;
2209 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2214 // Space between the TRD and RICH
2215 Double_t xRICH = 500.;
2216 xin = xout; xout = xRICH;
2217 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2218 for(Int_t i=0; i<steps; i++) {
2219 x = xin + i*dx + dx/2;
2220 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2230 //______________________________________________________
2232 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2235 // depending on the digitization parameters calculates "global"
2236 // time bin index for timebin <localTB> in plane <plane>
2239 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2240 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2241 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2243 Int_t tbAmp = fPar->GetTimeBefore();
2244 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2245 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2246 Int_t tbDrift = fPar->GetTimeMax();
2247 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2249 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2251 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2254 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2255 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2262 //______________________________________________________
2264 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2267 // For all sensitive time bins sets corresponding layer index
2268 // in the array fTimeBins
2273 for(Int_t i = 0; i < fN; i++) {
2274 index = fLayers[i]->GetTimeBinIndex();
2276 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2278 if(index < 0) continue;
2279 if(index >= (Int_t) kMaxTimeBinIndex) {
2280 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2281 printf(" index %d exceeds allowed maximum of %d!\n",
2282 index, kMaxTimeBinIndex-1);
2285 fTimeBinIndex[index] = i;
2288 Double_t x1, dx1, x2, dx2, gap;
2290 for(Int_t i = 0; i < fN-1; i++) {
2291 x1 = fLayers[i]->GetX();
2292 dx1 = fLayers[i]->GetdX();
2293 x2 = fLayers[i+1]->GetX();
2294 dx2 = fLayers[i+1]->GetdX();
2295 gap = (x2 - dx2/2) - (x1 + dx1/2);
2297 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2298 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2301 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2302 printf(" (%f - %f) - (%f + %f) = %f\n",
2303 x2, dx2/2, x1, dx1, gap);
2309 //______________________________________________________
2312 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2315 // Returns the number of time bin which in radial position is closest to <x>
2318 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2319 if(x <= fLayers[0]->GetX()) return 0;
2321 Int_t b=0, e=fN-1, m=(b+e)/2;
2322 for (; b<e; m=(b+e)/2) {
2323 if (x > fLayers[m]->GetX()) b=m+1;
2326 if(TMath::Abs(x - fLayers[m]->GetX()) >
2327 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2332 //______________________________________________________
2334 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2337 // Returns number of the innermost SENSITIVE propagation layer
2340 return GetLayerNumber(0);
2343 //______________________________________________________
2345 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2348 // Returns number of the outermost SENSITIVE time bin
2351 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2354 //______________________________________________________
2356 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2359 // Returns number of SENSITIVE time bins
2363 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2364 layer = GetLayerNumber(tb);
2370 //______________________________________________________
2372 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2375 // Insert layer <pl> in fLayers array.
2376 // Layers are sorted according to X coordinate.
2378 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2379 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2382 if (fN==0) {fLayers[fN++] = pl; return;}
2383 Int_t i=Find(pl->GetX());
2385 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2386 fLayers[i]=pl; fN++;
2390 //______________________________________________________
2392 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2395 // Returns index of the propagation layer nearest to X
2398 if (x <= fLayers[0]->GetX()) return 0;
2399 if (x > fLayers[fN-1]->GetX()) return fN;
2400 Int_t b=0, e=fN-1, m=(b+e)/2;
2401 for (; b<e; m=(b+e)/2) {
2402 if (x > fLayers[m]->GetX()) b=m+1;
2408 //______________________________________________________
2409 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2412 // set centers and the width of sectors
2413 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2414 fZc[icham] = center[icham];
2415 fZmax[icham] = w[icham];
2416 fZmaxSensitive[icham] = wsensitive[icham];
2417 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2420 //______________________________________________________
2422 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2425 // set centers and the width of sectors
2427 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2428 fIsHole[icham] = holes[icham];
2429 if (holes[icham]) fHole = kTRUE;
2435 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2436 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
2437 Bool_t &lookForCluster) const
2440 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2441 // and sensitivity <lookForCluster> in point <y,z>
2447 lookForCluster = kFALSE;
2449 // check dead regions in sensitive volume
2450 if(fTimeBinIndex >= 0) {
2453 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2454 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2456 lookForCluster = !(fIsHole[zone]);
2457 if(TMath::Abs(y) > fYmaxSensitive){
2458 lookForCluster = kFALSE;
2460 if (fIsHole[zone]) {
2472 if (fHole==kFALSE) return;
2474 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2475 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2486 Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2490 if (fTimeBinIndex < 0) return -20; //unknown
2491 Int_t zone=-10; // dead zone
2492 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2493 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2500 //______________________________________________________
2502 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2505 // Insert cluster in cluster array.
2506 // Clusters are sorted according to Y coordinate.
2508 if(fTimeBinIndex < 0) {
2509 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2513 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2514 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2517 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2518 Int_t i=Find(c->GetY());
2519 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2520 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2521 fIndex[i]=index; fClusters[i]=c; fN++;
2524 //______________________________________________________
2526 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2528 // Returns index of the cluster nearest in Y
2530 if (y <= fClusters[0]->GetY()) return 0;
2531 if (y > fClusters[fN-1]->GetY()) return fN;
2532 Int_t b=0, e=fN-1, m=(b+e)/2;
2533 for (; b<e; m=(b+e)/2) {
2534 if (y > fClusters[m]->GetY()) b=m+1;
2540 //---------------------------------------------------------
2542 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2544 // Returns correction factor for tilted pads geometry
2547 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2548 Int_t det = c->GetDetector();
2549 Int_t plane = fGeom->GetPlane(det);
2551 if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2553 if(fNoTilt) h01 = 0;