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 "AliTRDpadPlane.h"
33 #include "AliTRDgeometryFull.h"
34 #include "AliTRDcluster.h"
35 #include "AliTRDtrack.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 AliTRDgeometryFull();
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();
144 // fGeom->SetT0(fTzero);
147 fClusters = new TObjArray(2000);
149 fSeeds = new TObjArray(2000);
151 fTracks = new TObjArray(1000);
153 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
154 Int_t trS = CookSectorIndex(geomS);
155 fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS, fPar);
156 for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
157 fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
161 AliTRDpadPlane *padPlane = fPar->GetPadPlane(0,0);
162 Float_t tiltAngle = TMath::Abs(padPlane->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->GetDriftVelocity()
179 / fPar->GetSamplingFrequency();
180 Int_t tbAmp = fPar->GetTimeBefore();
181 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
182 if(kTRUE) maxAmp = 0; // intentional until we change the parameter class
183 Int_t tbDrift = fPar->GetTimeMax();
184 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
186 tbDrift = TMath::Min(tbDrift,maxDrift);
187 tbAmp = TMath::Min(tbAmp,maxAmp);
189 fTimeBinsPerPlane = tbAmp + tbDrift;
190 fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fgkSkipDepth);
197 //___________________________________________________________________
198 AliTRDtracker::~AliTRDtracker()
201 // Destructor of AliTRDtracker
219 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
220 delete fTrSec[geomS];
225 //_____________________________________________________________________
227 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track) {
229 // Rotates the track when necessary
232 Double_t alpha = AliTRDgeometry::GetAlpha();
233 Double_t y = track->GetY();
234 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
236 //Int_t ns = AliTRDgeometry::kNsect;
237 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
241 if (!track->Rotate(alpha)) return kFALSE;
242 } else if (y <-ymax) {
244 if (!track->Rotate(-alpha)) return kFALSE;
250 //_____________________________________________________________________
251 inline Double_t f1trd(Double_t x1,Double_t y1,
252 Double_t x2,Double_t y2,
253 Double_t x3,Double_t y3)
256 // Initial approximation of the track curvature
258 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
259 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
260 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
261 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
262 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
264 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
266 return -xr*yr/sqrt(xr*xr+yr*yr);
269 //_____________________________________________________________________
270 inline Double_t f2trd(Double_t x1,Double_t y1,
271 Double_t x2,Double_t y2,
272 Double_t x3,Double_t y3)
275 // Initial approximation of the track curvature times X coordinate
276 // of the center of curvature
279 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
280 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
281 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
282 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
283 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
285 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
287 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
290 //_____________________________________________________________________
291 inline Double_t f3trd(Double_t x1,Double_t y1,
292 Double_t x2,Double_t y2,
293 Double_t z1,Double_t z2)
296 // Initial approximation of the tangent of the track dip angle
299 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
303 AliTRDcluster * AliTRDtracker::GetCluster(AliTRDtrack * track, Int_t plane, Int_t timebin){
305 //try to find cluster in the backup list
307 AliTRDcluster * cl =0;
308 UInt_t *indexes = track->GetBackupIndexes();
309 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
310 if (indexes[i]==0) break;
311 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
313 if (cli->GetLocalTimeBin()!=timebin) continue;
314 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
324 Int_t AliTRDtracker::GetLastPlane(AliTRDtrack * track){
326 //return last updated plane
328 UInt_t *indexes = track->GetBackupIndexes();
329 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
330 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
332 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
333 if (iplane>lastplane) {
339 //___________________________________________________________________
340 Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
343 // Finds tracks within the TRD. The ESD event is expected to contain seeds
344 // at the outer part of the TRD. The seeds
345 // are found within the TRD if fAddTRDseeds is TRUE.
346 // The tracks are propagated to the innermost time bin
347 // of the TRD and the ESD event is updated
350 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
351 Float_t foundMin = fgkMinClustersInTrack * timeBins;
354 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
356 Int_t n = event->GetNumberOfTracks();
357 for (Int_t i=0; i<n; i++) {
358 AliESDtrack* seed=event->GetTrack(i);
359 ULong_t status=seed->GetStatus();
360 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
361 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
364 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
365 //seed2->ResetCovariance();
366 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
368 FollowProlongation(t, innerTB);
369 if (t.GetNumberOfClusters() >= foundMin) {
371 CookLabel(pt, 1-fgkLabelFraction);
375 // cout<<found<<'\r';
377 if(PropagateToTPC(t)) {
378 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
384 cout<<"Number of loaded seeds: "<<nseed<<endl;
385 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
387 // after tracks from loaded seeds are found and the corresponding
388 // clusters are used, look for additional seeds from TRD
391 // Find tracks for the seeds in the TRD
392 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
394 Int_t nSteps = (Int_t) (fgkSeedDepth / fgkSeedStep);
395 Int_t gap = (Int_t) (timeBins * fgkSeedGap);
396 Int_t step = (Int_t) (timeBins * fgkSeedStep);
398 // make a first turn with tight cut on initial curvature
399 for(Int_t turn = 1; turn <= 2; turn++) {
401 nSteps = (Int_t) (fgkSeedDepth / (3*fgkSeedStep));
402 step = (Int_t) (timeBins * (3*fgkSeedStep));
404 for(Int_t i=0; i<nSteps; i++) {
405 Int_t outer=timeBins-1-i*step;
406 Int_t inner=outer-gap;
408 nseed=fSeeds->GetEntriesFast();
410 MakeSeeds(inner, outer, turn);
412 nseed=fSeeds->GetEntriesFast();
413 // printf("\n turn %d, step %d: number of seeds for TRD inward %d\n",
416 for (Int_t i=0; i<nseed; i++) {
417 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
418 FollowProlongation(t,innerTB);
419 if (t.GetNumberOfClusters() >= foundMin) {
421 CookLabel(pt, 1-fgkLabelFraction);
424 // cout<<found<<'\r';
425 if(PropagateToTPC(t)) {
427 track.UpdateTrackParams(pt,AliESDtrack::kTRDin);
428 event->AddTrack(&track);
429 // track.SetTRDtrack(new AliTRDtrack(*pt));
432 delete fSeeds->RemoveAt(i);
439 cout<<"Total number of found tracks: "<<found<<endl;
446 //_____________________________________________________________________________
447 Int_t AliTRDtracker::PropagateBack(AliESD* event) {
449 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
450 // backpropagated by the TPC tracker. Each seed is first propagated
451 // to the TRD, and then its prolongation is searched in the TRD.
452 // If sufficiently long continuation of the track is found in the TRD
453 // the track is updated, otherwise it's stored as originaly defined
454 // by the TPC tracker.
458 Float_t foundMin = 20;
459 Int_t n = event->GetNumberOfTracks();
462 Float_t *quality =new Float_t[n];
463 Int_t *index =new Int_t[n];
464 for (Int_t i=0; i<n; i++) {
465 AliESDtrack* seed=event->GetTrack(i);
466 Double_t covariance[15];
467 seed->GetExternalCovariance(covariance);
468 quality[i] = covariance[0]+covariance[2];
470 TMath::Sort(n,quality,index,kFALSE);
472 for (Int_t i=0; i<n; i++) {
473 // AliESDtrack* seed=event->GetTrack(i);
474 AliESDtrack* seed=event->GetTrack(index[i]);
476 ULong_t status=seed->GetStatus();
477 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
478 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
480 Int_t lbl = seed->GetLabel();
481 AliTRDtrack *track = new AliTRDtrack(*seed);
482 track->SetSeedLabel(lbl);
483 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup); //make backup
485 Float_t p4 = track->GetC();
487 Int_t expectedClr = FollowBackProlongation(*track);
489 // only debug purpose
490 if (track->GetNumberOfClusters()<expectedClr/3){
491 AliTRDtrack *track1 = new AliTRDtrack(*seed);
492 track1->SetSeedLabel(lbl);
493 FollowBackProlongation(*track1);
494 AliTRDtrack *track2= new AliTRDtrack(*seed);
495 track->SetSeedLabel(lbl);
496 FollowBackProlongation(*track2);
501 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)<0.2 || TMath::Abs(track->GetPt())>0.8 ) {
503 //make backup for back propagation
505 Int_t foundClr = track->GetNumberOfClusters();
506 if (foundClr >= foundMin) {
508 CookLabel(track, 1-fgkLabelFraction);
509 if(track->GetChi2()/track->GetNumberOfClusters()<4) { // sign only gold tracks
510 if (seed->GetKinkIndex(0)==0&&TMath::Abs(track->GetPt())<1.5 ) UseClusters(track);
512 Bool_t isGold = kFALSE;
514 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
515 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
518 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
519 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
522 if (!isGold && track->GetBackupTrack()){
523 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
524 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
525 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
532 //Propagation to the TOF (I.Belikov)
533 CookdEdxTimBin(*track);
534 if (track->GetStop()==kFALSE){
537 Double_t c2=track->GetC()*xtof - track->GetEta();
538 if (TMath::Abs(c2)>=0.99) {
542 Double_t xTOF0 = 365. ;
543 PropagateToOuterPlane(*track,xTOF0);
545 //energy losses taken to the account - check one more time
546 c2=track->GetC()*xtof - track->GetEta();
547 if (TMath::Abs(c2)>=0.99) {
553 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
554 Double_t y=track->GetYat(xtof);
556 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
560 } else if (y <-ymax) {
561 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
567 if (track->PropagateTo(xtof)) {
568 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
569 for (Int_t i=0;i<kNPlane;i++) {
570 seed->SetTRDsignals(track->GetPIDsignals(i),i);
571 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
573 seed->SetTRDtrack(new AliTRDtrack(*track));
574 if (track->GetNumberOfClusters()>foundMin) found++;
577 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
578 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
579 //seed->SetStatus(AliESDtrack::kTRDStop);
580 for (Int_t i=0;i<kNPlane;i++) {
581 seed->SetTRDsignals(track->GetPIDsignals(i),i);
582 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
584 seed->SetTRDtrack(new AliTRDtrack(*track));
588 seed->SetTRDQuality(track->StatusForTOF());
591 //End of propagation to the TOF
592 //if (foundClr>foundMin)
593 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
598 cerr<<"Number of seeds: "<<fNseeds<<endl;
599 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
601 fSeeds->Clear(); fNseeds=0;
609 //_____________________________________________________________________________
610 Int_t AliTRDtracker::RefitInward(AliESD* event)
613 // Refits tracks within the TRD. The ESD event is expected to contain seeds
614 // at the outer part of the TRD.
615 // The tracks are propagated to the innermost time bin
616 // of the TRD and the ESD event is updated
617 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
620 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
621 Float_t foundMin = fgkMinClustersInTrack * timeBins;
624 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
627 Int_t n = event->GetNumberOfTracks();
628 for (Int_t i=0; i<n; i++) {
629 AliESDtrack* seed=event->GetTrack(i);
630 new(&seed2) AliTRDtrack(*seed);
631 if (seed2.GetX()<270){
632 seed->UpdateTrackParams(&seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
636 ULong_t status=seed->GetStatus();
637 if ( (status & AliESDtrack::kTRDout ) == 0 ) {
640 if ( (status & AliESDtrack::kTRDin) != 0 ) {
644 if (1/seed2.Get1Pt()>5.&& seed2.GetX()>260.) {
645 Double_t oldx = seed2.GetX();
646 seed2.PropagateTo(500.);
647 seed2.ResetCovariance(1.);
648 seed2.PropagateTo(oldx);
651 seed2.ResetCovariance(5.);
654 AliTRDtrack *pt = new AliTRDtrack(seed2,seed2.GetAlpha());
655 UInt_t * indexes2 = seed2.GetIndexes();
656 // for (Int_t i=0;i<kNPlane;i++) {
657 // pt->SetPIDsignals(seed2.GetPIDsignals(i),i);
658 // pt->SetPIDTimBin(seed2.GetPIDTimBin(i),i);
661 UInt_t * indexes3 = pt->GetBackupIndexes();
662 for (Int_t i=0;i<200;i++) {
663 if (indexes2[i]==0) break;
664 indexes3[i] = indexes2[i];
666 //AliTRDtrack *pt = seed2;
668 FollowProlongation(t, innerTB);
669 if (t.GetNumberOfClusters() >= foundMin) {
671 //CookLabel(pt, 1-fgkLabelFraction);
675 // cout<<found<<'\r';
677 if(PropagateToTPC(t)) {
678 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
679 // for (Int_t i=0;i<kNPlane;i++) {
680 // seed->SetTRDsignals(pt->GetPIDsignals(i),i);
681 // seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
684 //if not prolongation to TPC - propagate without update
685 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
686 seed2->ResetCovariance(5.);
687 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
689 if (PropagateToTPC(*pt2)) {
690 pt2->CookdEdx(0.,1.);
691 CookdEdxTimBin(*pt2);
692 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
693 // for (Int_t i=0;i<kNPlane;i++) {
694 // seed->SetTRDsignals(pt2->GetPIDsignals(i),i);
695 // seed->SetTRDTimBin(pt2->GetPIDTimBin(i),i);
703 cout<<"Number of loaded seeds: "<<nseed<<endl;
704 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
711 //---------------------------------------------------------------------------
712 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
714 // Starting from current position on track=t this function tries
715 // to extrapolate the track up to timeBin=0 and to confirm prolongation
716 // if a close cluster is found. Returns the number of clusters
717 // expected to be found in sensitive layers
719 Float_t wIndex, wTB, wChi2;
720 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
721 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
722 Float_t wPx, wPy, wPz, wC;
724 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
725 Int_t lastplane = GetLastPlane(&t);
727 Int_t trackIndex = t.GetLabel();
729 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
731 Int_t tryAgain=fMaxGap;
733 Double_t alpha=t.GetAlpha();
734 alpha = TVector2::Phi_0_2pi(alpha);
736 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
737 Double_t radLength, rho, x, dx, y, ymax, z;
739 Int_t expectedNumberOfClusters = 0;
740 Bool_t lookForCluster;
742 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
745 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
747 y = t.GetY(); z = t.GetZ();
749 // first propagate to the inner surface of the current time bin
750 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
751 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
752 if(!t.PropagateTo(x,radLength,rho)) break;
754 ymax = x*TMath::Tan(0.5*alpha);
757 if (!t.Rotate(alpha)) break;
758 if(!t.PropagateTo(x,radLength,rho)) break;
759 } else if (y <-ymax) {
761 if (!t.Rotate(-alpha)) break;
762 if(!t.PropagateTo(x,radLength,rho)) break;
765 y = t.GetY(); z = t.GetZ();
767 // now propagate to the middle plane of the next time bin
768 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
769 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
770 if(!t.PropagateTo(x,radLength,rho)) break;
772 ymax = x*TMath::Tan(0.5*alpha);
775 if (!t.Rotate(alpha)) break;
776 if(!t.PropagateTo(x,radLength,rho)) break;
777 } else if (y <-ymax) {
779 if (!t.Rotate(-alpha)) break;
780 if(!t.PropagateTo(x,radLength,rho)) break;
786 expectedNumberOfClusters++;
787 wIndex = (Float_t) t.GetLabel();
790 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
792 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
793 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
796 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
797 else return expectedNumberOfClusters;
801 wYwindow = (Float_t) road;
802 t.GetPxPyPz(px,py,pz);
806 wC = (Float_t) t.GetC();
807 wSigmaC2 = (Float_t) t.GetSigmaC2();
808 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
809 wSigmaY2 = (Float_t) t.GetSigmaY2();
810 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
817 Double_t maxChi2=fgkMaxChi2;
819 wYclosest = 12345678;
820 wYcorrect = 12345678;
821 wZclosest = 12345678;
822 wZcorrect = 12345678;
823 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
825 // Find the closest correct cluster for debugging purposes
826 if (timeBin&&fVocal) {
827 Float_t minDY = 1000000;
828 for (Int_t i=0; i<timeBin; i++) {
829 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
830 if((c->GetLabel(0) != trackIndex) &&
831 (c->GetLabel(1) != trackIndex) &&
832 (c->GetLabel(2) != trackIndex)) continue;
833 if(TMath::Abs(c->GetY() - y) > minDY) continue;
834 minDY = TMath::Abs(c->GetY() - y);
835 wYcorrect = c->GetY();
836 wZcorrect = c->GetZ();
838 Double_t h01 = GetTiltFactor(c);
839 wChi2 = t.GetPredictedChi2(c, h01);
843 // Now go for the real cluster search
847 //find cluster in history
850 AliTRDcluster * cl0 = timeBin[0];
854 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
855 if (plane>lastplane) continue;
856 Int_t timebin = cl0->GetLocalTimeBin();
857 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
860 Double_t h01 = GetTiltFactor(cl);
861 maxChi2=t.GetPredictedChi2(cl,h01);
863 if ((!cl) && road>fgkWideRoad) {
864 //if (t.GetNumberOfClusters()>4)
865 // cerr<<t.GetNumberOfClusters()
866 // <<"FindProlongation warning: Too broad road !\n";
872 Int_t maxn = timeBin;
873 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
874 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
875 if (c->GetY() > y+road) break;
876 if (c->IsUsed() > 0) continue;
877 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
879 Double_t h01 = GetTiltFactor(c);
880 Double_t chi2=t.GetPredictedChi2(c,h01);
882 if (chi2 > maxChi2) continue;
885 index=timeBin.GetIndex(i);
890 Int_t maxn = timeBin;
891 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
892 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
894 if (c->GetY() > y+road) break;
895 if (c->IsUsed() > 0) continue;
896 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
898 Double_t h01 = GetTiltFactor(c);
899 Double_t chi2=t.GetPredictedChi2(c, h01);
901 if (chi2 > maxChi2) continue;
904 index=timeBin.GetIndex(i);
910 wYclosest = cl->GetY();
911 wZclosest = cl->GetZ();
912 Double_t h01 = GetTiltFactor(cl);
914 if (cl->GetNPads()<5)
915 t.SetSampledEdx(cl->GetQ()/dx);
916 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
917 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
918 Int_t det = cl->GetDetector();
919 Int_t plane = fGeom->GetPlane(det);
921 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
922 //if(!t.Update(cl,maxChi2,index,h01)) {
923 //if(!tryAgain--) return 0;
925 else tryAgain=fMaxGap;
928 //if (tryAgain==0) break;
934 return expectedNumberOfClusters;
939 //___________________________________________________________________
941 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
943 // Starting from current radial position of track <t> this function
944 // extrapolates the track up to outer timebin and in the sensitive
945 // layers confirms prolongation if a close cluster is found.
946 // Returns the number of clusters expected to be found in sensitive layers
949 Float_t wIndex, wTB, wChi2;
950 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
951 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
952 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
954 Int_t trackIndex = t.GetLabel();
955 Int_t tryAgain=fMaxGap;
957 Double_t alpha=t.GetAlpha();
958 TVector2::Phi_0_2pi(alpha);
962 Int_t clusters[1000];
963 for (Int_t i=0;i<1000;i++) clusters[i]=-1;
965 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
966 Double_t radLength, rho, x, dx, y, ymax = 0, z;
967 Bool_t lookForCluster;
969 Int_t expectedNumberOfClusters = 0;
972 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
974 Int_t nRefPlane = fgkFirstPlane;
975 Bool_t isNewLayer = kFALSE;
981 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
986 // first propagate to the outer surface of the current time bin
989 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
990 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
994 if(!t.PropagateTo(x,radLength,rho)) break;
995 // if (!AdjustSector(&t)) break;
997 // MI -fix untill correct material desription will be implemented
999 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material
1000 if (!AdjustSector(&t)) break;
1001 Int_t cross = kFALSE;
1002 Int_t crosz = kFALSE;
1003 if (TMath::Abs(angle - t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
1004 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
1005 if (currentzone==-10) {cross = kTRUE,crosz=kTRUE;} // we are in the frame
1006 if (currentzone>-10){ // layer knows where we are
1007 if (zone==-10) zone = currentzone;
1008 if (zone!=currentzone) {
1013 if (TMath::Abs(t.GetSnp())>0.8 && t.GetBackupTrack()==0) t.MakeBackupTrack();
1015 if (t.GetNCross()==0 && t.GetBackupTrack()==0) t.MakeBackupTrack();
1017 if (t.GetNCross()>4) break;
1022 if (!t.PropagateTo(x,radLength,rho)) break;
1027 // Barrel Tracks [SR, 04.04.2003]
1030 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1031 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
1033 // if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
1036 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1037 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1039 } else {isNewLayer = kFALSE;}
1044 // now propagate to the middle plane of the next time bin
1045 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1047 rho = 1000*2.7; radLength = 24.01; //TEMPORARY - aluminium in between z - will be detected using GeoModeler in future versions
1049 x = fTrSec[s]->GetLayer(nr+1)->GetX();
1050 if(!t.PropagateTo(x,radLength,rho)) break;
1051 if (!AdjustSector(&t)) break;
1053 if(!t.PropagateTo(x,radLength,rho)) break;
1055 if (TMath::Abs(t.GetSnp())>0.95) break;
1060 // if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
1061 // printf("label %d, pl %d, lookForCluster %d \n",
1062 // trackIndex, nr+1, lookForCluster);
1064 if(lookForCluster) {
1065 // if (clusters[nr]==-1) {
1066 // FindClusters(s,nr,nr+30,&t,clusters);
1069 expectedNumberOfClusters++;
1071 if (t.fX>345) t.fNExpectedLast++;
1072 wIndex = (Float_t) t.GetLabel();
1073 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1075 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
1076 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1077 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1078 if((t.GetSigmaY2() + sy2) < 0) break;
1079 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1080 Double_t y=t.GetY(), z=t.GetZ();
1084 wYwindow = (Float_t) road;
1085 wSigmaC2 = (Float_t) t.GetSigmaC2();
1086 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1087 wSigmaY2 = (Float_t) t.GetSigmaY2();
1088 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1091 if (road>fgkWideRoad) {
1092 if (t.GetNumberOfClusters()>4)
1093 cerr<<t.GetNumberOfClusters()
1094 <<"FindProlongation warning: Too broad road !\n";
1098 AliTRDcluster *cl=0;
1101 Double_t maxChi2=fgkMaxChi2;
1106 maxChi2 = 10 * fgkMaxChi2;
1109 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1110 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1111 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
1114 wYclosest = 12345678;
1115 wYcorrect = 12345678;
1116 wZclosest = 12345678;
1117 wZcorrect = 12345678;
1118 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1120 // Find the closest correct cluster for debugging purposes
1121 if (timeBin&&fVocal) {
1123 for (Int_t i=0; i<timeBin; i++) {
1124 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1125 if((c->GetLabel(0) != trackIndex) &&
1126 (c->GetLabel(1) != trackIndex) &&
1127 (c->GetLabel(2) != trackIndex)) continue;
1128 if(TMath::Abs(c->GetY() - y) > minDY) continue;
1129 //minDY = TMath::Abs(c->GetY() - y);
1130 minDY = c->GetY() - y;
1131 wYcorrect = c->GetY();
1132 wZcorrect = c->GetZ();
1134 Double_t h01 = GetTiltFactor(c);
1135 wChi2 = t.GetPredictedChi2(c, h01);
1139 // Now go for the real cluster search
1143 if (clusters[nr+1]>0) {
1144 index = clusters[nr+1];
1145 cl = (AliTRDcluster*)GetCluster(index);
1146 Double_t h01 = GetTiltFactor(cl);
1147 maxChi2=t.GetPredictedChi2(cl,h01);
1152 Int_t maxn = timeBin;
1153 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
1154 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1155 if (c->GetY() > y+road) break;
1156 if (c->IsUsed() > 0) continue;
1157 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1159 Double_t h01 = GetTiltFactor(c);
1160 chi2=t.GetPredictedChi2(c,h01);
1162 if (chi2 > maxChi2) continue;
1165 index=timeBin.GetIndex(i);
1168 if((c->GetLabel(0) != trackIndex) &&
1169 (c->GetLabel(1) != trackIndex) &&
1170 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
1174 Int_t maxn = timeBin;
1175 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
1176 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
1178 if (c->GetY() > y+road) break;
1179 if (c->IsUsed() > 0) continue;
1180 // if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1181 if((c->GetZ()-z)*(c->GetZ()-z) > 12. * sz2) continue;
1184 Double_t h01 = GetTiltFactor(c);
1185 chi2=t.GetPredictedChi2(c,h01);
1187 if (chi2 > maxChi2) continue;
1190 index=timeBin.GetIndex(i);
1195 wYclosest = cl->GetY();
1196 wZclosest = cl->GetZ();
1197 if (cl->GetNPads()<5)
1198 t.SetSampledEdx(cl->GetQ()/dx);
1199 Double_t h01 = GetTiltFactor(cl);
1200 Int_t det = cl->GetDetector();
1201 Int_t plane = fGeom->GetPlane(det);
1204 t.fChi2Last+=maxChi2;
1206 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1207 if(!t.Update(cl,maxChi2,index,h01)) {
1208 if(!tryAgain--) return 0;
1211 else tryAgain=fMaxGap;
1214 if (tryAgain==0) break;
1218 isNewLayer = kFALSE;
1227 return expectedNumberOfClusters;
1232 //---------------------------------------------------------------------------
1233 Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1235 // Starting from current position on track=t this function tries
1236 // to extrapolate the track up to timeBin=0 and to reuse already
1237 // assigned clusters. Returns the number of clusters
1238 // expected to be found in sensitive layers
1239 // get indices of assigned clusters for each layer
1240 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1243 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1244 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1245 Int_t index = t.GetClusterIndex(i);
1246 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1248 Int_t detector=cl->GetDetector();
1249 Int_t localTimeBin=cl->GetLocalTimeBin();
1250 Int_t sector=fGeom->GetSector(detector);
1251 Int_t plane=fGeom->GetPlane(detector);
1253 Int_t trackingSector = CookSectorIndex(sector);
1255 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1256 if(gtb < 0) continue;
1257 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1258 iCluster[layer] = index;
1262 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1264 Double_t alpha=t.GetAlpha();
1265 alpha = TVector2::Phi_0_2pi(alpha);
1267 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1268 Double_t radLength, rho, x, dx, y, ymax, z;
1270 Int_t expectedNumberOfClusters = 0;
1271 Bool_t lookForCluster;
1273 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1276 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1278 y = t.GetY(); z = t.GetZ();
1280 // first propagate to the inner surface of the current time bin
1281 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1282 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1283 if(!t.PropagateTo(x,radLength,rho)) break;
1285 ymax = x*TMath::Tan(0.5*alpha);
1288 if (!t.Rotate(alpha)) break;
1289 if(!t.PropagateTo(x,radLength,rho)) break;
1290 } else if (y <-ymax) {
1292 if (!t.Rotate(-alpha)) break;
1293 if(!t.PropagateTo(x,radLength,rho)) break;
1296 y = t.GetY(); z = t.GetZ();
1298 // now propagate to the middle plane of the next time bin
1299 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1300 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1301 if(!t.PropagateTo(x,radLength,rho)) break;
1303 ymax = x*TMath::Tan(0.5*alpha);
1306 if (!t.Rotate(alpha)) break;
1307 if(!t.PropagateTo(x,radLength,rho)) break;
1308 } else if (y <-ymax) {
1310 if (!t.Rotate(-alpha)) break;
1311 if(!t.PropagateTo(x,radLength,rho)) break;
1314 if(lookForCluster) expectedNumberOfClusters++;
1316 // use assigned cluster
1317 if (!iCluster[nr-1]) continue;
1318 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1319 Double_t h01 = GetTiltFactor(cl);
1320 Double_t chi2=t.GetPredictedChi2(cl, h01);
1321 if (cl->GetNPads()<5) t.SetSampledEdx(cl->GetQ()/dx);
1323 //t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1324 t.Update(cl,chi2,iCluster[nr-1],h01);
1327 return expectedNumberOfClusters;
1330 //___________________________________________________________________
1332 Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1334 // Starting from current radial position of track <t> this function
1335 // extrapolates the track up to radial position <xToGo>.
1336 // Returns 1 if track reaches the plane, and 0 otherwise
1338 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1340 Double_t alpha=t.GetAlpha();
1342 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1343 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1345 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1347 Bool_t lookForCluster;
1348 Double_t radLength, rho, x, dx, y, ymax, z;
1352 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1354 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1356 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1358 y = t.GetY(); z = t.GetZ();
1360 // first propagate to the outer surface of the current time bin
1361 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1362 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
1363 if(!t.PropagateTo(x,radLength,rho)) return 0;
1365 ymax = x*TMath::Tan(0.5*alpha);
1368 if (!t.Rotate(alpha)) return 0;
1369 } else if (y <-ymax) {
1371 if (!t.Rotate(-alpha)) return 0;
1373 if(!t.PropagateTo(x,radLength,rho)) return 0;
1375 y = t.GetY(); z = t.GetZ();
1377 // now propagate to the middle plane of the next time bin
1378 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1379 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
1380 if(!t.PropagateTo(x,radLength,rho)) return 0;
1382 ymax = x*TMath::Tan(0.5*alpha);
1385 if (!t.Rotate(alpha)) return 0;
1386 } else if (y <-ymax) {
1388 if (!t.Rotate(-alpha)) return 0;
1390 if(!t.PropagateTo(x,radLength,rho)) return 0;
1395 //___________________________________________________________________
1397 Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1399 // Starting from current radial position of track <t> this function
1400 // extrapolates the track up to radial position of the outermost
1401 // padrow of the TPC.
1402 // Returns 1 if track reaches the TPC, and 0 otherwise
1404 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1406 Double_t alpha=t.GetAlpha();
1407 alpha = TVector2::Phi_0_2pi(alpha);
1409 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1411 Bool_t lookForCluster;
1412 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
1416 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1417 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1419 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1424 // first propagate to the outer surface of the current time bin
1425 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1426 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1428 if(!t.PropagateTo(x,radLength,rho)) return 0;
1430 if(!t.PropagateTo(x,radLength,rho)) return 0;
1435 // now propagate to the middle plane of the next time bin
1436 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1437 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1439 if(!t.PropagateTo(x,radLength,rho)) return 0;
1441 if(!t.PropagateTo(x,radLength,rho)) return 0;
1446 //_____________________________________________________________________________
1447 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1449 // Fills clusters into TRD tracking_sectors
1450 // Note that the numbering scheme for the TRD tracking_sectors
1451 // differs from that of TRD sectors
1452 cout<<"\n Read Sectors clusters"<<endl;
1453 if (ReadClusters(fClusters,cTree)) {
1454 Error("LoadClusters","Problem with reading the clusters !");
1457 Int_t ncl=fClusters->GetEntriesFast();
1459 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1462 for (Int_t ichamber=0;ichamber<5;ichamber++)
1463 for (Int_t isector=0;isector<18;isector++){
1464 fHoles[ichamber][isector]=kTRUE;
1469 // printf("\r %d left ",ncl);
1470 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1471 Int_t detector=c->GetDetector();
1472 Int_t localTimeBin=c->GetLocalTimeBin();
1473 Int_t sector=fGeom->GetSector(detector);
1474 Int_t plane=fGeom->GetPlane(detector);
1476 Int_t trackingSector = CookSectorIndex(sector);
1477 if (c->GetLabel(0)>0){
1478 Int_t chamber = fGeom->GetChamber(detector);
1479 fHoles[chamber][trackingSector]=kFALSE;
1482 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1483 if(gtb < 0) continue;
1484 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1487 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1493 for (Int_t isector=0;isector<18;isector++){
1494 for (Int_t ichamber=0;ichamber<5;ichamber++)
1495 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1496 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1497 fGeom->IsHole(0,ichamber,17-isector));
1503 //_____________________________________________________________________________
1504 void AliTRDtracker::UnloadClusters()
1507 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1512 nentr = fClusters->GetEntriesFast();
1513 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1516 nentr = fSeeds->GetEntriesFast();
1517 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1519 nentr = fTracks->GetEntriesFast();
1520 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1522 Int_t nsec = AliTRDgeometry::kNsect;
1524 for (i = 0; i < nsec; i++) {
1525 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1526 fTrSec[i]->GetLayer(pl)->Clear();
1532 //__________________________________________________________________________
1533 void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1535 // Creates track seeds using clusters in timeBins=i1,i2
1538 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1542 Double_t x[5], c[15];
1543 Int_t maxSec=AliTRDgeometry::kNsect;
1545 Double_t alpha=AliTRDgeometry::GetAlpha();
1546 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1547 Double_t cs=cos(alpha), sn=sin(alpha);
1548 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1551 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1552 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1554 Double_t x1 =fTrSec[0]->GetX(i1);
1555 Double_t xx2=fTrSec[0]->GetX(i2);
1557 for (Int_t ns=0; ns<maxSec; ns++) {
1559 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1560 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1561 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
1562 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1563 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1565 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1567 for (Int_t is=0; is < r1; is++) {
1568 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1570 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
1572 const AliTRDcluster *cl;
1573 Double_t x2, y2, z2;
1574 Double_t x3=0., y3=0.;
1577 if(turn != 2) continue;
1578 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1580 y2=cl->GetY(); z2=cl->GetZ();
1585 else if (js<nl2+nl) {
1586 if(turn != 1) continue;
1587 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
1589 y2=cl->GetY(); z2=cl->GetZ();
1594 else if (js<nl2+nl+nm) {
1595 if(turn != 1) continue;
1596 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1598 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1600 else if (js<nl2+nl+nm+nu) {
1601 if(turn != 1) continue;
1602 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1603 cl=r2[js-nl2-nl-nm];
1604 y2=cl->GetY(); z2=cl->GetZ();
1610 if(turn != 2) continue;
1611 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
1612 cl=r2[js-nl2-nl-nm-nu];
1613 y2=cl->GetY(); z2=cl->GetZ();
1619 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
1621 Double_t zz=z1 - z1/x1*(x1-x2);
1623 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
1625 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1626 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1630 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1632 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
1634 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1636 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1638 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1640 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
1642 Double_t a=asin(x[2]);
1643 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1645 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
1647 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1648 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
1649 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
1652 Double_t h01 = GetTiltFactor(r1[is]);
1653 Double_t xuFactor = 100.;
1659 sy1=sy1+sz1*h01*h01;
1660 Double_t syz=sz1*(-h01);
1661 // end of tilt changes
1663 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1664 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1665 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1666 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1667 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1668 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1669 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1670 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1671 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1672 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1676 // c[1]=0.; c[2]=sz1;
1677 c[1]=syz; c[2]=sz1*xuFactor;
1678 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1679 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1680 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1681 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1682 c[13]=f30*sy1*f40+f32*sy2*f42;
1683 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
1685 UInt_t index=r1.GetIndex(is);
1687 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1689 Int_t rc=FollowProlongation(*track, i2);
1692 (track->GetNumberOfClusters() <
1693 (outer-inner)*fgkMinClustersInSeed)) delete track;
1695 fSeeds->AddLast(track); fNseeds++;
1696 // cerr<<"\r found seed "<<fNseeds;
1703 //_____________________________________________________________________________
1704 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
1707 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1708 // from the file. The names of the cluster tree and branches
1709 // should match the ones used in AliTRDclusterizer::WriteClusters()
1711 Int_t nsize = Int_t(ClusterTree->GetTotBytes()/(sizeof(AliTRDcluster)));
1712 TObjArray *clusterArray = new TObjArray(nsize+1000);
1714 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1716 Error("ReadClusters","Can't get the branch !");
1719 branch->SetAddress(&clusterArray);
1721 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
1722 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
1724 // Loop through all entries in the tree
1726 AliTRDcluster *c = 0;
1728 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1731 nbytes += ClusterTree->GetEvent(iEntry);
1733 // Get the number of points in the detector
1734 Int_t nCluster = clusterArray->GetEntriesFast();
1735 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
1737 // Loop through all TRD digits
1738 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1739 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
1740 if (c->GetNPads()>3&&(iCluster%3>0)) {
1741 delete clusterArray->RemoveAt(iCluster);
1744 // AliTRDcluster *co = new AliTRDcluster(*c); //remove unnecesary coping - + clusters are together in memory
1745 AliTRDcluster *co = c;
1746 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
1747 Int_t ltb = co->GetLocalTimeBin();
1748 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1749 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
1751 // delete clusterArray->RemoveAt(iCluster);
1752 clusterArray->RemoveAt(iCluster);
1755 // cout<<"Allocated"<<nsize<<"\tLoaded"<<array->GetEntriesFast()<<"\n";
1757 delete clusterArray;
1762 //__________________________________________________________________
1763 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1766 // This cooks a label. Mmmmh, smells good...
1769 Int_t label=123456789, index, i, j;
1770 Int_t ncl=pt->GetNumberOfClusters();
1771 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
1775 // Int_t s[kRange][2];
1776 Int_t **s = new Int_t* [kRange];
1777 for (i=0; i<kRange; i++) {
1778 s[i] = new Int_t[2];
1780 for (i=0; i<kRange; i++) {
1786 for (i=0; i<ncl; i++) {
1787 index=pt->GetClusterIndex(i);
1788 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1794 for (i=0; i<ncl; i++) {
1795 index=pt->GetClusterIndex(i);
1796 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1797 for (Int_t k=0; k<3; k++) {
1798 label=c->GetLabel(k);
1799 labelAdded=kFALSE; j=0;
1801 while ( (!labelAdded) && ( j < kRange ) ) {
1802 if (s[j][0]==label || s[j][1]==0) {
1816 for (i=0; i<kRange; i++) {
1818 max=s[i][1]; label=s[i][0];
1822 for (i=0; i<kRange; i++) {
1828 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1830 pt->SetLabel(label);
1835 //__________________________________________________________________
1836 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1839 // Use clusters, but don't abuse them!
1842 Int_t ncl=t->GetNumberOfClusters();
1843 for (Int_t i=from; i<ncl; i++) {
1844 Int_t index = t->GetClusterIndex(i);
1845 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1851 //_____________________________________________________________________
1852 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
1854 // Parametrised "expected" error of the cluster reconstruction in Y
1856 Double_t s = 0.08 * 0.08;
1860 //_____________________________________________________________________
1861 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
1863 // Parametrised "expected" error of the cluster reconstruction in Z
1865 Double_t s = 9 * 9 /12.;
1869 //_____________________________________________________________________
1870 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
1873 // Returns radial position which corresponds to time bin <localTB>
1874 // in tracking sector <sector> and plane <plane>
1877 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
1878 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1879 return fTrSec[sector]->GetLayer(pl)->GetX();
1884 //_______________________________________________________
1885 AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
1886 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
1889 // AliTRDpropagationLayer constructor
1892 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1893 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
1896 for(Int_t i=0; i < (Int_t) kZones; i++) {
1897 fZc[i]=0; fZmax[i] = 0;
1902 if(fTimeBinIndex >= 0) {
1903 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1904 fIndex = new UInt_t[kMaxClusterPerTimeBin];
1907 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
1918 //_______________________________________________________
1919 void AliTRDtracker::AliTRDpropagationLayer::SetHole(
1920 Double_t Zmax, Double_t Ymax, Double_t rho,
1921 Double_t radLength, Double_t Yc, Double_t Zc)
1924 // Sets hole in the layer
1932 fHoleX0 = radLength;
1936 //_______________________________________________________
1937 AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
1940 // AliTRDtrackingSector Constructor
1943 AliTRDpadPlane *padPlane = 0;
1951 // get holes description from geometry
1952 Bool_t holes[AliTRDgeometry::kNcham];
1953 //printf("sector\t%d\t",gs);
1954 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
1955 holes[icham] = fGeom->IsHole(0,icham,gs);
1956 //printf("%d",holes[icham]);
1960 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
1963 AliTRDpropagationLayer* ppl;
1965 Double_t x, xin, xout, dx, rho, radLength;
1968 // set time bins in the gas of the TPC
1970 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
1971 rho = 0.9e-3; radLength = 28.94;
1973 for(Int_t i=0; i<steps; i++) {
1974 x = xin + i*dx + dx/2;
1975 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
1979 // set time bins in the outer field cage vessel
1981 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
1982 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
1985 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
1986 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
1989 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
1990 steps = 5; dx = (xout - xin)/steps;
1991 for(Int_t i=0; i<steps; i++) {
1992 x = xin + i*dx + dx/2;
1993 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
1997 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
1998 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2001 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2002 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2006 // set time bins in CO2
2008 xin = xout; xout = 275.0;
2009 steps = 50; dx = (xout - xin)/steps;
2010 rho = 1.977e-3; radLength = 36.2;
2012 for(Int_t i=0; i<steps; i++) {
2013 x = xin + i*dx + dx/2;
2014 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2018 // set time bins in the outer containment vessel
2020 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2021 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2024 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2025 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2028 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2029 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2032 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
2033 steps = 10; dx = (xout - xin)/steps;
2034 for(Int_t i=0; i<steps; i++) {
2035 x = xin + i*dx + dx/2;
2036 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2040 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2041 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2044 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2045 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2048 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2049 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2052 Double_t xtrd = (Double_t) fGeom->Rmin();
2054 // add layers between TPC and TRD (Air temporarily)
2055 xin = xout; xout = xtrd;
2056 steps = 50; dx = (xout - xin)/steps;
2057 rho = 1.2e-3; radLength = 36.66;
2059 for(Int_t i=0; i<steps; i++) {
2060 x = xin + i*dx + dx/2;
2061 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2066 // Double_t alpha=AliTRDgeometry::GetAlpha();
2068 // add layers for each of the planes
2070 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2071 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2072 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2073 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2074 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2075 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2076 Double_t dxPlane = dxTEC + dxSpace;
2079 const Int_t kNchambers = AliTRDgeometry::Ncham();
2082 Double_t ymaxsensitive=0;
2083 Double_t *zc = new Double_t[kNchambers];
2084 Double_t *zmax = new Double_t[kNchambers];
2085 Double_t *zmaxsensitive = new Double_t[kNchambers];
2086 // Double_t holeZmax = 1000.; // the whole sector is missing
2088 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2091 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
2092 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2093 for(Int_t i=0; i<steps; i++) {
2094 x = xin + i*dx + dx/2;
2095 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2099 ymax = fGeom->GetChamberWidth(plane)/2.;
2101 // Modidified for new pad plane class, 22.04.05 (C.B.)
2102 // ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2103 padPlane = fPar->GetPadPlane(plane,0);
2104 ymaxsensitive = (padPlane->GetColSize(1)*padPlane->GetNcols()-4)/2.;
2106 for(Int_t ch = 0; ch < kNchambers; ch++) {
2107 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2109 // Modidified for new pad plane class, 22.04.05 (C.B.)
2110 //Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2111 Float_t pad = padPlane->GetRowSize(1);
2112 Float_t row0 = fPar->GetRow0(plane,ch,0);
2113 Int_t nPads = fPar->GetRowMax(plane,ch,0);
2114 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2115 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2116 zc[ch] = (pad * nPads)/2 + row0;
2117 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2121 dx = fPar->GetDriftVelocity()
2122 / fPar->GetSamplingFrequency();
2123 rho = 0.00295 * 0.85; radLength = 11.0;
2125 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2126 Double_t xbottom = x0 - dxDrift;
2127 Double_t xtop = x0 + dxAmp;
2129 // Amplification region
2130 steps = (Int_t) (dxAmp/dx);
2132 for(tb = 0; tb < steps; tb++) {
2133 x = x0 + tb * dx + dx/2;
2134 tbIndex = CookTimeBinIndex(plane, -tb-1);
2135 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2136 ppl->SetYmax(ymax,ymaxsensitive);
2137 ppl->SetZ(zc, zmax, zmaxsensitive);
2138 ppl->SetHoles(holes);
2141 tbIndex = CookTimeBinIndex(plane, -steps);
2142 x = (x + dx/2 + xtop)/2;
2144 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2145 ppl->SetYmax(ymax,ymaxsensitive);
2146 ppl->SetZ(zc, zmax,zmaxsensitive);
2147 ppl->SetHoles(holes);
2151 dx = fPar->GetDriftVelocity()
2152 / fPar->GetSamplingFrequency();
2153 steps = (Int_t) (dxDrift/dx);
2155 for(tb = 0; tb < steps; tb++) {
2156 x = x0 - tb * dx - dx/2;
2157 tbIndex = CookTimeBinIndex(plane, tb);
2159 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2160 ppl->SetYmax(ymax,ymaxsensitive);
2161 ppl->SetZ(zc, zmax, zmaxsensitive);
2162 ppl->SetHoles(holes);
2165 tbIndex = CookTimeBinIndex(plane, steps);
2166 x = (x - dx/2 + xbottom)/2;
2168 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
2169 ppl->SetYmax(ymax,ymaxsensitive);
2170 ppl->SetZ(zc, zmax, zmaxsensitive);
2171 ppl->SetHoles(holes);
2175 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2176 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
2177 ppl->SetYmax(ymax,ymaxsensitive);
2178 ppl->SetZ(zc, zmax,zmax);
2179 ppl->SetHoles(holes);
2183 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
2184 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
2185 for(Int_t i=0; i<steps; i++) {
2186 x = xin + i*dx + dx/2;
2187 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2188 ppl->SetYmax(ymax,ymaxsensitive);
2189 ppl->SetZ(zc, zmax,zmax);
2190 ppl->SetHoles(holes);
2194 // Space between the chambers, air
2195 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
2196 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2197 for(Int_t i=0; i<steps; i++) {
2198 x = xin + i*dx + dx/2;
2199 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2204 // Space between the TRD and RICH
2205 Double_t xRICH = 500.;
2206 xin = xout; xout = xRICH;
2207 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
2208 for(Int_t i=0; i<steps; i++) {
2209 x = xin + i*dx + dx/2;
2210 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
2217 delete [] zmaxsensitive;
2221 //______________________________________________________
2223 Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2226 // depending on the digitization parameters calculates "global"
2227 // time bin index for timebin <localTB> in plane <plane>
2230 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2231 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2232 Double_t dx = (Double_t) fPar->GetDriftVelocity()
2233 / fPar->GetSamplingFrequency();
2235 Int_t tbAmp = fPar->GetTimeBefore();
2236 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
2237 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
2238 Int_t tbDrift = fPar->GetTimeMax();
2239 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2241 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
2243 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
2246 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2247 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
2254 //______________________________________________________
2256 void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2259 // For all sensitive time bins sets corresponding layer index
2260 // in the array fTimeBins
2265 for(Int_t i = 0; i < fN; i++) {
2266 index = fLayers[i]->GetTimeBinIndex();
2268 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2270 if(index < 0) continue;
2271 if(index >= (Int_t) kMaxTimeBinIndex) {
2272 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2273 printf(" index %d exceeds allowed maximum of %d!\n",
2274 index, kMaxTimeBinIndex-1);
2277 fTimeBinIndex[index] = i;
2280 Double_t x1, dx1, x2, dx2, gap;
2282 for(Int_t i = 0; i < fN-1; i++) {
2283 x1 = fLayers[i]->GetX();
2284 dx1 = fLayers[i]->GetdX();
2285 x2 = fLayers[i+1]->GetX();
2286 dx2 = fLayers[i+1]->GetdX();
2287 gap = (x2 - dx2/2) - (x1 + dx1/2);
2288 // if(gap < -0.01) {
2289 // printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2290 // printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2293 // printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2294 // printf(" (%f - %f) - (%f + %f) = %f\n",
2295 // x2, dx2/2, x1, dx1, gap);
2301 //______________________________________________________
2304 Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2307 // Returns the number of time bin which in radial position is closest to <x>
2310 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2311 if(x <= fLayers[0]->GetX()) return 0;
2313 Int_t b=0, e=fN-1, m=(b+e)/2;
2314 for (; b<e; m=(b+e)/2) {
2315 if (x > fLayers[m]->GetX()) b=m+1;
2318 if(TMath::Abs(x - fLayers[m]->GetX()) >
2319 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2324 //______________________________________________________
2326 Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2329 // Returns number of the innermost SENSITIVE propagation layer
2332 return GetLayerNumber(0);
2335 //______________________________________________________
2337 Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2340 // Returns number of the outermost SENSITIVE time bin
2343 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2346 //______________________________________________________
2348 Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2351 // Returns number of SENSITIVE time bins
2355 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2356 layer = GetLayerNumber(tb);
2362 //______________________________________________________
2364 void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2367 // Insert layer <pl> in fLayers array.
2368 // Layers are sorted according to X coordinate.
2370 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2371 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2374 if (fN==0) {fLayers[fN++] = pl; return;}
2375 Int_t i=Find(pl->GetX());
2377 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2378 fLayers[i]=pl; fN++;
2382 //______________________________________________________
2384 Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2387 // Returns index of the propagation layer nearest to X
2390 if (x <= fLayers[0]->GetX()) return 0;
2391 if (x > fLayers[fN-1]->GetX()) return fN;
2392 Int_t b=0, e=fN-1, m=(b+e)/2;
2393 for (; b<e; m=(b+e)/2) {
2394 if (x > fLayers[m]->GetX()) b=m+1;
2400 //______________________________________________________
2401 void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2404 // set centers and the width of sectors
2405 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2406 fZc[icham] = center[icham];
2407 fZmax[icham] = w[icham];
2408 fZmaxSensitive[icham] = wsensitive[icham];
2409 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2412 //______________________________________________________
2414 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2417 // set centers and the width of sectors
2419 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2420 fIsHole[icham] = holes[icham];
2421 if (holes[icham]) fHole = kTRUE;
2427 void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
2428 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
2429 Bool_t &lookForCluster) const
2432 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
2433 // and sensitivity <lookForCluster> in point <y,z>
2436 Double_t alpha = AliTRDgeometry::GetAlpha();
2437 Double_t ymax = fX*TMath::Tan(0.5*alpha);
2443 lookForCluster = kFALSE;
2444 Bool_t cross =kFALSE;
2447 if ( (ymax-TMath::Abs(y))<3.){ //cross material
2453 // check dead regions in sensitive volume
2456 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2457 if (TMath::Abs(z - fZc[ch]) > fZmax[ch]) continue; //not in given zone
2459 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2460 if (fTimeBinIndex>=0) lookForCluster = !(fIsHole[zone]);
2461 if(TMath::Abs(y) > fYmaxSensitive){
2462 lookForCluster = kFALSE;
2464 if (fIsHole[zone]) {
2470 rho = 2.7; radLength = 24.01; //aluminium in between
2474 if (fTimeBinIndex>=0) return;
2478 if (fHole==kFALSE) return;
2480 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2481 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2492 Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2496 if (fTimeBinIndex < 0) return -20; //unknown
2497 Int_t zone=-10; // dead zone
2498 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2499 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2506 //______________________________________________________
2508 void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
2511 // Insert cluster in cluster array.
2512 // Clusters are sorted according to Y coordinate.
2514 if(fTimeBinIndex < 0) {
2515 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2519 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2520 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2523 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2524 Int_t i=Find(c->GetY());
2525 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2526 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2527 fIndex[i]=index; fClusters[i]=c; fN++;
2530 //______________________________________________________
2532 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2534 // Returns index of the cluster nearest in Y
2536 if (y <= fClusters[0]->GetY()) return 0;
2537 if (y > fClusters[fN-1]->GetY()) return fN;
2538 Int_t b=0, e=fN-1, m=(b+e)/2;
2539 for (; b<e; m=(b+e)/2) {
2540 if (y > fClusters[m]->GetY()) b=m+1;
2546 //---------------------------------------------------------
2548 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2550 // Returns correction factor for tilted pads geometry
2553 AliTRDpadPlane *padPlane = fPar->GetPadPlane(0,0);
2554 Double_t h01 = sin(TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
2555 Int_t det = c->GetDetector();
2556 Int_t plane = fGeom->GetPlane(det);
2558 //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2559 if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
2561 if(fNoTilt) h01 = 0;
2567 void AliTRDtracker::CookdEdxTimBin(AliTRDtrack& TRDtrack)
2569 // *** ADDED TO GET MORE INFORMATION FOR TRD PID ---- PS
2570 // This is setting fdEdxPlane and fTimBinPlane
2571 // Sums up the charge in each plane for track TRDtrack and also get the
2572 // Time bin for Max. Cluster
2573 // Prashant Shukla (shukla@physi.uni-heidelberg.de)
2575 // const Int_t kNPlane = AliTRDgeometry::Nplan();
2576 // const Int_t kNPlane = 6;
2577 Double_t clscharge[kNPlane], maxclscharge[kNPlane];
2578 Int_t nCluster[kNPlane], timebin[kNPlane];
2580 //Initialization of cluster charge per plane.
2581 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2582 clscharge[iPlane] = 0.0;
2583 nCluster[iPlane] = 0;
2584 timebin[iPlane] = -1;
2585 maxclscharge[iPlane] = 0.0;
2588 // Loop through all clusters associated to track TRDtrack
2589 Int_t nClus = TRDtrack.GetNumberOfClusters(); // from Kalmantrack
2590 for (Int_t iClus = 0; iClus < nClus; iClus++) {
2591 Double_t charge = TRDtrack.GetClusterdQdl(iClus);
2592 Int_t index = TRDtrack.GetClusterIndex(iClus);
2593 AliTRDcluster *TRDcluster = (AliTRDcluster *) GetCluster(index);
2594 if (!TRDcluster) continue;
2595 Int_t tb = TRDcluster->GetLocalTimeBin();
2597 Int_t detector = TRDcluster->GetDetector();
2598 Int_t iPlane = fGeom->GetPlane(detector);
2599 clscharge[iPlane] = clscharge[iPlane]+charge;
2600 if(charge > maxclscharge[iPlane]) {
2601 maxclscharge[iPlane] = charge;
2602 timebin[iPlane] = tb;
2605 } // end of loop over cluster
2607 // Setting the fdEdxPlane and fTimBinPlane variabales
2608 Double_t Total_ch = 0;
2609 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2610 // Quality control of TRD track.
2611 if (nCluster[iPlane]<= 5) {
2612 clscharge[iPlane]=0.0;
2615 if (nCluster[iPlane]) clscharge[iPlane] /= nCluster[iPlane];
2616 TRDtrack.SetPIDsignals(clscharge[iPlane], iPlane);
2617 TRDtrack.SetPIDTimBin(timebin[iPlane], iPlane);
2618 Total_ch= Total_ch+clscharge[iPlane];
2621 // Int_t nc=TRDtrack.GetNumberOfClusters();
2623 // for (i=0; i<nc; i++) dedx += TRDtrack.GetClusterdQdl(i);
2625 // for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2626 // TRDtrack.SetPIDsignals(dedx, iPlane);
2627 // TRDtrack.SetPIDTimBin(timbin[iPlane], iPlane);
2630 } // end of function
2633 Int_t AliTRDtracker::FindClusters(Int_t sector, Int_t t0, Int_t t1, AliTRDtrack * track, Int_t *clusters)
2637 // try to find nearest clusters to the track in timebins from t0 to t1
2640 Double_t x[1000],yt[1000],zt[10000];
2641 Double_t dz[1000],dy[10000];
2642 Int_t indexes[2][10000];
2643 AliTRDcluster *cl[2][10000];
2645 for (Int_t it=t0;it<t1; it++){
2647 indexes[0][it]=-2; //reset indexes1
2648 indexes[1][it]=-2; //reset indexes1
2658 Double_t x0 = track->GetX();
2659 Double_t sy2=ExpectedSigmaY2(x0,track->GetTgl(),track->GetPt());
2660 Double_t sz2=ExpectedSigmaZ2(x0,track->GetTgl());
2661 Double_t road = 10.*sqrt(track->GetSigmaY2() + sy2);
2665 for (Int_t it=t0;it<t1;it++){
2666 Double_t maxChi2=fgkMaxChi2;
2667 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(it));
2668 if (timeBin==0) continue; // no indexes1
2669 Int_t maxn = timeBin;
2670 x[it] = timeBin.GetX();
2672 if (!track->GetProlongation(x[it],y,z)) continue;
2675 Double_t chi2 =1000000;
2678 // find nearest cluster at given pad
2679 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
2680 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
2681 Double_t h01 = GetTiltFactor(c);
2682 if (c->GetY() > y+road) break;
2683 if (c->IsUsed() > 0) continue;
2684 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
2685 chi2=track->GetPredictedChi2(c,h01);
2686 if (chi2 > maxChi2) continue;
2689 indexes[0][it] =timeBin.GetIndex(i);
2692 // find nearest cluster also in adjacent 2 pads
2694 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
2695 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
2696 Double_t h01 = GetTiltFactor(c);
2697 if (c->GetY() > y+road) break;
2698 if (c->IsUsed() > 0) continue;
2699 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
2700 chi2=track->GetPredictedChi2(c,h01);
2701 if (chi2 > maxChi2) continue;
2704 indexes[1][it]=timeBin.GetIndex(i);
2707 indexes[0][it]=timeBin.GetIndex(i);
2710 if (cl[0][it]||cl[1][it]) nfound++;
2713 if (nfound<5) return 0;
2715 // choose one of the variants
2717 Int_t changes[2]={0,0};
2718 Float_t sigmay[2]={1000,1000};
2719 Float_t meany[2] ={1000,1000};
2720 Float_t meanz[2] ={1000,1000};
2721 Int_t sumall[2] ={0,0};
2722 Int_t ngood[2] ={0,0};
2723 Int_t nbad[2] ={0,0};
2726 for (Int_t ih=0; ih<2;ih++){
2727 Float_t lastz =-10000;
2733 // how many changes ++ mean z ++mean y ++ sigma y
2734 for (Int_t it=t0;it<t1;it++){
2735 if (!cl[ih][it]) continue;
2737 if (lastz<-9999) lastz = cl[ih][it]->GetZ();
2738 if (TMath::Abs(lastz-cl[ih][it]->GetZ())>1) {
2739 lastz = cl[ih][it]->GetZ();
2742 sumz = cl[ih][it]->GetZ();
2744 Double_t h01 = GetTiltFactor(cl[ih][it]);
2745 dz[it] = cl[ih][it]->GetZ()- zt[it];
2746 dy[it] = cl[ih][it]->GetY()+ dz[it]*h01 -yt[it];
2749 sumdy2+= dy[it]*dy[it];
2750 Int_t label = TMath::Abs(track->GetLabel());
2751 if (cl[ih][it]->GetLabel(0)==label || cl[ih][it]->GetLabel(1)==label||cl[ih][it]->GetLabel(2)==label){
2759 meanz[ih] = sumz/sum;
2760 meany[ih] = sumdy/sum;
2761 sigmay[ih] = TMath::Sqrt(sumdy2/sum-meany[ih]*meany[ih]);
2766 if (sumall[0]<sumall[1]-2&&sigmay[1]<0.1){
2767 if (sigmay[1]<sigmay[0]) best = 1; // if sigma is better
2770 Float_t quality0 = (sigmay[0]+TMath::Abs(meany[0]))*(1.+Float_t(changes[0])/Float_t(sumall[0]));
2771 Float_t quality1 = (sigmay[1]+TMath::Abs(meany[1]))*(1.+Float_t(changes[1])/Float_t(sumall[1]));
2773 if (quality0>quality1){
2779 for (Int_t it=t0;it<t1;it++){
2780 if (!cl[best][it]) continue;
2781 Double_t h01 = GetTiltFactor(cl[best][it]);
2782 dz[it] = cl[best][it]->GetZ()- zt[it];
2783 dy[it] = cl[best][it]->GetY()+ dz[it]*h01 -yt[it];
2785 if (TMath::Abs(dy[it])<2.5*sigmay[best])
2786 clusters[it] = indexes[best][it];