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) {
576 track->CookdEdx(0.,1.);
577 CookdEdxTimBin(*track);
579 CookLabel(track, 1-fgkLabelFraction);
580 if(track->GetChi2()/track->GetNumberOfClusters()<6) { // sign only gold tracks
583 Bool_t isGold = kFALSE;
585 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
586 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
589 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
590 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
593 if (!isGold && track->GetBackupTrack()){
594 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
595 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
596 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
606 if (track->GetStop()==kFALSE){
609 Double_t c2=track->GetC()*xtof - track->GetEta();
610 if (TMath::Abs(c2)>=0.85) {
614 Double_t xTOF0 = 371. ;
615 PropagateToOuterPlane(*track,xTOF0);
617 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
618 Double_t y=track->GetYat(xtof);
620 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
624 } else if (y <-ymax) {
625 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
631 if (track->PropagateTo(xtof)) {
632 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
633 for (Int_t i=0;i<kNPlane;i++) {
634 seed->SetTRDsignals(track->GetPIDsignals(i),i);
635 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
637 seed->SetTRDtrack(new AliTRDtrack(*track));
638 if (track->GetNumberOfClusters()>foundMin) found++;
641 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
642 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
643 //seed->SetStatus(AliESDtrack::kTRDStop);
644 for (Int_t i=0;i<kNPlane;i++) {
645 seed->SetTRDsignals(track->GetPIDsignals(i),i);
646 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
648 seed->SetTRDtrack(new AliTRDtrack(*track));
655 //End of propagation to the TOF
656 //if (foundClr>foundMin)
657 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
662 cerr<<"Number of seeds: "<<fNseeds<<endl;
663 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
665 fSeeds->Clear(); fNseeds=0;
671 //_____________________________________________________________________________
672 Int_t AliTRDtracker::RefitInward(AliESD* event)
675 // Refits tracks within the TRD. The ESD event is expected to contain seeds
676 // at the outer part of the TRD.
677 // The tracks are propagated to the innermost time bin
678 // of the TRD and the ESD event is updated
679 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
682 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
683 Float_t foundMin = fgkMinClustersInTrack * timeBins;
686 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
688 Int_t n = event->GetNumberOfTracks();
689 for (Int_t i=0; i<n; i++) {
690 AliESDtrack* seed=event->GetTrack(i);
691 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
692 if (seed2->GetX()<270){
693 seed->UpdateTrackParams(seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
697 ULong_t status=seed->GetStatus();
698 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
699 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
701 seed2->ResetCovariance(5.);
702 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
703 for (Int_t i=0;i<kNPlane;i++) {
704 pt->SetPIDsignals(seed2->GetPIDsignals(i),i);
705 pt->SetPIDTimBin(seed2->GetPIDTimBin(i),i);
708 UInt_t * indexes2 = seed2->GetIndexes();
709 UInt_t * indexes3 = pt->GetBackupIndexes();
710 for (Int_t i=0;i<200;i++) {
711 if (indexes2[i]==0) break;
712 indexes3[i] = indexes2[i];
714 //AliTRDtrack *pt = seed2;
716 FollowProlongation(t, innerTB);
718 if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
719 // debug - why we dont go back?
720 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
721 UInt_t * indexes2 = seed2->GetIndexes();
722 UInt_t * indexes3 = pt2->GetBackupIndexes();
723 for (Int_t i=0;i<200;i++) {
724 if (indexes2[i]==0) break;
725 indexes3[i] = indexes2[i];
727 FollowProlongation(*pt2, innerTB);
731 if (t.GetNumberOfClusters() >= foundMin) {
733 //CookLabel(pt, 1-fgkLabelFraction);
737 // cout<<found<<'\r';
739 if(PropagateToTPC(t)) {
740 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
741 for (Int_t i=0;i<kNPlane;i++) {
742 seed->SetTRDsignals(pt->GetPIDsignals(i),i);
743 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
746 //if not prolongation to TPC - propagate without update
747 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
748 seed2->ResetCovariance(5.);
749 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
751 if (PropagateToTPC(*pt2)) {
752 pt2->CookdEdx(0.,1.);
753 CookdEdxTimBin(*pt2);
754 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
755 for (Int_t i=0;i<kNPlane;i++) {
756 seed->SetTRDsignals(seed2->GetPIDsignals(i),i);
757 seed->SetTRDTimBin(seed2->GetPIDTimBin(i),i);
767 cout<<"Number of loaded seeds: "<<nseed<<endl;
768 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
775 //---------------------------------------------------------------------------
776 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
778 // Starting from current position on track=t this function tries
779 // to extrapolate the track up to timeBin=0 and to confirm prolongation
780 // if a close cluster is found. Returns the number of clusters
781 // expected to be found in sensitive layers
783 Float_t wIndex, wTB, wChi2;
784 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
785 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
786 Float_t wPx, wPy, wPz, wC;
788 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
789 Int_t lastplane = GetLastPlane(&t);
791 Int_t trackIndex = t.GetLabel();
793 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
795 Int_t tryAgain=fMaxGap;
797 Double_t alpha=t.GetAlpha();
798 alpha = TVector2::Phi_0_2pi(alpha);
800 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
801 Double_t radLength, rho, x, dx, y, ymax, z;
803 Int_t expectedNumberOfClusters = 0;
804 Bool_t lookForCluster;
806 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
809 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
811 y = t.GetY(); z = t.GetZ();
813 // first propagate to the inner surface of the current time bin
814 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
815 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
816 if(!t.PropagateTo(x,radLength,rho)) break;
818 ymax = x*TMath::Tan(0.5*alpha);
821 if (!t.Rotate(alpha)) break;
822 if(!t.PropagateTo(x,radLength,rho)) break;
823 } else if (y <-ymax) {
825 if (!t.Rotate(-alpha)) break;
826 if(!t.PropagateTo(x,radLength,rho)) break;
829 y = t.GetY(); z = t.GetZ();
831 // now propagate to the middle plane of the next time bin
832 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
833 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
834 if(!t.PropagateTo(x,radLength,rho)) break;
836 ymax = x*TMath::Tan(0.5*alpha);
839 if (!t.Rotate(alpha)) break;
840 if(!t.PropagateTo(x,radLength,rho)) break;
841 } else if (y <-ymax) {
843 if (!t.Rotate(-alpha)) break;
844 if(!t.PropagateTo(x,radLength,rho)) break;
850 expectedNumberOfClusters++;
851 wIndex = (Float_t) t.GetLabel();
854 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
856 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
857 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
860 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
861 else return expectedNumberOfClusters;
865 wYwindow = (Float_t) road;
866 t.GetPxPyPz(px,py,pz);
870 wC = (Float_t) t.GetC();
871 wSigmaC2 = (Float_t) t.GetSigmaC2();
872 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
873 wSigmaY2 = (Float_t) t.GetSigmaY2();
874 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
881 Double_t maxChi2=fgkMaxChi2;
883 wYclosest = 12345678;
884 wYcorrect = 12345678;
885 wZclosest = 12345678;
886 wZcorrect = 12345678;
887 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
889 // Find the closest correct cluster for debugging purposes
891 Float_t minDY = 1000000;
892 for (Int_t i=0; i<timeBin; i++) {
893 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
894 if((c->GetLabel(0) != trackIndex) &&
895 (c->GetLabel(1) != trackIndex) &&
896 (c->GetLabel(2) != trackIndex)) continue;
897 if(TMath::Abs(c->GetY() - y) > minDY) continue;
898 minDY = TMath::Abs(c->GetY() - y);
899 wYcorrect = c->GetY();
900 wZcorrect = c->GetZ();
902 Double_t h01 = GetTiltFactor(c);
903 wChi2 = t.GetPredictedChi2(c, h01);
907 // Now go for the real cluster search
911 //find cluster in history
914 AliTRDcluster * cl0 = timeBin[0];
918 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
919 if (plane>lastplane) continue;
920 Int_t timebin = cl0->GetLocalTimeBin();
921 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
924 Double_t h01 = GetTiltFactor(cl);
925 maxChi2=t.GetPredictedChi2(cl,h01);
927 if ((!cl) && road>fgkWideRoad) {
928 //if (t.GetNumberOfClusters()>4)
929 // cerr<<t.GetNumberOfClusters()
930 // <<"FindProlongation warning: Too broad road !\n";
937 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
938 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
939 if (c->GetY() > y+road) break;
940 if (c->IsUsed() > 0) continue;
941 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
943 Double_t h01 = GetTiltFactor(c);
944 Double_t chi2=t.GetPredictedChi2(c,h01);
946 if (chi2 > maxChi2) continue;
949 index=timeBin.GetIndex(i);
955 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
956 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
958 if (c->GetY() > y+road) break;
959 if (c->IsUsed() > 0) continue;
960 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
962 Double_t h01 = GetTiltFactor(c);
963 Double_t chi2=t.GetPredictedChi2(c, h01);
965 if (chi2 > maxChi2) continue;
968 index=timeBin.GetIndex(i);
972 wYclosest = cl->GetY();
973 wZclosest = cl->GetZ();
974 Double_t h01 = GetTiltFactor(cl);
976 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
977 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
978 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
979 Int_t det = cl->GetDetector();
980 Int_t plane = fGeom->GetPlane(det);
982 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
983 //if(!t.Update(cl,maxChi2,index,h01)) {
984 //if(!tryAgain--) return 0;
986 else tryAgain=fMaxGap;
989 //if (tryAgain==0) break;
994 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
996 printf(" %f", wIndex); //1
997 printf(" %f", wTB); //2
998 printf(" %f", wYrt); //3
999 printf(" %f", wYclosest); //4
1000 printf(" %f", wYcorrect); //5
1001 printf(" %f", wYwindow); //6
1002 printf(" %f", wZrt); //7
1003 printf(" %f", wZclosest); //8
1004 printf(" %f", wZcorrect); //9
1005 printf(" %f", wZwindow); //10
1006 printf(" %f", wPx); //11
1007 printf(" %f", wPy); //12
1008 printf(" %f", wPz); //13
1009 printf(" %f", wSigmaC2*1000000); //14
1010 printf(" %f", wSigmaTgl2*1000); //15
1011 printf(" %f", wSigmaY2); //16
1012 // printf(" %f", wSigmaZ2); //17
1013 printf(" %f", wChi2); //17
1014 printf(" %f", wC); //18
1021 return expectedNumberOfClusters;
1026 //___________________________________________________________________
1028 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
1030 // Starting from current radial position of track <t> this function
1031 // extrapolates the track up to outer timebin and in the sensitive
1032 // layers confirms prolongation if a close cluster is found.
1033 // Returns the number of clusters expected to be found in sensitive layers
1036 Float_t wIndex, wTB, wChi2;
1037 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
1038 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
1039 Float_t wPx, wPy, wPz, wC;
1040 Double_t px, py, pz;
1041 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
1043 Int_t trackIndex = t.GetLabel();
1044 Int_t tryAgain=fMaxGap;
1046 Double_t alpha=t.GetAlpha();
1047 TVector2::Phi_0_2pi(alpha);
1051 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
1052 Double_t radLength, rho, x, dx, y, ymax = 0, z;
1053 Bool_t lookForCluster;
1055 Int_t expectedNumberOfClusters = 0;
1058 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1060 Int_t nRefPlane = fgkFirstPlane;
1061 Bool_t isNewLayer = kFALSE;
1067 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
1072 // first propagate to the outer surface of the current time bin
1075 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1076 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
1080 if(!t.PropagateTo(x,radLength,rho)) break;
1081 // if (!AdjustSector(&t)) break;
1083 // MI -fix untill correct material desription will be implemented
1085 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material
1086 if (!AdjustSector(&t)) break;
1087 Int_t cross = kFALSE;
1089 if (TMath::Abs(angle - t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
1090 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1091 if (currentzone==-10) cross = kTRUE; // we are in the frame
1092 if (currentzone>-10){ // layer knows where we are
1093 if (zone==-10) zone = currentzone;
1094 if (zone!=currentzone) cross=kTRUE;
1098 if (t.GetNCross()==1) t.MakeBackupTrack();
1099 if (t.GetNCross()>2) break;
1105 if (!t.PropagateTo(x,radLength,rho)) break;
1110 // Barrel Tracks [SR, 04.04.2003]
1113 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1114 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1116 // if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1119 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1120 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1122 } else {isNewLayer = kFALSE;}
1127 // now propagate to the middle plane of the next time bin
1128 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1130 x = fTrSec[s]->GetLayer(nr+1)->GetX();
1131 if(!t.PropagateTo(x,radLength,rho)) break;
1132 if (!AdjustSector(&t)) break;
1134 if(!t.PropagateTo(x,radLength,rho)) break;
1139 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1140 // printf("label %d, pl %d, lookForCluster %d \n",
1141 // trackIndex, nr+1, lookForCluster);
1143 if(lookForCluster) {
1144 expectedNumberOfClusters++;
1146 wIndex = (Float_t) t.GetLabel();
1147 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1149 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1150 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1151 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1152 if((t.GetSigmaY2() + sy2) < 0) break;
1153 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1154 Double_t y=t.GetY(), z=t.GetZ();
1158 wYwindow = (Float_t) road;
1159 t.GetPxPyPz(px,py,pz);
1163 wC = (Float_t) t.GetC();
1164 wSigmaC2 = (Float_t) t.GetSigmaC2();
1165 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1166 wSigmaY2 = (Float_t) t.GetSigmaY2();
1167 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1170 if (road>fgkWideRoad) {
1171 if (t.GetNumberOfClusters()>4)
1172 cerr<<t.GetNumberOfClusters()
1173 <<"FindProlongation warning: Too broad road !\n";
1177 AliTRDcluster *cl=0;
1180 Double_t maxChi2=fgkMaxChi2;
1185 maxChi2 = 10 * fgkMaxChi2;
1188 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1189 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1190 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1193 wYclosest = 12345678;
1194 wYcorrect = 12345678;
1195 wZclosest = 12345678;
1196 wZcorrect = 12345678;
1197 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1199 // Find the closest correct cluster for debugging purposes
1202 for (Int_t i=0; i<timeBin; i++) {
1203 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1204 if((c->GetLabel(0) != trackIndex) &&
1205 (c->GetLabel(1) != trackIndex) &&
1206 (c->GetLabel(2) != trackIndex)) continue;
1207 if(TMath::Abs(c->GetY() - y) > minDY) continue;
1208 //minDY = TMath::Abs(c->GetY() - y);
1209 minDY = c->GetY() - y;
1210 wYcorrect = c->GetY();
1211 wZcorrect = c->GetZ();
1213 Double_t h01 = GetTiltFactor(c);
1214 wChi2 = t.GetPredictedChi2(c, h01);
1218 // Now go for the real cluster search
1222 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1223 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1224 if (c->GetY() > y+road) break;
1225 if (c->IsUsed() > 0) continue;
1226 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1228 Double_t h01 = GetTiltFactor(c);
1229 chi2=t.GetPredictedChi2(c,h01);
1231 if (chi2 > maxChi2) continue;
1234 index=timeBin.GetIndex(i);
1237 if((c->GetLabel(0) != trackIndex) &&
1238 (c->GetLabel(1) != trackIndex) &&
1239 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1244 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1245 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1247 if (c->GetY() > y+road) break;
1248 if (c->IsUsed() > 0) continue;
1249 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1251 Double_t h01 = GetTiltFactor(c);
1252 chi2=t.GetPredictedChi2(c,h01);
1254 if (chi2 > maxChi2) continue;
1257 index=timeBin.GetIndex(i);
1262 wYclosest = cl->GetY();
1263 wZclosest = cl->GetZ();
1265 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1266 Double_t h01 = GetTiltFactor(cl);
1267 Int_t det = cl->GetDetector();
1268 Int_t plane = fGeom->GetPlane(det);
1270 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1271 //if(!t.Update(cl,maxChi2,index,h01)) {
1272 if(!tryAgain--) return 0;
1274 else tryAgain=fMaxGap;
1277 if (tryAgain==0) break;
1280 //if (minDY < 1000000 && isNewLayer)
1281 //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() << "\t" <<
1282 // road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
1286 isNewLayer = kFALSE;
1289 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1291 printf(" %f", wIndex); //1
1292 printf(" %f", wTB); //2
1293 printf(" %f", wYrt); //3
1294 printf(" %f", wYclosest); //4
1295 printf(" %f", wYcorrect); //5
1296 printf(" %f", wYwindow); //6
1297 printf(" %f", wZrt); //7
1298 printf(" %f", wZclosest); //8
1299 printf(" %f", wZcorrect); //9
1300 printf(" %f", wZwindow); //10
1301 printf(" %f", wPx); //11
1302 printf(" %f", wPy); //12
1303 printf(" %f", wPz); //13
1304 printf(" %f", wSigmaC2*1000000); //14
1305 printf(" %f", wSigmaTgl2*1000); //15
1306 printf(" %f", wSigmaY2); //16
1307 // printf(" %f", wSigmaZ2); //17
1308 printf(" %f", wChi2); //17
1309 printf(" %f", wC); //18
1320 return expectedNumberOfClusters;
1325 //---------------------------------------------------------------------------
1326 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1328 // Starting from current position on track=t this function tries
1329 // to extrapolate the track up to timeBin=0 and to reuse already
1330 // assigned clusters. Returns the number of clusters
1331 // expected to be found in sensitive layers
1332 // get indices of assigned clusters for each layer
1333 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1336 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1337 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1338 Int_t index = t.GetClusterIndex(i);
1339 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1341 Int_t detector=cl->GetDetector();
1342 Int_t localTimeBin=cl->GetLocalTimeBin();
1343 Int_t sector=fGeom->GetSector(detector);
1344 Int_t plane=fGeom->GetPlane(detector);
1346 Int_t trackingSector = CookSectorIndex(sector);
1348 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1349 if(gtb < 0) continue;
1350 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1351 iCluster[layer] = index;
1355 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1357 Double_t alpha=t.GetAlpha();
1358 alpha = TVector2::Phi_0_2pi(alpha);
1360 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1361 Double_t radLength, rho, x, dx, y, ymax, z;
1363 Int_t expectedNumberOfClusters = 0;
1364 Bool_t lookForCluster;
1366 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1369 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1371 y = t.GetY(); z = t.GetZ();
1373 // first propagate to the inner surface of the current time bin
1374 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1375 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1376 if(!t.PropagateTo(x,radLength,rho)) break;
1378 ymax = x*TMath::Tan(0.5*alpha);
1381 if (!t.Rotate(alpha)) break;
1382 if(!t.PropagateTo(x,radLength,rho)) break;
1383 } else if (y <-ymax) {
1385 if (!t.Rotate(-alpha)) break;
1386 if(!t.PropagateTo(x,radLength,rho)) break;
1389 y = t.GetY(); z = t.GetZ();
1391 // now propagate to the middle plane of the next time bin
1392 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1393 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1394 if(!t.PropagateTo(x,radLength,rho)) break;
1396 ymax = x*TMath::Tan(0.5*alpha);
1399 if (!t.Rotate(alpha)) break;
1400 if(!t.PropagateTo(x,radLength,rho)) break;
1401 } else if (y <-ymax) {
1403 if (!t.Rotate(-alpha)) break;
1404 if(!t.PropagateTo(x,radLength,rho)) break;
1407 if(lookForCluster) expectedNumberOfClusters++;
1409 // use assigned cluster
1410 if (!iCluster[nr-1]) continue;
1411 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1412 Double_t h01 = GetTiltFactor(cl);
1413 Double_t chi2=t.GetPredictedChi2(cl, h01);
1414 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1415 t.Update(cl,chi2,iCluster[nr-1],h01);
1418 return expectedNumberOfClusters;
1421 //___________________________________________________________________
1423 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1425 // Starting from current radial position of track <t> this function
1426 // extrapolates the track up to radial position <xToGo>.
1427 // Returns 1 if track reaches the plane, and 0 otherwise
1429 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1431 Double_t alpha=t.GetAlpha();
1433 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1434 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1436 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1438 Bool_t lookForCluster;
1439 Double_t radLength, rho, x, dx, y, ymax, z;
1443 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1445 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1447 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1449 y = t.GetY(); z = t.GetZ();
1451 // first propagate to the outer surface of the current time bin
1452 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1453 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1454 if(!t.PropagateTo(x,radLength,rho)) return 0;
1456 ymax = x*TMath::Tan(0.5*alpha);
1459 if (!t.Rotate(alpha)) return 0;
1460 } else if (y <-ymax) {
1462 if (!t.Rotate(-alpha)) return 0;
1464 if(!t.PropagateTo(x,radLength,rho)) return 0;
1466 y = t.GetY(); z = t.GetZ();
1468 // now propagate to the middle plane of the next time bin
1469 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1470 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1471 if(!t.PropagateTo(x,radLength,rho)) return 0;
1473 ymax = x*TMath::Tan(0.5*alpha);
1476 if (!t.Rotate(alpha)) return 0;
1477 } else if (y <-ymax) {
1479 if (!t.Rotate(-alpha)) return 0;
1481 if(!t.PropagateTo(x,radLength,rho)) return 0;
1486 //___________________________________________________________________
1488 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1490 // Starting from current radial position of track <t> this function
1491 // extrapolates the track up to radial position of the outermost
1492 // padrow of the TPC.
1493 // Returns 1 if track reaches the TPC, and 0 otherwise
1495 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1497 Double_t alpha=t.GetAlpha();
1498 alpha = TVector2::Phi_0_2pi(alpha);
1500 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1502 Bool_t lookForCluster;
1503 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1507 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1508 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1510 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1515 // first propagate to the outer surface of the current time bin
1516 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1517 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1519 if(!t.PropagateTo(x,radLength,rho)) return 0;
1521 if(!t.PropagateTo(x,radLength,rho)) return 0;
1526 // now propagate to the middle plane of the next time bin
1527 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1528 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1530 if(!t.PropagateTo(x,radLength,rho)) return 0;
1532 if(!t.PropagateTo(x,radLength,rho)) return 0;
1537 //_____________________________________________________________________________
1538 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1540 // Fills clusters into TRD tracking_sectors
1541 // Note that the numbering scheme for the TRD tracking_sectors
1542 // differs from that of TRD sectors
1544 if (ReadClusters(fClusters,cTree)) {
1545 Error("LoadClusters","Problem with reading the clusters !");
1548 Int_t ncl=fClusters->GetEntriesFast();
1550 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1553 for (Int_t ichamber=0;ichamber<5;ichamber++)
1554 for (Int_t isector=0;isector<18;isector++){
1555 fHoles[ichamber][isector]=kTRUE;
1560 // printf("\r %d left ",ncl);
1561 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1562 Int_t detector=c->GetDetector();
1563 Int_t localTimeBin=c->GetLocalTimeBin();
1564 Int_t sector=fGeom->GetSector(detector);
1565 Int_t plane=fGeom->GetPlane(detector);
1567 Int_t trackingSector = CookSectorIndex(sector);
1568 if (c->GetLabel(0)>0){
1569 Int_t chamber = fGeom->GetChamber(detector);
1570 fHoles[chamber][trackingSector]=kFALSE;
1573 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1574 if(gtb < 0) continue;
1575 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1578 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1584 for (Int_t isector=0;isector<18;isector++){
1585 for (Int_t ichamber=0;ichamber<5;ichamber++)
1586 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1587 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1588 fGeom->IsHole(0,ichamber,17-isector));
1594 //_____________________________________________________________________________
1595 void AliTRDtracker::UnloadClusters()
1598 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1603 nentr = fClusters->GetEntriesFast();
1604 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1607 nentr = fSeeds->GetEntriesFast();
1608 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1610 nentr = fTracks->GetEntriesFast();
1611 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1613 Int_t nsec = AliTRDgeometry::kNsect;
1615 for (i = 0; i < nsec; i++) {
1616 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1617 fTrSec[i]->GetLayer(pl)->Clear();
1623 //__________________________________________________________________________
1624 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1626 // Creates track seeds using clusters in timeBins=i1,i2
1629 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1633 Double_t x[5], c[15];
1634 Int_t maxSec=AliTRDgeometry::kNsect;
1636 Double_t alpha=AliTRDgeometry::GetAlpha();
1637 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1638 Double_t cs=cos(alpha), sn=sin(alpha);
1639 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1642 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1643 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1645 Double_t x1 =fTrSec[0]->GetX(i1);
1646 Double_t xx2=fTrSec[0]->GetX(i2);
1648 for (Int_t ns=0; ns<maxSec; ns++) {
1650 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1651 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1652 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1653 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1654 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1656 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1658 for (Int_t is=0; is < r1; is++) {
1659 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1661 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1663 const AliTRDcluster *cl;
1664 Double_t x2, y2, z2;
1665 Double_t x3=0., y3=0.;
1668 if(turn != 2) continue;
1669 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1671 y2=cl->GetY(); z2=cl->GetZ();
1676 else if (js<nl2+nl) {
1677 if(turn != 1) continue;
1678 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1680 y2=cl->GetY(); z2=cl->GetZ();
1685 else if (js<nl2+nl+nm) {
1686 if(turn != 1) continue;
1687 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1689 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1691 else if (js<nl2+nl+nm+nu) {
1692 if(turn != 1) continue;
1693 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1694 cl=r2[js-nl2-nl-nm];
1695 y2=cl->GetY(); z2=cl->GetZ();
1701 if(turn != 2) continue;
1702 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1703 cl=r2[js-nl2-nl-nm-nu];
1704 y2=cl->GetY(); z2=cl->GetZ();
1710 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1712 Double_t zz=z1 - z1/x1*(x1-x2);
1714 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1716 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1717 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1721 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1723 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
1725 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1727 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1729 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1731 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1733 Double_t a=asin(x[2]);
1734 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1736 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1738 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1739 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1740 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
1743 Double_t h01 = GetTiltFactor(r1[is]);
1744 Double_t xuFactor = 100.;
1750 sy1=sy1+sz1*h01*h01;
1751 Double_t syz=sz1*(-h01);
1752 // end of tilt changes
1754 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1755 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1756 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1757 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1758 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1759 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1760 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1761 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1762 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1763 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1767 // c[1]=0.; c[2]=sz1;
1768 c[1]=syz; c[2]=sz1*xuFactor;
1769 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1770 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1771 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1772 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1773 c[13]=f30*sy1*f40+f32*sy2*f42;
1774 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1776 UInt_t index=r1.GetIndex(is);
1778 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1780 Int_t rc=FollowProlongation(*track, i2);
1783 (track->GetNumberOfClusters() <
1784 (outer-inner)*fgkMinClustersInSeed)) delete track;
1786 fSeeds->AddLast(track); fNseeds++;
1787 // cerr<<"\r found seed "<<fNseeds;
1794 //_____________________________________________________________________________
1795 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1798 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1799 // from the file. The names of the cluster tree and branches
1800 // should match the ones used in AliTRDclusterizer::WriteClusters()
1802 TObjArray *clusterArray = new TObjArray(400);
1804 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1806 Error("ReadClusters","Can't get the branch !");
1809 branch->SetAddress(&clusterArray);
1811 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1812 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1814 // Loop through all entries in the tree
1816 AliTRDcluster *c = 0;
1819 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1822 nbytes += ClusterTree->GetEvent(iEntry);
1824 // Get the number of points in the detector
1825 Int_t nCluster = clusterArray->GetEntriesFast();
1826 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1828 // Loop through all TRD digits
1829 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1830 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1831 AliTRDcluster *co = new AliTRDcluster(*c);
1832 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1833 Int_t ltb = co->GetLocalTimeBin();
1834 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1835 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1837 delete clusterArray->RemoveAt(iCluster);
1841 delete clusterArray;
1846 //__________________________________________________________________
1847 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1850 // This cooks a label. Mmmmh, smells good...
1853 Int_t label=123456789, index, i, j;
1854 Int_t ncl=pt->GetNumberOfClusters();
1855 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1859 // Int_t s[kRange][2];
1860 Int_t **s = new Int_t* [kRange];
1861 for (i=0; i<kRange; i++) {
1862 s[i] = new Int_t[2];
1864 for (i=0; i<kRange; i++) {
1870 for (i=0; i<ncl; i++) {
1871 index=pt->GetClusterIndex(i);
1872 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1878 for (i=0; i<ncl; i++) {
1879 index=pt->GetClusterIndex(i);
1880 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1881 for (Int_t k=0; k<3; k++) {
1882 label=c->GetLabel(k);
1883 labelAdded=kFALSE; j=0;
1885 while ( (!labelAdded) && ( j < kRange ) ) {
1886 if (s[j][0]==label || s[j][1]==0) {
1900 for (i=0; i<kRange; i++) {
1902 max=s[i][1]; label=s[i][0];
1906 for (i=0; i<kRange; i++) {
1912 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1914 pt->SetLabel(label);
1919 //__________________________________________________________________
1920 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1923 // Use clusters, but don't abuse them!
1926 Int_t ncl=t->GetNumberOfClusters();
1927 for (Int_t i=from; i<ncl; i++) {
1928 Int_t index = t->GetClusterIndex(i);
1929 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1935 //_____________________________________________________________________
1936 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1938 // Parametrised "expected" error of the cluster reconstruction in Y
1940 Double_t s = 0.08 * 0.08;
1944 //_____________________________________________________________________
1945 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1947 // Parametrised "expected" error of the cluster reconstruction in Z
1949 Double_t s = 9 * 9 /12.;
1953 //_____________________________________________________________________
1954 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
1957 // Returns radial position which corresponds to time bin <localTB>
1958 // in tracking sector <sector> and plane <plane>
1961 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
1962 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1963 return fTrSec[sector]->GetLayer(pl)->GetX();
1968 //_______________________________________________________
1969 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1970 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1973 // AliTRDpropagationLayer constructor
1976 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1977 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1980 for(Int_t i=0; i < (Int_t) kZones; i++) {
1981 fZc[i]=0; fZmax[i] = 0;
1986 if(fTimeBinIndex >= 0) {
1987 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1988 fIndex = new UInt_t[kMaxClusterPerTimeBin];
1991 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
2002 //_______________________________________________________
2003 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
2004 Double_t Zmax, Double_t Ymax, Double_t rho,
2005 Double_t radLength, Double_t Yc, Double_t Zc)
2008 // Sets hole in the layer
2016 fHoleX0 = radLength;
2020 //_______________________________________________________
2021 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
2024 // AliTRDtrackingSector Constructor
2033 // get holes description from geometry
2034 Bool_t holes[AliTRDgeometry::kNcham];
2035 //printf("sector\t%d\t",gs);
2036 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
2037 holes[icham] = fGeom->IsHole(0,icham,gs);
2038 //printf("%d",holes[icham]);
2042 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
2045 AliTRDpropagationLayer* ppl;
2047 Double_t x, xin, xout, dx, rho, radLength;
2050 // set time bins in the gas of the TPC
2052 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
2053 rho = 0.9e-3; radLength = 28.94;
2055 for(Int_t i=0; i<steps; i++) {
2056 x = xin + i*dx + dx/2;
2057 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2061 // set time bins in the outer field cage vessel
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 = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2068 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2071 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2072 steps = 5; dx = (xout - xin)/steps;
2073 for(Int_t i=0; i<steps; i++) {
2074 x = xin + i*dx + dx/2;
2075 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2079 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2080 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2083 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2084 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2088 // set time bins in CO2
2090 xin = xout; xout = 275.0;
2091 steps = 50; dx = (xout - xin)/steps;
2092 rho = 1.977e-3; radLength = 36.2;
2094 for(Int_t i=0; i<steps; i++) {
2095 x = xin + i*dx + dx/2;
2096 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2100 // set time bins in the outer containment vessel
2102 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2103 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2106 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2107 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2110 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2111 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2114 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2115 steps = 10; dx = (xout - xin)/steps;
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);
2122 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2123 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2126 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2127 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2130 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2131 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2134 Double_t xtrd = (Double_t) fGeom->Rmin();
2136 // add layers between TPC and TRD (Air temporarily)
2137 xin = xout; xout = xtrd;
2138 steps = 50; dx = (xout - xin)/steps;
2139 rho = 1.2e-3; radLength = 36.66;
2141 for(Int_t i=0; i<steps; i++) {
2142 x = xin + i*dx + dx/2;
2143 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2148 // Double_t alpha=AliTRDgeometry::GetAlpha();
2150 // add layers for each of the planes
2152 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2153 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2154 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2155 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2156 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2157 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2158 Double_t dxPlane = dxTEC + dxSpace;
2161 const Int_t kNchambers = AliTRDgeometry::Ncham();
2164 Double_t ymaxsensitive=0;
2165 Double_t *zc = new Double_t[kNchambers];
2166 Double_t *zmax = new Double_t[kNchambers];
2167 Double_t *zmaxsensitive = new Double_t[kNchambers];
2168 // Double_t holeZmax = 1000.; // the whole sector is missing
2170 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2173 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2174 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2175 for(Int_t i=0; i<steps; i++) {
2176 x = xin + i*dx + dx/2;
2177 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2181 ymax = fGeom->GetChamberWidth(plane)/2.;
2182 ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2184 for(Int_t ch = 0; ch < kNchambers; ch++) {
2185 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2186 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2187 Float_t row0 = fPar->GetRow0(plane,ch,0);
2188 Int_t nPads = fPar->GetRowMax(plane,ch,0);
2189 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2190 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2191 zc[ch] = (pad * nPads)/2 + row0;
2192 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2196 dx = fPar->GetTimeBinSize();
2197 rho = 0.00295 * 0.85; radLength = 11.0;
2199 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2200 Double_t xbottom = x0 - dxDrift;
2201 Double_t xtop = x0 + dxAmp;
2203 // Amplification region
2204 steps = (Int_t) (dxAmp/dx);
2206 for(tb = 0; tb < steps; tb++) {
2207 x = x0 + tb * dx + dx/2;
2208 tbIndex = CookTimeBinIndex(plane, -tb-1);
2209 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2210 ppl->SetYmax(ymax,ymaxsensitive);
2211 ppl->SetZ(zc, zmax, zmaxsensitive);
2212 ppl->SetHoles(holes);
2215 tbIndex = CookTimeBinIndex(plane, -steps);
2216 x = (x + dx/2 + xtop)/2;
2218 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2219 ppl->SetYmax(ymax,ymaxsensitive);
2220 ppl->SetZ(zc, zmax,zmaxsensitive);
2221 ppl->SetHoles(holes);
2225 dx = fPar->GetTimeBinSize();
2226 steps = (Int_t) (dxDrift/dx);
2228 for(tb = 0; tb < steps; tb++) {
2229 x = x0 - tb * dx - dx/2;
2230 tbIndex = CookTimeBinIndex(plane, tb);
2232 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2233 ppl->SetYmax(ymax,ymaxsensitive);
2234 ppl->SetZ(zc, zmax, zmaxsensitive);
2235 ppl->SetHoles(holes);
2238 tbIndex = CookTimeBinIndex(plane, steps);
2239 x = (x - dx/2 + xbottom)/2;
2241 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2242 ppl->SetYmax(ymax,ymaxsensitive);
2243 ppl->SetZ(zc, zmax, zmaxsensitive);
2244 ppl->SetHoles(holes);
2248 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2249 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2250 ppl->SetYmax(ymax,ymaxsensitive);
2251 ppl->SetZ(zc, zmax,zmax);
2252 ppl->SetHoles(holes);
2256 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2257 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2258 for(Int_t i=0; i<steps; i++) {
2259 x = xin + i*dx + dx/2;
2260 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2261 ppl->SetYmax(ymax,ymaxsensitive);
2262 ppl->SetZ(zc, zmax,zmax);
2263 ppl->SetHoles(holes);
2267 // Space between the chambers, air
2268 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2269 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2270 for(Int_t i=0; i<steps; i++) {
2271 x = xin + i*dx + dx/2;
2272 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2277 // Space between the TRD and RICH
2278 Double_t xRICH = 500.;
2279 xin = xout; xout = xRICH;
2280 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2281 for(Int_t i=0; i<steps; i++) {
2282 x = xin + i*dx + dx/2;
2283 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2293 //______________________________________________________
2295 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2298 // depending on the digitization parameters calculates "global"
2299 // time bin index for timebin <localTB> in plane <plane>
2302 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2303 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2304 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2306 Int_t tbAmp = fPar->GetTimeBefore();
2307 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2308 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2309 Int_t tbDrift = fPar->GetTimeMax();
2310 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2312 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2314 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2317 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2318 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2325 //______________________________________________________
2327 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2330 // For all sensitive time bins sets corresponding layer index
2331 // in the array fTimeBins
2336 for(Int_t i = 0; i < fN; i++) {
2337 index = fLayers[i]->GetTimeBinIndex();
2339 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2341 if(index < 0) continue;
2342 if(index >= (Int_t) kMaxTimeBinIndex) {
2343 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2344 printf(" index %d exceeds allowed maximum of %d!\n",
2345 index, kMaxTimeBinIndex-1);
2348 fTimeBinIndex[index] = i;
2351 Double_t x1, dx1, x2, dx2, gap;
2353 for(Int_t i = 0; i < fN-1; i++) {
2354 x1 = fLayers[i]->GetX();
2355 dx1 = fLayers[i]->GetdX();
2356 x2 = fLayers[i+1]->GetX();
2357 dx2 = fLayers[i+1]->GetdX();
2358 gap = (x2 - dx2/2) - (x1 + dx1/2);
2360 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2361 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2364 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2365 printf(" (%f - %f) - (%f + %f) = %f\n",
2366 x2, dx2/2, x1, dx1, gap);
2372 //______________________________________________________
2375 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2378 // Returns the number of time bin which in radial position is closest to <x>
2381 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2382 if(x <= fLayers[0]->GetX()) return 0;
2384 Int_t b=0, e=fN-1, m=(b+e)/2;
2385 for (; b<e; m=(b+e)/2) {
2386 if (x > fLayers[m]->GetX()) b=m+1;
2389 if(TMath::Abs(x - fLayers[m]->GetX()) >
2390 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2395 //______________________________________________________
2397 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2400 // Returns number of the innermost SENSITIVE propagation layer
2403 return GetLayerNumber(0);
2406 //______________________________________________________
2408 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2411 // Returns number of the outermost SENSITIVE time bin
2414 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2417 //______________________________________________________
2419 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2422 // Returns number of SENSITIVE time bins
2426 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2427 layer = GetLayerNumber(tb);
2433 //______________________________________________________
2435 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2438 // Insert layer <pl> in fLayers array.
2439 // Layers are sorted according to X coordinate.
2441 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2442 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2445 if (fN==0) {fLayers[fN++] = pl; return;}
2446 Int_t i=Find(pl->GetX());
2448 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2449 fLayers[i]=pl; fN++;
2453 //______________________________________________________
2455 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2458 // Returns index of the propagation layer nearest to X
2461 if (x <= fLayers[0]->GetX()) return 0;
2462 if (x > fLayers[fN-1]->GetX()) return fN;
2463 Int_t b=0, e=fN-1, m=(b+e)/2;
2464 for (; b<e; m=(b+e)/2) {
2465 if (x > fLayers[m]->GetX()) b=m+1;
2471 //______________________________________________________
2472 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2475 // set centers and the width of sectors
2476 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2477 fZc[icham] = center[icham];
2478 fZmax[icham] = w[icham];
2479 fZmaxSensitive[icham] = wsensitive[icham];
2480 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2483 //______________________________________________________
2485 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2488 // set centers and the width of sectors
2490 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2491 fIsHole[icham] = holes[icham];
2492 if (holes[icham]) fHole = kTRUE;
2498 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2499 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
2500 Bool_t &lookForCluster) const
2503 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2504 // and sensitivity <lookForCluster> in point <y,z>
2510 lookForCluster = kFALSE;
2512 // check dead regions in sensitive volume
2513 if(fTimeBinIndex >= 0) {
2516 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2517 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2519 lookForCluster = !(fIsHole[zone]);
2520 if(TMath::Abs(y) > fYmaxSensitive){
2521 lookForCluster = kFALSE;
2523 if (fIsHole[zone]) {
2535 if (fHole==kFALSE) return;
2537 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2538 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2549 Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2553 if (fTimeBinIndex < 0) return -20; //unknown
2554 Int_t zone=-10; // dead zone
2555 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2556 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2563 //______________________________________________________
2565 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2568 // Insert cluster in cluster array.
2569 // Clusters are sorted according to Y coordinate.
2571 if(fTimeBinIndex < 0) {
2572 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2576 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2577 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2580 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2581 Int_t i=Find(c->GetY());
2582 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2583 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2584 fIndex[i]=index; fClusters[i]=c; fN++;
2587 //______________________________________________________
2589 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2591 // Returns index of the cluster nearest in Y
2593 if (y <= fClusters[0]->GetY()) return 0;
2594 if (y > fClusters[fN-1]->GetY()) return fN;
2595 Int_t b=0, e=fN-1, m=(b+e)/2;
2596 for (; b<e; m=(b+e)/2) {
2597 if (y > fClusters[m]->GetY()) b=m+1;
2603 //---------------------------------------------------------
2605 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2607 // Returns correction factor for tilted pads geometry
2610 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2611 Int_t det = c->GetDetector();
2612 Int_t plane = fGeom->GetPlane(det);
2614 //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2615 if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
2617 if(fNoTilt) h01 = 0;
2623 void AliTRDtracker::CookdEdxTimBin(AliTRDtrack& TRDtrack)
2625 // *** ADDED TO GET MORE INFORMATION FOR TRD PID ---- PS
2626 // This is setting fdEdxPlane and fTimBinPlane
2627 // Sums up the charge in each plane for track TRDtrack and also get the
2628 // Time bin for Max. Cluster
2629 // Prashant Shukla (shukla@physi.uni-heidelberg.de)
2631 // const Int_t kNPlane = AliTRDgeometry::Nplan();
2632 // const Int_t kNPlane = 6;
2633 Double_t clscharge[kNPlane], maxclscharge[kNPlane];
2634 Int_t nCluster[kNPlane], timebin[kNPlane];
2636 //Initialization of cluster charge per plane.
2637 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2638 clscharge[iPlane] = 0.0;
2639 nCluster[iPlane] = 0;
2640 timebin[iPlane] = -1;
2641 maxclscharge[iPlane] = 0.0;
2644 // Loop through all clusters associated to track TRDtrack
2645 Int_t nClus = TRDtrack.GetNumberOfClusters(); // from Kalmantrack
2646 for (Int_t iClus = 0; iClus < nClus; iClus++) {
2647 Double_t charge = TRDtrack.GetClusterdQdl(iClus);
2648 Int_t index = TRDtrack.GetClusterIndex(iClus);
2649 AliTRDcluster *TRDcluster = (AliTRDcluster *) GetCluster(index);
2650 if (!TRDcluster) continue;
2651 Int_t tb = TRDcluster->GetLocalTimeBin();
2653 Int_t detector = TRDcluster->GetDetector();
2654 Int_t iPlane = fGeom->GetPlane(detector);
2655 clscharge[iPlane] = clscharge[iPlane]+charge;
2656 if(charge > maxclscharge[iPlane]) {
2657 maxclscharge[iPlane] = charge;
2658 timebin[iPlane] = tb;
2661 } // end of loop over cluster
2663 // Setting the fdEdxPlane and fTimBinPlane variabales
2664 Double_t Total_ch = 0;
2665 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2666 if (nCluster[iPlane]) clscharge[iPlane] /= nCluster[iPlane];
2667 TRDtrack.SetPIDsignals(clscharge[iPlane], iPlane);
2668 TRDtrack.SetPIDTimBin(timebin[iPlane], iPlane);
2669 Total_ch= Total_ch+clscharge[iPlane];
2672 // Int_t nc=TRDtrack.GetNumberOfClusters();
2674 // for (i=0; i<nc; i++) dedx += TRDtrack.GetClusterdQdl(i);
2676 // for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2677 // TRDtrack.SetPIDsignals(dedx, iPlane);
2678 // TRDtrack.SetPIDTimBin(timbin[iPlane], iPlane);
2681 } // end of function