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 //
21 // Based on Kalman filltering approach //
23 ///////////////////////////////////////////////////////////////////////////////
25 #include <Riostream.h>
29 #include <TObjArray.h>
31 #include "AliTRDgeometry.h"
32 #include "AliTRDpadPlane.h"
33 #include "AliTRDgeometryFull.h"
34 #include "AliTRDcluster.h"
35 #include "AliTRDtrack.h"
36 #include "AliTRDseed.h"
39 #include "AliTRDcalibDB.h"
40 #include "AliTRDCommonParam.h"
42 #include "TTreeStream.h"
44 #include "AliTRDtracker.h"
45 #include "TLinearFitter.h"
46 #include "AliRieman.h"
47 #include "AliTrackPointArray.h"
48 #include "AliAlignObj.h"
49 #include "AliTRDReconstructor.h"
51 ClassImp(AliTRDtracker)
53 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5;
54 const Float_t AliTRDtracker::fgkLabelFraction = 0.8;
55 const Double_t AliTRDtracker::fgkMaxChi2 = 12.;
56 const Double_t AliTRDtracker::fgkMaxSnp = 0.95; // correspond to tan = 3
57 const Double_t AliTRDtracker::fgkMaxStep = 2.; // maximal step size in propagation
59 //_____________________________________________________________________________
60 AliTRDtracker::AliTRDtracker():AliTracker(),
73 // Default constructor
76 for(Int_t i=0;i<kTrackingSectors;i++) fTrSec[i]=0;
77 for(Int_t j=0;j<5;j++)
78 for(Int_t k=0;k<18;k++) fHoles[j][k]=kFALSE;
83 //_____________________________________________________________________________
84 AliTRDtracker::AliTRDtracker(const TFile *geomfile):AliTracker()
90 fAddTRDseeds = kFALSE;
94 TDirectory *savedir=gDirectory;
95 TFile *in=(TFile*)geomfile;
97 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
98 printf(" FULL TRD geometry and DEFAULT TRD parameter will be used\n");
102 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
106 // printf("Found geometry version %d on file \n", fGeom->IsVersion());
109 printf("AliTRDtracker::AliTRDtracker(): can't find TRD geometry!\n");
110 fGeom = new AliTRDgeometryFull();
112 fGeom->ReadGeoMatrices();
117 fClusters = new TObjArray(2000);
119 fSeeds = new TObjArray(2000);
121 fTracks = new TObjArray(1000);
123 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
124 Int_t trS = CookSectorIndex(geomS);
125 fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS);
126 for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
127 fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
130 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(0,0);
131 Float_t tiltAngle = TMath::Abs(padPlane->GetTiltingAngle());
132 if(tiltAngle < 0.1) {
136 fTimeBinsPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
138 fDebugStreamer = new TTreeSRedirector("TRDdebug.root");
144 //_____________________________________________________________________________
145 AliTRDtracker::~AliTRDtracker()
148 // Destructor of AliTRDtracker
165 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
166 delete fTrSec[geomS];
168 if (fDebugStreamer) {
169 //fDebugStreamer->Close();
170 delete fDebugStreamer;
175 //_____________________________________________________________________________
176 Int_t AliTRDtracker::LocalToGlobalID(Int_t lid)
179 // Transform internal TRD ID to global detector ID
182 Int_t isector = fGeom->GetSector(lid);
183 Int_t ichamber= fGeom->GetChamber(lid);
184 Int_t iplan = fGeom->GetPlane(lid);
186 AliAlignObj::ELayerID iLayer = AliAlignObj::kTRD1;
189 iLayer = AliAlignObj::kTRD1;
192 iLayer = AliAlignObj::kTRD2;
195 iLayer = AliAlignObj::kTRD3;
198 iLayer = AliAlignObj::kTRD4;
201 iLayer = AliAlignObj::kTRD5;
204 iLayer = AliAlignObj::kTRD6;
207 Int_t modId = isector*fGeom->Ncham()+ichamber;
208 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,modId);
214 //_____________________________________________________________________________
215 Int_t AliTRDtracker::GlobalToLocalID(Int_t gid)
218 // Transform global detector ID to local detector ID
222 AliAlignObj::ELayerID layerId = AliAlignObj::VolUIDToLayer(gid, modId);
223 Int_t isector = modId/fGeom->Ncham();
224 Int_t ichamber = modId%fGeom->Ncham();
227 case AliAlignObj::kTRD1:
230 case AliAlignObj::kTRD2:
233 case AliAlignObj::kTRD3:
236 case AliAlignObj::kTRD4:
239 case AliAlignObj::kTRD5:
242 case AliAlignObj::kTRD6:
248 if (iLayer<0) return -1;
249 Int_t lid = fGeom->GetDetector(iLayer,ichamber,isector);
255 //_____________________________________________________________________________
256 Bool_t AliTRDtracker::Transform(AliTRDcluster * cluster)
259 // Transform something ... whatever ...
262 const Double_t kX0shift = 2.52; // magic constants for geo manager transformation
263 const Double_t kX0shift5 = 3.05; //
266 // apply alignment and calibration to transform cluster
269 Int_t detector = cluster->GetDetector();
270 Int_t plane = fGeom->GetPlane(cluster->GetDetector());
271 Int_t chamber = fGeom->GetChamber(cluster->GetDetector());
272 Int_t sector = fGeom->GetSector(cluster->GetDetector());
274 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
275 Double_t driftX = TMath::Max(cluster->GetX()-dxAmp*0.5,0.); // drift distance
279 Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
280 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift);
282 AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
283 AliTRDpadPlane * padPlane = commonParam->GetPadPlane(plane,chamber);
284 Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
285 Double_t localPos[3], localPosTracker[3];
286 localPos[0] = -cluster->GetX();
287 localPos[1] = cluster->GetY() - driftX*exB;
288 localPos[2] = cluster->GetZ() -zshiftIdeal;
290 cluster->SetY(cluster->GetY() - driftX*exB);
291 Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
292 cluster->SetX(xplane- cluster->GetX());
294 TGeoHMatrix * matrix = fGeom->GetCorrectionMatrix(cluster->GetDetector());
296 // no matrix found - if somebody used geometry with holes
297 AliError("Invalid Geometry - Default Geometry used\n");
300 matrix->LocalToMaster(localPos, localPosTracker);
304 if (AliTRDReconstructor::StreamLevel()>1){
305 (*fDebugStreamer)<<"Transform"<<
308 "Detector="<<detector<<
311 "Chamber="<<chamber<<
312 "lx0="<<localPosTracker[0]<<
313 "ly0="<<localPosTracker[1]<<
314 "lz0="<<localPosTracker[2]<<
319 cluster->SetX(localPosTracker[0]+kX0shift5);
321 cluster->SetX(localPosTracker[0]+kX0shift);
323 cluster->SetY(localPosTracker[1]);
324 cluster->SetZ(localPosTracker[2]);
330 //_____________________________________________________________________________
331 // Bool_t AliTRDtracker::Transform(AliTRDcluster * cluster)
335 // const Double_t kDriftCorrection = 1.01; // drift coeficient correction
336 // const Double_t kTime0Cor = 0.32; // time0 correction
338 // const Double_t kX0shift = 2.52;
339 // const Double_t kX0shift5 = 3.05;
342 // // apply alignment and calibration to transform cluster
345 // Int_t detector = cluster->GetDetector();
346 // Int_t plane = fGeom->GetPlane(cluster->GetDetector());
347 // Int_t chamber = fGeom->GetChamber(cluster->GetDetector());
348 // Int_t sector = fGeom->GetSector(cluster->GetDetector());
350 // Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
351 // Double_t driftX = TMath::Max(cluster->GetX()-dxAmp*0.5,0.); // drift distance
355 // Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
356 // Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift);
359 // AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
360 // AliTRDpadPlane * padPlane = commonParam->GetPadPlane(plane,chamber);
361 // Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
362 // Double_t localPos[3], globalPos[3], localPosTracker[3], localPosTracker2[3];
363 // localPos[2] = -cluster->GetX();
364 // localPos[0] = cluster->GetY() - driftX*exB;
365 // localPos[1] = cluster->GetZ() -zshiftIdeal;
366 // TGeoHMatrix * matrix = fGeom->GetGeoMatrix(cluster->GetDetector());
367 // matrix->LocalToMaster(localPos, globalPos);
369 // Double_t sectorAngle = 20.*(sector%18)+10;
370 // TGeoHMatrix rotSector;
371 // rotSector.RotateZ(sectorAngle);
372 // rotSector.LocalToMaster(globalPos, localPosTracker);
375 // TGeoHMatrix matrix2(*matrix);
376 // matrix2.MultiplyLeft(&rotSector);
377 // matrix2.LocalToMaster(localPos,localPosTracker2);
381 // cluster->SetY(cluster->GetY() - driftX*exB);
382 // Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
383 // cluster->SetX(xplane- kDriftCorrection*(cluster->GetX()-kTime0Cor));
384 // (*fDebugStreamer)<<"Transform"<<
386 // "matrix.="<<matrix<<
387 // "matrix2.="<<&matrix2<<
388 // "Detector="<<detector<<
389 // "Sector="<<sector<<
391 // "Chamber="<<chamber<<
392 // "lx0="<<localPosTracker[0]<<
393 // "ly0="<<localPosTracker[1]<<
394 // "lz0="<<localPosTracker[2]<<
395 // "lx2="<<localPosTracker2[0]<<
396 // "ly2="<<localPosTracker2[1]<<
397 // "lz2="<<localPosTracker2[2]<<
401 // cluster->SetX(localPosTracker[0]+kX0shift5);
403 // cluster->SetX(localPosTracker[0]+kX0shift);
405 // cluster->SetY(localPosTracker[1]);
406 // cluster->SetZ(localPosTracker[2]);
410 //_____________________________________________________________________________
411 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track)
414 // Rotates the track when necessary
417 Double_t alpha = AliTRDgeometry::GetAlpha();
418 Double_t y = track->GetY();
419 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
421 //Int_t ns = AliTRDgeometry::kNsect;
422 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
426 if (!track->Rotate(alpha)) return kFALSE;
427 } else if (y <-ymax) {
429 if (!track->Rotate(-alpha)) return kFALSE;
436 //_____________________________________________________________________________
437 AliTRDcluster *AliTRDtracker::GetCluster(AliTRDtrack *track, Int_t plane
438 , Int_t timebin, UInt_t &index)
441 // Try to find cluster in the backup list
444 AliTRDcluster * cl =0;
445 Int_t *indexes = track->GetBackupIndexes();
446 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
447 if (indexes[i]==0) break;
448 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
450 if (cli->GetLocalTimeBin()!=timebin) continue;
451 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
463 //_____________________________________________________________________________
464 Int_t AliTRDtracker::GetLastPlane(AliTRDtrack * track)
467 // Return last updated plane
471 Int_t *indexes = track->GetBackupIndexes();
472 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
473 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
475 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
476 if (iplane>lastplane) {
485 //_____________________________________________________________________________
486 Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
489 // Finds tracks within the TRD. The ESD event is expected to contain seeds
490 // at the outer part of the TRD. The seeds
491 // are found within the TRD if fAddTRDseeds is TRUE.
492 // The tracks are propagated to the innermost time bin
493 // of the TRD and the ESD event is updated
496 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
497 Float_t foundMin = fgkMinClustersInTrack * timeBins;
500 // Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
502 Int_t n = event->GetNumberOfTracks();
503 for (Int_t i=0; i<n; i++) {
504 AliESDtrack* seed=event->GetTrack(i);
505 ULong_t status=seed->GetStatus();
506 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
507 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
510 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
511 //seed2->ResetCovariance();
512 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
514 FollowProlongation(t);
515 if (t.GetNumberOfClusters() >= foundMin) {
517 CookLabel(pt, 1-fgkLabelFraction);
521 // cout<<found<<'\r';
524 if (PropagateToX(t,xTPC,fgkMaxStep)) {
525 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
531 cout<<"Number of loaded seeds: "<<nseed<<endl;
532 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
534 // after tracks from loaded seeds are found and the corresponding
535 // clusters are used, look for additional seeds from TRD
538 cout<<"Total number of found tracks: "<<found<<endl;
544 //_____________________________________________________________________________
545 Int_t AliTRDtracker::PropagateBack(AliESD* event)
548 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
549 // backpropagated by the TPC tracker. Each seed is first propagated
550 // to the TRD, and then its prolongation is searched in the TRD.
551 // If sufficiently long continuation of the track is found in the TRD
552 // the track is updated, otherwise it's stored as originaly defined
553 // by the TPC tracker.
557 Float_t foundMin = 20;
558 Int_t n = event->GetNumberOfTracks();
561 Float_t *quality =new Float_t[n];
562 Int_t *index =new Int_t[n];
563 for (Int_t i=0; i<n; i++) {
564 AliESDtrack* seed=event->GetTrack(i);
565 Double_t covariance[15];
566 seed->GetExternalCovariance(covariance);
567 quality[i] = covariance[0]+covariance[2];
569 TMath::Sort(n,quality,index,kFALSE);
571 for (Int_t i=0; i<n; i++) {
572 // AliESDtrack* seed=event->GetTrack(i);
573 AliESDtrack* seed=event->GetTrack(index[i]);
575 ULong_t status=seed->GetStatus();
576 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
577 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
579 Int_t lbl = seed->GetLabel();
580 AliTRDtrack *track = new AliTRDtrack(*seed);
581 track->SetSeedLabel(lbl);
582 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup); //make backup
584 Float_t p4 = track->GetC();
586 Int_t expectedClr = FollowBackProlongation(*track);
587 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)<0.2 || TMath::Abs(track->GetPt())>0.8 ) {
589 //make backup for back propagation
591 Int_t foundClr = track->GetNumberOfClusters();
592 if (foundClr >= foundMin) {
594 CookdEdxTimBin(*track);
595 CookLabel(track, 1-fgkLabelFraction);
596 if (track->GetBackupTrack()) UseClusters(track->GetBackupTrack());
597 if(track->GetChi2()/track->GetNumberOfClusters()<4) { // sign only gold tracks
598 if (seed->GetKinkIndex(0)==0&&TMath::Abs(track->GetPt())<1.5 ) UseClusters(track);
600 Bool_t isGold = kFALSE;
602 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
603 // seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
604 if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
607 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
608 // seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
609 if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
612 if (!isGold && track->GetBackupTrack()){
613 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
614 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
615 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
619 if (track->StatusForTOF()>0 &&track->fNCross==0 && Float_t(track->fN)/Float_t(track->fNExpected)>0.4){
620 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
624 // Debug part of tracking
625 TTreeSRedirector& cstream = *fDebugStreamer;
626 Int_t eventNr = event->GetEventNumber();
627 if (AliTRDReconstructor::StreamLevel()>0){
628 if (track->GetBackupTrack()){
630 "EventNr="<<eventNr<<
633 "trdback.="<<track->GetBackupTrack()<<
637 "EventNr="<<eventNr<<
645 //Propagation to the TOF (I.Belikov)
646 if (track->GetStop()==kFALSE){
649 Double_t c2=track->GetC()*xtof - track->GetEta();
650 if (TMath::Abs(c2)>=0.99) {
654 Double_t xTOF0 = 370. ;
655 PropagateToX(*track,xTOF0,fgkMaxStep);
657 //energy losses taken to the account - check one more time
658 c2=track->GetC()*xtof - track->GetEta();
659 if (TMath::Abs(c2)>=0.99) {
665 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
666 Double_t y=track->GetYat(xtof);
668 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
672 } else if (y <-ymax) {
673 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
679 if (track->PropagateTo(xtof)) {
680 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
681 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
682 seed->SetTRDsignals(track->GetPIDsignals(i),i);
683 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
685 // seed->SetTRDtrack(new AliTRDtrack(*track));
686 if (track->GetNumberOfClusters()>foundMin) found++;
689 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
690 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
691 //seed->SetStatus(AliESDtrack::kTRDStop);
692 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
693 seed->SetTRDsignals(track->GetPIDsignals(i),i);
694 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
696 //seed->SetTRDtrack(new AliTRDtrack(*track));
700 seed->SetTRDQuality(track->StatusForTOF());
701 seed->SetTRDBudget(track->fBudget[0]);
705 //End of propagation to the TOF
706 //if (foundClr>foundMin)
707 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
712 cerr<<"Number of seeds: "<<fNseeds<<endl;
713 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
715 if (AliTRDReconstructor::SeedingOn()) MakeSeedsMI(3,5,event); //new seeding
717 fSeeds->Clear(); fNseeds=0;
725 //_____________________________________________________________________________
726 Int_t AliTRDtracker::RefitInward(AliESD* event)
729 // Refits tracks within the TRD. The ESD event is expected to contain seeds
730 // at the outer part of the TRD.
731 // The tracks are propagated to the innermost time bin
732 // of the TRD and the ESD event is updated
733 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
736 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
737 Float_t foundMin = fgkMinClustersInTrack * timeBins;
740 // Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
743 Int_t n = event->GetNumberOfTracks();
744 for (Int_t i=0; i<n; i++) {
745 AliESDtrack* seed=event->GetTrack(i);
746 new(&seed2) AliTRDtrack(*seed);
747 if (seed2.GetX()<270){
748 seed->UpdateTrackParams(&seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
752 ULong_t status=seed->GetStatus();
753 if ( (status & AliESDtrack::kTRDout ) == 0 ) {
756 if ( (status & AliESDtrack::kTRDin) != 0 ) {
760 // if (1/seed2.Get1Pt()>1.5&& seed2.GetX()>260.) {
761 // Double_t oldx = seed2.GetX();
762 // seed2.PropagateTo(500.);
763 // seed2.ResetCovariance(1.);
764 // seed2.PropagateTo(oldx);
767 // seed2.ResetCovariance(5.);
770 AliTRDtrack *pt = new AliTRDtrack(seed2,seed2.GetAlpha());
771 Int_t * indexes2 = seed2.GetIndexes();
772 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
773 pt->SetPIDsignals(seed2.GetPIDsignals(i),i);
774 pt->SetPIDTimBin(seed2.GetPIDTimBin(i),i);
777 Int_t * indexes3 = pt->GetBackupIndexes();
778 for (Int_t i=0;i<200;i++) {
779 if (indexes2[i]==0) break;
780 indexes3[i] = indexes2[i];
782 //AliTRDtrack *pt = seed2;
784 FollowProlongation(t);
785 if (t.GetNumberOfClusters() >= foundMin) {
787 //CookLabel(pt, 1-fgkLabelFraction);
792 // cout<<found<<'\r';
794 if(PropagateToX(t,xTPC,fgkMaxStep)) {
795 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
796 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
797 seed->SetTRDsignals(pt->GetPIDsignals(i),i);
798 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
801 //if not prolongation to TPC - propagate without update
802 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
803 seed2->ResetCovariance(5.);
804 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
806 if (PropagateToX(*pt2,xTPC,fgkMaxStep)) {
807 //pt2->CookdEdx(0.,1.);
808 pt2->CookdEdx( ); // Modification by PS
809 CookdEdxTimBin(*pt2);
810 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
811 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
812 seed->SetTRDsignals(pt2->GetPIDsignals(i),i);
813 seed->SetTRDTimBin(pt2->GetPIDTimBin(i),i);
821 cout<<"Number of loaded seeds: "<<nseed<<endl;
822 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
828 //_____________________________________________________________________________
829 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t)
832 // Starting from current position on track=t this function tries
833 // to extrapolate the track up to timeBin=0 and to confirm prolongation
834 // if a close cluster is found. Returns the number of clusters
835 // expected to be found in sensitive layers
836 // GeoManager used to estimate mean density
840 Int_t lastplane = GetLastPlane(&t);
841 Double_t radLength = 0.0;
843 Int_t expectedNumberOfClusters = 0;
847 for (Int_t iplane = lastplane; iplane>=0; iplane--){
849 Int_t row0 = GetGlobalTimeBin(0, iplane,GetTimeBinsPerPlane()-1);
850 Int_t rowlast = GetGlobalTimeBin(0, iplane,0);
852 // propagate track close to the plane if neccessary
854 Double_t currentx = fTrSec[0]->GetLayer(rowlast)->GetX();
855 if (currentx < -fgkMaxStep +t.GetX()){
856 //propagate closer to chamber - safety space fgkMaxStep
857 if (!PropagateToX(t, currentx+fgkMaxStep, fgkMaxStep)) break;
859 if (!AdjustSector(&t)) break;
861 // get material budget
863 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
864 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
865 // end global position
866 x = fTrSec[0]->GetLayer(row0)->GetX();
867 if (!t.GetProlongation(x,y,z)) break;
868 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
869 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
871 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
873 radLength = param[1]; // get mean propagation parameters
875 // propagate nad update
877 sector = t.GetSector();
878 // for (Int_t itime=GetTimeBinsPerPlane()-1;itime>=0;itime--) {
879 for (Int_t itime=0 ;itime<GetTimeBinsPerPlane();itime++) {
880 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
881 expectedNumberOfClusters++;
883 if (t.fX>345) t.fNExpectedLast++;
884 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(ilayer));
887 Double_t maxChi2=fgkMaxChi2;
890 AliTRDcluster * cl0 = timeBin[0];
891 if (!cl0) continue; // no clusters in given time bin
892 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
893 if (plane>lastplane) continue;
894 Int_t timebin = cl0->GetLocalTimeBin();
895 AliTRDcluster * cl2= GetCluster(&t,plane, timebin,index);
899 Double_t h01 = GetTiltFactor(cl);
900 maxChi2=t.GetPredictedChi2(cl,h01);
903 // if (cl->GetNPads()<5)
904 Double_t dxsample = timeBin.GetdX();
905 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
906 Double_t h01 = GetTiltFactor(cl);
907 Int_t det = cl->GetDetector();
908 Int_t plane = fGeom->GetPlane(det);
911 t.fChi2Last+=maxChi2;
913 Double_t xcluster = cl->GetX();
914 t.PropagateTo(xcluster,radLength,rho);
915 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
922 return expectedNumberOfClusters;
926 //_____________________________________________________________________________
927 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
930 // Starting from current radial position of track <t> this function
931 // extrapolates the track up to outer timebin and in the sensitive
932 // layers confirms prolongation if a close cluster is found.
933 // Returns the number of clusters expected to be found in sensitive layers
934 // Use GEO manager for material Description
938 Int_t clusters[1000];
939 for (Int_t i=0;i<1000;i++) clusters[i]=-1;
940 Double_t radLength = 0.0;
942 Int_t expectedNumberOfClusters = 0;
944 AliTRDtracklet tracklet;
947 for (Int_t iplane = 0; iplane<AliESDtrack::kNPlane; iplane++){
948 Int_t row0 = GetGlobalTimeBin(0, iplane,GetTimeBinsPerPlane()-1);
949 Int_t rowlast = GetGlobalTimeBin(0, iplane,0);
951 Double_t currentx = fTrSec[0]->GetLayer(row0)->GetX();
952 if (currentx<t.GetX()) continue;
954 // propagate closer to chamber if neccessary
956 if (currentx > fgkMaxStep +t.GetX()){
957 if (!PropagateToX(t, currentx-fgkMaxStep, fgkMaxStep)) break;
959 if (!AdjustSector(&t)) break;
960 if (TMath::Abs(t.GetSnp())>fgkMaxSnp) break;
962 // get material budget inside of chamber
964 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
965 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
966 // end global position
967 x = fTrSec[0]->GetLayer(rowlast)->GetX();
968 if (!t.GetProlongation(x,y,z)) break;
969 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
970 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
972 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
974 radLength = param[1]; // get mean propagation parameters
978 sector = t.GetSector();
979 Float_t ncl = FindClusters(sector,row0,rowlast,&t,clusters,tracklet);
980 if (tracklet.GetN()<GetTimeBinsPerPlane()/3) continue;
982 // Propagate and update track
984 for (Int_t itime= GetTimeBinsPerPlane()-1;itime>=0;itime--) {
985 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
986 expectedNumberOfClusters++;
988 if (t.fX>345) t.fNExpectedLast++;
989 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(ilayer));
992 Double_t maxChi2=fgkMaxChi2;
996 if (clusters[ilayer]>0) {
997 index = clusters[ilayer];
998 cl = (AliTRDcluster*)GetCluster(index);
999 Double_t h01 = GetTiltFactor(cl);
1000 maxChi2=t.GetPredictedChi2(cl,h01);
1004 // if (cl->GetNPads()<5)
1005 Double_t dxsample = timeBin.GetdX();
1006 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
1007 Double_t h01 = GetTiltFactor(cl);
1008 Int_t det = cl->GetDetector();
1009 Int_t plane = fGeom->GetPlane(det);
1012 t.fChi2Last+=maxChi2;
1014 Double_t xcluster = cl->GetX();
1015 t.PropagateTo(xcluster,radLength,rho);
1016 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1017 if(!t.Update(cl,maxChi2,index,h01)) {
1021 // reset material budget if 2 consecutive gold
1023 if (t.fTracklets[plane].GetN()+t.fTracklets[plane-1].GetN()>20){
1029 ratio0 = ncl/Float_t(fTimeBinsPerPlane);
1030 Float_t ratio1 = Float_t(t.fN+1)/Float_t(t.fNExpected+1.);
1031 if (tracklet.GetChi2()<18.&&ratio0>0.8 && ratio1>0.6 && ratio0+ratio1>1.5 && t.GetNCross()==0 && TMath::Abs(t.GetSnp())<0.85&&t.fN>20){
1032 t.MakeBackupTrack(); // make backup of the track until is gold
1037 return expectedNumberOfClusters;
1041 //_____________________________________________________________________________
1042 Int_t AliTRDtracker::PropagateToX(AliTRDtrack& t, Double_t xToGo, Double_t maxStep)
1045 // Starting from current radial position of track <t> this function
1046 // extrapolates the track up to radial position <xToGo>.
1047 // Returns 1 if track reaches the plane, and 0 otherwise
1050 const Double_t kEpsilon = 0.00001;
1051 // Double_t tanmax = TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
1052 Double_t xpos = t.GetX();
1053 Double_t dir = (xpos<xToGo) ? 1.:-1.;
1055 while ( (xToGo-xpos)*dir > kEpsilon){
1056 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1058 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
1059 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
1062 if (!t.GetProlongation(x,y,z)) return 0; // no prolongation
1064 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
1065 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
1068 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
1069 if (!t.PropagateTo(x,param[1],param[0])) return 0;
1078 //_____________________________________________________________________________
1079 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1082 // Fills clusters into TRD tracking_sectors
1083 // Note that the numbering scheme for the TRD tracking_sectors
1084 // differs from that of TRD sectors
1087 cout<<"\n Read Sectors clusters"<<endl;
1088 if (ReadClusters(fClusters,cTree)) {
1089 Error("LoadClusters","Problem with reading the clusters !");
1092 Int_t ncl=fClusters->GetEntriesFast();
1094 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1097 for (Int_t ichamber=0;ichamber<5;ichamber++)
1098 for (Int_t isector=0;isector<18;isector++){
1099 fHoles[ichamber][isector]=kTRUE;
1104 // printf("\r %d left ",ncl);
1105 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
1106 Int_t detector=c->GetDetector();
1107 Int_t localTimeBin=c->GetLocalTimeBin();
1108 Int_t sector=fGeom->GetSector(detector);
1109 Int_t plane=fGeom->GetPlane(detector);
1111 Int_t trackingSector = CookSectorIndex(sector);
1112 if (c->GetLabel(0)>0){
1113 Int_t chamber = fGeom->GetChamber(detector);
1114 fHoles[chamber][trackingSector]=kFALSE;
1117 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1118 if(gtb < 0) continue;
1119 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1123 // apply pos correction
1125 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1132 //_____________________________________________________________________________
1133 void AliTRDtracker::UnloadClusters()
1136 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1141 nentr = fClusters->GetEntriesFast();
1142 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
1145 nentr = fSeeds->GetEntriesFast();
1146 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1148 nentr = fTracks->GetEntriesFast();
1149 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1151 Int_t nsec = AliTRDgeometry::kNsect;
1153 for (i = 0; i < nsec; i++) {
1154 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1155 fTrSec[i]->GetLayer(pl)->Clear();
1161 //_____________________________________________________________________________
1162 void AliTRDtracker::MakeSeedsMI(Int_t /*inner*/, Int_t /*outer*/, AliESD * esd)
1165 // Creates seeds using clusters between position inner plane and outer plane
1168 const Double_t kMaxTheta = 1;
1169 const Double_t kMaxPhi = 2.0;
1171 const Double_t kRoad0y = 6; // road for middle cluster
1172 const Double_t kRoad0z = 8.5; // road for middle cluster
1174 const Double_t kRoad1y = 2; // road in y for seeded cluster
1175 const Double_t kRoad1z = 20; // road in z for seeded cluster
1177 const Double_t kRoad2y = 3; // road in y for extrapolated cluster
1178 const Double_t kRoad2z = 20; // road in z for extrapolated cluster
1179 const Int_t kMaxSeed = 3000;
1180 Int_t maxSec=AliTRDgeometry::kNsect;
1183 // linear fitters in planes
1184 TLinearFitter fitterTC(2,"hyp2"); // fitting with tilting pads - kz fixed - kz= Z/x, + vertex const
1185 TLinearFitter fitterT2(4,"hyp4"); // fitting with tilting pads - kz not fixed
1186 fitterTC.StoreData(kTRUE);
1187 fitterT2.StoreData(kTRUE);
1188 AliRieman rieman(1000); // rieman fitter
1189 AliRieman rieman2(1000); // rieman fitter
1191 // find the maximal and minimal layer for the planes
1194 AliTRDpropagationLayer* reflayers[6];
1195 for (Int_t i=0;i<6;i++){layers[i][0]=10000; layers[i][1]=0;}
1196 for (Int_t ns=0;ns<maxSec;ns++){
1197 for (Int_t ilayer=0;ilayer<fTrSec[ns]->GetNumberOfLayers();ilayer++){
1198 AliTRDpropagationLayer& layer=*(fTrSec[ns]->GetLayer(ilayer));
1199 if (layer==0) continue;
1200 Int_t det = layer[0]->GetDetector();
1201 Int_t plane = fGeom->GetPlane(det);
1202 if (ilayer<layers[plane][0]) layers[plane][0] = ilayer;
1203 if (ilayer>layers[plane][1]) layers[plane][1] = ilayer;
1207 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(0,0);
1208 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
1209 Double_t hL[6]; // tilting angle
1210 Double_t xcl[6]; // x - position of reference cluster
1211 Double_t ycl[6]; // y - position of reference cluster
1212 Double_t zcl[6]; // z - position of reference cluster
1213 AliTRDcluster *cl[6]={0,0,0,0,0,0}; // seeding clusters
1214 Float_t padlength[6]={10,10,10,10,10,10}; //current pad-length
1215 Double_t chi2R =0, chi2Z=0;
1216 Double_t chi2RF =0, chi2ZF=0;
1218 Int_t nclusters; // total number of clusters
1219 for (Int_t i=0;i<6;i++) {hL[i]=h01; if (i%2==1) hL[i]*=-1.;}
1223 AliTRDseed *pseed = new AliTRDseed[kMaxSeed*6];
1224 AliTRDseed *seed[kMaxSeed];
1225 for (Int_t iseed=0;iseed<kMaxSeed;iseed++) seed[iseed]= &pseed[iseed*6];
1226 AliTRDseed *cseed = seed[0];
1228 Double_t seedquality[kMaxSeed];
1229 Double_t seedquality2[kMaxSeed];
1230 Double_t seedparams[kMaxSeed][7];
1231 Int_t seedlayer[kMaxSeed];
1232 Int_t registered =0;
1233 Int_t sort[kMaxSeed];
1237 for (Int_t ns = 0; ns<maxSec; ns++){ //loop over sectors
1238 //for (Int_t ns = 0; ns<5; ns++){ //loop over sectors
1239 registered = 0; // reset registerd seed counter
1240 cseed = seed[registered];
1242 for (Int_t sLayer=2; sLayer>=0;sLayer--){
1243 //for (Int_t dseed=5;dseed<15; dseed+=3){ //loop over central seeding time bins
1245 Int_t dseed = 5+Int_t(iter)*3;
1246 // Initialize seeding layers
1247 for (Int_t ilayer=0;ilayer<6;ilayer++){
1248 reflayers[ilayer] = fTrSec[ns]->GetLayer(layers[ilayer][1]-dseed);
1249 xcl[ilayer] = reflayers[ilayer]->GetX();
1252 Double_t xref = (xcl[sLayer+1] + xcl[sLayer+2])*0.5;
1253 AliTRDpropagationLayer& layer0=*reflayers[sLayer+0];
1254 AliTRDpropagationLayer& layer1=*reflayers[sLayer+1];
1255 AliTRDpropagationLayer& layer2=*reflayers[sLayer+2];
1256 AliTRDpropagationLayer& layer3=*reflayers[sLayer+3];
1258 Int_t maxn3 = layer3;
1259 for (Int_t icl3=0;icl3<maxn3;icl3++){
1260 AliTRDcluster *cl3 = layer3[icl3];
1262 padlength[sLayer+3] = TMath::Sqrt(cl3->GetSigmaZ2()*12.);
1263 ycl[sLayer+3] = cl3->GetY();
1264 zcl[sLayer+3] = cl3->GetZ();
1265 Float_t yymin0 = ycl[sLayer+3] - 1- kMaxPhi *(xcl[sLayer+3]-xcl[sLayer+0]);
1266 Float_t yymax0 = ycl[sLayer+3] + 1+ kMaxPhi *(xcl[sLayer+3]-xcl[sLayer+0]);
1267 Int_t maxn0 = layer0; //
1268 for (Int_t icl0=layer0.Find(yymin0);icl0<maxn0;icl0++){
1269 AliTRDcluster *cl0 = layer0[icl0];
1271 if (cl3->IsUsed()&&cl0->IsUsed()) continue;
1272 ycl[sLayer+0] = cl0->GetY();
1273 zcl[sLayer+0] = cl0->GetZ();
1274 if ( ycl[sLayer+0]>yymax0) break;
1275 Double_t tanphi = (ycl[sLayer+3]-ycl[sLayer+0])/(xcl[sLayer+3]-xcl[sLayer+0]);
1276 if (TMath::Abs(tanphi)>kMaxPhi) continue;
1277 Double_t tantheta = (zcl[sLayer+3]-zcl[sLayer+0])/(xcl[sLayer+3]-xcl[sLayer+0]);
1278 if (TMath::Abs(tantheta)>kMaxTheta) continue;
1279 padlength[sLayer+0] = TMath::Sqrt(cl0->GetSigmaZ2()*12.);
1281 // expected position in 1 layer
1282 Double_t y1exp = ycl[sLayer+0]+(tanphi) *(xcl[sLayer+1]-xcl[sLayer+0]);
1283 Double_t z1exp = zcl[sLayer+0]+(tantheta)*(xcl[sLayer+1]-xcl[sLayer+0]);
1284 Float_t yymin1 = y1exp - kRoad0y-tanphi;
1285 Float_t yymax1 = y1exp + kRoad0y+tanphi;
1286 Int_t maxn1 = layer1; //
1288 for (Int_t icl1=layer1.Find(yymin1);icl1<maxn1;icl1++){
1289 AliTRDcluster *cl1 = layer1[icl1];
1292 if (cl3->IsUsed()) nusedCl++;
1293 if (cl0->IsUsed()) nusedCl++;
1294 if (cl1->IsUsed()) nusedCl++;
1295 if (nusedCl>1) continue;
1296 ycl[sLayer+1] = cl1->GetY();
1297 zcl[sLayer+1] = cl1->GetZ();
1298 if ( ycl[sLayer+1]>yymax1) break;
1299 if (TMath::Abs(ycl[sLayer+1]-y1exp)>kRoad0y+tanphi) continue;
1300 if (TMath::Abs(zcl[sLayer+1]-z1exp)>kRoad0z) continue;
1301 padlength[sLayer+1] = TMath::Sqrt(cl1->GetSigmaZ2()*12.);
1303 Double_t y2exp = ycl[sLayer+0]+(tanphi) *(xcl[sLayer+2]-xcl[sLayer+0])+(ycl[sLayer+1]-y1exp);
1304 Double_t z2exp = zcl[sLayer+0]+(tantheta)*(xcl[sLayer+2]-xcl[sLayer+0]);
1305 Int_t index2 = layer2.FindNearestCluster(y2exp,z2exp,kRoad1y, kRoad1z);
1306 if (index2<=0) continue;
1307 AliTRDcluster *cl2 = (AliTRDcluster*)GetCluster(index2);
1308 padlength[sLayer+2] = TMath::Sqrt(cl2->GetSigmaZ2()*12.);
1309 ycl[sLayer+2] = cl2->GetY();
1310 zcl[sLayer+2] = cl2->GetZ();
1311 if (TMath::Abs(cl2->GetZ()-z2exp)>kRoad0z) continue;
1314 rieman.AddPoint(xcl[sLayer+0],ycl[sLayer+0],zcl[sLayer+0],1,10);
1315 rieman.AddPoint(xcl[sLayer+1],ycl[sLayer+1],zcl[sLayer+1],1,10);
1316 rieman.AddPoint(xcl[sLayer+3],ycl[sLayer+3],zcl[sLayer+3],1,10);
1317 rieman.AddPoint(xcl[sLayer+2],ycl[sLayer+2],zcl[sLayer+2],1,10);
1321 for (Int_t iLayer=0;iLayer<6;iLayer++){
1322 cseed[iLayer].Reset();
1324 chi2Z =0.; chi2R=0.;
1325 for (Int_t iLayer=0;iLayer<4;iLayer++){
1326 cseed[sLayer+iLayer].fZref[0] = rieman.GetZat(xcl[sLayer+iLayer]);
1327 chi2Z += (cseed[sLayer+iLayer].fZref[0]- zcl[sLayer+iLayer])*
1328 (cseed[sLayer+iLayer].fZref[0]- zcl[sLayer+iLayer]);
1329 cseed[sLayer+iLayer].fZref[1] = rieman.GetDZat(xcl[sLayer+iLayer]);
1330 cseed[sLayer+iLayer].fYref[0] = rieman.GetYat(xcl[sLayer+iLayer]);
1331 chi2R += (cseed[sLayer+iLayer].fYref[0]- ycl[sLayer+iLayer])*
1332 (cseed[sLayer+iLayer].fYref[0]- ycl[sLayer+iLayer]);
1333 cseed[sLayer+iLayer].fYref[1] = rieman.GetDYat(xcl[sLayer+iLayer]);
1335 if (TMath::Sqrt(chi2R)>1./iter) continue;
1336 if (TMath::Sqrt(chi2Z)>7./iter) continue;
1340 Float_t minmax[2]={-100,100};
1341 for (Int_t iLayer=0;iLayer<4;iLayer++){
1342 Float_t max = zcl[sLayer+iLayer]+padlength[sLayer+iLayer]*0.5+1 -cseed[sLayer+iLayer].fZref[0];
1343 if (max<minmax[1]) minmax[1]=max;
1344 Float_t min = zcl[sLayer+iLayer]-padlength[sLayer+iLayer]*0.5-1 -cseed[sLayer+iLayer].fZref[0];
1345 if (min>minmax[0]) minmax[0]=min;
1347 Bool_t isFake = kFALSE;
1348 if (cl0->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
1349 if (cl1->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
1350 if (cl2->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
1351 if (AliTRDReconstructor::StreamLevel()>0){
1352 if ((!isFake) || (icl3%10)==0 ){ //debugging print
1353 TTreeSRedirector& cstream = *fDebugStreamer;
1361 "X0="<<xcl[sLayer+0]<<
1362 "X1="<<xcl[sLayer+1]<<
1363 "X2="<<xcl[sLayer+2]<<
1364 "X3="<<xcl[sLayer+3]<<
1369 "Seed0.="<<&cseed[sLayer+0]<<
1370 "Seed1.="<<&cseed[sLayer+1]<<
1371 "Seed2.="<<&cseed[sLayer+2]<<
1372 "Seed3.="<<&cseed[sLayer+3]<<
1373 "Zmin="<<minmax[0]<<
1374 "Zmax="<<minmax[1]<<
1379 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1380 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1381 //<<<<<<<<<<<<<<<<<< FIT SEEDING PART <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1382 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1388 for (Int_t jLayer=0;jLayer<4;jLayer++){
1389 cseed[sLayer+jLayer].fTilt = hL[sLayer+jLayer];
1390 cseed[sLayer+jLayer].fPadLength = padlength[sLayer+jLayer];
1391 cseed[sLayer+jLayer].fX0 = xcl[sLayer+jLayer];
1392 for (Int_t iter=0; iter<2; iter++){
1394 // in iteration 0 we try only one pad-row
1395 // if quality not sufficient we try 2 pad-rows - about 5% of tracks cross 2 pad-rows
1397 AliTRDseed tseed = cseed[sLayer+jLayer];
1398 Float_t roadz = padlength[sLayer+jLayer]*0.5;
1399 if (iter>0) roadz = padlength[sLayer+jLayer];
1401 Float_t quality =10000;
1402 for (Int_t iTime=2;iTime<20;iTime++){
1403 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[sLayer+jLayer][1]-iTime));
1404 Double_t dxlayer= layer.GetX()-xcl[sLayer+jLayer];
1405 Double_t zexp = cl[sLayer+jLayer]->GetZ() ;
1407 // try 2 pad-rows in second iteration
1408 zexp = tseed.fZref[0]+ tseed.fZref[1]*dxlayer;
1409 if (zexp>cl[sLayer+jLayer]->GetZ()) zexp = cl[sLayer+jLayer]->GetZ()+padlength[sLayer+jLayer]*0.5;
1410 if (zexp<cl[sLayer+jLayer]->GetZ()) zexp = cl[sLayer+jLayer]->GetZ()-padlength[sLayer+jLayer]*0.5;
1413 Double_t yexp = tseed.fYref[0]+
1414 tseed.fYref[1]*dxlayer;
1415 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y, roadz);
1416 if (index<=0) continue;
1417 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1419 tseed.fIndexes[iTime] = index;
1420 tseed.fClusters[iTime] = cl; // register cluster
1421 tseed.fX[iTime] = dxlayer; // register cluster
1422 tseed.fY[iTime] = cl->GetY(); // register cluster
1423 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1426 //count the number of clusters and distortions into quality
1427 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1428 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1429 TMath::Abs(tseed.fYfit[0]-tseed.fYref[0])/0.2+
1430 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1431 if (iter==0 && tseed.IsOK()) {
1432 cseed[sLayer+jLayer] = tseed;
1434 if (tquality<5) break;
1436 if (tseed.IsOK() && tquality<quality)
1437 cseed[sLayer+jLayer] = tseed;
1439 if (!cseed[sLayer+jLayer].IsOK()){
1443 cseed[sLayer+jLayer].CookLabels();
1444 cseed[sLayer+jLayer].UpdateUsed();
1445 nusedCl+= cseed[sLayer+jLayer].fNUsed;
1452 if (!isOK) continue;
1454 for (Int_t iLayer=0;iLayer<4;iLayer++){
1455 if (cseed[sLayer+iLayer].IsOK()){
1456 nclusters+=cseed[sLayer+iLayer].fN2;
1462 for (Int_t iLayer=0;iLayer<4;iLayer++){
1463 rieman.AddPoint(xcl[sLayer+iLayer],cseed[sLayer+iLayer].fYfitR[0],
1464 cseed[sLayer+iLayer].fZProb,1,10);
1470 for (Int_t iLayer=0;iLayer<4;iLayer++){
1471 cseed[sLayer+iLayer].fYref[0] = rieman.GetYat(xcl[sLayer+iLayer]);
1472 chi2R += (cseed[sLayer+iLayer].fYref[0]-cseed[sLayer+iLayer].fYfitR[0])*
1473 (cseed[sLayer+iLayer].fYref[0]-cseed[sLayer+iLayer].fYfitR[0]);
1474 cseed[sLayer+iLayer].fYref[1] = rieman.GetDYat(xcl[sLayer+iLayer]);
1475 cseed[sLayer+iLayer].fZref[0] = rieman.GetZat(xcl[sLayer+iLayer]);
1476 chi2Z += (cseed[sLayer+iLayer].fZref[0]- cseed[sLayer+iLayer].fMeanz)*
1477 (cseed[sLayer+iLayer].fZref[0]- cseed[sLayer+iLayer].fMeanz);
1478 cseed[sLayer+iLayer].fZref[1] = rieman.GetDZat(xcl[sLayer+iLayer]);
1480 Double_t curv = rieman.GetC();
1485 TMath::Abs(cseed[sLayer+0].fYfitR[1]- cseed[sLayer+0].fYref[1])+
1486 TMath::Abs(cseed[sLayer+1].fYfitR[1]- cseed[sLayer+1].fYref[1])+
1487 TMath::Abs(cseed[sLayer+2].fYfitR[1]- cseed[sLayer+2].fYref[1])+
1488 TMath::Abs(cseed[sLayer+3].fYfitR[1]- cseed[sLayer+3].fYref[1]);
1489 Double_t likea = TMath::Exp(-sumda*10.6);
1490 Double_t likechi2 = 0.0000000001;
1491 if (chi2R<0.5) likechi2+=TMath::Exp(-TMath::Sqrt(chi2R)*7.73);
1492 Double_t likechi2z = TMath::Exp(-chi2Z*0.088)/TMath::Exp(-chi2Z*0.019);
1493 Double_t likeN = TMath::Exp(-(72-nclusters)*0.19);
1494 Double_t like = likea*likechi2*likechi2z*likeN;
1496 Double_t likePrimY = TMath::Exp(-TMath::Abs(cseed[sLayer+0].fYref[1]-130*curv)*1.9);
1497 Double_t likePrimZ = TMath::Exp(-TMath::Abs(cseed[sLayer+0].fZref[1]-
1498 cseed[sLayer+0].fZref[0]/xcl[sLayer+0])*5.9);
1499 Double_t likePrim = TMath::Max(likePrimY*likePrimZ,0.0005);
1501 seedquality[registered] = like;
1502 seedlayer[registered] = sLayer;
1503 if (TMath::Log(0.000000000000001+like)<-15) continue;
1504 AliTRDseed seedb[6];
1505 for (Int_t iLayer=0;iLayer<6;iLayer++){
1506 seedb[iLayer] = cseed[iLayer];
1509 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1510 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1511 //<<<<<<<<<<<<<<< FULL TRACK FIT PART <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1512 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1518 // add new layers - avoid long extrapolation
1520 Int_t tLayer[2]={0,0};
1521 if (sLayer==2) {tLayer[0]=1; tLayer[1]=0;}
1522 if (sLayer==1) {tLayer[0]=5; tLayer[1]=0;}
1523 if (sLayer==0) {tLayer[0]=4; tLayer[1]=5;}
1525 for (Int_t iLayer=0;iLayer<2;iLayer++){
1526 Int_t jLayer = tLayer[iLayer]; // set tracking layer
1527 cseed[jLayer].Reset();
1528 cseed[jLayer].fTilt = hL[jLayer];
1529 cseed[jLayer].fPadLength = padlength[jLayer];
1530 cseed[jLayer].fX0 = xcl[jLayer];
1531 // get pad length and rough cluster
1532 Int_t indexdummy = reflayers[jLayer]->FindNearestCluster(cseed[jLayer].fYref[0],
1533 cseed[jLayer].fZref[0],kRoad2y,kRoad2z);
1534 if (indexdummy<=0) continue;
1535 AliTRDcluster *cldummy = (AliTRDcluster*)GetCluster(indexdummy);
1536 padlength[jLayer] = TMath::Sqrt(cldummy->GetSigmaZ2()*12.);
1538 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1540 for (Int_t iLayer=0;iLayer<2;iLayer++){
1541 Int_t jLayer = tLayer[iLayer]; // set tracking layer
1542 if ( (jLayer==0) && !(cseed[1].IsOK())) continue; // break not allowed
1543 if ( (jLayer==5) && !(cseed[4].IsOK())) continue; // break not allowed
1544 Float_t zexp = cseed[jLayer].fZref[0];
1545 Double_t zroad = padlength[jLayer]*0.5+1.;
1548 for (Int_t iter=0;iter<2;iter++){
1549 AliTRDseed tseed = cseed[jLayer];
1550 Float_t quality = 10000;
1551 for (Int_t iTime=2;iTime<20;iTime++){
1552 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[jLayer][1]-iTime));
1553 Double_t dxlayer = layer.GetX()-xcl[jLayer];
1554 Double_t yexp = tseed.fYref[0]+tseed.fYref[1]*dxlayer;
1555 Float_t yroad = kRoad1y;
1556 Int_t index = layer.FindNearestCluster(yexp,zexp, yroad, zroad);
1557 if (index<=0) continue;
1558 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1560 tseed.fIndexes[iTime] = index;
1561 tseed.fClusters[iTime] = cl; // register cluster
1562 tseed.fX[iTime] = dxlayer; // register cluster
1563 tseed.fY[iTime] = cl->GetY(); // register cluster
1564 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1568 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1569 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1570 TMath::Abs(tseed.fYfit[0]-tseed.fYref[0])/0.2+
1571 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1573 if (tquality<quality){
1574 cseed[jLayer]=tseed;
1580 if ( cseed[jLayer].IsOK()){
1581 cseed[jLayer].CookLabels();
1582 cseed[jLayer].UpdateUsed();
1583 nusedf+= cseed[jLayer].fNUsed;
1584 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1590 AliTRDseed bseed[6];
1591 for (Int_t jLayer=0;jLayer<6;jLayer++){
1592 bseed[jLayer] = cseed[jLayer];
1594 Float_t lastquality = 10000;
1595 Float_t lastchi2 = 10000;
1596 Float_t chi2 = 1000;
1599 for (Int_t iter =0; iter<4;iter++){
1601 // sort tracklets according "quality", try to "improve" 4 worst
1603 Float_t sumquality = 0;
1604 Float_t squality[6];
1605 Int_t sortindexes[6];
1606 for (Int_t jLayer=0;jLayer<6;jLayer++){
1607 if (bseed[jLayer].IsOK()){
1608 AliTRDseed &tseed = bseed[jLayer];
1609 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1610 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1611 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1612 TMath::Abs(tseed.fYfit[0]-(tseed.fYref[0]-zcor))/0.2+
1613 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1614 squality[jLayer] = tquality;
1616 else squality[jLayer]=-1;
1617 sumquality +=squality[jLayer];
1620 if (sumquality>=lastquality || chi2>lastchi2) break;
1621 lastquality = sumquality;
1624 for (Int_t jLayer=0;jLayer<6;jLayer++){
1625 cseed[jLayer] = bseed[jLayer];
1628 TMath::Sort(6,squality,sortindexes,kFALSE);
1631 for (Int_t jLayer=5;jLayer>1;jLayer--){
1632 Int_t bLayer = sortindexes[jLayer];
1633 AliTRDseed tseed = bseed[bLayer];
1634 for (Int_t iTime=2;iTime<20;iTime++){
1635 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[bLayer][1]-iTime));
1636 Double_t dxlayer= layer.GetX()-xcl[bLayer];
1638 Double_t zexp = tseed.fZref[0];
1639 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1641 Float_t roadz = padlength[bLayer]+1;
1642 if (TMath::Abs(tseed.fZProb-zexp)> padlength[bLayer]*0.5) {roadz = padlength[bLayer]*0.5;}
1643 if (tseed.fZfit[1]*tseed.fZref[1]<0) {roadz = padlength[bLayer]*0.5;}
1644 if (TMath::Abs(tseed.fZProb-zexp)<0.1*padlength[bLayer]) {
1645 zexp = tseed.fZProb;
1646 roadz = padlength[bLayer]*0.5;
1649 Double_t yexp = tseed.fYref[0]+
1650 tseed.fYref[1]*dxlayer-zcor;
1651 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y, roadz);
1652 if (index<=0) continue;
1653 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1655 tseed.fIndexes[iTime] = index;
1656 tseed.fClusters[iTime] = cl; // register cluster
1657 tseed.fX[iTime] = dxlayer; // register cluster
1658 tseed.fY[iTime] = cl->GetY(); // register cluster
1659 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1663 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1664 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1666 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1667 TMath::Abs(tseed.fYfit[0]-(tseed.fYref[0]-zcor))/0.2+
1668 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1670 if (tquality<squality[bLayer])
1671 bseed[bLayer] = tseed;
1674 chi2 = AliTRDseed::FitRiemanTilt(bseed, kTRUE);
1682 for (Int_t iLayer=0;iLayer<6;iLayer++) {
1683 if (TMath::Abs(cseed[iLayer].fYref[0]/cseed[iLayer].fX0)<0.15)
1685 if (cseed[iLayer].IsOK()){
1686 nclusters+=cseed[iLayer].fN2;
1690 if (nlayers<3) continue;
1692 for (Int_t iLayer=0;iLayer<6;iLayer++){
1693 if (cseed[iLayer].IsOK()) rieman.AddPoint(xcl[iLayer],cseed[iLayer].fYfitR[0],
1694 cseed[iLayer].fZProb,1,10);
1700 for (Int_t iLayer=0;iLayer<6;iLayer++){
1701 if (cseed[iLayer].IsOK()){
1702 cseed[iLayer].fYref[0] = rieman.GetYat(xcl[iLayer]);
1703 chi2RF += (cseed[iLayer].fYref[0]-cseed[iLayer].fYfitR[0])*
1704 (cseed[iLayer].fYref[0]-cseed[iLayer].fYfitR[0]);
1705 cseed[iLayer].fYref[1] = rieman.GetDYat(xcl[iLayer]);
1706 cseed[iLayer].fZref[0] = rieman.GetZat(xcl[iLayer]);
1707 chi2ZF += (cseed[iLayer].fZref[0]- cseed[iLayer].fMeanz)*
1708 (cseed[iLayer].fZref[0]- cseed[iLayer].fMeanz);
1709 cseed[iLayer].fZref[1] = rieman.GetDZat(xcl[iLayer]);
1712 chi2RF/=TMath::Max((nlayers-3.),1.);
1713 chi2ZF/=TMath::Max((nlayers-3.),1.);
1714 curv = rieman.GetC();
1718 Double_t xref2 = (xcl[2]+xcl[3])*0.5; // middle of the chamber
1719 Double_t dzmf = rieman.GetDZat(xref2);
1720 Double_t zmf = rieman.GetZat(xref2);
1725 fitterTC.ClearPoints();
1726 fitterT2.ClearPoints();
1728 for (Int_t iLayer=0; iLayer<6;iLayer++){
1729 if (!cseed[iLayer].IsOK()) continue;
1730 for (Int_t itime=0;itime<25;itime++){
1731 if (!cseed[iLayer].fUsable[itime]) continue;
1732 Double_t x = cseed[iLayer].fX[itime]+cseed[iLayer].fX0-xref2; // x relative to the midle chamber
1733 Double_t y = cseed[iLayer].fY[itime];
1734 Double_t z = cseed[iLayer].fZ[itime];
1735 // ExB correction to the correction
1739 Double_t x2 = cseed[iLayer].fX[itime]+cseed[iLayer].fX0; // global x
1741 Double_t t = 1./(x2*x2+y*y);
1743 uvt[0] = 2.*x2*uvt[1]; // u
1745 uvt[2] = 2.0*hL[iLayer]*uvt[1];
1746 uvt[3] = 2.0*hL[iLayer]*x*uvt[1];
1747 uvt[4] = 2.0*(y+hL[iLayer]*z)*uvt[1];
1749 Double_t error = 2*0.2*uvt[1];
1750 fitterT2.AddPoint(uvt,uvt[4],error);
1752 // constrained rieman
1754 z =cseed[iLayer].fZ[itime];
1755 uvt[0] = 2.*x2*t; // u
1756 uvt[1] = 2*hL[iLayer]*x2*uvt[1];
1757 uvt[2] = 2*(y+hL[iLayer]*(z-GetZ()))*t;
1758 fitterTC.AddPoint(uvt,uvt[2],error);
1760 rieman2.AddPoint(x2,y,z,1,10);
1767 Double_t rpolz0 = fitterT2.GetParameter(3);
1768 Double_t rpolz1 = fitterT2.GetParameter(4);
1770 // linear fitter - not possible to make boundaries
1771 // non accept non possible z and dzdx combination
1773 Bool_t acceptablez =kTRUE;
1774 for (Int_t iLayer=0; iLayer<6;iLayer++){
1775 if (cseed[iLayer].IsOK()){
1776 Double_t zT2 = rpolz0+rpolz1*(xcl[iLayer] - xref2);
1777 if (TMath::Abs(cseed[iLayer].fZProb-zT2)>padlength[iLayer]*0.5+1)
1778 acceptablez = kFALSE;
1782 fitterT2.FixParameter(3,zmf);
1783 fitterT2.FixParameter(4,dzmf);
1785 fitterT2.ReleaseParameter(3);
1786 fitterT2.ReleaseParameter(4);
1787 rpolz0 = fitterT2.GetParameter(3);
1788 rpolz1 = fitterT2.GetParameter(4);
1791 Double_t chi2TR = fitterT2.GetChisquare()/Float_t(npointsT);
1792 Double_t chi2TC = fitterTC.GetChisquare()/Float_t(npointsT);
1794 Double_t polz1c = fitterTC.GetParameter(2);
1795 Double_t polz0c = polz1c*xref2;
1797 Double_t aC = fitterTC.GetParameter(0);
1798 Double_t bC = fitterTC.GetParameter(1);
1799 Double_t cC = aC/TMath::Sqrt(bC*bC+1.); // curvature
1801 Double_t aR = fitterT2.GetParameter(0);
1802 Double_t bR = fitterT2.GetParameter(1);
1803 Double_t dR = fitterT2.GetParameter(2);
1804 Double_t cR = 1+bR*bR-dR*aR;
1807 dca = -dR/(TMath::Sqrt(1+bR*bR-dR*aR)+TMath::Sqrt(1+bR*bR));
1808 cR = aR/TMath::Sqrt(cR);
1811 Double_t chi2ZT2=0, chi2ZTC=0;
1812 for (Int_t iLayer=0; iLayer<6;iLayer++){
1813 if (cseed[iLayer].IsOK()){
1814 Double_t zT2 = rpolz0+rpolz1*(xcl[iLayer] - xref2);
1815 Double_t zTC = polz0c+polz1c*(xcl[iLayer] - xref2);
1816 chi2ZT2 += TMath::Abs(cseed[iLayer].fMeanz-zT2);
1817 chi2ZTC += TMath::Abs(cseed[iLayer].fMeanz-zTC);
1820 chi2ZT2/=TMath::Max((nlayers-3.),1.);
1821 chi2ZTC/=TMath::Max((nlayers-3.),1.);
1825 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1827 for (Int_t iLayer=0;iLayer<6;iLayer++){
1828 if (cseed[iLayer].IsOK())
1829 sumdaf += TMath::Abs((cseed[iLayer].fYfit[1]-cseed[iLayer].fYref[1])/cseed[iLayer].fSigmaY2);
1831 sumdaf /= Float_t (nlayers-2.);
1833 // likelihoods for full track
1835 Double_t likezf = TMath::Exp(-chi2ZF*0.14);
1836 Double_t likechi2C = TMath::Exp(-chi2TC*0.677);
1837 Double_t likechi2TR = TMath::Exp(-chi2TR*0.78);
1838 Double_t likeaf = TMath::Exp(-sumdaf*3.23);
1839 seedquality2[registered] = likezf*likechi2TR*likeaf;
1840 // Bool_t isGold = kFALSE;
1842 // if (nlayers == 6 && TMath::Log(0.000000001+seedquality2[index])<-5.) isGold =kTRUE; // gold
1843 // if (nlayers == findable && TMath::Log(0.000000001+seedquality2[index])<-4.) isGold =kTRUE; // gold
1844 // if (isGold &&nusedf<10){
1845 // for (Int_t jLayer=0;jLayer<6;jLayer++){
1846 // if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.1)
1847 // seed[index][jLayer].UseClusters(); //sign gold
1854 if (!cseed[0].IsOK()){
1856 if (!cseed[1].IsOK()) index0 = 2;
1858 seedparams[registered][0] = cseed[index0].fX0;
1859 seedparams[registered][1] = cseed[index0].fYref[0];
1860 seedparams[registered][2] = cseed[index0].fZref[0];
1861 seedparams[registered][5] = cR;
1862 seedparams[registered][3] = cseed[index0].fX0*cR - TMath::Sin(TMath::ATan(cseed[0].fYref[1]));
1863 seedparams[registered][4] = cseed[index0].fZref[1]/
1864 TMath::Sqrt(1+cseed[index0].fYref[1]*cseed[index0].fYref[1]);
1865 seedparams[registered][6] = ns;
1868 Int_t labels[12], outlab[24];
1870 for (Int_t iLayer=0;iLayer<6;iLayer++){
1871 if (!cseed[iLayer].IsOK()) continue;
1872 if (cseed[iLayer].fLabels[0]>=0) {
1873 labels[nlab] = cseed[iLayer].fLabels[0];
1876 if (cseed[iLayer].fLabels[1]>=0) {
1877 labels[nlab] = cseed[iLayer].fLabels[1];
1881 Freq(nlab,labels,outlab,kFALSE);
1882 Int_t label = outlab[0];
1883 Int_t frequency = outlab[1];
1884 for (Int_t iLayer=0;iLayer<6;iLayer++){
1885 cseed[iLayer].fFreq = frequency;
1886 cseed[iLayer].fC = cR;
1887 cseed[iLayer].fCC = cC;
1888 cseed[iLayer].fChi2 = chi2TR;
1889 cseed[iLayer].fChi2Z = chi2ZF;
1892 if (1||(!isFake)){ //debugging print
1893 Float_t zvertex = GetZ();
1894 TTreeSRedirector& cstream = *fDebugStreamer;
1895 if (AliTRDReconstructor::StreamLevel()>0)
1898 "Vertex="<<zvertex<<
1899 "Rieman2.="<<&rieman2<<
1900 "Rieman.="<<&rieman<<
1910 "Chi2RF="<<chi2RF<< //chi2 of trackletes on full track
1911 "Chi2ZF="<<chi2ZF<< //chi2 z on tracklets on full track
1912 "Chi2ZT2="<<chi2ZT2<< //chi2 z on tracklets on full track - rieman tilt
1913 "Chi2ZTC="<<chi2ZTC<< //chi2 z on tracklets on full track - rieman tilt const
1915 "Chi2TR="<<chi2TR<< //chi2 without vertex constrain
1916 "Chi2TC="<<chi2TC<< //chi2 with vertex constrain
1917 "C="<<curv<< // non constrained - no tilt correction
1918 "DR="<<dR<< // DR parameter - tilt correction
1919 "DCA="<<dca<< // DCA - tilt correction
1920 "CR="<<cR<< // non constrained curvature - tilt correction
1921 "CC="<<cC<< // constrained curvature
1927 "Nlayers="<<nlayers<<
1928 "NUsedS="<<nusedCl<<
1930 "Findable="<<findable<<
1932 "LikePrim="<<likePrim<<
1933 "Likechi2C="<<likechi2C<<
1934 "Likechi2TR="<<likechi2TR<<
1936 "LikeF="<<seedquality2[registered]<<
1943 "SB0.="<<&seedb[0]<<
1944 "SB1.="<<&seedb[1]<<
1945 "SB2.="<<&seedb[2]<<
1946 "SB3.="<<&seedb[3]<<
1947 "SB4.="<<&seedb[4]<<
1948 "SB5.="<<&seedb[5]<<
1950 "Freq="<<frequency<<
1954 if (registered<kMaxSeed-1) {
1956 cseed = seed[registered];
1958 }// end of loop over layer 1
1959 } // end of loop over layer 0
1960 } // end of loop over layer 3
1961 } // end of loop over seeding time bins
1965 TMath::Sort(registered,seedquality2,sort,kTRUE);
1966 Bool_t signedseed[kMaxSeed];
1967 for (Int_t i=0;i<registered;i++){
1968 signedseed[i]= kFALSE;
1970 for (Int_t iter=0; iter<5; iter++){
1971 for (Int_t iseed=0;iseed<registered;iseed++){
1972 Int_t index = sort[iseed];
1973 if (signedseed[index]) continue;
1974 Int_t labelsall[1000];
1977 Int_t sLayer = seedlayer[index];
1982 for (Int_t jLayer=0;jLayer<6;jLayer++){
1983 if (TMath::Abs(seed[index][jLayer].fYref[0]/xcl[jLayer])<0.15)
1985 if (seed[index][jLayer].IsOK()){
1986 seed[index][jLayer].UpdateUsed();
1987 ncl +=seed[index][jLayer].fN2;
1988 nused +=seed[index][jLayer].fNUsed;
1991 for (Int_t itime=0;itime<25;itime++){
1992 if (seed[index][jLayer].fUsable[itime]){
1994 for (Int_t ilab=0;ilab<3;ilab++){
1995 Int_t tindex = seed[index][jLayer].fClusters[itime]->GetLabel(ilab);
1997 labelsall[nlabelsall] = tindex;
2006 if (nused>30) continue;
2009 if (nlayers<6) continue;
2010 if (TMath::Log(0.000000001+seedquality2[index])<-5.) continue; // gold
2014 if (nlayers<findable) continue;
2015 if (TMath::Log(0.000000001+seedquality2[index])<-4.) continue; //
2020 if (nlayers==findable || nlayers==6) continue;
2021 if (TMath::Log(0.000000001+seedquality2[index])<-6.) continue;
2025 if (TMath::Log(0.000000001+seedquality2[index])<-5.) continue;
2029 if (TMath::Log(0.000000001+seedquality2[index])-nused/(nlayers-3.)<-15.) continue;
2032 signedseed[index] = kTRUE;
2034 Int_t labels[1000], outlab[1000];
2036 for (Int_t iLayer=0;iLayer<6;iLayer++){
2037 if (seed[index][iLayer].IsOK()){
2038 if (seed[index][iLayer].fLabels[0]>=0) {
2039 labels[nlab] = seed[index][iLayer].fLabels[0];
2042 if (seed[index][iLayer].fLabels[1]>=0) {
2043 labels[nlab] = seed[index][iLayer].fLabels[1];
2048 Freq(nlab,labels,outlab,kFALSE);
2049 Int_t label = outlab[0];
2050 Int_t frequency = outlab[1];
2051 Freq(nlabelsall,labelsall,outlab,kFALSE);
2052 Int_t label1 = outlab[0];
2053 Int_t label2 = outlab[2];
2054 Float_t fakeratio = (naccepted-outlab[1])/Float_t(naccepted);
2055 Float_t ratio = Float_t(nused)/Float_t(ncl);
2057 for (Int_t jLayer=0;jLayer<6;jLayer++){
2058 if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.2 )
2059 seed[index][jLayer].UseClusters(); //sign gold
2063 Int_t eventNr = esd->GetEventNumber();
2064 TTreeSRedirector& cstream = *fDebugStreamer;
2068 AliTRDtrack * track = RegisterSeed(seed[index],seedparams[index]);
2070 if (!track) track=&dummy;
2072 AliESDtrack esdtrack;
2073 esdtrack.UpdateTrackParams(track, AliESDtrack::kTRDout);
2074 esdtrack.SetLabel(label);
2075 esd->AddTrack(&esdtrack);
2076 TTreeSRedirector& cstream = *fDebugStreamer;
2077 if (AliTRDReconstructor::StreamLevel()>0)
2079 "EventNr="<<eventNr<<
2080 "ESD.="<<&esdtrack<<
2082 "trdback.="<<track<<
2085 if (AliTRDReconstructor::StreamLevel()>0)
2089 "Like="<<seedquality[index]<<
2090 "LikeF="<<seedquality2[index]<<
2091 "S0.="<<&seed[index][0]<<
2092 "S1.="<<&seed[index][1]<<
2093 "S2.="<<&seed[index][2]<<
2094 "S3.="<<&seed[index][3]<<
2095 "S4.="<<&seed[index][4]<<
2096 "S5.="<<&seed[index][5]<<
2100 "FakeRatio="<<fakeratio<<
2101 "Freq="<<frequency<<
2103 "Nlayers="<<nlayers<<
2104 "Findable="<<findable<<
2107 "EventNr="<<eventNr<<
2111 } // end of loop over sectors
2117 //_____________________________________________________________________________
2118 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
2121 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
2122 // from the file. The names of the cluster tree and branches
2123 // should match the ones used in AliTRDclusterizer::WriteClusters()
2126 Int_t nsize = Int_t(ClusterTree->GetTotBytes()/(sizeof(AliTRDcluster)));
2127 TObjArray *clusterArray = new TObjArray(nsize+1000);
2129 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
2131 Error("ReadClusters","Can't get the branch !");
2134 branch->SetAddress(&clusterArray);
2136 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
2137 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
2139 // Loop through all entries in the tree
2141 AliTRDcluster *c = 0;
2143 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
2146 nbytes += ClusterTree->GetEvent(iEntry);
2148 // Get the number of points in the detector
2149 Int_t nCluster = clusterArray->GetEntriesFast();
2150 // printf("\r Read %d clusters from entry %d", nCluster, iEntry);
2152 // Loop through all TRD digits
2153 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
2154 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
2155 AliTRDcluster *co = c;
2157 // delete clusterArray->RemoveAt(iCluster);
2158 clusterArray->RemoveAt(iCluster);
2161 // cout<<"Allocated"<<nsize<<"\tLoaded"<<array->GetEntriesFast()<<"\n";
2163 delete clusterArray;
2169 //_____________________________________________________________________________
2170 Bool_t AliTRDtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
2173 // Get track space point with index i
2174 // Origin: C.Cheshkov
2177 AliTRDcluster *cl = (AliTRDcluster*)fClusters->UncheckedAt(index);
2178 Int_t idet = cl->GetDetector();
2179 Int_t isector = fGeom->GetSector(idet);
2180 Int_t ichamber= fGeom->GetChamber(idet);
2181 Int_t iplan = fGeom->GetPlane(idet);
2183 local[0]=GetX(isector,iplan,cl->GetLocalTimeBin());
2184 local[1]=cl->GetY();
2185 local[2]=cl->GetZ();
2187 fGeom->RotateBack(idet,local,global);
2188 p.SetXYZ(global[0],global[1],global[2]);
2189 AliAlignObj::ELayerID iLayer = AliAlignObj::kTRD1;
2192 iLayer = AliAlignObj::kTRD1;
2195 iLayer = AliAlignObj::kTRD2;
2198 iLayer = AliAlignObj::kTRD3;
2201 iLayer = AliAlignObj::kTRD4;
2204 iLayer = AliAlignObj::kTRD5;
2207 iLayer = AliAlignObj::kTRD6;
2210 Int_t modId = isector*fGeom->Ncham()+ichamber;
2211 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,modId);
2212 p.SetVolumeID(volid);
2218 //_____________________________________________________________________________
2219 void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
2222 // This cooks a label. Mmmmh, smells good...
2225 Int_t label=123456789, index, i, j;
2226 Int_t ncl=pt->GetNumberOfClusters();
2227 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
2231 // Int_t s[kRange][2];
2232 Int_t **s = new Int_t* [kRange];
2233 for (i=0; i<kRange; i++) {
2234 s[i] = new Int_t[2];
2236 for (i=0; i<kRange; i++) {
2242 for (i=0; i<ncl; i++) {
2243 index=pt->GetClusterIndex(i);
2244 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
2250 for (i=0; i<ncl; i++) {
2251 index=pt->GetClusterIndex(i);
2252 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
2253 for (Int_t k=0; k<3; k++) {
2254 label=c->GetLabel(k);
2255 labelAdded=kFALSE; j=0;
2257 while ( (!labelAdded) && ( j < kRange ) ) {
2258 if (s[j][0]==label || s[j][1]==0) {
2272 for (i=0; i<kRange; i++) {
2274 max=s[i][1]; label=s[i][0];
2278 for (i=0; i<kRange; i++) {
2284 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
2286 pt->SetLabel(label);
2290 //_____________________________________________________________________________
2291 void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
2294 // Use clusters, but don't abuse them!
2297 const Float_t kmaxchi2 =18;
2298 const Float_t kmincl =10;
2299 AliTRDtrack * track = (AliTRDtrack*)t;
2301 Int_t ncl=t->GetNumberOfClusters();
2302 for (Int_t i=from; i<ncl; i++) {
2303 Int_t index = t->GetClusterIndex(i);
2304 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
2306 Int_t iplane = fGeom->GetPlane(c->GetDetector());
2307 if (track->fTracklets[iplane].GetChi2()>kmaxchi2) continue;
2308 if (track->fTracklets[iplane].GetN()<kmincl) continue;
2309 if (!(c->IsUsed())) c->Use();
2314 //_____________________________________________________________________________
2315 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
2318 // Parametrised "expected" error of the cluster reconstruction in Y
2321 Double_t s = 0.08 * 0.08;
2326 //_____________________________________________________________________________
2327 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
2330 // Parametrised "expected" error of the cluster reconstruction in Z
2333 Double_t s = 9 * 9 /12.;
2338 //_____________________________________________________________________________
2339 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
2342 // Returns radial position which corresponds to time bin <localTB>
2343 // in tracking sector <sector> and plane <plane>
2346 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
2347 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
2348 return fTrSec[sector]->GetLayer(pl)->GetX();
2352 //_____________________________________________________________________________
2353 AliTRDtracker::AliTRDpropagationLayer
2354 ::AliTRDpropagationLayer(Double_t x, Double_t dx, Double_t rho
2355 , Double_t radLength, Int_t tbIndex, Int_t plane)
2358 // AliTRDpropagationLayer constructor
2361 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
2362 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
2365 for(Int_t i=0; i < (Int_t) kZones; i++) {
2366 fZc[i]=0; fZmax[i] = 0;
2371 if(fTimeBinIndex >= 0) {
2372 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
2373 fIndex = new UInt_t[kMaxClusterPerTimeBin];
2376 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
2387 //_____________________________________________________________________________
2388 void AliTRDtracker::AliTRDpropagationLayer
2389 ::SetHole(Double_t Zmax, Double_t Ymax, Double_t rho
2390 , Double_t radLength, Double_t Yc, Double_t Zc)
2393 // Sets hole in the layer
2402 fHoleX0 = radLength;
2406 //_____________________________________________________________________________
2407 AliTRDtracker::AliTRDtrackingSector
2408 ::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs)
2411 // AliTRDtrackingSector Constructor
2414 AliTRDpadPlane *padPlane = 0;
2420 // get holes description from geometry
2421 Bool_t holes[AliTRDgeometry::kNcham];
2422 //printf("sector\t%d\t",gs);
2423 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
2424 holes[icham] = fGeom->IsHole(0,icham,gs);
2425 //printf("%d",holes[icham]);
2429 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
2432 AliTRDpropagationLayer* ppl;
2434 Double_t x, dx, rho, radLength;
2437 // add layers for each of the planes
2438 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2439 //Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2442 const Int_t kNchambers = AliTRDgeometry::Ncham();
2444 Double_t ymaxsensitive=0;
2445 Double_t *zc = new Double_t[kNchambers];
2446 Double_t *zmax = new Double_t[kNchambers];
2447 Double_t *zmaxsensitive = new Double_t[kNchambers];
2449 AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
2452 printf("<AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector> ");
2453 printf("Could not get common params\n");
2457 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2459 ymax = fGeom->GetChamberWidth(plane)/2.;
2460 // Modidified for new pad plane class, 22.04.05 (C.B.)
2461 padPlane = commonParam->GetPadPlane(plane,0);
2462 ymaxsensitive = (padPlane->GetColSize(1)*padPlane->GetNcols()-4)/2.;
2463 for(Int_t ch = 0; ch < kNchambers; ch++) {
2464 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
2466 // Modidified for new pad plane class, 22.04.05 (C.B.)
2467 Float_t pad = padPlane->GetRowSize(1);
2468 Float_t row0 = commonParam->GetRow0(plane,ch,0);
2469 Int_t nPads = commonParam->GetRowMax(plane,ch,0);
2470 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
2471 zc[ch] = -(pad * nPads)/2 + row0;
2474 dx = AliTRDcalibDB::Instance()->GetVdrift(0,0,0)
2475 / AliTRDcalibDB::Instance()->GetSamplingFrequency();
2476 rho = 0.00295 * 0.85; radLength = 11.0;
2478 Double_t x0 = (Double_t) AliTRDgeometry::GetTime0(plane);
2479 //Double_t xbottom = x0 - dxDrift;
2480 //Double_t xtop = x0 + dxAmp;
2482 Int_t nTimeBins = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
2483 for (Int_t iTime = 0; iTime<nTimeBins; iTime++){
2484 Double_t xlayer = iTime*dx - dxAmp;
2485 //if (xlayer<0) xlayer=dxAmp/2.;
2488 tbIndex = CookTimeBinIndex(plane, iTime);
2489 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex, plane);
2490 ppl->SetYmax(ymax,ymaxsensitive);
2491 ppl->SetZ(zc, zmax, zmaxsensitive);
2492 ppl->SetHoles(holes);
2500 delete [] zmaxsensitive;
2504 //_____________________________________________________________________________
2505 Int_t AliTRDtracker::AliTRDtrackingSector
2506 ::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2509 // depending on the digitization parameters calculates "global"
2510 // time bin index for timebin <localTB> in plane <plane>
2514 Int_t tbPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
2515 Int_t gtb = (plane+1) * tbPerPlane - localTB -1;
2516 if (localTB<0) return -1;
2517 if (gtb<0) return -1;
2523 //_____________________________________________________________________________
2524 void AliTRDtracker::AliTRDtrackingSector
2525 ::MapTimeBinLayers()
2528 // For all sensitive time bins sets corresponding layer index
2529 // in the array fTimeBins
2534 for(Int_t i = 0; i < fN; i++) {
2535 index = fLayers[i]->GetTimeBinIndex();
2537 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2539 if(index < 0) continue;
2540 if(index >= (Int_t) kMaxTimeBinIndex) {
2541 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2542 printf(" index %d exceeds allowed maximum of %d!\n",
2543 index, kMaxTimeBinIndex-1);
2546 fTimeBinIndex[index] = i;
2551 //_____________________________________________________________________________
2552 Int_t AliTRDtracker::AliTRDtrackingSector
2553 ::GetLayerNumber(Double_t x) const
2556 // Returns the number of time bin which in radial position is closest to <x>
2559 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2560 if(x <= fLayers[0]->GetX()) return 0;
2562 Int_t b=0, e=fN-1, m=(b+e)/2;
2563 for (; b<e; m=(b+e)/2) {
2564 if (x > fLayers[m]->GetX()) b=m+1;
2567 if(TMath::Abs(x - fLayers[m]->GetX()) >
2568 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2574 //_____________________________________________________________________________
2575 Int_t AliTRDtracker::AliTRDtrackingSector
2576 ::GetInnerTimeBin() const
2579 // Returns number of the innermost SENSITIVE propagation layer
2582 return GetLayerNumber(0);
2586 //_____________________________________________________________________________
2587 Int_t AliTRDtracker::AliTRDtrackingSector
2588 ::GetOuterTimeBin() const
2591 // Returns number of the outermost SENSITIVE time bin
2594 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2598 //_____________________________________________________________________________
2599 Int_t AliTRDtracker::AliTRDtrackingSector
2600 ::GetNumberOfTimeBins() const
2603 // Returns number of SENSITIVE time bins
2607 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
2608 layer = GetLayerNumber(tb);
2616 //_____________________________________________________________________________
2617 void AliTRDtracker::AliTRDtrackingSector
2618 ::InsertLayer(AliTRDpropagationLayer* pl)
2621 // Insert layer <pl> in fLayers array.
2622 // Layers are sorted according to X coordinate.
2625 if ( fN == ((Int_t) kMaxLayersPerSector)) {
2626 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2629 if (fN==0) {fLayers[fN++] = pl; return;}
2630 Int_t i=Find(pl->GetX());
2632 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2633 fLayers[i]=pl; fN++;
2637 //_____________________________________________________________________________
2638 Int_t AliTRDtracker::AliTRDtrackingSector
2639 ::Find(Double_t x) const
2642 // Returns index of the propagation layer nearest to X
2645 if (x <= fLayers[0]->GetX()) return 0;
2646 if (x > fLayers[fN-1]->GetX()) return fN;
2647 Int_t b=0, e=fN-1, m=(b+e)/2;
2648 for (; b<e; m=(b+e)/2) {
2649 if (x > fLayers[m]->GetX()) b=m+1;
2657 //_____________________________________________________________________________
2658 void AliTRDtracker::AliTRDpropagationLayer
2659 ::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2662 // set centers and the width of sectors
2665 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2666 fZc[icham] = center[icham];
2667 fZmax[icham] = w[icham];
2668 fZmaxSensitive[icham] = wsensitive[icham];
2669 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2674 //_____________________________________________________________________________
2675 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2678 // set centers and the width of sectors
2682 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2683 fIsHole[icham] = holes[icham];
2684 if (holes[icham]) fHole = kTRUE;
2689 //_____________________________________________________________________________
2690 void AliTRDtracker::AliTRDpropagationLayer
2691 ::InsertCluster(AliTRDcluster* c, UInt_t index)
2694 // Insert cluster in cluster array.
2695 // Clusters are sorted according to Y coordinate.
2698 if(fTimeBinIndex < 0) {
2699 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2703 if (fN== (Int_t) kMaxClusterPerTimeBin) {
2704 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2707 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2708 Int_t i=Find(c->GetY());
2709 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2710 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2711 fIndex[i]=index; fClusters[i]=c; fN++;
2715 //_____________________________________________________________________________
2716 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Float_t y) const
2719 // Returns index of the cluster nearest in Y
2722 if (fN<=0) return 0;
2723 if (y <= fClusters[0]->GetY()) return 0;
2724 if (y > fClusters[fN-1]->GetY()) return fN;
2725 Int_t b=0, e=fN-1, m=(b+e)/2;
2726 for (; b<e; m=(b+e)/2) {
2727 if (y > fClusters[m]->GetY()) b=m+1;
2735 //_____________________________________________________________________________
2736 Int_t AliTRDtracker::AliTRDpropagationLayer
2737 ::FindNearestCluster(Float_t y, Float_t z, Float_t maxroad
2738 , Float_t maxroadz) const
2741 // Returns index of the cluster nearest to the given y,z
2746 Float_t mindist = maxroad;
2748 for (Int_t i=Find(y-maxroad); i<maxn; i++) {
2749 AliTRDcluster* c=(AliTRDcluster*)(fClusters[i]);
2750 Float_t ycl = c->GetY();
2752 if (ycl > y+maxroad) break;
2753 if (TMath::Abs(c->GetZ()-z) > maxroadz) continue;
2754 if (TMath::Abs(ycl-y)<mindist){
2755 mindist = TMath::Abs(ycl-y);
2764 //_____________________________________________________________________________
2765 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c)
2768 // Returns correction factor for tilted pads geometry
2771 Int_t det = c->GetDetector();
2772 Int_t plane = fGeom->GetPlane(det);
2773 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(plane,0);
2774 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
2776 if(fNoTilt) h01 = 0;
2782 //_____________________________________________________________________________
2783 void AliTRDtracker::CookdEdxTimBin(AliTRDtrack& TRDtrack)
2786 // *** ADDED TO GET MORE INFORMATION FOR TRD PID ---- PS
2787 // This is setting fdEdxPlane and fTimBinPlane
2788 // Sums up the charge in each plane for track TRDtrack and also get the
2789 // Time bin for Max. Cluster
2790 // Prashant Shukla (shukla@physi.uni-heidelberg.de)
2793 Double_t clscharge[AliESDtrack::kNPlane], maxclscharge[AliESDtrack::kNPlane];
2794 Int_t nCluster[AliESDtrack::kNPlane], timebin[AliESDtrack::kNPlane];
2796 //Initialization of cluster charge per plane.
2797 for (Int_t iPlane = 0; iPlane < AliESDtrack::kNPlane; iPlane++) {
2798 clscharge[iPlane] = 0.0;
2799 nCluster[iPlane] = 0;
2800 timebin[iPlane] = -1;
2801 maxclscharge[iPlane] = 0.0;
2804 // Loop through all clusters associated to track TRDtrack
2805 Int_t nClus = TRDtrack.GetNumberOfClusters(); // from Kalmantrack
2806 for (Int_t iClus = 0; iClus < nClus; iClus++) {
2807 Double_t charge = TRDtrack.GetClusterdQdl(iClus);
2808 Int_t index = TRDtrack.GetClusterIndex(iClus);
2809 AliTRDcluster *pTRDcluster = (AliTRDcluster *) GetCluster(index);
2810 if (!pTRDcluster) continue;
2811 Int_t tb = pTRDcluster->GetLocalTimeBin();
2813 Int_t detector = pTRDcluster->GetDetector();
2814 Int_t iPlane = fGeom->GetPlane(detector);
2815 clscharge[iPlane] = clscharge[iPlane]+charge;
2816 if(charge > maxclscharge[iPlane]) {
2817 maxclscharge[iPlane] = charge;
2818 timebin[iPlane] = tb;
2821 } // end of loop over cluster
2823 // Setting the fdEdxPlane and fTimBinPlane variabales
2824 Double_t totalCharge = 0;
2825 for (Int_t iPlane = 0; iPlane < AliESDtrack::kNPlane; iPlane++) {
2826 // Quality control of TRD track.
2827 if (nCluster[iPlane]<= 5) {
2828 clscharge[iPlane]=0.0;
2831 if (nCluster[iPlane]) clscharge[iPlane] /= nCluster[iPlane];
2832 TRDtrack.SetPIDsignals(clscharge[iPlane], iPlane);
2833 TRDtrack.SetPIDTimBin(timebin[iPlane], iPlane);
2834 totalCharge= totalCharge+clscharge[iPlane];
2837 // Int_t nc=TRDtrack.GetNumberOfClusters();
2839 // for (i=0; i<nc; i++) dedx += TRDtrack.GetClusterdQdl(i);
2841 // for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2842 // TRDtrack.SetPIDsignals(dedx, iPlane);
2843 // TRDtrack.SetPIDTimBin(timbin[iPlane], iPlane);
2848 //_____________________________________________________________________________
2849 Int_t AliTRDtracker::FindClusters(Int_t sector, Int_t t0, Int_t t1
2850 , AliTRDtrack * track
2851 , Int_t *clusters,AliTRDtracklet&tracklet)
2855 // Try to find nearest clusters to the track in timebins from t0 to t1
2859 // correction coeficients - depends on TRD parameters - to be changed according it
2862 Double_t x[100],yt[100],zt[100];
2863 Double_t xmean=0; //reference x
2864 Double_t dz[10][100],dy[10][100];
2865 Float_t zmean[100], nmean[100];
2867 Int_t indexes[10][100]; // indexes of the clusters in the road
2868 AliTRDcluster *cl[10][100]; // pointers to the clusters in the road
2869 Int_t best[10][100]; // index of best matching cluster
2873 for (Int_t it=0;it<100; it++){
2881 for (Int_t ih=0;ih<10;ih++){
2882 indexes[ih][it]=-2; //reset indexes1
2890 Double_t x0 = track->GetX();
2891 Double_t sigmaz = TMath::Sqrt(TMath::Abs(track->GetSigmaZ2()));
2897 Float_t padlength=0;
2898 AliTRDtrack track2(*track);
2899 Float_t snpy = track->GetSnp();
2900 Float_t tany = TMath::Sqrt(snpy*snpy/(1.-snpy*snpy));
2901 if (snpy<0) tany*=-1;
2903 Double_t sy2=ExpectedSigmaY2(x0,track->GetTgl(),track->GetPt());
2904 Double_t sz2=ExpectedSigmaZ2(x0,track->GetTgl());
2905 Double_t road = 15.*sqrt(track->GetSigmaY2() + sy2);
2906 if (road>6.) road=6.;
2909 for (Int_t it=0;it<t1-t0;it++){
2910 Double_t maxChi2[2]={fgkMaxChi2,fgkMaxChi2};
2911 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(it+t0));
2912 if (timeBin==0) continue; // no indexes1
2913 Int_t maxn = timeBin;
2914 x[it] = timeBin.GetX();
2915 track2.PropagateTo(x[it]);
2916 yt[it] = track2.GetY();
2917 zt[it] = track2.GetZ();
2919 Double_t y=yt[it],z=zt[it];
2920 Double_t chi2 =1000000;
2923 // find 2 nearest cluster at given time bin
2926 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
2927 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
2928 h01 = GetTiltFactor(c);
2930 Int_t det = c->GetDetector();
2931 plane = fGeom->GetPlane(det);
2932 padlength = TMath::Sqrt(c->GetSigmaZ2()*12.);
2934 // if (c->GetLocalTimeBin()==0) continue;
2935 if (c->GetY() > y+road) break;
2936 if((c->GetZ()-z)*(c->GetZ()-z) > 12. * sz2) continue;
2938 Double_t dist = TMath::Abs(c->GetZ()-z);
2939 if (dist> (0.5*padlength+6.*sigmaz)) continue; // 6 sigma boundary cut
2942 if (dist> (0.5*padlength-sigmaz)){ // sigma boundary cost function
2943 cost = (dist-0.5*padlength)/(2.*sigmaz);
2944 if (cost>-1) cost= (cost+1.)*(cost+1.);
2947 // Int_t label = TMath::Abs(track->GetLabel());
2948 // if (c->GetLabel(0)!=label && c->GetLabel(1)!=label&&c->GetLabel(2)!=label) continue;
2949 chi2=track2.GetPredictedChi2(c,h01)+cost;
2952 if (chi2 > maxChi2[1]) continue;
2953 detector = c->GetDetector();
2955 for (Int_t ih=2;ih<9; ih++){ //store the clusters in the road
2958 indexes[ih][it] =timeBin.GetIndex(i); // index - 9 - reserved for outliers
2963 if (chi2 <maxChi2[0]){
2964 maxChi2[1] = maxChi2[0];
2966 indexes[1][it] = indexes[0][it];
2967 cl[1][it] = cl[0][it];
2968 indexes[0][it] = timeBin.GetIndex(i);
2974 indexes[1][it] =timeBin.GetIndex(i);
2982 if (nfound<4) return 0;
2983 xmean /=Float_t(nfound); // middle x
2984 track2.PropagateTo(xmean); // propagate track to the center
2986 // choose one of the variants
2992 Double_t sumdy2 = 0;
3002 Double_t moffset[10]; // mean offset
3003 Double_t mean[10]; // mean value
3004 Double_t angle[10]; // angle
3006 Double_t smoffset[10]; // sigma of mean offset
3007 Double_t smean[10]; // sigma of mean value
3008 Double_t sangle[10]; // sigma of angle
3009 Double_t smeanangle[10]; // correlation
3011 Double_t sigmas[10];
3012 Double_t tchi2s[10]; // chi2s for tracklet
3014 for (Int_t it=0;it<10;it++) {
3020 moffset[it] = 0; // mean offset
3021 mean[it] = 0; // mean value
3022 angle[it] = 0; // angle
3024 smoffset[it] = 1e10; // sigma of mean offset
3025 smean[it] = 1e10; // sigma of mean value
3026 sangle[it] = 1e10; // sigma of angle
3027 smeanangle[it] = 0; // correlation
3030 tchi2s[it] = 1e10; // chi2s for tracklet
3037 for (Int_t it=0;it<t1-t0;it++){
3038 if (!cl[0][it]) continue;
3039 for (Int_t dt=-3;dt<=3;dt++){
3040 if (it+dt<0) continue;
3041 if (it+dt>t1-t0) continue;
3042 if (!cl[0][it+dt]) continue;
3043 zmean[it]+=cl[0][it+dt]->GetZ();
3046 zmean[it]/=nmean[it];
3049 for (Int_t it=0; it<t1-t0;it++){
3051 for (Int_t ih=0;ih<10;ih++){
3054 if (!cl[ih][it]) continue;
3055 Double_t xcluster = cl[ih][it]->GetX();
3056 Double_t ytrack,ztrack;
3057 track2.GetProlongation(xcluster, ytrack, ztrack );
3058 dz[ih][it] = cl[ih][it]->GetZ()- ztrack; // calculate distance from track in z
3059 dy[ih][it] = cl[ih][it]->GetY()+ dz[ih][it]*h01 -ytrack; // in y
3062 if (!cl[0][it]) continue;
3063 if (TMath::Abs(cl[0][it]->GetZ()-zmean[it])> padlength*0.8 &&cl[1][it])
3064 if (TMath::Abs(cl[1][it]->GetZ()-zmean[it])< padlength*0.5){
3069 // iterative choosing of "best path"
3072 Int_t label = TMath::Abs(track->GetLabel());
3075 for (Int_t iter=0;iter<9;iter++){
3078 sumz = 0; sum=0; sumdy=0;sumdy2=0;sumx=0;sumx2=0;sumxy=0;mpads=0; ngood[iter]=0; nbad[iter]=0;
3080 for (Int_t it=0;it<t1-t0;it++){
3081 if (!cl[best[iter][it]][it]) continue;
3082 //calculates pad-row changes
3083 Double_t zbefore= cl[best[iter][it]][it]->GetZ();
3084 Double_t zafter = cl[best[iter][it]][it]->GetZ();
3085 for (Int_t itd = it-1; itd>=0;itd--) {
3086 if (cl[best[iter][itd]][itd]) {
3087 zbefore= cl[best[iter][itd]][itd]->GetZ();
3091 for (Int_t itd = it+1; itd<t1-t0;itd++) {
3092 if (cl[best[iter][itd]][itd]) {
3093 zafter= cl[best[iter][itd]][itd]->GetZ();
3097 if (TMath::Abs(cl[best[iter][it]][it]->GetZ()-zbefore)>0.1&&TMath::Abs(cl[best[iter][it]][it]->GetZ()-zafter)>0.1) changes[iter]++;
3099 Double_t dx = x[it]-xmean; // distance to reference x
3100 sumz += cl[best[iter][it]][it]->GetZ();
3102 sumdy += dy[best[iter][it]][it];
3103 sumdy2+= dy[best[iter][it]][it]*dy[best[iter][it]][it];
3106 sumxy += dx*dy[best[iter][it]][it];
3107 mpads += cl[best[iter][it]][it]->GetNPads();
3108 if (cl[best[iter][it]][it]->GetLabel(0)==label || cl[best[iter][it]][it]->GetLabel(1)==label||cl[best[iter][it]][it]->GetLabel(2)==label){
3116 // calculates line parameters
3118 Double_t det = sum*sumx2-sumx*sumx;
3119 angle[iter] = (sum*sumxy-sumx*sumdy)/det;
3120 mean[iter] = (sumx2*sumdy-sumx*sumxy)/det;
3121 meanz[iter] = sumz/sum;
3122 moffset[iter] = sumdy/sum;
3123 mpads /= sum; // mean number of pads
3126 Double_t sigma2 = 0; // normalized residuals - for line fit
3127 Double_t sigma1 = 0; // normalized residuals - constant fit
3129 for (Int_t it=0;it<t1-t0;it++){
3130 if (!cl[best[iter][it]][it]) continue;
3131 Double_t dx = x[it]-xmean;
3132 Double_t ytr = mean[iter]+angle[iter]*dx;
3133 sigma2 += (dy[best[iter][it]][it]-ytr)*(dy[best[iter][it]][it]-ytr);
3134 sigma1 += (dy[best[iter][it]][it]-moffset[iter])*(dy[best[iter][it]][it]-moffset[iter]);
3137 sigma2 /=(sum-2); // normalized residuals
3138 sigma1 /=(sum-1); // normalized residuals
3140 smean[iter] = sigma2*(sumx2/det); // estimated error2 of mean
3141 sangle[iter] = sigma2*(sum/det); // estimated error2 of angle
3142 smeanangle[iter] = sigma2*(-sumx/det); // correlation
3145 sigmas[iter] = TMath::Sqrt(sigma1); //
3146 smoffset[iter]= (sigma1/sum)+0.01*0.01; // sigma of mean offset + unisochronity sigma
3148 // iterative choosing of "better path"
3150 for (Int_t it=0;it<t1-t0;it++){
3151 if (!cl[best[iter][it]][it]) continue;
3153 Double_t sigmatr2 = smoffset[iter]+0.5*tany*tany; //add unisochronity + angular effect contribution
3154 Double_t sweight = 1./sigmatr2+1./track->GetSigmaY2();
3155 Double_t weighty = (moffset[iter]/sigmatr2)/sweight; // weighted mean
3156 Double_t sigmacl = TMath::Sqrt(sigma1*sigma1+track->GetSigmaY2()); //
3157 Double_t mindist=100000;
3159 for (Int_t ih=0;ih<10;ih++){
3160 if (!cl[ih][it]) break;
3161 Double_t dist2 = (dy[ih][it]-weighty)/sigmacl;
3162 dist2*=dist2; //chi2 distance
3168 best[iter+1][it]=ihbest;
3171 // update best hypothesy if better chi2 according tracklet position and angle
3173 Double_t sy2 = smean[iter] + track->GetSigmaY2();
3174 Double_t sa2 = sangle[iter] + track->fCee;
3175 Double_t say = track->fCey;
3176 // Double_t chi20 = mean[bestiter]*mean[bestiter]/sy2+angle[bestiter]*angle[bestiter]/sa2;
3177 // Double_t chi21 = mean[iter]*mean[iter]/sy2+angle[iter]*angle[iter]/sa2;
3179 Double_t detchi = sy2*sa2-say*say;
3180 Double_t invers[3] = {sa2/detchi, sy2/detchi, -say/detchi}; //inverse value of covariance matrix
3182 Double_t chi20 = mean[bestiter]*mean[bestiter]*invers[0]+angle[bestiter]*angle[bestiter]*invers[1]+
3183 2.*mean[bestiter]*angle[bestiter]*invers[2];
3184 Double_t chi21 = mean[iter]*mean[iter]*invers[0]+angle[iter]*angle[iter]*invers[1]+
3185 2*mean[iter]*angle[iter]*invers[2];
3186 tchi2s[iter] =chi21;
3188 if (changes[iter]<=changes[bestiter] && chi21<chi20) {
3195 Double_t sigma2 = sigmas[0]; // choose as sigma from 0 iteration
3196 Short_t maxpos = -1;
3197 Float_t maxcharge = 0;
3198 Short_t maxpos4 = -1;
3199 Float_t maxcharge4 = 0;
3200 Short_t maxpos5 = -1;
3201 Float_t maxcharge5 = 0;
3203 //if (tchi2s[bestiter]>25.) sigma2*=tchi2s[bestiter]/25.;
3204 //if (tchi2s[bestiter]>25.) sigma2=1000.; // dont'accept
3206 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(AliTRDcalibDB::Instance()->GetVdrift(0,0,0));
3207 Double_t expectederr = sigma2*sigma2+0.01*0.01;
3208 if (mpads>3.5) expectederr += (mpads-3.5)*0.04;
3209 if (changes[bestiter]>1) expectederr+= changes[bestiter]*0.01;
3210 expectederr+=(0.03*(tany-exB)*(tany-exB))*15;
3211 // if (tchi2s[bestiter]>18.) expectederr*= tchi2s[bestiter]/18.;
3212 //expectederr+=10000;
3213 for (Int_t it=0;it<t1-t0;it++){
3214 if (!cl[best[bestiter][it]][it]) continue;
3215 cl[best[bestiter][it]][it]->SetSigmaY2(expectederr); // set cluster error
3216 if (!cl[best[bestiter][it]][it]->IsUsed()){
3217 cl[best[bestiter][it]][it]->SetY( cl[best[bestiter][it]][it]->GetY());
3218 // cl[best[bestiter][it]][it]->Use();
3221 // time bins with maximal charge
3222 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge){
3223 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3224 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3227 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge4){
3228 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=4){
3229 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3230 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3233 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge5){
3234 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=5){
3235 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3236 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3240 // time bins with maximal charge
3241 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge){
3242 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3243 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3246 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge4){
3247 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=4){
3248 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3249 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3252 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge5){
3253 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=5){
3254 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3255 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3258 clusters[it+t0] = indexes[best[bestiter][it]][it];
3259 //if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>4 && cl[best[bestiter][it]][it]->GetLocalTimeBin()<18) clusters[it+t0] = indexes[best[bestiter][it]][it]; //Test
3262 // set tracklet parameters
3264 Double_t trackleterr2 = smoffset[bestiter]+0.01*0.01;
3265 if (mpads>3.5) trackleterr2 += (mpads-3.5)*0.04;
3266 trackleterr2+= changes[bestiter]*0.01;
3267 trackleterr2*= TMath::Max(14.-nfound,1.);
3268 trackleterr2+= 0.2*(tany-exB)*(tany-exB);
3270 tracklet.Set(xmean, track2.GetY()+moffset[bestiter], meanz[bestiter], track2.GetAlpha(), trackleterr2); //set tracklet parameters
3271 tracklet.SetTilt(h01);
3272 tracklet.SetP0(mean[bestiter]);
3273 tracklet.SetP1(angle[bestiter]);
3274 tracklet.SetN(nfound);
3275 tracklet.SetNCross(changes[bestiter]);
3276 tracklet.SetPlane(plane);
3277 tracklet.SetSigma2(expectederr);
3278 tracklet.SetChi2(tchi2s[bestiter]);
3279 tracklet.SetMaxPos(maxpos,maxpos4,maxpos5);
3280 track->fTracklets[plane] = tracklet;
3281 track->fNWrong+=nbad[0];
3285 TClonesArray array0("AliTRDcluster");
3286 TClonesArray array1("AliTRDcluster");
3287 array0.ExpandCreateFast(t1-t0+1);
3288 array1.ExpandCreateFast(t1-t0+1);
3289 TTreeSRedirector& cstream = *fDebugStreamer;
3290 AliTRDcluster dummy;
3294 for (Int_t it=0;it<t1-t0;it++){
3295 dy0[it] = dy[0][it];
3296 dyb[it] = dy[best[bestiter][it]][it];
3298 new(array0[it]) AliTRDcluster(*cl[0][it]);
3301 new(array0[it]) AliTRDcluster(dummy);
3303 if(cl[best[bestiter][it]][it]) {
3304 new(array1[it]) AliTRDcluster(*cl[best[bestiter][it]][it]);
3307 new(array1[it]) AliTRDcluster(dummy);
3310 TGraph graph0(t1-t0,x,dy0);
3311 TGraph graph1(t1-t0,x,dyb);
3312 TGraph graphy(t1-t0,x,yt);
3313 TGraph graphz(t1-t0,x,zt);
3316 if (AliTRDReconstructor::StreamLevel()>0)
3317 cstream<<"tracklet"<<
3318 "track.="<<track<< // track parameters
3319 "tany="<<tany<< // tangent of the local track angle
3320 "xmean="<<xmean<< // xmean - reference x of tracklet
3321 "tilt="<<h01<< // tilt angle
3322 "nall="<<nall<< // number of foundable clusters
3323 "nfound="<<nfound<< // number of found clusters
3324 "clfound="<<clfound<< // total number of found clusters in road
3325 "mpads="<<mpads<< // mean number of pads per cluster
3326 "plane="<<plane<< // plane number
3327 "detector="<<detector<< // detector number
3328 "road="<<road<< // the width of the used road
3329 "graph0.="<<&graph0<< // x - y = dy for closest cluster
3330 "graph1.="<<&graph1<< // x - y = dy for second closest cluster
3331 "graphy.="<<&graphy<< // y position of the track
3332 "graphz.="<<&graphz<< // z position of the track
3333 // "fCl.="<<&array0<< // closest cluster
3334 //"fCl2.="<<&array1<< // second closest cluster
3335 "maxpos="<<maxpos<< // maximal charge postion
3336 "maxcharge="<<maxcharge<< // maximal charge
3337 "maxpos4="<<maxpos4<< // maximal charge postion - after bin 4
3338 "maxcharge4="<<maxcharge4<< // maximal charge - after bin 4
3339 "maxpos5="<<maxpos5<< // maximal charge postion - after bin 5
3340 "maxcharge5="<<maxcharge5<< // maximal charge - after bin 5
3342 "bestiter="<<bestiter<< // best iteration number
3343 "tracklet.="<<&tracklet<< // corrspond to the best iteration
3344 "tchi20="<<tchi2s[0]<< // chi2 of cluster in the 0 iteration
3345 "tchi2b="<<tchi2s[bestiter]<< // chi2 of cluster in the best iteration
3346 "sigmas0="<<sigmas[0]<< // residuals sigma
3347 "sigmasb="<<sigmas[bestiter]<< // residulas sigma
3349 "ngood0="<<ngood[0]<< // number of good clusters in 0 iteration
3350 "nbad0="<<nbad[0]<< // number of bad clusters in 0 iteration
3351 "ngoodb="<<ngood[bestiter]<< // in best iteration
3352 "nbadb="<<nbad[bestiter]<< // in best iteration
3354 "changes0="<<changes[0]<< // changes of pardrows in iteration number 0
3355 "changesb="<<changes[bestiter]<< // changes of pardrows in best iteration
3357 "moffset0="<<moffset[0]<< // offset fixing angle in iter=0
3358 "smoffset0="<<smoffset[0]<< // sigma of offset fixing angle in iter=0
3359 "moffsetb="<<moffset[bestiter]<< // offset fixing angle in iter=best
3360 "smoffsetb="<<smoffset[bestiter]<< // sigma of offset fixing angle in iter=best
3362 "mean0="<<mean[0]<< // mean dy in iter=0;
3363 "smean0="<<smean[0]<< // sigma of mean dy in iter=0
3364 "meanb="<<mean[bestiter]<< // mean dy in iter=best
3365 "smeanb="<<smean[bestiter]<< // sigma of mean dy in iter=best
3367 "angle0="<<angle[0]<< // angle deviation in the iteration number 0
3368 "sangle0="<<sangle[0]<< // sigma of angular deviation in iteration number 0
3369 "angleb="<<angle[bestiter]<< // angle deviation in the best iteration
3370 "sangleb="<<sangle[bestiter]<< // sigma of angle deviation in the best iteration
3372 "expectederr="<<expectederr<< // expected error of cluster position
3380 //_____________________________________________________________________________
3381 Int_t AliTRDtracker::Freq(Int_t n, const Int_t *inlist
3382 , Int_t *outlist, Bool_t down)
3385 // Sort eleements according occurancy
3386 // The size of output array has is 2*n
3389 Int_t * sindexS = new Int_t[n]; // temp array for sorting
3390 Int_t * sindexF = new Int_t[2*n];
3391 for (Int_t i=0;i<n;i++) sindexF[i]=0;
3393 TMath::Sort(n,inlist, sindexS, down);
3394 Int_t last = inlist[sindexS[0]];
3397 sindexF[0+n] = last;
3401 for(Int_t i=1;i<n; i++){
3402 val = inlist[sindexS[i]];
3403 if (last == val) sindexF[countPos]++;
3406 sindexF[countPos+n] = val;
3407 sindexF[countPos]++;
3411 if (last==val) countPos++;
3412 // sort according frequency
3413 TMath::Sort(countPos, sindexF, sindexS, kTRUE);
3414 for (Int_t i=0;i<countPos;i++){
3415 outlist[2*i ] = sindexF[sindexS[i]+n];
3416 outlist[2*i+1] = sindexF[sindexS[i]];
3425 //_____________________________________________________________________________
3426 AliTRDtrack *AliTRDtracker::RegisterSeed(AliTRDseed * seeds, Double_t * params)
3431 Double_t alpha=AliTRDgeometry::GetAlpha();
3432 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
3435 c[1] = 0 ; c[2] = 2;
3436 c[3] = 0 ; c[4] = 0; c[5] = 0.02;
3437 c[6] = 0 ; c[7] = 0; c[8] = 0; c[9] = 0.1;
3438 c[10] = 0 ; c[11] = 0; c[12] = 0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
3441 AliTRDcluster *cl =0;
3442 for (Int_t ilayer=0;ilayer<6;ilayer++){
3443 if (seeds[ilayer].IsOK()){
3444 for (Int_t itime=22;itime>0;itime--){
3445 if (seeds[ilayer].fIndexes[itime]>0){
3446 index = seeds[ilayer].fIndexes[itime];
3447 cl = seeds[ilayer].fClusters[itime];
3454 if (cl==0) return 0;
3455 AliTRDtrack * track = new AliTRDtrack(cl,index,¶ms[1],c, params[0],params[6]*alpha+shift);
3456 track->PropagateTo(params[0]-5.);
3457 track->ResetCovariance(1);
3459 Int_t rc=FollowBackProlongation(*track);
3465 CookdEdxTimBin(*track);
3466 CookLabel(track, 0.9);