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);
551 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup); //make backup
553 Float_t p4 = track->GetC();
555 Int_t expectedClr = FollowBackProlongation(*track);
557 // only debug purpose
558 if (track->GetNumberOfClusters()<expectedClr/3){
559 AliTRDtrack *track1 = new AliTRDtrack(*seed);
560 track1->SetSeedLabel(lbl);
561 FollowBackProlongation(*track1);
562 AliTRDtrack *track2= new AliTRDtrack(*seed);
563 track->SetSeedLabel(lbl);
564 FollowBackProlongation(*track2);
569 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)>0.2) {
571 continue; //too big change of curvature - to be checked
574 Int_t foundClr = track->GetNumberOfClusters();
575 if (foundClr >= foundMin) {
577 CookLabel(track, 1-fgkLabelFraction);
578 if(track->GetChi2()/track->GetNumberOfClusters()<6) { // sign only gold tracks
581 Bool_t isGold = kFALSE;
583 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
584 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
587 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
588 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
591 if (!isGold && track->GetBackupTrack()){
592 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
593 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
594 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
603 //Propagation to the TOF (I.Belikov)
605 if (track->GetStop()==kFALSE){
608 Double_t c2=track->GetC()*xtof - track->GetEta();
609 if (TMath::Abs(c2)>=0.85) {
613 Double_t xTOF0 = 371. ;
614 PropagateToOuterPlane(*track,xTOF0);
616 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
617 Double_t y=track->GetYat(xtof);
619 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
623 } else if (y <-ymax) {
624 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
630 if (track->PropagateTo(xtof)) {
631 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
632 seed->SetTRDtrack(new AliTRDtrack(*track));
633 if (track->GetNumberOfClusters()>foundMin) found++;
636 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
637 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
638 //seed->SetStatus(AliESDtrack::kTRDStop);
639 seed->SetTRDtrack(new AliTRDtrack(*track));
646 //End of propagation to the TOF
647 //if (foundClr>foundMin)
648 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
653 cerr<<"Number of seeds: "<<fNseeds<<endl;
654 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
656 fSeeds->Clear(); fNseeds=0;
662 //_____________________________________________________________________________
663 Int_t AliTRDtracker::RefitInward(AliESD* event)
666 // Refits tracks within the TRD. The ESD event is expected to contain seeds
667 // at the outer part of the TRD.
668 // The tracks are propagated to the innermost time bin
669 // of the TRD and the ESD event is updated
670 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
673 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
674 Float_t foundMin = fgkMinClustersInTrack * timeBins;
677 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
679 Int_t n = event->GetNumberOfTracks();
680 for (Int_t i=0; i<n; i++) {
681 AliESDtrack* seed=event->GetTrack(i);
682 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
683 if (seed2->GetX()<270){
684 seed->UpdateTrackParams(seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
688 ULong_t status=seed->GetStatus();
689 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
690 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
692 seed2->ResetCovariance(5.);
693 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
694 UInt_t * indexes2 = seed2->GetIndexes();
695 UInt_t * indexes3 = pt->GetBackupIndexes();
696 for (Int_t i=0;i<200;i++) {
697 if (indexes2[i]==0) break;
698 indexes3[i] = indexes2[i];
700 //AliTRDtrack *pt = seed2;
702 FollowProlongation(t, innerTB);
704 if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
705 // debug - why we dont go back?
706 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
707 UInt_t * indexes2 = seed2->GetIndexes();
708 UInt_t * indexes3 = pt2->GetBackupIndexes();
709 for (Int_t i=0;i<200;i++) {
710 if (indexes2[i]==0) break;
711 indexes3[i] = indexes2[i];
713 FollowProlongation(*pt2, innerTB);
717 if (t.GetNumberOfClusters() >= foundMin) {
719 //CookLabel(pt, 1-fgkLabelFraction);
723 // cout<<found<<'\r';
725 if(PropagateToTPC(t)) {
726 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
728 //if not prolongation to TPC - propagate without update
729 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
730 seed2->ResetCovariance(5.);
731 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
733 if (PropagateToTPC(*pt2))
734 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
742 cout<<"Number of loaded seeds: "<<nseed<<endl;
743 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
750 //---------------------------------------------------------------------------
751 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
753 // Starting from current position on track=t this function tries
754 // to extrapolate the track up to timeBin=0 and to confirm prolongation
755 // if a close cluster is found. Returns the number of clusters
756 // expected to be found in sensitive layers
758 Float_t wIndex, wTB, wChi2;
759 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
760 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
761 Float_t wPx, wPy, wPz, wC;
763 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
764 Int_t lastplane = GetLastPlane(&t);
766 Int_t trackIndex = t.GetLabel();
768 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
770 Int_t tryAgain=fMaxGap;
772 Double_t alpha=t.GetAlpha();
773 alpha = TVector2::Phi_0_2pi(alpha);
775 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
776 Double_t radLength, rho, x, dx, y, ymax, z;
778 Int_t expectedNumberOfClusters = 0;
779 Bool_t lookForCluster;
781 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
784 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
786 y = t.GetY(); z = t.GetZ();
788 // first propagate to the inner surface of the current time bin
789 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
790 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
791 if(!t.PropagateTo(x,radLength,rho)) break;
793 ymax = x*TMath::Tan(0.5*alpha);
796 if (!t.Rotate(alpha)) break;
797 if(!t.PropagateTo(x,radLength,rho)) break;
798 } else if (y <-ymax) {
800 if (!t.Rotate(-alpha)) break;
801 if(!t.PropagateTo(x,radLength,rho)) break;
804 y = t.GetY(); z = t.GetZ();
806 // now propagate to the middle plane of the next time bin
807 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
808 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
809 if(!t.PropagateTo(x,radLength,rho)) break;
811 ymax = x*TMath::Tan(0.5*alpha);
814 if (!t.Rotate(alpha)) break;
815 if(!t.PropagateTo(x,radLength,rho)) break;
816 } else if (y <-ymax) {
818 if (!t.Rotate(-alpha)) break;
819 if(!t.PropagateTo(x,radLength,rho)) break;
825 expectedNumberOfClusters++;
826 wIndex = (Float_t) t.GetLabel();
829 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
831 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
832 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
835 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
836 else return expectedNumberOfClusters;
840 wYwindow = (Float_t) road;
841 t.GetPxPyPz(px,py,pz);
845 wC = (Float_t) t.GetC();
846 wSigmaC2 = (Float_t) t.GetSigmaC2();
847 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
848 wSigmaY2 = (Float_t) t.GetSigmaY2();
849 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
856 Double_t maxChi2=fgkMaxChi2;
858 wYclosest = 12345678;
859 wYcorrect = 12345678;
860 wZclosest = 12345678;
861 wZcorrect = 12345678;
862 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
864 // Find the closest correct cluster for debugging purposes
866 Float_t minDY = 1000000;
867 for (Int_t i=0; i<timeBin; i++) {
868 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
869 if((c->GetLabel(0) != trackIndex) &&
870 (c->GetLabel(1) != trackIndex) &&
871 (c->GetLabel(2) != trackIndex)) continue;
872 if(TMath::Abs(c->GetY() - y) > minDY) continue;
873 minDY = TMath::Abs(c->GetY() - y);
874 wYcorrect = c->GetY();
875 wZcorrect = c->GetZ();
877 Double_t h01 = GetTiltFactor(c);
878 wChi2 = t.GetPredictedChi2(c, h01);
882 // Now go for the real cluster search
886 //find cluster in history
889 AliTRDcluster * cl0 = timeBin[0];
893 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
894 if (plane>lastplane) continue;
895 Int_t timebin = cl0->GetLocalTimeBin();
896 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
899 Double_t h01 = GetTiltFactor(cl);
900 maxChi2=t.GetPredictedChi2(cl,h01);
902 if ((!cl) && road>fgkWideRoad) {
903 //if (t.GetNumberOfClusters()>4)
904 // cerr<<t.GetNumberOfClusters()
905 // <<"FindProlongation warning: Too broad road !\n";
912 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
913 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) > 3 * 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);
930 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
931 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
933 if (c->GetY() > y+road) break;
934 if (c->IsUsed() > 0) continue;
935 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
937 Double_t h01 = GetTiltFactor(c);
938 Double_t chi2=t.GetPredictedChi2(c, h01);
940 if (chi2 > maxChi2) continue;
943 index=timeBin.GetIndex(i);
947 wYclosest = cl->GetY();
948 wZclosest = cl->GetZ();
949 Double_t h01 = GetTiltFactor(cl);
951 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
952 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
953 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
954 Int_t det = cl->GetDetector();
955 Int_t plane = fGeom->GetPlane(det);
957 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
958 //if(!t.Update(cl,maxChi2,index,h01)) {
959 //if(!tryAgain--) return 0;
961 else tryAgain=fMaxGap;
964 //if (tryAgain==0) break;
969 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
971 printf(" %f", wIndex); //1
972 printf(" %f", wTB); //2
973 printf(" %f", wYrt); //3
974 printf(" %f", wYclosest); //4
975 printf(" %f", wYcorrect); //5
976 printf(" %f", wYwindow); //6
977 printf(" %f", wZrt); //7
978 printf(" %f", wZclosest); //8
979 printf(" %f", wZcorrect); //9
980 printf(" %f", wZwindow); //10
981 printf(" %f", wPx); //11
982 printf(" %f", wPy); //12
983 printf(" %f", wPz); //13
984 printf(" %f", wSigmaC2*1000000); //14
985 printf(" %f", wSigmaTgl2*1000); //15
986 printf(" %f", wSigmaY2); //16
987 // printf(" %f", wSigmaZ2); //17
988 printf(" %f", wChi2); //17
989 printf(" %f", wC); //18
996 return expectedNumberOfClusters;
1001 //___________________________________________________________________
1003 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
1005 // Starting from current radial position of track <t> this function
1006 // extrapolates the track up to outer timebin and in the sensitive
1007 // layers confirms prolongation if a close cluster is found.
1008 // Returns the number of clusters expected to be found in sensitive layers
1011 Float_t wIndex, wTB, wChi2;
1012 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
1013 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
1014 Float_t wPx, wPy, wPz, wC;
1015 Double_t px, py, pz;
1016 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
1018 Int_t trackIndex = t.GetLabel();
1019 Int_t tryAgain=fMaxGap;
1021 Double_t alpha=t.GetAlpha();
1022 TVector2::Phi_0_2pi(alpha);
1026 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
1027 Double_t radLength, rho, x, dx, y, ymax = 0, z;
1028 Bool_t lookForCluster;
1030 Int_t expectedNumberOfClusters = 0;
1033 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1035 Int_t nRefPlane = fgkFirstPlane;
1036 Bool_t isNewLayer = kFALSE;
1042 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
1047 // first propagate to the outer surface of the current time bin
1050 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1051 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
1055 if(!t.PropagateTo(x,radLength,rho)) break;
1056 // if (!AdjustSector(&t)) break;
1058 // MI -fix untill correct material desription will be implemented
1060 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material
1061 if (!AdjustSector(&t)) break;
1062 Int_t cross = kFALSE;
1064 if (TMath::Abs(angle - t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
1065 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1066 if (currentzone==-10) cross = kTRUE; // we are in the frame
1067 if (currentzone>-10){ // layer knows where we are
1068 if (zone==-10) zone = currentzone;
1069 if (zone!=currentzone) cross=kTRUE;
1073 if (t.GetNCross()==1) t.MakeBackupTrack();
1074 if (t.GetNCross()>2) break;
1080 if (!t.PropagateTo(x,radLength,rho)) break;
1085 // Barrel Tracks [SR, 04.04.2003]
1088 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1089 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1091 // if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1094 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1095 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1097 } else {isNewLayer = kFALSE;}
1102 // now propagate to the middle plane of the next time bin
1103 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1105 x = fTrSec[s]->GetLayer(nr+1)->GetX();
1106 if(!t.PropagateTo(x,radLength,rho)) break;
1107 if (!AdjustSector(&t)) break;
1109 if(!t.PropagateTo(x,radLength,rho)) break;
1114 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1115 // printf("label %d, pl %d, lookForCluster %d \n",
1116 // trackIndex, nr+1, lookForCluster);
1118 if(lookForCluster) {
1119 expectedNumberOfClusters++;
1121 wIndex = (Float_t) t.GetLabel();
1122 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1124 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1125 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1126 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1127 if((t.GetSigmaY2() + sy2) < 0) break;
1128 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1129 Double_t y=t.GetY(), z=t.GetZ();
1133 wYwindow = (Float_t) road;
1134 t.GetPxPyPz(px,py,pz);
1138 wC = (Float_t) t.GetC();
1139 wSigmaC2 = (Float_t) t.GetSigmaC2();
1140 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1141 wSigmaY2 = (Float_t) t.GetSigmaY2();
1142 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1145 if (road>fgkWideRoad) {
1146 if (t.GetNumberOfClusters()>4)
1147 cerr<<t.GetNumberOfClusters()
1148 <<"FindProlongation warning: Too broad road !\n";
1152 AliTRDcluster *cl=0;
1155 Double_t maxChi2=fgkMaxChi2;
1160 maxChi2 = 10 * fgkMaxChi2;
1163 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1164 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1165 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1168 wYclosest = 12345678;
1169 wYcorrect = 12345678;
1170 wZclosest = 12345678;
1171 wZcorrect = 12345678;
1172 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1174 // Find the closest correct cluster for debugging purposes
1177 for (Int_t i=0; i<timeBin; i++) {
1178 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1179 if((c->GetLabel(0) != trackIndex) &&
1180 (c->GetLabel(1) != trackIndex) &&
1181 (c->GetLabel(2) != trackIndex)) continue;
1182 if(TMath::Abs(c->GetY() - y) > minDY) continue;
1183 //minDY = TMath::Abs(c->GetY() - y);
1184 minDY = c->GetY() - y;
1185 wYcorrect = c->GetY();
1186 wZcorrect = c->GetZ();
1188 Double_t h01 = GetTiltFactor(c);
1189 wChi2 = t.GetPredictedChi2(c, h01);
1193 // Now go for the real cluster search
1197 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1198 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1199 if (c->GetY() > y+road) break;
1200 if (c->IsUsed() > 0) continue;
1201 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1203 Double_t h01 = GetTiltFactor(c);
1204 chi2=t.GetPredictedChi2(c,h01);
1206 if (chi2 > maxChi2) continue;
1209 index=timeBin.GetIndex(i);
1212 if((c->GetLabel(0) != trackIndex) &&
1213 (c->GetLabel(1) != trackIndex) &&
1214 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1219 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1220 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1222 if (c->GetY() > y+road) break;
1223 if (c->IsUsed() > 0) continue;
1224 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1226 Double_t h01 = GetTiltFactor(c);
1227 chi2=t.GetPredictedChi2(c,h01);
1229 if (chi2 > maxChi2) continue;
1232 index=timeBin.GetIndex(i);
1237 wYclosest = cl->GetY();
1238 wZclosest = cl->GetZ();
1240 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1241 Double_t h01 = GetTiltFactor(cl);
1242 Int_t det = cl->GetDetector();
1243 Int_t plane = fGeom->GetPlane(det);
1245 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1246 //if(!t.Update(cl,maxChi2,index,h01)) {
1247 if(!tryAgain--) return 0;
1249 else tryAgain=fMaxGap;
1252 if (tryAgain==0) break;
1255 //if (minDY < 1000000 && isNewLayer)
1256 //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() << "\t" <<
1257 // road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
1261 isNewLayer = kFALSE;
1264 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1266 printf(" %f", wIndex); //1
1267 printf(" %f", wTB); //2
1268 printf(" %f", wYrt); //3
1269 printf(" %f", wYclosest); //4
1270 printf(" %f", wYcorrect); //5
1271 printf(" %f", wYwindow); //6
1272 printf(" %f", wZrt); //7
1273 printf(" %f", wZclosest); //8
1274 printf(" %f", wZcorrect); //9
1275 printf(" %f", wZwindow); //10
1276 printf(" %f", wPx); //11
1277 printf(" %f", wPy); //12
1278 printf(" %f", wPz); //13
1279 printf(" %f", wSigmaC2*1000000); //14
1280 printf(" %f", wSigmaTgl2*1000); //15
1281 printf(" %f", wSigmaY2); //16
1282 // printf(" %f", wSigmaZ2); //17
1283 printf(" %f", wChi2); //17
1284 printf(" %f", wC); //18
1295 return expectedNumberOfClusters;
1300 //---------------------------------------------------------------------------
1301 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1303 // Starting from current position on track=t this function tries
1304 // to extrapolate the track up to timeBin=0 and to reuse already
1305 // assigned clusters. Returns the number of clusters
1306 // expected to be found in sensitive layers
1307 // get indices of assigned clusters for each layer
1308 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1311 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1312 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1313 Int_t index = t.GetClusterIndex(i);
1314 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1316 Int_t detector=cl->GetDetector();
1317 Int_t localTimeBin=cl->GetLocalTimeBin();
1318 Int_t sector=fGeom->GetSector(detector);
1319 Int_t plane=fGeom->GetPlane(detector);
1321 Int_t trackingSector = CookSectorIndex(sector);
1323 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1324 if(gtb < 0) continue;
1325 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1326 iCluster[layer] = index;
1330 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1332 Double_t alpha=t.GetAlpha();
1333 alpha = TVector2::Phi_0_2pi(alpha);
1335 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1336 Double_t radLength, rho, x, dx, y, ymax, z;
1338 Int_t expectedNumberOfClusters = 0;
1339 Bool_t lookForCluster;
1341 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1344 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1346 y = t.GetY(); z = t.GetZ();
1348 // first propagate to the inner surface of the current time bin
1349 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1350 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1351 if(!t.PropagateTo(x,radLength,rho)) break;
1353 ymax = x*TMath::Tan(0.5*alpha);
1356 if (!t.Rotate(alpha)) break;
1357 if(!t.PropagateTo(x,radLength,rho)) break;
1358 } else if (y <-ymax) {
1360 if (!t.Rotate(-alpha)) break;
1361 if(!t.PropagateTo(x,radLength,rho)) break;
1364 y = t.GetY(); z = t.GetZ();
1366 // now propagate to the middle plane of the next time bin
1367 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1368 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1369 if(!t.PropagateTo(x,radLength,rho)) break;
1371 ymax = x*TMath::Tan(0.5*alpha);
1374 if (!t.Rotate(alpha)) break;
1375 if(!t.PropagateTo(x,radLength,rho)) break;
1376 } else if (y <-ymax) {
1378 if (!t.Rotate(-alpha)) break;
1379 if(!t.PropagateTo(x,radLength,rho)) break;
1382 if(lookForCluster) expectedNumberOfClusters++;
1384 // use assigned cluster
1385 if (!iCluster[nr-1]) continue;
1386 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1387 Double_t h01 = GetTiltFactor(cl);
1388 Double_t chi2=t.GetPredictedChi2(cl, h01);
1389 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1390 t.Update(cl,chi2,iCluster[nr-1],h01);
1393 return expectedNumberOfClusters;
1396 //___________________________________________________________________
1398 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1400 // Starting from current radial position of track <t> this function
1401 // extrapolates the track up to radial position <xToGo>.
1402 // Returns 1 if track reaches the plane, and 0 otherwise
1404 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1406 Double_t alpha=t.GetAlpha();
1408 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1409 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1411 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1413 Bool_t lookForCluster;
1414 Double_t radLength, rho, x, dx, y, ymax, z;
1418 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1420 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1422 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1424 y = t.GetY(); z = t.GetZ();
1426 // first propagate to the outer surface of the current time bin
1427 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1428 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1429 if(!t.PropagateTo(x,radLength,rho)) return 0;
1431 ymax = x*TMath::Tan(0.5*alpha);
1434 if (!t.Rotate(alpha)) return 0;
1435 } else if (y <-ymax) {
1437 if (!t.Rotate(-alpha)) return 0;
1439 if(!t.PropagateTo(x,radLength,rho)) return 0;
1441 y = t.GetY(); z = t.GetZ();
1443 // now propagate to the middle plane of the next time bin
1444 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1445 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1446 if(!t.PropagateTo(x,radLength,rho)) return 0;
1448 ymax = x*TMath::Tan(0.5*alpha);
1451 if (!t.Rotate(alpha)) return 0;
1452 } else if (y <-ymax) {
1454 if (!t.Rotate(-alpha)) return 0;
1456 if(!t.PropagateTo(x,radLength,rho)) return 0;
1461 //___________________________________________________________________
1463 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1465 // Starting from current radial position of track <t> this function
1466 // extrapolates the track up to radial position of the outermost
1467 // padrow of the TPC.
1468 // Returns 1 if track reaches the TPC, and 0 otherwise
1470 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1472 Double_t alpha=t.GetAlpha();
1473 alpha = TVector2::Phi_0_2pi(alpha);
1475 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1477 Bool_t lookForCluster;
1478 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1482 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1483 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1485 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1490 // first propagate to the outer surface of the current time bin
1491 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1492 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1494 if(!t.PropagateTo(x,radLength,rho)) return 0;
1496 if(!t.PropagateTo(x,radLength,rho)) return 0;
1501 // now propagate to the middle plane of the next time bin
1502 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1503 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1505 if(!t.PropagateTo(x,radLength,rho)) return 0;
1507 if(!t.PropagateTo(x,radLength,rho)) return 0;
1512 //_____________________________________________________________________________
1513 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1515 // Fills clusters into TRD tracking_sectors
1516 // Note that the numbering scheme for the TRD tracking_sectors
1517 // differs from that of TRD sectors
1519 if (ReadClusters(fClusters,cTree)) {
1520 Error("LoadClusters","Problem with reading the clusters !");
1523 Int_t ncl=fClusters->GetEntriesFast();
1525 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1528 for (Int_t ichamber=0;ichamber<5;ichamber++)
1529 for (Int_t isector=0;isector<18;isector++){
1530 fHoles[ichamber][isector]=kTRUE;
1535 // printf("\r %d left ",ncl);
1536 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1537 Int_t detector=c->GetDetector();
1538 Int_t localTimeBin=c->GetLocalTimeBin();
1539 Int_t sector=fGeom->GetSector(detector);
1540 Int_t plane=fGeom->GetPlane(detector);
1542 Int_t trackingSector = CookSectorIndex(sector);
1543 if (c->GetLabel(0)>0){
1544 Int_t chamber = fGeom->GetChamber(detector);
1545 fHoles[chamber][trackingSector]=kFALSE;
1548 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1549 if(gtb < 0) continue;
1550 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1553 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1559 for (Int_t isector=0;isector<18;isector++){
1560 for (Int_t ichamber=0;ichamber<5;ichamber++)
1561 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1562 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1563 fGeom->IsHole(0,ichamber,17-isector));
1569 //_____________________________________________________________________________
1570 void AliTRDtracker::UnloadClusters()
1573 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1578 nentr = fClusters->GetEntriesFast();
1579 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1582 nentr = fSeeds->GetEntriesFast();
1583 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1585 nentr = fTracks->GetEntriesFast();
1586 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1588 Int_t nsec = AliTRDgeometry::kNsect;
1590 for (i = 0; i < nsec; i++) {
1591 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1592 fTrSec[i]->GetLayer(pl)->Clear();
1598 //__________________________________________________________________________
1599 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1601 // Creates track seeds using clusters in timeBins=i1,i2
1604 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1608 Double_t x[5], c[15];
1609 Int_t maxSec=AliTRDgeometry::kNsect;
1611 Double_t alpha=AliTRDgeometry::GetAlpha();
1612 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1613 Double_t cs=cos(alpha), sn=sin(alpha);
1614 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1617 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1618 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1620 Double_t x1 =fTrSec[0]->GetX(i1);
1621 Double_t xx2=fTrSec[0]->GetX(i2);
1623 for (Int_t ns=0; ns<maxSec; ns++) {
1625 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1626 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1627 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1628 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1629 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1631 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1633 for (Int_t is=0; is < r1; is++) {
1634 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1636 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1638 const AliTRDcluster *cl;
1639 Double_t x2, y2, z2;
1640 Double_t x3=0., y3=0.;
1643 if(turn != 2) continue;
1644 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1646 y2=cl->GetY(); z2=cl->GetZ();
1651 else if (js<nl2+nl) {
1652 if(turn != 1) continue;
1653 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1655 y2=cl->GetY(); z2=cl->GetZ();
1660 else if (js<nl2+nl+nm) {
1661 if(turn != 1) continue;
1662 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1664 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1666 else if (js<nl2+nl+nm+nu) {
1667 if(turn != 1) continue;
1668 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1669 cl=r2[js-nl2-nl-nm];
1670 y2=cl->GetY(); z2=cl->GetZ();
1676 if(turn != 2) continue;
1677 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1678 cl=r2[js-nl2-nl-nm-nu];
1679 y2=cl->GetY(); z2=cl->GetZ();
1685 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1687 Double_t zz=z1 - z1/x1*(x1-x2);
1689 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1691 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1692 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1696 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1698 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
1700 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1702 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1704 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1706 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1708 Double_t a=asin(x[2]);
1709 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1711 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1713 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1714 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1715 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
1718 Double_t h01 = GetTiltFactor(r1[is]);
1719 Double_t xuFactor = 100.;
1725 sy1=sy1+sz1*h01*h01;
1726 Double_t syz=sz1*(-h01);
1727 // end of tilt changes
1729 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1730 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1731 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1732 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1733 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1734 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1735 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1736 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1737 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1738 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1742 // c[1]=0.; c[2]=sz1;
1743 c[1]=syz; c[2]=sz1*xuFactor;
1744 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1745 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1746 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1747 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1748 c[13]=f30*sy1*f40+f32*sy2*f42;
1749 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1751 UInt_t index=r1.GetIndex(is);
1753 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1755 Int_t rc=FollowProlongation(*track, i2);
1758 (track->GetNumberOfClusters() <
1759 (outer-inner)*fgkMinClustersInSeed)) delete track;
1761 fSeeds->AddLast(track); fNseeds++;
1762 // cerr<<"\r found seed "<<fNseeds;
1769 //_____________________________________________________________________________
1770 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1773 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1774 // from the file. The names of the cluster tree and branches
1775 // should match the ones used in AliTRDclusterizer::WriteClusters()
1777 TObjArray *clusterArray = new TObjArray(400);
1779 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1781 Error("ReadClusters","Can't get the branch !");
1784 branch->SetAddress(&clusterArray);
1786 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1787 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1789 // Loop through all entries in the tree
1791 AliTRDcluster *c = 0;
1794 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1797 nbytes += ClusterTree->GetEvent(iEntry);
1799 // Get the number of points in the detector
1800 Int_t nCluster = clusterArray->GetEntriesFast();
1801 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1803 // Loop through all TRD digits
1804 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1805 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1806 AliTRDcluster *co = new AliTRDcluster(*c);
1807 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1808 Int_t ltb = co->GetLocalTimeBin();
1809 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1810 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1812 delete clusterArray->RemoveAt(iCluster);
1816 delete clusterArray;
1821 //__________________________________________________________________
1822 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1825 // This cooks a label. Mmmmh, smells good...
1828 Int_t label=123456789, index, i, j;
1829 Int_t ncl=pt->GetNumberOfClusters();
1830 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1834 // Int_t s[kRange][2];
1835 Int_t **s = new Int_t* [kRange];
1836 for (i=0; i<kRange; i++) {
1837 s[i] = new Int_t[2];
1839 for (i=0; i<kRange; i++) {
1845 for (i=0; i<ncl; i++) {
1846 index=pt->GetClusterIndex(i);
1847 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1853 for (i=0; i<ncl; i++) {
1854 index=pt->GetClusterIndex(i);
1855 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1856 for (Int_t k=0; k<3; k++) {
1857 label=c->GetLabel(k);
1858 labelAdded=kFALSE; j=0;
1860 while ( (!labelAdded) && ( j < kRange ) ) {
1861 if (s[j][0]==label || s[j][1]==0) {
1875 for (i=0; i<kRange; i++) {
1877 max=s[i][1]; label=s[i][0];
1881 for (i=0; i<kRange; i++) {
1887 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1889 pt->SetLabel(label);
1894 //__________________________________________________________________
1895 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1898 // Use clusters, but don't abuse them!
1901 Int_t ncl=t->GetNumberOfClusters();
1902 for (Int_t i=from; i<ncl; i++) {
1903 Int_t index = t->GetClusterIndex(i);
1904 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1910 //_____________________________________________________________________
1911 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1913 // Parametrised "expected" error of the cluster reconstruction in Y
1915 Double_t s = 0.08 * 0.08;
1919 //_____________________________________________________________________
1920 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1922 // Parametrised "expected" error of the cluster reconstruction in Z
1924 Double_t s = 9 * 9 /12.;
1928 //_____________________________________________________________________
1929 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
1932 // Returns radial position which corresponds to time bin <localTB>
1933 // in tracking sector <sector> and plane <plane>
1936 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
1937 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1938 return fTrSec[sector]->GetLayer(pl)->GetX();
1943 //_______________________________________________________
1944 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1945 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1948 // AliTRDpropagationLayer constructor
1951 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1952 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1955 for(Int_t i=0; i < (Int_t) kZones; i++) {
1956 fZc[i]=0; fZmax[i] = 0;
1961 if(fTimeBinIndex >= 0) {
1962 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1963 fIndex = new UInt_t[kMaxClusterPerTimeBin];
1966 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
1977 //_______________________________________________________
1978 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1979 Double_t Zmax, Double_t Ymax, Double_t rho,
1980 Double_t radLength, Double_t Yc, Double_t Zc)
1983 // Sets hole in the layer
1991 fHoleX0 = radLength;
1995 //_______________________________________________________
1996 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1999 // AliTRDtrackingSector Constructor
2008 // get holes description from geometry
2009 Bool_t holes[AliTRDgeometry::kNcham];
2010 //printf("sector\t%d\t",gs);
2011 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
2012 holes[icham] = fGeom->IsHole(0,icham,gs);
2013 //printf("%d",holes[icham]);
2017 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
2020 AliTRDpropagationLayer* ppl;
2022 Double_t x, xin, xout, dx, rho, radLength;
2025 // set time bins in the gas of the TPC
2027 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
2028 rho = 0.9e-3; radLength = 28.94;
2030 for(Int_t i=0; i<steps; i++) {
2031 x = xin + i*dx + dx/2;
2032 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2036 // set time bins in the outer field cage vessel
2038 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2039 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2042 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2043 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2046 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2047 steps = 5; dx = (xout - xin)/steps;
2048 for(Int_t i=0; i<steps; i++) {
2049 x = xin + i*dx + dx/2;
2050 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2054 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2055 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2058 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2059 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2063 // set time bins in CO2
2065 xin = xout; xout = 275.0;
2066 steps = 50; dx = (xout - xin)/steps;
2067 rho = 1.977e-3; radLength = 36.2;
2069 for(Int_t i=0; i<steps; i++) {
2070 x = xin + i*dx + dx/2;
2071 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2075 // set time bins in the outer containment vessel
2077 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2078 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2081 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2082 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2085 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2086 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2089 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2090 steps = 10; dx = (xout - xin)/steps;
2091 for(Int_t i=0; i<steps; i++) {
2092 x = xin + i*dx + dx/2;
2093 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2097 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2098 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2101 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2102 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2105 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2106 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2109 Double_t xtrd = (Double_t) fGeom->Rmin();
2111 // add layers between TPC and TRD (Air temporarily)
2112 xin = xout; xout = xtrd;
2113 steps = 50; dx = (xout - xin)/steps;
2114 rho = 1.2e-3; radLength = 36.66;
2116 for(Int_t i=0; i<steps; i++) {
2117 x = xin + i*dx + dx/2;
2118 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2123 // Double_t alpha=AliTRDgeometry::GetAlpha();
2125 // add layers for each of the planes
2127 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2128 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2129 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2130 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2131 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2132 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2133 Double_t dxPlane = dxTEC + dxSpace;
2136 const Int_t kNchambers = AliTRDgeometry::Ncham();
2139 Double_t ymaxsensitive=0;
2140 Double_t *zc = new Double_t[kNchambers];
2141 Double_t *zmax = new Double_t[kNchambers];
2142 Double_t *zmaxsensitive = new Double_t[kNchambers];
2143 // Double_t holeZmax = 1000.; // the whole sector is missing
2145 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2148 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2149 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2150 for(Int_t i=0; i<steps; i++) {
2151 x = xin + i*dx + dx/2;
2152 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2156 ymax = fGeom->GetChamberWidth(plane)/2.;
2157 ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2159 for(Int_t ch = 0; ch < kNchambers; ch++) {
2160 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2161 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2162 Float_t row0 = fPar->GetRow0(plane,ch,0);
2163 Int_t nPads = fPar->GetRowMax(plane,ch,0);
2164 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2165 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2166 zc[ch] = (pad * nPads)/2 + row0;
2167 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2171 dx = fPar->GetTimeBinSize();
2172 rho = 0.00295 * 0.85; radLength = 11.0;
2174 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2175 Double_t xbottom = x0 - dxDrift;
2176 Double_t xtop = x0 + dxAmp;
2178 // Amplification region
2179 steps = (Int_t) (dxAmp/dx);
2181 for(tb = 0; tb < steps; tb++) {
2182 x = x0 + tb * dx + dx/2;
2183 tbIndex = CookTimeBinIndex(plane, -tb-1);
2184 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2185 ppl->SetYmax(ymax,ymaxsensitive);
2186 ppl->SetZ(zc, zmax, zmaxsensitive);
2187 ppl->SetHoles(holes);
2190 tbIndex = CookTimeBinIndex(plane, -steps);
2191 x = (x + dx/2 + xtop)/2;
2193 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2194 ppl->SetYmax(ymax,ymaxsensitive);
2195 ppl->SetZ(zc, zmax,zmaxsensitive);
2196 ppl->SetHoles(holes);
2200 dx = fPar->GetTimeBinSize();
2201 steps = (Int_t) (dxDrift/dx);
2203 for(tb = 0; tb < steps; tb++) {
2204 x = x0 - tb * dx - dx/2;
2205 tbIndex = CookTimeBinIndex(plane, tb);
2207 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2208 ppl->SetYmax(ymax,ymaxsensitive);
2209 ppl->SetZ(zc, zmax, zmaxsensitive);
2210 ppl->SetHoles(holes);
2213 tbIndex = CookTimeBinIndex(plane, steps);
2214 x = (x - dx/2 + xbottom)/2;
2216 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2217 ppl->SetYmax(ymax,ymaxsensitive);
2218 ppl->SetZ(zc, zmax, zmaxsensitive);
2219 ppl->SetHoles(holes);
2223 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2224 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2225 ppl->SetYmax(ymax,ymaxsensitive);
2226 ppl->SetZ(zc, zmax,zmax);
2227 ppl->SetHoles(holes);
2231 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2232 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2233 for(Int_t i=0; i<steps; i++) {
2234 x = xin + i*dx + dx/2;
2235 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2236 ppl->SetYmax(ymax,ymaxsensitive);
2237 ppl->SetZ(zc, zmax,zmax);
2238 ppl->SetHoles(holes);
2242 // Space between the chambers, air
2243 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2244 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2245 for(Int_t i=0; i<steps; i++) {
2246 x = xin + i*dx + dx/2;
2247 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2252 // Space between the TRD and RICH
2253 Double_t xRICH = 500.;
2254 xin = xout; xout = xRICH;
2255 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2256 for(Int_t i=0; i<steps; i++) {
2257 x = xin + i*dx + dx/2;
2258 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2268 //______________________________________________________
2270 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2273 // depending on the digitization parameters calculates "global"
2274 // time bin index for timebin <localTB> in plane <plane>
2277 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2278 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2279 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2281 Int_t tbAmp = fPar->GetTimeBefore();
2282 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2283 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2284 Int_t tbDrift = fPar->GetTimeMax();
2285 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2287 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2289 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2292 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2293 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2300 //______________________________________________________
2302 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2305 // For all sensitive time bins sets corresponding layer index
2306 // in the array fTimeBins
2311 for(Int_t i = 0; i < fN; i++) {
2312 index = fLayers[i]->GetTimeBinIndex();
2314 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2316 if(index < 0) continue;
2317 if(index >= (Int_t) kMaxTimeBinIndex) {
2318 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2319 printf(" index %d exceeds allowed maximum of %d!\n",
2320 index, kMaxTimeBinIndex-1);
2323 fTimeBinIndex[index] = i;
2326 Double_t x1, dx1, x2, dx2, gap;
2328 for(Int_t i = 0; i < fN-1; i++) {
2329 x1 = fLayers[i]->GetX();
2330 dx1 = fLayers[i]->GetdX();
2331 x2 = fLayers[i+1]->GetX();
2332 dx2 = fLayers[i+1]->GetdX();
2333 gap = (x2 - dx2/2) - (x1 + dx1/2);
2335 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2336 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2339 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2340 printf(" (%f - %f) - (%f + %f) = %f\n",
2341 x2, dx2/2, x1, dx1, gap);
2347 //______________________________________________________
2350 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2353 // Returns the number of time bin which in radial position is closest to <x>
2356 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2357 if(x <= fLayers[0]->GetX()) return 0;
2359 Int_t b=0, e=fN-1, m=(b+e)/2;
2360 for (; b<e; m=(b+e)/2) {
2361 if (x > fLayers[m]->GetX()) b=m+1;
2364 if(TMath::Abs(x - fLayers[m]->GetX()) >
2365 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2370 //______________________________________________________
2372 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2375 // Returns number of the innermost SENSITIVE propagation layer
2378 return GetLayerNumber(0);
2381 //______________________________________________________
2383 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2386 // Returns number of the outermost SENSITIVE time bin
2389 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2392 //______________________________________________________
2394 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2397 // Returns number of SENSITIVE time bins
2401 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2402 layer = GetLayerNumber(tb);
2408 //______________________________________________________
2410 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2413 // Insert layer <pl> in fLayers array.
2414 // Layers are sorted according to X coordinate.
2416 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2417 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2420 if (fN==0) {fLayers[fN++] = pl; return;}
2421 Int_t i=Find(pl->GetX());
2423 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2424 fLayers[i]=pl; fN++;
2428 //______________________________________________________
2430 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2433 // Returns index of the propagation layer nearest to X
2436 if (x <= fLayers[0]->GetX()) return 0;
2437 if (x > fLayers[fN-1]->GetX()) return fN;
2438 Int_t b=0, e=fN-1, m=(b+e)/2;
2439 for (; b<e; m=(b+e)/2) {
2440 if (x > fLayers[m]->GetX()) b=m+1;
2446 //______________________________________________________
2447 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2450 // set centers and the width of sectors
2451 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2452 fZc[icham] = center[icham];
2453 fZmax[icham] = w[icham];
2454 fZmaxSensitive[icham] = wsensitive[icham];
2455 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2458 //______________________________________________________
2460 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2463 // set centers and the width of sectors
2465 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2466 fIsHole[icham] = holes[icham];
2467 if (holes[icham]) fHole = kTRUE;
2473 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2474 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
2475 Bool_t &lookForCluster) const
2478 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2479 // and sensitivity <lookForCluster> in point <y,z>
2485 lookForCluster = kFALSE;
2487 // check dead regions in sensitive volume
2488 if(fTimeBinIndex >= 0) {
2491 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2492 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2494 lookForCluster = !(fIsHole[zone]);
2495 if(TMath::Abs(y) > fYmaxSensitive){
2496 lookForCluster = kFALSE;
2498 if (fIsHole[zone]) {
2510 if (fHole==kFALSE) return;
2512 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2513 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2524 Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2528 if (fTimeBinIndex < 0) return -20; //unknown
2529 Int_t zone=-10; // dead zone
2530 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2531 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2538 //______________________________________________________
2540 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2543 // Insert cluster in cluster array.
2544 // Clusters are sorted according to Y coordinate.
2546 if(fTimeBinIndex < 0) {
2547 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2551 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2552 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2555 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2556 Int_t i=Find(c->GetY());
2557 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2558 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2559 fIndex[i]=index; fClusters[i]=c; fN++;
2562 //______________________________________________________
2564 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2566 // Returns index of the cluster nearest in Y
2568 if (y <= fClusters[0]->GetY()) return 0;
2569 if (y > fClusters[fN-1]->GetY()) return fN;
2570 Int_t b=0, e=fN-1, m=(b+e)/2;
2571 for (; b<e; m=(b+e)/2) {
2572 if (y > fClusters[m]->GetY()) b=m+1;
2578 //---------------------------------------------------------
2580 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2582 // Returns correction factor for tilted pads geometry
2585 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2586 Int_t det = c->GetDetector();
2587 Int_t plane = fGeom->GetPlane(det);
2589 //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2590 if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
2592 if(fNoTilt) h01 = 0;