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 //
24 // M. Ivanov (Marian.Ivanov@cern.ch) //
25 // Y. Belikov (Jouri.Belikov@cern.ch) //
27 ///////////////////////////////////////////////////////////////////////////////
30 #include <Riostream.h>
36 #include <TLinearFitter.h>
37 #include <TObjArray.h>
40 #include <TTreeStream.h>
42 #include "AliESDEvent.h"
43 #include "AliESDtrack.h"
44 #include "AliAlignObj.h"
45 #include "AliRieman.h"
46 #include "AliTrackPointArray.h"
48 #include "AliTRDgeometry.h"
49 #include "AliTRDpadPlane.h"
50 #include "AliTRDgeometry.h"
51 #include "AliTRDcluster.h"
52 #include "AliTRDtrack.h"
53 #include "AliTRDseed.h"
54 #include "AliTRDcalibDB.h"
55 #include "AliTRDCommonParam.h"
56 #include "AliTRDtracker.h"
57 #include "AliTRDReconstructor.h"
58 #include "AliTRDCalibraFillHisto.h"
60 ClassImp(AliTRDtracker)
62 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5; //
63 const Float_t AliTRDtracker::fgkLabelFraction = 0.8; //
64 const Double_t AliTRDtracker::fgkMaxChi2 = 12.0; //
65 const Double_t AliTRDtracker::fgkMaxSnp = 0.95; // Corresponds to tan = 3
66 const Double_t AliTRDtracker::fgkMaxStep = 2.0; // Maximal step size in propagation
68 //_____________________________________________________________________________
69 AliTRDtracker::AliTRDtracker()
96 // Default constructor
99 for (Int_t i = 0; i < kTrackingSectors; i++) {
102 for (Int_t j = 0; j < 5; j++) {
103 for (Int_t k = 0; k < 18; k++) {
104 fHoles[j][k] = kFALSE;
112 //_____________________________________________________________________________
113 AliTRDtracker::AliTRDtracker(const AliTRDtracker &t)
134 ,fTimeBinsPerPlane(0)
135 ,fAddTRDseeds(kFALSE)
145 //_____________________________________________________________________________
146 AliTRDtracker::AliTRDtracker(const TFile */*geomfile*/)
162 ,fClusters(new TObjArray(2000))
164 ,fSeeds(new TObjArray(2000))
166 ,fTracks(new TObjArray(1000))
167 ,fTimeBinsPerPlane(0)
168 ,fAddTRDseeds(kFALSE)
176 TDirectory *savedir = gDirectory;
178 fGeom = new AliTRDgeometry();
179 fGeom->ReadGeoMatrices();
181 for (Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
183 fTrSec[trS] = new AliTRDtrackingSector(fGeom,geomS);
184 for (Int_t icham = 0; icham < AliTRDgeometry::kNcham; icham++) {
186 fHoles[icham][trS] = fGeom->IsHole(0,icham,geomS);
190 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(0,0);
191 Float_t tiltAngle = TMath::Abs(padPlane->GetTiltingAngle());
192 if (tiltAngle < 0.1) {
196 if (!AliTRDcalibDB::Instance()) {
197 AliFatal("Could not get calibration object");
199 fTimeBinsPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
201 fDebugStreamer = new TTreeSRedirector("TRDdebug.root");
209 //_____________________________________________________________________________
210 AliTRDtracker::~AliTRDtracker()
213 // Destructor of AliTRDtracker
235 for (Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
236 delete fTrSec[geomS];
239 if (fDebugStreamer) {
240 delete fDebugStreamer;
245 //_____________________________________________________________________________
246 Int_t AliTRDtracker::LocalToGlobalID(Int_t lid)
249 // Transform internal TRD ID to global detector ID
252 Int_t isector = fGeom->GetSector(lid);
253 Int_t ichamber = fGeom->GetChamber(lid);
254 Int_t iplan = fGeom->GetPlane(lid);
256 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
259 iLayer = AliGeomManager::kTRD1;
262 iLayer = AliGeomManager::kTRD2;
265 iLayer = AliGeomManager::kTRD3;
268 iLayer = AliGeomManager::kTRD4;
271 iLayer = AliGeomManager::kTRD5;
274 iLayer = AliGeomManager::kTRD6;
278 Int_t modId = isector * fGeom->Ncham() + ichamber;
279 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer,modId);
285 //_____________________________________________________________________________
286 Int_t AliTRDtracker::GlobalToLocalID(Int_t gid)
289 // Transform global detector ID to local detector ID
293 AliGeomManager::ELayerID layerId = AliGeomManager::VolUIDToLayer(gid,modId);
295 Int_t isector = modId / fGeom->Ncham();
296 Int_t ichamber = modId % fGeom->Ncham();
300 case AliGeomManager::kTRD1:
303 case AliGeomManager::kTRD2:
306 case AliGeomManager::kTRD3:
309 case AliGeomManager::kTRD4:
312 case AliGeomManager::kTRD5:
315 case AliGeomManager::kTRD6:
326 Int_t lid = fGeom->GetDetector(iLayer,ichamber,isector);
332 //_____________________________________________________________________________
333 Bool_t AliTRDtracker::Transform(AliTRDcluster *cluster)
336 // Transform from cluster system to tracking system
339 // Magic constants for geo manager transformation
340 const Double_t kX0shift = 2.52;
343 // 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.0); // Drift distance
356 Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
357 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift,-AliTracker::GetBz()*0.1);
359 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(plane,chamber);
360 Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
361 Double_t localPos[3];
362 Double_t localPosTracker[3];
363 localPos[0] = -cluster->GetX();
364 localPos[1] = cluster->GetY() - driftX * exB;
365 localPos[2] = cluster->GetZ() - zshiftIdeal;
367 cluster->SetY(cluster->GetY() - driftX*exB);
368 Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
369 cluster->SetX(xplane - cluster->GetX());
371 TGeoHMatrix *matrix = fGeom->GetCorrectionMatrix(cluster->GetDetector());
373 // No matrix found - if somebody used geometry with holes
374 AliError("Invalid Geometry - Default Geometry used\n");
377 matrix->LocalToMaster(localPos,localPosTracker);
379 if (AliTRDReconstructor::StreamLevel() > 1) {
380 (* fDebugStreamer) << "Transform"
382 << "matrix.=" << matrix
383 << "Detector=" << detector
384 << "Sector=" << sector
386 << "Chamber=" << chamber
387 << "lx0=" << localPosTracker[0]
388 << "ly0=" << localPosTracker[1]
389 << "lz0=" << localPosTracker[2]
393 cluster->SetX(localPosTracker[0]+kX0shift);
394 cluster->SetY(localPosTracker[1]);
395 cluster->SetZ(localPosTracker[2]);
401 //_____________________________________________________________________________
402 // Bool_t AliTRDtracker::Transform(AliTRDcluster *cluster)
405 // // Is this still needed ????
407 // const Double_t kDriftCorrection = 1.01; // drift coeficient correction
408 // const Double_t kTime0Cor = 0.32; // time0 correction
410 // const Double_t kX0shift = 2.52;
411 // const Double_t kX0shift5 = 3.05;
414 // // apply alignment and calibration to transform cluster
417 // Int_t detector = cluster->GetDetector();
418 // Int_t plane = fGeom->GetPlane(cluster->GetDetector());
419 // Int_t chamber = fGeom->GetChamber(cluster->GetDetector());
420 // Int_t sector = fGeom->GetSector(cluster->GetDetector());
422 // Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
423 // Double_t driftX = TMath::Max(cluster->GetX()-dxAmp*0.5,0.); // drift distance
427 // Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
428 // Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift,-AliTracker::GetBz()*0.1);
431 // AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
432 // AliTRDpadPlane * padPlane = commonParam->GetPadPlane(plane,chamber);
433 // Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
434 // Double_t localPos[3], globalPos[3], localPosTracker[3], localPosTracker2[3];
435 // localPos[2] = -cluster->GetX();
436 // localPos[0] = cluster->GetY() - driftX*exB;
437 // localPos[1] = cluster->GetZ() -zshiftIdeal;
438 // TGeoHMatrix * matrix = fGeom->GetGeoMatrix(cluster->GetDetector());
439 // matrix->LocalToMaster(localPos, globalPos);
441 // Double_t sectorAngle = 20.*(sector%18)+10;
442 // TGeoHMatrix rotSector;
443 // rotSector.RotateZ(sectorAngle);
444 // rotSector.LocalToMaster(globalPos, localPosTracker);
447 // TGeoHMatrix matrix2(*matrix);
448 // matrix2.MultiplyLeft(&rotSector);
449 // matrix2.LocalToMaster(localPos,localPosTracker2);
453 // cluster->SetY(cluster->GetY() - driftX*exB);
454 // Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
455 // cluster->SetX(xplane- kDriftCorrection*(cluster->GetX()-kTime0Cor));
456 // (*fDebugStreamer)<<"Transform"<<
458 // "matrix.="<<matrix<<
459 // "matrix2.="<<&matrix2<<
460 // "Detector="<<detector<<
461 // "Sector="<<sector<<
463 // "Chamber="<<chamber<<
464 // "lx0="<<localPosTracker[0]<<
465 // "ly0="<<localPosTracker[1]<<
466 // "lz0="<<localPosTracker[2]<<
467 // "lx2="<<localPosTracker2[0]<<
468 // "ly2="<<localPosTracker2[1]<<
469 // "lz2="<<localPosTracker2[2]<<
473 // cluster->SetX(localPosTracker[0]+kX0shift5);
475 // cluster->SetX(localPosTracker[0]+kX0shift);
477 // cluster->SetY(localPosTracker[1]);
478 // cluster->SetZ(localPosTracker[2]);
482 //_____________________________________________________________________________
483 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track)
486 // Rotates the track when necessary
489 Double_t alpha = AliTRDgeometry::GetAlpha();
490 Double_t y = track->GetY();
491 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
493 // Is this still needed ????
494 //Int_t ns = AliTRDgeometry::kNsect;
495 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
499 if (!track->Rotate( alpha)) {
503 else if (y < -ymax) {
505 if (!track->Rotate(-alpha)) {
514 //_____________________________________________________________________________
515 AliTRDcluster *AliTRDtracker::GetCluster(AliTRDtrack *track, Int_t plane
516 , Int_t timebin, UInt_t &index)
519 // Try to find cluster in the backup list
522 AliTRDcluster *cl =0;
523 Int_t *indexes = track->GetBackupIndexes();
525 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
526 if (indexes[i] == 0) {
529 AliTRDcluster *cli = (AliTRDcluster *) fClusters->UncheckedAt(indexes[i]);
533 if (cli->GetLocalTimeBin() != timebin) {
536 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
537 if (iplane == plane) {
548 //_____________________________________________________________________________
549 Int_t AliTRDtracker::GetLastPlane(AliTRDtrack *track)
552 // Return last updated plane
556 Int_t *indexes = track->GetBackupIndexes();
558 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
559 AliTRDcluster *cli = (AliTRDcluster *) fClusters->UncheckedAt(indexes[i]);
563 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
564 if (iplane > lastplane) {
573 //_____________________________________________________________________________
574 Int_t AliTRDtracker::Clusters2Tracks(AliESDEvent *event)
577 // Finds tracks within the TRD. The ESD event is expected to contain seeds
578 // at the outer part of the TRD. The seeds
579 // are found within the TRD if fAddTRDseeds is TRUE.
580 // The tracks are propagated to the innermost time bin
581 // of the TRD and the ESD event is updated
584 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
585 Float_t foundMin = fgkMinClustersInTrack * timeBins;
588 //Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
590 Int_t n = event->GetNumberOfTracks();
591 for (Int_t i = 0; i < n; i++) {
593 AliESDtrack *seed = event->GetTrack(i);
594 ULong_t status = seed->GetStatus();
595 if ((status & AliESDtrack::kTRDout) == 0) {
598 if ((status & AliESDtrack::kTRDin) != 0) {
603 AliTRDtrack *seed2 = new AliTRDtrack(*seed);
604 //seed2->ResetCovariance();
605 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
606 AliTRDtrack &t = *pt;
607 FollowProlongation(t);
608 if (t.GetNumberOfClusters() >= foundMin) {
610 CookLabel(pt,1 - fgkLabelFraction);
615 Double_t xTPC = 250.0;
616 if (PropagateToX(t,xTPC,fgkMaxStep)) {
617 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
624 AliInfo(Form("Number of loaded seeds: %d",nseed));
625 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
626 AliInfo(Form("Total number of found tracks: %d",found));
632 //_____________________________________________________________________________
633 Int_t AliTRDtracker::PropagateBack(AliESDEvent *event)
636 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
637 // backpropagated by the TPC tracker. Each seed is first propagated
638 // to the TRD, and then its prolongation is searched in the TRD.
639 // If sufficiently long continuation of the track is found in the TRD
640 // the track is updated, otherwise it's stored as originaly defined
641 // by the TPC tracker.
644 Int_t found = 0; // number of tracks found
645 Float_t foundMin = 20.0;
646 Int_t n = event->GetNumberOfTracks();
649 Float_t *quality = new Float_t[n];
650 Int_t *index = new Int_t[n];
651 for (Int_t i = 0; i < n; i++) {
652 AliESDtrack *seed = event->GetTrack(i);
653 Double_t covariance[15];
654 seed->GetExternalCovariance(covariance);
655 quality[i] = covariance[0]+covariance[2];
656 //quality[i] = covariance[0];
658 TMath::Sort(n,quality,index,kFALSE);
660 for (Int_t i = 0; i < n; i++) {
662 //AliESDtrack *seed = event->GetTrack(i);
663 AliESDtrack *seed = event->GetTrack(index[i]);
666 ULong_t status = seed->GetStatus();
667 if ((status & AliESDtrack::kTPCout) == 0) {
672 if ((status & AliESDtrack::kTRDout) != 0) {
677 Int_t lbl = seed->GetLabel();
678 AliTRDtrack *track = new AliTRDtrack(*seed);
679 track->SetSeedLabel(lbl);
680 seed->UpdateTrackParams(track,AliESDtrack::kTRDbackup); // Make backup
682 Float_t p4 = track->GetC();
683 Int_t expectedClr = FollowBackProlongation(*track);
686 fHX->Fill(track->GetX());
689 // store the last measurement
691 fHNClTrack->Fill(track->GetNumberOfClusters());
692 if (track->GetNumberOfClusters() >= foundMin) {
696 CookdEdxTimBin(*track);
697 CookLabel(track,1 - fgkLabelFraction);
698 if (track->GetBackupTrack()) {
699 //fHBackfit->Fill(5);
700 UseClusters(track->GetBackupTrack());
701 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
707 // inter-tracks competition ???
708 if ((TMath::Abs(track->GetC() - p4) / TMath::Abs(p4) < 0.2) ||
709 (TMath::Abs(track->GetPt()) > 0.8)) {
714 // Make backup for back propagation
717 Int_t foundClr = track->GetNumberOfClusters();
718 if (foundClr >= foundMin) {
720 track->CookdEdxTimBin(); // A.Bercuci 25.07.07
721 CookLabel(track,1 - fgkLabelFraction);
722 if (track->GetBackupTrack()) {
723 UseClusters(track->GetBackupTrack());
726 // Sign only gold tracks
727 if (track->GetChi2() / track->GetNumberOfClusters() < 4) {
728 if ((seed->GetKinkIndex(0) == 0) &&
729 (TMath::Abs(track->GetPt()) < 1.5)) {
733 Bool_t isGold = kFALSE;
736 if (track->GetChi2() / track->GetNumberOfClusters() < 5) {
737 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
738 if (track->GetBackupTrack()) {
739 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
748 (track->GetNCross() == 0) &&
749 (track->GetChi2() / track->GetNumberOfClusters() < 7)) {
750 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
751 if (track->GetBackupTrack()) {
752 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
758 (track->GetBackupTrack())) {
759 if ((track->GetBackupTrack()->GetNumberOfClusters() > foundMin) &&
760 ((track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
761 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
766 if ((track->StatusForTOF() > 0) &&
767 (track->GetNCross() == 0) &&
768 (Float_t(track->GetNumberOfClusters()) / Float_t(track->GetNExpected()) > 0.4)) {
769 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
777 // Debug part of tracking
778 TTreeSRedirector &cstream = *fDebugStreamer;
779 Int_t eventNrInFile = event->GetEventNumberInFile(); // This is most likely NOT the event number you'd like to use. It has nothing to do with the 'real' event number.
780 if (AliTRDReconstructor::StreamLevel() > 0) {
781 if (track->GetBackupTrack()) {
783 << "EventNrInFile=" << eventNrInFile
786 << "trdback.=" << track->GetBackupTrack()
791 << "EventNrInFile=" << eventNrInFile
794 << "trdback.=" << track
800 // Propagation to the TOF (I.Belikov)
801 if (track->GetStop() == kFALSE) {
804 Double_t xtof = 371.0;
805 Double_t xTOF0 = 370.0;
807 Double_t c2 = track->GetSnp() + track->GetC() * (xtof - track->GetX());
808 if (TMath::Abs(c2) >= 0.99) {
815 PropagateToX(*track,xTOF0,fgkMaxStep);
817 // Energy losses taken to the account - check one more time
818 c2 = track->GetSnp() + track->GetC() * (xtof - track->GetX());
819 if (TMath::Abs(c2) >= 0.99) {
826 //if (!PropagateToX(*track,xTOF0,fgkMaxStep)) {
827 // fHBackfit->Fill(7);
832 Double_t ymax = xtof * TMath::Tan(0.5 * AliTRDgeometry::GetAlpha());
834 track->GetYAt(xtof,GetBz(),y);
836 if (!track->Rotate( AliTRDgeometry::GetAlpha())) {
842 else if (y < -ymax) {
843 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
850 if (track->PropagateTo(xtof)) {
851 seed->UpdateTrackParams(track,AliESDtrack::kTRDout);
854 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
855 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
856 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
858 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
860 //seed->SetTRDtrack(new AliTRDtrack(*track));
861 if (track->GetNumberOfClusters() > foundMin) {
872 if ((track->GetNumberOfClusters() > 15) &&
873 (track->GetNumberOfClusters() > 0.5*expectedClr)) {
875 seed->UpdateTrackParams(track,AliESDtrack::kTRDout);
878 //seed->SetStatus(AliESDtrack::kTRDStop);
879 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
880 for (Int_t j = 0; j <AliESDtrack::kNSlice; j++) {
881 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
883 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
885 //seed->SetTRDtrack(new AliTRDtrack(*track));
891 seed->SetTRDQuality(track->StatusForTOF());
892 seed->SetTRDBudget(track->GetBudget(0));
898 AliInfo(Form("Number of seeds: %d",fNseeds));
899 AliInfo(Form("Number of back propagated TRD tracks: %d",found));
902 if (AliTRDReconstructor::SeedingOn()) {
903 MakeSeedsMI(3,5,event);
917 //_____________________________________________________________________________
918 Int_t AliTRDtracker::RefitInward(AliESDEvent *event)
921 // Refits tracks within the TRD. The ESD event is expected to contain seeds
922 // at the outer part of the TRD.
923 // The tracks are propagated to the innermost time bin
924 // of the TRD and the ESD event is updated
925 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
928 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
929 Float_t foundMin = fgkMinClustersInTrack * timeBins;
932 //Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
935 Int_t n = event->GetNumberOfTracks();
936 for (Int_t i = 0; i < n; i++) {
938 AliESDtrack *seed = event->GetTrack(i);
939 new (&seed2) AliTRDtrack(*seed);
942 if (seed2.GetX() < 270.0) {
943 seed->UpdateTrackParams(&seed2,AliESDtrack::kTRDbackup); // Backup TPC track - only update
948 ULong_t status = seed->GetStatus();
949 if ((status & AliESDtrack::kTRDout) == 0) {
953 if ((status & AliESDtrack::kTRDin) != 0) {
961 seed2.ResetCovariance(50.0);
963 AliTRDtrack *pt = new AliTRDtrack(seed2,seed2.GetAlpha());
964 Int_t *indexes2 = seed2.GetIndexes();
965 for (Int_t i = 0; i < AliESDtrack::kNPlane;i++) {
966 for (Int_t j = 0; j < AliESDtrack::kNSlice;j++) {
967 pt->SetPIDsignals(seed2.GetPIDsignals(i,j),i,j);
969 pt->SetPIDTimBin(seed2.GetPIDTimBin(i),i);
972 Int_t *indexes3 = pt->GetBackupIndexes();
973 for (Int_t i = 0; i < 200;i++) {
974 if (indexes2[i] == 0) {
977 indexes3[i] = indexes2[i];
980 //AliTRDtrack *pt = seed2;
981 //AliTRDtrack &t = *pt;
982 FollowProlongation(*pt);
983 if (pt->GetNumberOfClusters() >= foundMin) {
985 //CookLabel(pt, 1-fgkLabelFraction);
987 pt->CookdEdxTimBin();
989 //pt->Calibrate(); // slot for calibration
993 Double_t xTPC = 250.0;
994 if (PropagateToX(*pt,xTPC,fgkMaxStep)) {
996 seed->UpdateTrackParams(pt,AliESDtrack::kTRDrefit);
999 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
1000 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
1001 seed->SetTRDsignals(pt->GetPIDsignals(i,j),i,j);
1003 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
1007 // If not prolongation to TPC - propagate without update
1009 AliTRDtrack *seed2 = new AliTRDtrack(*seed);
1010 seed2->ResetCovariance(5.0);
1011 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
1013 if (PropagateToX(*pt2,xTPC,fgkMaxStep)) {
1015 pt2->CookdEdxTimBin();
1016 seed->UpdateTrackParams(pt2,AliESDtrack::kTRDrefit);
1019 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
1020 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
1021 seed->SetTRDsignals(pt2->GetPIDsignals(i,j),i,j);
1023 seed->SetTRDTimBin(pt2->GetPIDTimBin(i),i);
1027 // Add TRD track to ESDfriendTrack - maybe this tracks are not useful for post-processing - TODO make decission
1028 if (AliTRDReconstructor::StreamLevel() > 0) {
1029 seed->AddCalibObject(new AliTRDtrack(*pt2/*, kTRUE*/));
1034 // Add TRD track to ESDfriendTrack
1035 if (AliTRDReconstructor::StreamLevel() > 0) {
1036 seed->AddCalibObject(new AliTRDtrack(*pt/*, kTRUE*/));
1042 AliInfo(Form("Number of loaded seeds: %d",nseed));
1043 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
1049 //_____________________________________________________________________________
1050 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack &t)
1053 // Starting from current position on track=t this function tries
1054 // to extrapolate the track up to timeBin=0 and to confirm prolongation
1055 // if a close cluster is found. Returns the number of clusters
1056 // expected to be found in sensitive layers
1057 // GeoManager used to estimate mean density
1061 Int_t lastplane = GetLastPlane(&t);
1062 Double_t radLength = 0.0;
1064 Int_t expectedNumberOfClusters = 0;
1066 for (Int_t iplane = lastplane; iplane >= 0; iplane--) {
1068 Int_t row0 = GetGlobalTimeBin(0,iplane,GetTimeBinsPerPlane()-1);
1069 Int_t rowlast = GetGlobalTimeBin(0,iplane,0);
1072 // Propagate track close to the plane if neccessary
1074 Double_t currentx = fTrSec[0]->GetLayer(rowlast)->GetX();
1075 if (currentx < (-fgkMaxStep + t.GetX())) {
1076 // Propagate closer to chamber - safety space fgkMaxStep
1077 if (!PropagateToX(t,currentx+fgkMaxStep,fgkMaxStep)) {
1082 if (!AdjustSector(&t)) {
1087 // Get material budget
1096 // Starting global position
1098 // End global position
1099 x = fTrSec[0]->GetLayer(row0)->GetX();
1100 if (!t.GetProlongation(x,y,z)) {
1103 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1104 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1106 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
1108 radLength = param[1]; // Get mean propagation parameters
1109 // A.Bercuci 25.07.07
1110 // Flags for marking the track momentum and direction. The track is
1111 // marked only if it has at least 1 cluster picked up in the current
1113 Bool_t kUPDATED = kFALSE, kMARKED = kFALSE;
1116 // Propagate and update
1118 sector = t.GetSector();
1119 //for (Int_t itime=GetTimeBinsPerPlane()-1;itime>=0;itime--) {
1120 for (Int_t itime = 0 ; itime < GetTimeBinsPerPlane(); itime++) {
1121 // A.Bercuci 25.07.07
1122 // Mark track kinematics
1123 if(itime > 10 && kUPDATED && !kMARKED){
1124 t.SetTrackSegmentDirMom(iplane);
1128 Int_t ilayer = GetGlobalTimeBin(0,iplane,itime);
1129 expectedNumberOfClusters++;
1130 t.SetNExpected(t.GetNExpected() + 1);
1131 if (t.GetX() > 345.0) {
1132 t.SetNExpectedLast(t.GetNExpectedLast() + 1);
1134 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(ilayer));
1135 AliTRDcluster *cl = 0;
1137 Double_t maxChi2 = fgkMaxChi2;
1142 AliTRDcluster *cl0 = timeBin[0];
1144 // No clusters in given time bin
1148 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
1149 if (plane > lastplane) {
1153 Int_t timebin = cl0->GetLocalTimeBin();
1154 AliTRDcluster *cl2 = GetCluster(&t,plane,timebin,index);
1158 //Double_t h01 = GetTiltFactor(cl); //I.B's fix
1159 //maxChi2=t.GetPredictedChi2(cl,h01);
1162 //if (cl->GetNPads()<5)
1163 Double_t dxsample = timeBin.GetdX();
1164 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
1165 Double_t h01 = GetTiltFactor(cl);
1166 Int_t det = cl->GetDetector();
1167 Int_t plane = fGeom->GetPlane(det);
1168 if (t.GetX() > 345.0) {
1169 t.SetNLast(t.GetNLast() + 1);
1170 t.SetChi2Last(t.GetChi2Last() + maxChi2);
1173 Double_t xcluster = cl->GetX();
1174 t.PropagateTo(xcluster,radLength,rho);
1176 if (!AdjustSector(&t)) {
1179 maxChi2 = t.GetPredictedChi2(cl,h01);
1182 if (!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1185 // A.Bercuci 25.07.07
1186 //SetCluster(cl, GetNumberOfClusters()-1);
1194 return expectedNumberOfClusters;
1197 //_____________________________________________________________________________
1198 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack &t)
1201 // Starting from current radial position of track <t> this function
1202 // extrapolates the track up to outer timebin and in the sensitive
1203 // layers confirms prolongation if a close cluster is found.
1204 // Returns the number of clusters expected to be found in sensitive layers
1205 // Use GEO manager for material Description
1207 // return number of assigned clusters ?
1212 Int_t clusters[1000];
1213 Double_t radLength = 0.0;
1215 Int_t expectedNumberOfClusters = 0;
1216 Float_t ratio0 = 0.0;
1217 AliTRDtracklet tracklet;
1219 // Calibration fill 2D
1220 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
1222 AliInfo("Could not get Calibra instance\n");
1224 if (calibra->GetMITracking()) {
1225 calibra->ResetTrack();
1228 for (Int_t i = 0; i < 1000; i++) {
1232 for (Int_t iplane = 0; iplane < AliESDtrack::kNPlane; iplane++) {
1234 int hb = iplane * 10;
1235 fHClSearch->Fill(hb);
1237 Int_t row0 = GetGlobalTimeBin(0,iplane,GetTimeBinsPerPlane()-1);
1238 Int_t rowlast = GetGlobalTimeBin(0,iplane,0);
1239 Double_t currentx = fTrSec[0]->GetLayer(row0)->GetX();
1240 if (currentx < t.GetX()) {
1241 fHClSearch->Fill(hb+1);
1246 // Propagate closer to chamber if neccessary
1248 if (currentx > (fgkMaxStep + t.GetX())) {
1249 if (!PropagateToX(t,currentx-fgkMaxStep,fgkMaxStep)) {
1250 fHClSearch->Fill(hb+2);
1254 if (!AdjustSector(&t)) {
1255 fHClSearch->Fill(hb+3);
1258 if (TMath::Abs(t.GetSnp()) > fgkMaxSnp) {
1259 fHClSearch->Fill(hb+4);
1264 // Get material budget inside of chamber
1272 // Starting global position
1274 // End global position
1275 x = fTrSec[0]->GetLayer(rowlast)->GetX();
1276 if (!t.GetProlongation(x,y,z)) {
1277 fHClSearch->Fill(hb+5);
1280 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1281 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1283 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
1285 radLength = param[1]; // Get mean propagation parameters
1290 sector = t.GetSector();
1291 Float_t ncl = FindClusters(sector,row0,rowlast,&t,clusters,tracklet);
1292 fHNCl->Fill(tracklet.GetN());
1294 if (tracklet.GetN() < GetTimeBinsPerPlane()/3) {
1295 fHClSearch->Fill(hb+6);
1300 // Propagate and update track
1302 for (Int_t itime = GetTimeBinsPerPlane()-1; itime >= 0; itime--) {
1304 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
1305 expectedNumberOfClusters++;
1306 t.SetNExpected(t.GetNExpected() + 1);
1307 if (t.GetX() > 345.0) {
1308 t.SetNExpectedLast(t.GetNExpectedLast() + 1);
1310 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(ilayer));
1311 AliTRDcluster *cl = 0;
1313 Double_t maxChi2 = fgkMaxChi2;
1318 if (clusters[ilayer] > 0) {
1319 index = clusters[ilayer];
1320 cl = (AliTRDcluster *)GetCluster(index);
1321 //Double_t h01 = GetTiltFactor(cl); // I.B's fix
1322 //maxChi2=t.GetPredictedChi2(cl,h01); //
1327 //if (cl->GetNPads() < 5)
1328 Double_t dxsample = timeBin.GetdX();
1329 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
1330 Double_t h01 = GetTiltFactor(cl);
1331 Int_t det = cl->GetDetector();
1332 Int_t plane = fGeom->GetPlane(det);
1333 if (t.GetX() > 345.0) {
1334 t.SetNLast(t.GetNLast() + 1);
1335 t.SetChi2Last(t.GetChi2Last() + maxChi2);
1337 Double_t xcluster = cl->GetX();
1338 t.PropagateTo(xcluster,radLength,rho);
1339 maxChi2 = t.GetPredictedChi2(cl,h01);
1342 if (!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1343 if (!t.Update(cl,maxChi2,index,h01)) {
1346 } // else SetCluster(cl, GetNumberOfClusters()-1); // A.Bercuci 25.07.07
1349 if (calibra->GetMITracking()) {
1350 calibra->UpdateHistograms(cl,&t);
1353 // Reset material budget if 2 consecutive gold
1355 if ((t.GetTracklets(plane).GetN() + t.GetTracklets(plane-1).GetN()) > 20) {
1366 ratio0 = ncl / Float_t(fTimeBinsPerPlane);
1367 Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
1368 if ((tracklet.GetChi2() < 18.0) &&
1371 (ratio0+ratio1 > 1.5) &&
1372 (t.GetNCross() == 0) &&
1373 (TMath::Abs(t.GetSnp()) < 0.85) &&
1374 (t.GetNumberOfClusters() > 20)){
1375 //if (ratio0 > 0.8) {
1376 t.MakeBackupTrack(); // Make backup of the track until is gold
1381 return expectedNumberOfClusters;
1384 //_____________________________________________________________________________
1385 Int_t AliTRDtracker::PropagateToX(AliTRDtrack &t, Double_t xToGo, Double_t maxStep)
1388 // Starting from current radial position of track <t> this function
1389 // extrapolates the track up to radial position <xToGo>.
1390 // Returns 1 if track reaches the plane, and 0 otherwise
1393 const Double_t kEpsilon = 0.00001;
1394 //Double_t tanmax = TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
1395 Double_t xpos = t.GetX();
1396 Double_t dir = (xpos<xToGo) ? 1.0 : -1.0;
1398 while (((xToGo-xpos)*dir) > kEpsilon) {
1400 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1408 // Starting global position
1412 if (!t.GetProlongation(x,y,z)) {
1413 return 0; // No prolongation
1416 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1417 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1420 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
1421 if (!t.PropagateTo(x,param[1],param[0])) {
1433 //_____________________________________________________________________________
1434 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1437 // Fills clusters into TRD tracking_sectors
1438 // Note that the numbering scheme for the TRD tracking_sectors
1439 // differs from that of TRD sectors
1442 if (ReadClusters(fClusters,cTree)) {
1443 AliError("Problem with reading the clusters !");
1446 Int_t ncl = fClusters->GetEntriesFast();
1448 AliInfo(Form("Sorting %d clusters",ncl));
1451 for (Int_t ichamber = 0; ichamber < 5; ichamber++) {
1452 for (Int_t isector = 0; isector < 18; isector++) {
1453 fHoles[ichamber][isector] = kTRUE;
1459 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(ncl);
1460 Int_t detector = c->GetDetector();
1461 Int_t localTimeBin = c->GetLocalTimeBin();
1462 Int_t sector = fGeom->GetSector(detector);
1463 Int_t plane = fGeom->GetPlane(detector);
1464 Int_t trackingSector = sector;
1466 //if (c->GetLabel(0) > 0) {
1467 if (c->GetQ() > 10) {
1468 Int_t chamber = fGeom->GetChamber(detector);
1469 fHoles[chamber][trackingSector] = kFALSE;
1472 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1476 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1480 // Apply pos correction
1482 fHXCl->Fill(c->GetX());
1484 fTrSec[trackingSector]->GetLayer(layer)->SetX(c->GetX());
1485 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1492 //_____________________________________________________________________________
1493 void AliTRDtracker::UnloadClusters()
1496 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1502 nentr = fClusters->GetEntriesFast();
1503 for (i = 0; i < nentr; i++) {
1504 delete fClusters->RemoveAt(i);
1508 nentr = fSeeds->GetEntriesFast();
1509 for (i = 0; i < nentr; i++) {
1510 delete fSeeds->RemoveAt(i);
1513 nentr = fTracks->GetEntriesFast();
1514 for (i = 0; i < nentr; i++) {
1515 delete fTracks->RemoveAt(i);
1518 Int_t nsec = AliTRDgeometry::kNsect;
1519 for (i = 0; i < nsec; i++) {
1520 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1521 fTrSec[i]->GetLayer(pl)->Clear();
1527 //_____________________________________________________________________________
1528 void AliTRDtracker::MakeSeedsMI(Int_t /*inner*/, Int_t /*outer*/, AliESDEvent *esd)
1531 // Creates seeds using clusters between position inner plane and outer plane
1534 const Double_t kMaxTheta = 1.0;
1535 const Double_t kMaxPhi = 2.0;
1537 const Double_t kRoad0y = 6.0; // Road for middle cluster
1538 const Double_t kRoad0z = 8.5; // Road for middle cluster
1540 const Double_t kRoad1y = 2.0; // Road in y for seeded cluster
1541 const Double_t kRoad1z = 20.0; // Road in z for seeded cluster
1543 const Double_t kRoad2y = 3.0; // Road in y for extrapolated cluster
1544 const Double_t kRoad2z = 20.0; // Road in z for extrapolated cluster
1545 const Int_t kMaxSeed = 3000;
1547 Int_t maxSec = AliTRDgeometry::kNsect;
1549 // Linear fitters in planes
1550 TLinearFitter fitterTC(2,"hyp2"); // Fitting with tilting pads - kz fixed - kz= Z/x, + vertex const
1551 TLinearFitter fitterT2(4,"hyp4"); // Fitting with tilting pads - kz not fixed
1552 fitterTC.StoreData(kTRUE);
1553 fitterT2.StoreData(kTRUE);
1554 AliRieman rieman(1000); // Rieman fitter
1555 AliRieman rieman2(1000); // Rieman fitter
1557 // Find the maximal and minimal layer for the planes
1559 AliTRDpropagationLayer *reflayers[6];
1560 for (Int_t i = 0; i < 6; i++) {
1561 layers[i][0] = 10000;
1564 for (Int_t ns = 0; ns < maxSec; ns++) {
1565 for (Int_t ilayer = 0; ilayer < fTrSec[ns]->GetNumberOfLayers(); ilayer++) {
1566 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(ilayer));
1570 Int_t det = layer[0]->GetDetector();
1571 Int_t plane = fGeom->GetPlane(det);
1572 if (ilayer < layers[plane][0]) {
1573 layers[plane][0] = ilayer;
1575 if (ilayer > layers[plane][1]) {
1576 layers[plane][1] = ilayer;
1581 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(0,0);
1582 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
1583 Double_t hL[6]; // Tilting angle
1584 Double_t xcl[6]; // X - position of reference cluster
1585 Double_t ycl[6]; // Y - position of reference cluster
1586 Double_t zcl[6]; // Z - position of reference cluster
1588 AliTRDcluster *cl[6] = { 0, 0, 0, 0, 0, 0 }; // Seeding clusters
1589 Float_t padlength[6] = { 10.0, 10.0, 10.0, 10.0, 10.0, 10.0 }; // Current pad-length
1591 Double_t chi2R = 0.0;
1592 Double_t chi2Z = 0.0;
1593 Double_t chi2RF = 0.0;
1594 Double_t chi2ZF = 0.0;
1596 Int_t nclusters; // Total number of clusters
1597 for (Int_t i = 0; i < 6; i++) {
1605 AliTRDseed *pseed = new AliTRDseed[kMaxSeed*6];
1606 AliTRDseed *seed[kMaxSeed];
1607 for (Int_t iseed = 0; iseed < kMaxSeed; iseed++) {
1608 seed[iseed]= &pseed[iseed*6];
1610 AliTRDseed *cseed = seed[0];
1612 Double_t seedquality[kMaxSeed];
1613 Double_t seedquality2[kMaxSeed];
1614 Double_t seedparams[kMaxSeed][7];
1615 Int_t seedlayer[kMaxSeed];
1616 Int_t registered = 0;
1617 Int_t sort[kMaxSeed];
1622 for (Int_t ns = 0; ns < maxSec; ns++) { // Loop over sectors
1623 //for (Int_t ns = 0; ns < 5; ns++) { // Loop over sectors
1625 registered = 0; // Reset registerd seed counter
1626 cseed = seed[registered];
1629 for (Int_t sLayer = 2; sLayer >= 0; sLayer--) {
1630 //for (Int_t dseed = 5; dseed < 15; dseed += 3) {
1633 Int_t dseed = 5 + Int_t(iter) * 3;
1635 // Initialize seeding layers
1636 for (Int_t ilayer = 0; ilayer < 6; ilayer++) {
1637 reflayers[ilayer] = fTrSec[ns]->GetLayer(layers[ilayer][1]-dseed);
1638 xcl[ilayer] = reflayers[ilayer]->GetX();
1641 Double_t xref = (xcl[sLayer+1] + xcl[sLayer+2]) * 0.5;
1642 AliTRDpropagationLayer &layer0 = *reflayers[sLayer+0];
1643 AliTRDpropagationLayer &layer1 = *reflayers[sLayer+1];
1644 AliTRDpropagationLayer &layer2 = *reflayers[sLayer+2];
1645 AliTRDpropagationLayer &layer3 = *reflayers[sLayer+3];
1647 Int_t maxn3 = layer3;
1648 for (Int_t icl3 = 0; icl3 < maxn3; icl3++) {
1650 AliTRDcluster *cl3 = layer3[icl3];
1654 padlength[sLayer+3] = TMath::Sqrt(cl3->GetSigmaZ2() * 12.0);
1655 ycl[sLayer+3] = cl3->GetY();
1656 zcl[sLayer+3] = cl3->GetZ();
1657 Float_t yymin0 = ycl[sLayer+3] - 1.0 - kMaxPhi * (xcl[sLayer+3]-xcl[sLayer+0]);
1658 Float_t yymax0 = ycl[sLayer+3] + 1.0 + kMaxPhi * (xcl[sLayer+3]-xcl[sLayer+0]);
1659 Int_t maxn0 = layer0;
1661 for (Int_t icl0 = layer0.Find(yymin0); icl0 < maxn0; icl0++) {
1663 AliTRDcluster *cl0 = layer0[icl0];
1667 if (cl3->IsUsed() && cl0->IsUsed()) {
1670 ycl[sLayer+0] = cl0->GetY();
1671 zcl[sLayer+0] = cl0->GetZ();
1672 if (ycl[sLayer+0] > yymax0) {
1675 Double_t tanphi = (ycl[sLayer+3]-ycl[sLayer+0]) / (xcl[sLayer+3]-xcl[sLayer+0]);
1676 if (TMath::Abs(tanphi) > kMaxPhi) {
1679 Double_t tantheta = (zcl[sLayer+3]-zcl[sLayer+0]) / (xcl[sLayer+3]-xcl[sLayer+0]);
1680 if (TMath::Abs(tantheta) > kMaxTheta) {
1683 padlength[sLayer+0] = TMath::Sqrt(cl0->GetSigmaZ2() * 12.0);
1685 // Expected position in 1 layer
1686 Double_t y1exp = ycl[sLayer+0] + (tanphi) * (xcl[sLayer+1]-xcl[sLayer+0]);
1687 Double_t z1exp = zcl[sLayer+0] + (tantheta) * (xcl[sLayer+1]-xcl[sLayer+0]);
1688 Float_t yymin1 = y1exp - kRoad0y - tanphi;
1689 Float_t yymax1 = y1exp + kRoad0y + tanphi;
1690 Int_t maxn1 = layer1;
1692 for (Int_t icl1 = layer1.Find(yymin1); icl1 < maxn1; icl1++) {
1694 AliTRDcluster *cl1 = layer1[icl1];
1699 if (cl3->IsUsed()) nusedCl++;
1700 if (cl0->IsUsed()) nusedCl++;
1701 if (cl1->IsUsed()) nusedCl++;
1705 ycl[sLayer+1] = cl1->GetY();
1706 zcl[sLayer+1] = cl1->GetZ();
1707 if (ycl[sLayer+1] > yymax1) {
1710 if (TMath::Abs(ycl[sLayer+1]-y1exp) > kRoad0y+tanphi) {
1713 if (TMath::Abs(zcl[sLayer+1]-z1exp) > kRoad0z) {
1716 padlength[sLayer+1] = TMath::Sqrt(cl1->GetSigmaZ2() * 12.0);
1718 Double_t y2exp = ycl[sLayer+0]+(tanphi) * (xcl[sLayer+2]-xcl[sLayer+0]) + (ycl[sLayer+1]-y1exp);
1719 Double_t z2exp = zcl[sLayer+0]+(tantheta) * (xcl[sLayer+2]-xcl[sLayer+0]);
1720 Int_t index2 = layer2.FindNearestCluster(y2exp,z2exp,kRoad1y,kRoad1z);
1724 AliTRDcluster *cl2 = (AliTRDcluster *) GetCluster(index2);
1725 padlength[sLayer+2] = TMath::Sqrt(cl2->GetSigmaZ2() * 12.0);
1726 ycl[sLayer+2] = cl2->GetY();
1727 zcl[sLayer+2] = cl2->GetZ();
1728 if (TMath::Abs(cl2->GetZ()-z2exp) > kRoad0z) {
1733 rieman.AddPoint(xcl[sLayer+0],ycl[sLayer+0],zcl[sLayer+0],1,10);
1734 rieman.AddPoint(xcl[sLayer+1],ycl[sLayer+1],zcl[sLayer+1],1,10);
1735 rieman.AddPoint(xcl[sLayer+3],ycl[sLayer+3],zcl[sLayer+3],1,10);
1736 rieman.AddPoint(xcl[sLayer+2],ycl[sLayer+2],zcl[sLayer+2],1,10);
1740 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1741 cseed[iLayer].Reset();
1746 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1747 cseed[sLayer+iLayer].SetZref(0,rieman.GetZat(xcl[sLayer+iLayer]));
1748 chi2Z += (cseed[sLayer+iLayer].GetZref(0)- zcl[sLayer+iLayer])
1749 * (cseed[sLayer+iLayer].GetZref(0)- zcl[sLayer+iLayer]);
1750 cseed[sLayer+iLayer].SetZref(1,rieman.GetDZat(xcl[sLayer+iLayer]));
1751 cseed[sLayer+iLayer].SetYref(0,rieman.GetYat(xcl[sLayer+iLayer]));
1752 chi2R += (cseed[sLayer+iLayer].GetYref(0)- ycl[sLayer+iLayer])
1753 * (cseed[sLayer+iLayer].GetYref(0)- ycl[sLayer+iLayer]);
1754 cseed[sLayer+iLayer].SetYref(1,rieman.GetDYat(xcl[sLayer+iLayer]));
1756 if (TMath::Sqrt(chi2R) > 1.0/iter) {
1759 if (TMath::Sqrt(chi2Z) > 7.0/iter) {
1763 Float_t minmax[2] = { -100.0, 100.0 };
1764 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1765 Float_t max = zcl[sLayer+iLayer]+padlength[sLayer+iLayer] * 0.5
1766 + 1.0 - cseed[sLayer+iLayer].GetZref(0);
1767 if (max < minmax[1]) {
1770 Float_t min = zcl[sLayer+iLayer]-padlength[sLayer+iLayer] * 0.5
1771 - 1.0 - cseed[sLayer+iLayer].GetZref(0);
1772 if (min > minmax[0]) {
1777 Bool_t isFake = kFALSE;
1778 if (cl0->GetLabel(0) != cl3->GetLabel(0)) {
1781 if (cl1->GetLabel(0) != cl3->GetLabel(0)) {
1784 if (cl2->GetLabel(0) != cl3->GetLabel(0)) {
1788 if (AliTRDReconstructor::StreamLevel() > 0) {
1789 if ((!isFake) || ((icl3%10) == 0)) { // Debugging print
1790 TTreeSRedirector &cstream = *fDebugStreamer;
1792 << "isFake=" << isFake
1798 << "X0=" << xcl[sLayer+0]
1799 << "X1=" << xcl[sLayer+1]
1800 << "X2=" << xcl[sLayer+2]
1801 << "X3=" << xcl[sLayer+3]
1802 << "Y2exp=" << y2exp
1803 << "Z2exp=" << z2exp
1804 << "Chi2R=" << chi2R
1805 << "Chi2Z=" << chi2Z
1806 << "Seed0.=" << &cseed[sLayer+0]
1807 << "Seed1.=" << &cseed[sLayer+1]
1808 << "Seed2.=" << &cseed[sLayer+2]
1809 << "Seed3.=" << &cseed[sLayer+3]
1810 << "Zmin=" << minmax[0]
1811 << "Zmax=" << minmax[1]
1816 ////////////////////////////////////////////////////////////////////////////////////
1820 ////////////////////////////////////////////////////////////////////////////////////
1826 Bool_t isOK = kTRUE;
1828 for (Int_t jLayer = 0; jLayer < 4; jLayer++) {
1830 cseed[sLayer+jLayer].SetTilt(hL[sLayer+jLayer]);
1831 cseed[sLayer+jLayer].SetPadLength(padlength[sLayer+jLayer]);
1832 cseed[sLayer+jLayer].SetX0(xcl[sLayer+jLayer]);
1834 for (Int_t iter = 0; iter < 2; iter++) {
1837 // In iteration 0 we try only one pad-row
1838 // If quality not sufficient we try 2 pad-rows - about 5% of tracks cross 2 pad-rows
1840 AliTRDseed tseed = cseed[sLayer+jLayer];
1841 Float_t roadz = padlength[sLayer+jLayer] * 0.5;
1843 roadz = padlength[sLayer+jLayer];
1846 Float_t quality = 10000.0;
1848 for (Int_t iTime = 2; iTime < 20; iTime++) {
1850 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[sLayer+jLayer][1]-iTime));
1851 Double_t dxlayer = layer.GetX() - xcl[sLayer+jLayer];
1852 Double_t zexp = cl[sLayer+jLayer]->GetZ();
1855 // Try 2 pad-rows in second iteration
1856 zexp = tseed.GetZref(0) + tseed.GetZref(1) * dxlayer;
1857 if (zexp > cl[sLayer+jLayer]->GetZ()) {
1858 zexp = cl[sLayer+jLayer]->GetZ() + padlength[sLayer+jLayer]*0.5;
1860 if (zexp < cl[sLayer+jLayer]->GetZ()) {
1861 zexp = cl[sLayer+jLayer]->GetZ() - padlength[sLayer+jLayer]*0.5;
1865 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer;
1866 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y,roadz);
1870 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
1872 tseed.SetIndexes(iTime,index);
1873 tseed.SetClusters(iTime,cl); // Register cluster
1874 tseed.SetX(iTime,dxlayer); // Register cluster
1875 tseed.SetY(iTime,cl->GetY()); // Register cluster
1876 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
1882 // Count the number of clusters and distortions into quality
1883 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
1884 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0 + TMath::Abs(dangle) / 0.1
1885 + TMath::Abs(tseed.GetYfit(0) - tseed.GetYref(0)) / 0.2
1886 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
1887 if ((iter == 0) && tseed.IsOK()) {
1888 cseed[sLayer+jLayer] = tseed;
1894 if (tseed.IsOK() && (tquality < quality)) {
1895 cseed[sLayer+jLayer] = tseed;
1900 if (!cseed[sLayer+jLayer].IsOK()) {
1905 cseed[sLayer+jLayer].CookLabels();
1906 cseed[sLayer+jLayer].UpdateUsed();
1907 nusedCl += cseed[sLayer+jLayer].GetNUsed();
1919 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1920 if (cseed[sLayer+iLayer].IsOK()) {
1921 nclusters += cseed[sLayer+iLayer].GetN2();
1927 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1928 rieman.AddPoint(xcl[sLayer+iLayer]
1929 ,cseed[sLayer+iLayer].GetYfitR(0)
1930 ,cseed[sLayer+iLayer].GetZProb()
1939 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1940 cseed[sLayer+iLayer].SetYref(0,rieman.GetYat(xcl[sLayer+iLayer]));
1941 chi2R += (cseed[sLayer+iLayer].GetYref(0) - cseed[sLayer+iLayer].GetYfitR(0))
1942 * (cseed[sLayer+iLayer].GetYref(0) - cseed[sLayer+iLayer].GetYfitR(0));
1943 cseed[sLayer+iLayer].SetYref(1,rieman.GetDYat(xcl[sLayer+iLayer]));
1944 cseed[sLayer+iLayer].SetZref(0,rieman.GetZat(xcl[sLayer+iLayer]));
1945 chi2Z += (cseed[sLayer+iLayer].GetZref(0) - cseed[sLayer+iLayer].GetMeanz())
1946 * (cseed[sLayer+iLayer].GetZref(0) - cseed[sLayer+iLayer].GetMeanz());
1947 cseed[sLayer+iLayer].SetZref(1,rieman.GetDZat(xcl[sLayer+iLayer]));
1949 Double_t curv = rieman.GetC();
1954 Double_t sumda = TMath::Abs(cseed[sLayer+0].GetYfitR(1) - cseed[sLayer+0].GetYref(1))
1955 + TMath::Abs(cseed[sLayer+1].GetYfitR(1) - cseed[sLayer+1].GetYref(1))
1956 + TMath::Abs(cseed[sLayer+2].GetYfitR(1) - cseed[sLayer+2].GetYref(1))
1957 + TMath::Abs(cseed[sLayer+3].GetYfitR(1) - cseed[sLayer+3].GetYref(1));
1958 Double_t likea = TMath::Exp(-sumda*10.6);
1959 Double_t likechi2 = 0.0000000001;
1961 likechi2 += TMath::Exp(-TMath::Sqrt(chi2R) * 7.73);
1963 Double_t likechi2z = TMath::Exp(-chi2Z * 0.088) / TMath::Exp(-chi2Z * 0.019);
1964 Double_t likeN = TMath::Exp(-(72 - nclusters) * 0.19);
1965 Double_t like = likea * likechi2 * likechi2z * likeN;
1966 Double_t likePrimY = TMath::Exp(-TMath::Abs(cseed[sLayer+0].GetYref(1) - 130.0*curv) * 1.9);
1967 Double_t likePrimZ = TMath::Exp(-TMath::Abs(cseed[sLayer+0].GetZref(1)
1968 - cseed[sLayer+0].GetZref(0) / xcl[sLayer+0]) * 5.9);
1969 Double_t likePrim = TMath::Max(likePrimY*likePrimZ,0.0005);
1971 seedquality[registered] = like;
1972 seedlayer[registered] = sLayer;
1973 if (TMath::Log(0.000000000000001 + like) < -15) {
1976 AliTRDseed seedb[6];
1977 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1978 seedb[iLayer] = cseed[iLayer];
1981 ////////////////////////////////////////////////////////////////////////////////////
1983 // Full track fit part
1985 ////////////////////////////////////////////////////////////////////////////////////
1992 // Add new layers - avoid long extrapolation
1994 Int_t tLayer[2] = { 0, 0 };
2008 for (Int_t iLayer = 0; iLayer < 2; iLayer++) {
2009 Int_t jLayer = tLayer[iLayer]; // Set tracking layer
2010 cseed[jLayer].Reset();
2011 cseed[jLayer].SetTilt(hL[jLayer]);
2012 cseed[jLayer].SetPadLength(padlength[jLayer]);
2013 cseed[jLayer].SetX0(xcl[jLayer]);
2014 // Get pad length and rough cluster
2015 Int_t indexdummy = reflayers[jLayer]->FindNearestCluster(cseed[jLayer].GetYref(0)
2016 ,cseed[jLayer].GetZref(0)
2019 if (indexdummy <= 0) {
2022 AliTRDcluster *cldummy = (AliTRDcluster *) GetCluster(indexdummy);
2023 padlength[jLayer] = TMath::Sqrt(cldummy->GetSigmaZ2() * 12.0);
2025 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
2027 for (Int_t iLayer = 0; iLayer < 2; iLayer++) {
2029 Int_t jLayer = tLayer[iLayer]; // set tracking layer
2030 if ((jLayer == 0) && !(cseed[1].IsOK())) {
2031 continue; // break not allowed
2033 if ((jLayer == 5) && !(cseed[4].IsOK())) {
2034 continue; // break not allowed
2036 Float_t zexp = cseed[jLayer].GetZref(0);
2037 Double_t zroad = padlength[jLayer] * 0.5 + 1.0;
2039 for (Int_t iter = 0; iter < 2; iter++) {
2041 AliTRDseed tseed = cseed[jLayer];
2042 Float_t quality = 10000.0;
2044 for (Int_t iTime = 2; iTime < 20; iTime++) {
2045 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[jLayer][1]-iTime));
2046 Double_t dxlayer = layer.GetX()-xcl[jLayer];
2047 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer;
2048 Float_t yroad = kRoad1y;
2049 Int_t index = layer.FindNearestCluster(yexp,zexp,yroad,zroad);
2053 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
2054 tseed.SetIndexes(iTime,index);
2055 tseed.SetClusters(iTime,cl); // Register cluster
2056 tseed.SetX(iTime,dxlayer); // Register cluster
2057 tseed.SetY(iTime,cl->GetY()); // Register cluster
2058 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
2063 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
2064 Float_t tquality = (18.0 - tseed.GetN2())/2.0 + TMath::Abs(dangle) / 0.1
2065 + TMath::Abs(tseed.GetYfit(0) - tseed.GetYref(0)) / 0.2
2066 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
2067 if (tquality < quality) {
2068 cseed[jLayer] = tseed;
2077 if ( cseed[jLayer].IsOK()) {
2078 cseed[jLayer].CookLabels();
2079 cseed[jLayer].UpdateUsed();
2080 nusedf += cseed[jLayer].GetNUsed();
2081 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
2087 AliTRDseed bseed[6];
2088 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2089 bseed[jLayer] = cseed[jLayer];
2091 Float_t lastquality = 10000.0;
2092 Float_t lastchi2 = 10000.0;
2093 Float_t chi2 = 1000.0;
2095 for (Int_t iter = 0; iter < 4; iter++) {
2097 // Sort tracklets according "quality", try to "improve" 4 worst
2098 Float_t sumquality = 0.0;
2099 Float_t squality[6];
2100 Int_t sortindexes[6];
2102 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2103 if (bseed[jLayer].IsOK()) {
2104 AliTRDseed &tseed = bseed[jLayer];
2105 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
2106 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
2107 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0 + TMath::Abs(dangle) / 0.1
2108 + TMath::Abs(tseed.GetYfit(0) - (tseed.GetYref(0) - zcor)) / 0.2
2109 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
2110 squality[jLayer] = tquality;
2113 squality[jLayer] = -1.0;
2115 sumquality +=squality[jLayer];
2118 if ((sumquality >= lastquality) ||
2119 (chi2 > lastchi2)) {
2122 lastquality = sumquality;
2125 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2126 cseed[jLayer] = bseed[jLayer];
2129 TMath::Sort(6,squality,sortindexes,kFALSE);
2131 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
2133 Int_t bLayer = sortindexes[jLayer];
2134 AliTRDseed tseed = bseed[bLayer];
2136 for (Int_t iTime = 2; iTime < 20; iTime++) {
2138 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[bLayer][1]-iTime));
2139 Double_t dxlayer = layer.GetX() - xcl[bLayer];
2140 Double_t zexp = tseed.GetZref(0);
2141 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
2142 Float_t roadz = padlength[bLayer] + 1;
2143 if (TMath::Abs(tseed.GetZProb() - zexp) > 0.5*padlength[bLayer]) {
2144 roadz = padlength[bLayer] * 0.5;
2146 if (tseed.GetZfit(1)*tseed.GetZref(1) < 0.0) {
2147 roadz = padlength[bLayer] * 0.5;
2149 if (TMath::Abs(tseed.GetZProb() - zexp) < 0.1*padlength[bLayer]) {
2150 zexp = tseed.GetZProb();
2151 roadz = padlength[bLayer] * 0.5;
2154 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer - zcor;
2155 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y,roadz);
2159 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
2161 tseed.SetIndexes(iTime,index);
2162 tseed.SetClusters(iTime,cl); // Register cluster
2163 tseed.SetX(iTime,dxlayer); // Register cluster
2164 tseed.SetY(iTime,cl->GetY()); // Register cluster
2165 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
2171 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
2172 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
2173 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0
2174 + TMath::Abs(dangle) / 0.1
2175 + TMath::Abs(tseed.GetYfit(0) - (tseed.GetYref(0) - zcor)) / 0.2
2176 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
2177 if (tquality<squality[bLayer]) {
2178 bseed[bLayer] = tseed;
2184 chi2 = AliTRDseed::FitRiemanTilt(bseed,kTRUE);
2191 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2192 if (TMath::Abs(cseed[iLayer].GetYref(0) / cseed[iLayer].GetX0()) < 0.15) {
2195 if (cseed[iLayer].IsOK()) {
2196 nclusters += cseed[iLayer].GetN2();
2204 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2205 if (cseed[iLayer].IsOK()) {
2206 rieman.AddPoint(xcl[iLayer]
2207 ,cseed[iLayer].GetYfitR(0)
2208 ,cseed[iLayer].GetZProb()
2217 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2218 if (cseed[iLayer].IsOK()) {
2219 cseed[iLayer].SetYref(0,rieman.GetYat(xcl[iLayer]));
2220 chi2RF += (cseed[iLayer].GetYref(0) - cseed[iLayer].GetYfitR(0))
2221 * (cseed[iLayer].GetYref(0) - cseed[iLayer].GetYfitR(0));
2222 cseed[iLayer].SetYref(1,rieman.GetDYat(xcl[iLayer]));
2223 cseed[iLayer].SetZref(0,rieman.GetZat(xcl[iLayer]));
2224 chi2ZF += (cseed[iLayer].GetZref(0) - cseed[iLayer].GetMeanz())
2225 * (cseed[iLayer].GetZref(0) - cseed[iLayer].GetMeanz());
2226 cseed[iLayer].SetZref(1,rieman.GetDZat(xcl[iLayer]));
2229 chi2RF /= TMath::Max((nlayers - 3.0),1.0);
2230 chi2ZF /= TMath::Max((nlayers - 3.0),1.0);
2231 curv = rieman.GetC();
2233 Double_t xref2 = (xcl[2] + xcl[3]) * 0.5; // Middle of the chamber
2234 Double_t dzmf = rieman.GetDZat(xref2);
2235 Double_t zmf = rieman.GetZat(xref2);
2241 fitterTC.ClearPoints();
2242 fitterT2.ClearPoints();
2245 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2247 if (!cseed[iLayer].IsOK()) {
2251 for (Int_t itime = 0; itime < 25; itime++) {
2253 if (!cseed[iLayer].IsUsable(itime)) {
2256 // X relative to the middle chamber
2257 Double_t x = cseed[iLayer].GetX(itime) + cseed[iLayer].GetX0() - xref2;
2258 Double_t y = cseed[iLayer].GetY(itime);
2259 Double_t z = cseed[iLayer].GetZ(itime);
2260 // ExB correction to the correction
2264 Double_t x2 = cseed[iLayer].GetX(itime) + cseed[iLayer].GetX0();
2265 Double_t t = 1.0 / (x2*x2 + y*y);
2267 uvt[0] = 2.0 * x2 * uvt[1]; // u
2268 uvt[2] = 2.0 * hL[iLayer] * uvt[1];
2269 uvt[3] = 2.0 * hL[iLayer] * x * uvt[1];
2270 uvt[4] = 2.0 * (y + hL[iLayer]*z) * uvt[1];
2271 Double_t error = 2.0 * 0.2 * uvt[1];
2272 fitterT2.AddPoint(uvt,uvt[4],error);
2275 // Constrained rieman
2277 z = cseed[iLayer].GetZ(itime);
2278 uvt[0] = 2.0 * x2 * t; // u
2279 uvt[1] = 2.0 * hL[iLayer] * x2 * uvt[1];
2280 uvt[2] = 2.0 * (y + hL[iLayer] * (z - GetZ())) * t;
2281 fitterTC.AddPoint(uvt,uvt[2],error);
2282 rieman2.AddPoint(x2,y,z,1,10);
2292 Double_t rpolz0 = fitterT2.GetParameter(3);
2293 Double_t rpolz1 = fitterT2.GetParameter(4);
2296 // Linear fitter - not possible to make boundaries
2297 // Do not accept non possible z and dzdx combinations
2299 Bool_t acceptablez = kTRUE;
2300 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2301 if (cseed[iLayer].IsOK()) {
2302 Double_t zT2 = rpolz0 + rpolz1 * (xcl[iLayer] - xref2);
2303 if (TMath::Abs(cseed[iLayer].GetZProb() - zT2) > padlength[iLayer] * 0.5 + 1.0) {
2304 acceptablez = kFALSE;
2309 fitterT2.FixParameter(3,zmf);
2310 fitterT2.FixParameter(4,dzmf);
2312 fitterT2.ReleaseParameter(3);
2313 fitterT2.ReleaseParameter(4);
2314 rpolz0 = fitterT2.GetParameter(3);
2315 rpolz1 = fitterT2.GetParameter(4);
2318 Double_t chi2TR = fitterT2.GetChisquare() / Float_t(npointsT);
2319 Double_t chi2TC = fitterTC.GetChisquare() / Float_t(npointsT);
2320 Double_t polz1c = fitterTC.GetParameter(2);
2321 Double_t polz0c = polz1c * xref2;
2322 Double_t aC = fitterTC.GetParameter(0);
2323 Double_t bC = fitterTC.GetParameter(1);
2324 Double_t cC = aC / TMath::Sqrt(bC * bC + 1.0); // Curvature
2325 Double_t aR = fitterT2.GetParameter(0);
2326 Double_t bR = fitterT2.GetParameter(1);
2327 Double_t dR = fitterT2.GetParameter(2);
2328 Double_t cR = 1.0 + bR*bR - dR*aR;
2331 dca = -dR / (TMath::Sqrt(1.0 + bR*bR - dR*aR) + TMath::Sqrt(1.0 + bR*bR));
2332 cR = aR / TMath::Sqrt(cR);
2335 Double_t chi2ZT2 = 0.0;
2336 Double_t chi2ZTC = 0.0;
2337 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2338 if (cseed[iLayer].IsOK()) {
2339 Double_t zT2 = rpolz0 + rpolz1 * (xcl[iLayer] - xref2);
2340 Double_t zTC = polz0c + polz1c * (xcl[iLayer] - xref2);
2341 chi2ZT2 += TMath::Abs(cseed[iLayer].GetMeanz() - zT2);
2342 chi2ZTC += TMath::Abs(cseed[iLayer].GetMeanz() - zTC);
2345 chi2ZT2 /= TMath::Max((nlayers - 3.0),1.0);
2346 chi2ZTC /= TMath::Max((nlayers - 3.0),1.0);
2348 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
2349 Float_t sumdaf = 0.0;
2350 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2351 if (cseed[iLayer].IsOK()) {
2352 sumdaf += TMath::Abs((cseed[iLayer].GetYfit(1) - cseed[iLayer].GetYref(1))
2353 / cseed[iLayer].GetSigmaY2());
2356 sumdaf /= Float_t (nlayers - 2.0);
2359 // Likelihoods for full track
2361 Double_t likezf = TMath::Exp(-chi2ZF * 0.14);
2362 Double_t likechi2C = TMath::Exp(-chi2TC * 0.677);
2363 Double_t likechi2TR = TMath::Exp(-chi2TR * 0.78);
2364 Double_t likeaf = TMath::Exp(-sumdaf * 3.23);
2365 seedquality2[registered] = likezf * likechi2TR * likeaf;
2367 // Still needed ????
2368 // Bool_t isGold = kFALSE;
2370 // if (nlayers == 6 && TMath::Log(0.000000001+seedquality2[index])<-5.) isGold =kTRUE; // gold
2371 // if (nlayers == findable && TMath::Log(0.000000001+seedquality2[index])<-4.) isGold =kTRUE; // gold
2372 // if (isGold &&nusedf<10){
2373 // for (Int_t jLayer=0;jLayer<6;jLayer++){
2374 // if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.1)
2375 // seed[index][jLayer].UseClusters(); //sign gold
2380 if (!cseed[0].IsOK()) {
2382 if (!cseed[1].IsOK()) {
2386 seedparams[registered][0] = cseed[index0].GetX0();
2387 seedparams[registered][1] = cseed[index0].GetYref(0);
2388 seedparams[registered][2] = cseed[index0].GetZref(0);
2389 seedparams[registered][5] = cR;
2390 seedparams[registered][3] = cseed[index0].GetX0() * cR - TMath::Sin(TMath::ATan(cseed[0].GetYref(1)));
2391 seedparams[registered][4] = cseed[index0].GetZref(1)
2392 / TMath::Sqrt(1.0 + cseed[index0].GetYref(1) * cseed[index0].GetYref(1));
2393 seedparams[registered][6] = ns;
2398 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2399 if (!cseed[iLayer].IsOK()) {
2402 if (cseed[iLayer].GetLabels(0) >= 0) {
2403 labels[nlab] = cseed[iLayer].GetLabels(0);
2406 if (cseed[iLayer].GetLabels(1) >= 0) {
2407 labels[nlab] = cseed[iLayer].GetLabels(1);
2411 Freq(nlab,labels,outlab,kFALSE);
2412 Int_t label = outlab[0];
2413 Int_t frequency = outlab[1];
2414 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2415 cseed[iLayer].SetFreq(frequency);
2416 cseed[iLayer].SetC(cR);
2417 cseed[iLayer].SetCC(cC);
2418 cseed[iLayer].SetChi2(chi2TR);
2419 cseed[iLayer].SetChi2Z(chi2ZF);
2423 if (1 || (!isFake)) {
2424 Float_t zvertex = GetZ();
2425 TTreeSRedirector &cstream = *fDebugStreamer;
2426 if (AliTRDReconstructor::StreamLevel() > 0) {
2428 << "isFake=" << isFake
2429 << "Vertex=" << zvertex
2430 << "Rieman2.=" << &rieman2
2431 << "Rieman.=" << &rieman
2439 << "Chi2R=" << chi2R
2440 << "Chi2Z=" << chi2Z
2441 << "Chi2RF=" << chi2RF // Chi2 of trackletes on full track
2442 << "Chi2ZF=" << chi2ZF // Chi2 z on tracklets on full track
2443 << "Chi2ZT2=" << chi2ZT2 // Chi2 z on tracklets on full track - rieman tilt
2444 << "Chi2ZTC=" << chi2ZTC // Chi2 z on tracklets on full track - rieman tilt const
2445 << "Chi2TR=" << chi2TR // Chi2 without vertex constrain
2446 << "Chi2TC=" << chi2TC // Chi2 with vertex constrain
2447 << "C=" << curv // Non constrained - no tilt correction
2448 << "DR=" << dR // DR parameter - tilt correction
2449 << "DCA=" << dca // DCA - tilt correction
2450 << "CR=" << cR // Non constrained curvature - tilt correction
2451 << "CC=" << cC // Constrained curvature
2452 << "Polz0=" << polz0c
2453 << "Polz1=" << polz1c
2454 << "RPolz0=" << rpolz0
2455 << "RPolz1=" << rpolz1
2456 << "Ncl=" << nclusters
2457 << "Nlayers=" << nlayers
2458 << "NUsedS=" << nusedCl
2459 << "NUsed=" << nusedf
2460 << "Findable=" << findable
2462 << "LikePrim=" << likePrim
2463 << "Likechi2C=" << likechi2C
2464 << "Likechi2TR=" << likechi2TR
2465 << "Likezf=" << likezf
2466 << "LikeF=" << seedquality2[registered]
2467 << "S0.=" << &cseed[0]
2468 << "S1.=" << &cseed[1]
2469 << "S2.=" << &cseed[2]
2470 << "S3.=" << &cseed[3]
2471 << "S4.=" << &cseed[4]
2472 << "S5.=" << &cseed[5]
2473 << "SB0.=" << &seedb[0]
2474 << "SB1.=" << &seedb[1]
2475 << "SB2.=" << &seedb[2]
2476 << "SB3.=" << &seedb[3]
2477 << "SB4.=" << &seedb[4]
2478 << "SB5.=" << &seedb[5]
2479 << "Label=" << label
2480 << "Freq=" << frequency
2481 << "sLayer=" << sLayer
2486 if (registered<kMaxSeed - 1) {
2488 cseed = seed[registered];
2491 } // End of loop over layer 1
2493 } // End of loop over layer 0
2495 } // End of loop over layer 3
2497 } // End of loop over seeding time bins
2503 TMath::Sort(registered,seedquality2,sort,kTRUE);
2504 Bool_t signedseed[kMaxSeed];
2505 for (Int_t i = 0; i < registered; i++) {
2506 signedseed[i] = kFALSE;
2509 for (Int_t iter = 0; iter < 5; iter++) {
2511 for (Int_t iseed = 0; iseed < registered; iseed++) {
2513 Int_t index = sort[iseed];
2514 if (signedseed[index]) {
2517 Int_t labelsall[1000];
2518 Int_t nlabelsall = 0;
2519 Int_t naccepted = 0;;
2520 Int_t sLayer = seedlayer[index];
2526 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2528 if (TMath::Abs(seed[index][jLayer].GetYref(0) / xcl[jLayer]) < 0.15) {
2531 if (seed[index][jLayer].IsOK()) {
2532 seed[index][jLayer].UpdateUsed();
2533 ncl +=seed[index][jLayer].GetN2();
2534 nused +=seed[index][jLayer].GetNUsed();
2537 for (Int_t itime = 0; itime < 25; itime++) {
2538 if (seed[index][jLayer].IsUsable(itime)) {
2540 for (Int_t ilab = 0; ilab < 3; ilab++) {
2541 Int_t tindex = seed[index][jLayer].GetClusters(itime)->GetLabel(ilab);
2543 labelsall[nlabelsall] = tindex;
2561 if (TMath::Log(0.000000001+seedquality2[index]) < -5.0) {
2567 if (nlayers < findable) {
2570 if (TMath::Log(0.000000001+seedquality2[index]) < -4.0) {
2576 if ((nlayers == findable) ||
2580 if (TMath::Log(0.000000001+seedquality2[index]) < -6.0) {
2586 if (TMath::Log(0.000000001+seedquality2[index]) < -5.0) {
2592 if (TMath::Log(0.000000001+seedquality2[index]) - nused/(nlayers-3.0) < -15.0) {
2597 signedseed[index] = kTRUE;
2602 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2603 if (seed[index][iLayer].IsOK()) {
2604 if (seed[index][iLayer].GetLabels(0) >= 0) {
2605 labels[nlab] = seed[index][iLayer].GetLabels(0);
2608 if (seed[index][iLayer].GetLabels(1) >= 0) {
2609 labels[nlab] = seed[index][iLayer].GetLabels(1);
2614 Freq(nlab,labels,outlab,kFALSE);
2615 Int_t label = outlab[0];
2616 Int_t frequency = outlab[1];
2617 Freq(nlabelsall,labelsall,outlab,kFALSE);
2618 Int_t label1 = outlab[0];
2619 Int_t label2 = outlab[2];
2620 Float_t fakeratio = (naccepted - outlab[1]) / Float_t(naccepted);
2621 Float_t ratio = Float_t(nused) / Float_t(ncl);
2623 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2624 if ((seed[index][jLayer].IsOK()) &&
2625 (TMath::Abs(seed[index][jLayer].GetYfit(1) - seed[index][jLayer].GetYfit(1)) < 0.2)) {
2626 seed[index][jLayer].UseClusters(); // Sign gold
2631 Int_t eventNrInFile = esd->GetEventNumberInFile(); // This is most likely NOT the event number you'd like to use. It has nothing to do with the 'real' event number.
2632 TTreeSRedirector &cstream = *fDebugStreamer;
2637 AliTRDtrack *track = RegisterSeed(seed[index],seedparams[index]);
2643 AliESDtrack esdtrack;
2644 esdtrack.UpdateTrackParams(track,AliESDtrack::kTRDout);
2645 esdtrack.SetLabel(label);
2646 esd->AddTrack(&esdtrack);
2647 TTreeSRedirector &cstream = *fDebugStreamer;
2648 if (AliTRDReconstructor::StreamLevel() > 0) {
2650 << "EventNrInFile=" << eventNrInFile
2651 << "ESD.=" << &esdtrack
2653 << "trdback.=" << track
2658 if (AliTRDReconstructor::StreamLevel() > 0) {
2661 << "Track.=" << track
2662 << "Like=" << seedquality[index]
2663 << "LikeF=" << seedquality2[index]
2664 << "S0.=" << &seed[index][0]
2665 << "S1.=" << &seed[index][1]
2666 << "S2.=" << &seed[index][2]
2667 << "S3.=" << &seed[index][3]
2668 << "S4.=" << &seed[index][4]
2669 << "S5.=" << &seed[index][5]
2670 << "Label=" << label
2671 << "Label1=" << label1
2672 << "Label2=" << label2
2673 << "FakeRatio=" << fakeratio
2674 << "Freq=" << frequency
2676 << "Nlayers=" << nlayers
2677 << "Findable=" << findable
2678 << "NUsed=" << nused
2679 << "sLayer=" << sLayer
2680 << "EventNrInFile=" << eventNrInFile
2688 } // End of loop over sectors
2694 //_____________________________________________________________________________
2695 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *clusterTree) const
2698 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
2699 // from the file. The names of the cluster tree and branches
2700 // should match the ones used in AliTRDclusterizer::WriteClusters()
2703 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
2704 TObjArray *clusterArray = new TObjArray(nsize+1000);
2706 TBranch *branch = clusterTree->GetBranch("TRDcluster");
2708 AliError("Can't get the branch !");
2711 branch->SetAddress(&clusterArray);
2713 // Loop through all entries in the tree
2714 Int_t nEntries = (Int_t) clusterTree->GetEntries();
2716 AliTRDcluster *c = 0;
2717 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
2720 nbytes += clusterTree->GetEvent(iEntry);
2722 // Get the number of points in the detector
2723 Int_t nCluster = clusterArray->GetEntriesFast();
2725 // Loop through all TRD digits
2726 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
2727 c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster);
2728 AliTRDcluster *co = c;
2730 clusterArray->RemoveAt(iCluster);
2735 delete clusterArray;
2741 //_____________________________________________________________________________
2742 Bool_t AliTRDtracker::GetTrackPoint(Int_t index, AliTrackPoint &p) const
2745 // Get track space point with index i
2746 // Origin: C.Cheshkov
2749 AliTRDcluster *cl = (AliTRDcluster *) fClusters->UncheckedAt(index);
2750 Int_t idet = cl->GetDetector();
2751 Int_t isector = fGeom->GetSector(idet);
2752 Int_t ichamber = fGeom->GetChamber(idet);
2753 Int_t iplan = fGeom->GetPlane(idet);
2755 local[0] = GetX(isector,iplan,cl->GetLocalTimeBin());
2756 local[1] = cl->GetY();
2757 local[2] = cl->GetZ();
2759 fGeom->RotateBack(idet,local,global);
2760 p.SetXYZ(global[0],global[1],global[2]);
2761 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
2764 iLayer = AliGeomManager::kTRD1;
2767 iLayer = AliGeomManager::kTRD2;
2770 iLayer = AliGeomManager::kTRD3;
2773 iLayer = AliGeomManager::kTRD4;
2776 iLayer = AliGeomManager::kTRD5;
2779 iLayer = AliGeomManager::kTRD6;
2782 Int_t modId = isector * fGeom->Ncham() + ichamber;
2783 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer,modId);
2784 p.SetVolumeID(volid);
2790 //_____________________________________________________________________________
2791 void AliTRDtracker::CookLabel(AliKalmanTrack *pt, Float_t wrong) const
2794 // This cooks a label. Mmmmh, smells good...
2797 Int_t label = 123456789;
2801 Int_t ncl = pt->GetNumberOfClusters();
2803 const Int_t kRange = fTrSec[0]->GetOuterTimeBin() + 1;
2807 Int_t **s = new Int_t* [kRange];
2808 for (i = 0; i < kRange; i++) {
2809 s[i] = new Int_t[2];
2811 for (i = 0; i < kRange; i++) {
2820 for (i = 0; i < ncl; i++) {
2821 index = pt->GetClusterIndex(i);
2822 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2828 for (i = 0; i < ncl; i++) {
2829 index = pt->GetClusterIndex(i);
2830 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2831 for (Int_t k = 0; k < 3; k++) {
2832 label = c->GetLabel(k);
2833 labelAdded = kFALSE;
2836 while ((!labelAdded) && (j < kRange)) {
2837 if ((s[j][0] == label) ||
2840 s[j][1] = s[j][1] + 1;
2852 for (i = 0; i < kRange; i++) {
2853 if (s[i][1] > max) {
2859 for (i = 0; i < kRange; i++) {
2865 if ((1.0 - Float_t(max)/ncl) > wrong) {
2869 pt->SetLabel(label);
2873 //_____________________________________________________________________________
2874 void AliTRDtracker::UseClusters(const AliKalmanTrack *t, Int_t from) const
2877 // Use clusters, but don't abuse them!
2880 const Float_t kmaxchi2 = 18;
2881 const Float_t kmincl = 10;
2883 AliTRDtrack *track = (AliTRDtrack *) t;
2885 Int_t ncl = t->GetNumberOfClusters();
2886 for (Int_t i = from; i < ncl; i++) {
2887 Int_t index = t->GetClusterIndex(i);
2888 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2889 Int_t iplane = fGeom->GetPlane(c->GetDetector());
2890 if (track->GetTracklets(iplane).GetChi2() > kmaxchi2) {
2893 if (track->GetTracklets(iplane).GetN() < kmincl) {
2896 if (!(c->IsUsed())) {
2903 //_____________________________________________________________________________
2904 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
2907 // Parametrised "expected" error of the cluster reconstruction in Y
2910 Double_t s = 0.08 * 0.08;
2915 //_____________________________________________________________________________
2916 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
2919 // Parametrised "expected" error of the cluster reconstruction in Z
2922 Double_t s = 9.0 * 9.0 / 12.0;
2927 //_____________________________________________________________________________
2928 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
2931 // Returns radial position which corresponds to time bin <localTB>
2932 // in tracking sector <sector> and plane <plane>
2935 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
2936 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
2938 return fTrSec[sector]->GetLayer(pl)->GetX();
2942 //_____________________________________________________________________________
2943 AliTRDtracker::AliTRDpropagationLayer
2944 ::AliTRDpropagationLayer(Double_t x, Double_t dx, Double_t rho
2945 , Double_t radLength, Int_t tbIndex, Int_t plane)
2954 ,fTimeBinIndex(tbIndex)
2967 // AliTRDpropagationLayer constructor
2970 for (Int_t i = 0; i < (Int_t) kZones; i++) {
2975 if (fTimeBinIndex >= 0) {
2976 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
2977 fIndex = new UInt_t[kMaxClusterPerTimeBin];
2980 for (Int_t i = 0; i < 5; i++) {
2981 fIsHole[i] = kFALSE;
2986 //_____________________________________________________________________________
2987 void AliTRDtracker::AliTRDpropagationLayer
2988 ::SetHole(Double_t Zmax, Double_t Ymax, Double_t rho
2989 , Double_t radLength, Double_t Yc, Double_t Zc)
2992 // Sets hole in the layer
3001 fHoleX0 = radLength;
3005 //_____________________________________________________________________________
3006 AliTRDtracker::AliTRDtrackingSector
3007 ::AliTRDtrackingSector(AliTRDgeometry *geo, Int_t gs)
3013 // AliTRDtrackingSector Constructor
3016 AliTRDpadPlane *padPlane = 0;
3017 AliTRDpropagationLayer *ppl = 0;
3019 // Get holes description from geometry
3020 Bool_t holes[AliTRDgeometry::kNcham];
3021 for (Int_t icham = 0; icham < AliTRDgeometry::kNcham; icham++) {
3022 holes[icham] = fGeom->IsHole(0,icham,gs);
3025 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
3026 fTimeBinIndex[i] = -1;
3034 // Add layers for each of the planes
3035 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
3036 //Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
3038 const Int_t kNchambers = AliTRDgeometry::Ncham();
3041 Double_t ymaxsensitive = 0;
3042 Double_t *zc = new Double_t[kNchambers];
3043 Double_t *zmax = new Double_t[kNchambers];
3044 Double_t *zmaxsensitive = new Double_t[kNchambers];
3046 for (Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
3048 ymax = fGeom->GetChamberWidth(plane) / 2.0;
3049 padPlane = fGeom->GetPadPlane(plane,0);
3050 ymaxsensitive = (padPlane->GetColSize(1) * padPlane->GetNcols() - 4.0) / 2.0;
3052 for (Int_t ch = 0; ch < kNchambers; ch++) {
3053 zmax[ch] = fGeom->GetChamberLength(plane,ch) / 2.0;
3054 Float_t pad = padPlane->GetRowSize(1);
3055 Float_t row0 = fGeom->GetRow0(plane,ch,0);
3056 Int_t nPads = fGeom->GetRowMax(plane,ch,0);
3057 zmaxsensitive[ch] = Float_t(nPads) * pad / 2.0;
3058 zc[ch] = -(pad * nPads) / 2.0 + row0;
3061 dx = AliTRDcalibDB::Instance()->GetVdrift(0,0,0)
3062 / AliTRDCommonParam::Instance()->GetSamplingFrequency();
3063 rho = 0.00295 * 0.85; //????
3066 Double_t x0 = (Double_t) AliTRDgeometry::GetTime0(plane);
3067 //Double_t xbottom = x0 - dxDrift;
3068 //Double_t xtop = x0 + dxAmp;
3070 Int_t nTimeBins = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
3071 for (Int_t iTime = 0; iTime < nTimeBins; iTime++) {
3073 Double_t xlayer = iTime * dx - dxAmp;
3074 //if (xlayer<0) xlayer = dxAmp / 2.0;
3077 tbIndex = CookTimeBinIndex(plane,iTime);
3078 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex,plane);
3079 ppl->SetYmax(ymax,ymaxsensitive);
3080 ppl->SetZ(zc,zmax,zmaxsensitive);
3081 ppl->SetHoles(holes);
3092 delete [] zmaxsensitive;
3096 //_____________________________________________________________________________
3097 AliTRDtracker::AliTRDtrackingSector
3098 ::AliTRDtrackingSector(const AliTRDtrackingSector &/*t*/)
3109 //_____________________________________________________________________________
3110 Int_t AliTRDtracker::AliTRDtrackingSector
3111 ::CookTimeBinIndex(Int_t plane, Int_t localTB) const
3114 // depending on the digitization parameters calculates "global"
3115 // time bin index for timebin <localTB> in plane <plane>
3119 Int_t tbPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
3120 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1;
3132 //_____________________________________________________________________________
3133 void AliTRDtracker::AliTRDtrackingSector
3134 ::MapTimeBinLayers()
3137 // For all sensitive time bins sets corresponding layer index
3138 // in the array fTimeBins
3143 for (Int_t i = 0; i < fN; i++) {
3145 index = fLayers[i]->GetTimeBinIndex();
3150 if (index >= (Int_t) kMaxTimeBinIndex) {
3151 //AliWarning(Form("Index %d exceeds allowed maximum of %d!\n"
3152 // ,index,kMaxTimeBinIndex-1));
3156 fTimeBinIndex[index] = i;
3162 //_____________________________________________________________________________
3163 Int_t AliTRDtracker::AliTRDtrackingSector
3164 ::GetLayerNumber(Double_t x) const
3167 // Returns the number of time bin which in radial position is closest to <x>
3170 if (x >= fLayers[fN-1]->GetX()) {
3173 if (x <= fLayers[ 0]->GetX()) {
3179 Int_t m = (b + e) / 2;
3181 for ( ; b < e; m = (b + e) / 2) {
3182 if (x > fLayers[m]->GetX()) {
3190 if (TMath::Abs(x - fLayers[m]->GetX()) > TMath::Abs(x - fLayers[m+1]->GetX())) {
3199 //_____________________________________________________________________________
3200 Int_t AliTRDtracker::AliTRDtrackingSector
3201 ::GetInnerTimeBin() const
3204 // Returns number of the innermost SENSITIVE propagation layer
3207 return GetLayerNumber(0);
3211 //_____________________________________________________________________________
3212 Int_t AliTRDtracker::AliTRDtrackingSector
3213 ::GetOuterTimeBin() const
3216 // Returns number of the outermost SENSITIVE time bin
3219 return GetLayerNumber(GetNumberOfTimeBins() - 1);
3223 //_____________________________________________________________________________
3224 Int_t AliTRDtracker::AliTRDtrackingSector
3225 ::GetNumberOfTimeBins() const
3228 // Returns number of SENSITIVE time bins
3234 for (tb = kMaxTimeBinIndex - 1; tb >= 0; tb--) {
3235 layer = GetLayerNumber(tb);
3245 //_____________________________________________________________________________
3246 void AliTRDtracker::AliTRDtrackingSector
3247 ::InsertLayer(AliTRDpropagationLayer *pl)
3250 // Insert layer <pl> in fLayers array.
3251 // Layers are sorted according to X coordinate.
3254 if (fN == ((Int_t) kMaxLayersPerSector)) {
3255 //AliWarning("Too many layers !\n");
3264 Int_t i = Find(pl->GetX());
3266 memmove(fLayers+i+1,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
3273 //_____________________________________________________________________________
3274 Int_t AliTRDtracker::AliTRDtrackingSector
3275 ::Find(Double_t x) const
3278 // Returns index of the propagation layer nearest to X
3281 if (x <= fLayers[0]->GetX()) {
3285 if (x > fLayers[fN-1]->GetX()) {
3291 Int_t m = (b + e) / 2;
3293 for (; b < e; m = (b + e) / 2) {
3294 if (x > fLayers[m]->GetX()) {
3306 //_____________________________________________________________________________
3307 void AliTRDtracker::AliTRDpropagationLayer
3308 ::SetZ(Double_t *center, Double_t *w, Double_t *wsensitive )
3311 // set centers and the width of sectors
3314 for (Int_t icham = 0; icham < AliTRDgeometry::kNcham; icham++) {
3315 fZc[icham] = center[icham];
3316 fZmax[icham] = w[icham];
3317 fZmaxSensitive[icham] = wsensitive[icham];
3322 //_____________________________________________________________________________
3323 void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
3326 // set centers and the width of sectors
3331 for (Int_t icham = 0; icham < AliTRDgeometry::kNcham; icham++) {
3332 fIsHole[icham] = holes[icham];
3340 //_____________________________________________________________________________
3341 void AliTRDtracker::AliTRDpropagationLayer
3342 ::InsertCluster(AliTRDcluster *c, UInt_t index)
3345 // Insert cluster in cluster array.
3346 // Clusters are sorted according to Y coordinate.
3349 if (fTimeBinIndex < 0) {
3350 //AliWarning("Attempt to insert cluster into non-sensitive time bin!\n");
3354 if (fN == (Int_t) kMaxClusterPerTimeBin) {
3355 //AliWarning("Too many clusters !\n");
3361 fClusters[fN++] = c;
3365 Int_t i = Find(c->GetY());
3366 memmove(fClusters+i+1,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
3367 memmove(fIndex +i+1,fIndex +i,(fN-i)*sizeof(UInt_t));
3374 //_____________________________________________________________________________
3375 Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Float_t y) const
3378 // Returns index of the cluster nearest in Y
3384 if (y <= fClusters[0]->GetY()) {
3387 if (y > fClusters[fN-1]->GetY()) {
3393 Int_t m = (b + e) / 2;
3395 for ( ; b < e; m = (b + e) / 2) {
3396 if (y > fClusters[m]->GetY()) {
3408 //_____________________________________________________________________________
3409 Int_t AliTRDtracker::AliTRDpropagationLayer
3410 ::FindNearestCluster(Float_t y, Float_t z, Float_t maxroad
3411 , Float_t maxroadz) const
3414 // Returns index of the cluster nearest to the given y,z
3419 Float_t mindist = maxroad;
3421 for (Int_t i = Find(y-maxroad); i < maxn; i++) {
3422 AliTRDcluster *c = (AliTRDcluster *) (fClusters[i]);
3423 Float_t ycl = c->GetY();
3424 if (ycl > (y + maxroad)) {
3427 if (TMath::Abs(c->GetZ() - z) > maxroadz) {
3430 if (TMath::Abs(ycl - y) < mindist) {
3431 mindist = TMath::Abs(ycl - y);
3440 //_____________________________________________________________________________
3441 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster *c)
3444 // Returns correction factor for tilted pads geometry
3447 Int_t det = c->GetDetector();
3448 Int_t plane = fGeom->GetPlane(det);
3449 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(plane,0);
3450 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
3461 //_____________________________________________________________________________
3462 Int_t AliTRDtracker::FindClusters(Int_t sector, Int_t t0, Int_t t1
3463 , AliTRDtrack *track
3464 , Int_t *clusters, AliTRDtracklet &tracklet)
3468 // Try to find nearest clusters to the track in timebins from t0 to t1
3470 // Correction coeficients - depend on TRD parameters - to be changed accordingly
3476 Double_t xmean = 0.0; // Reference x
3477 Double_t dz[10][100];
3478 Double_t dy[10][100];
3482 Int_t indexes[10][100]; // Indexes of the clusters in the road
3483 Int_t best[10][100]; // Index of best matching cluster
3484 AliTRDcluster *cl[10][100]; // Pointers to the clusters in the road
3486 for (Int_t it = 0; it < 100; it++) {
3493 for (Int_t ih = 0; ih < 10;ih++) {
3494 indexes[ih][it] = -2; // Reset indexes1
3496 dz[ih][it] = -100.0;
3497 dy[ih][it] = -100.0;
3502 Double_t x0 = track->GetX();
3503 Double_t sigmaz = TMath::Sqrt(TMath::Abs(track->GetSigmaZ2()));
3508 Int_t detector = -1;
3509 Float_t padlength = 0.0;
3510 AliTRDtrack track2(* track);
3511 Float_t snpy = track->GetSnp();
3512 Float_t tany = TMath::Sqrt(snpy*snpy / (1.0 - snpy*snpy));
3517 Double_t sy2 = ExpectedSigmaY2(x0,track->GetTgl(),track->GetPt());
3518 Double_t sz2 = ExpectedSigmaZ2(x0,track->GetTgl());
3519 Double_t road = 15.0 * TMath::Sqrt(track->GetSigmaY2() + sy2);
3525 for (Int_t it = 0; it < t1-t0; it++) {
3527 Double_t maxChi2[2] = { fgkMaxChi2, fgkMaxChi2 };
3528 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(it+t0));
3530 continue; // No indexes1
3533 Int_t maxn = timeBin;
3534 x[it] = timeBin.GetX();
3535 track2.PropagateTo(x[it]);
3536 yt[it] = track2.GetY();
3537 zt[it] = track2.GetZ();
3539 Double_t y = yt[it];
3540 Double_t z = zt[it];
3541 Double_t chi2 = 1000000.0;
3545 // Find 2 nearest cluster at given time bin
3547 int checkPoint[4] = {0,0,0,0};
3548 double minY = 123456789;
3549 double minD[2] = {1,1};
3551 for (Int_t i = timeBin.Find(y - road); i < maxn; i++) {
3552 //for (Int_t i = 0; i < maxn; i++) {
3554 AliTRDcluster *c = (AliTRDcluster *) (timeBin[i]);
3555 h01 = GetTiltFactor(c);
3557 Int_t det = c->GetDetector();
3558 plane = fGeom->GetPlane(det);
3559 padlength = TMath::Sqrt(c->GetSigmaZ2() * 12.0);
3562 //if (c->GetLocalTimeBin()==0) continue;
3563 if (c->GetY() > (y + road)) {
3567 fHDeltaX->Fill(c->GetX() - x[it]);
3568 //printf("%f\t%f\t%f \n", c->GetX(), x[it], c->GetX()-x[it]);
3570 if (TMath::Abs(c->GetY()-y) < TMath::Abs(minY)) {
3572 minD[0] = c->GetY()-y;
3573 minD[1] = c->GetZ()-z;
3578 fHMinZ->Fill(c->GetZ() - z);
3579 if ((c->GetZ() - z) * (c->GetZ() - z) > 2 * (12.0 * sz2)) {
3584 Double_t dist = TMath::Abs(c->GetZ() - z);
3585 if (dist > (0.5 * padlength + 6.0 * sigmaz)) { // 0.5
3586 continue; // 6 sigma boundary cut
3590 Double_t cost = 0.0;
3591 // Sigma boundary cost function
3592 if (dist> (0.5 * padlength - sigmaz)){
3593 cost = (dist - 0.5*padlength) / (2.0 * sigmaz);
3595 cost = (cost + 1.0) * (cost + 1.0);
3601 //Int_t label = TMath::Abs(track->GetLabel());
3602 //if (c->GetLabel(0)!=label && c->GetLabel(1)!=label&&c->GetLabel(2)!=label) continue;
3603 chi2 = track2.GetPredictedChi2(c,h01) + cost;
3606 if (chi2 > maxChi2[1]) {
3611 detector = c->GetDetector();
3612 // Store the clusters in the road
3613 for (Int_t ih = 2; ih < 9; ih++) {
3614 if (cl[ih][it] == 0) {
3616 indexes[ih][it] = timeBin.GetIndex(i); // Index - 9 - reserved for outliers
3621 if (chi2 < maxChi2[0]) {
3622 maxChi2[1] = maxChi2[0];
3624 indexes[1][it] = indexes[0][it];
3625 cl[1][it] = cl[0][it];
3626 indexes[0][it] = timeBin.GetIndex(i);
3632 indexes[1][it] = timeBin.GetIndex(i);
3636 for(int iCheckPoint = 0; iCheckPoint<4; iCheckPoint++)
3637 fHFindCl[iCheckPoint]->Fill(checkPoint[iCheckPoint]);
3639 if (checkPoint[3]) {
3640 if (track->GetPt() > 0) fHMinYPos->Fill(minY);
3641 else fHMinYNeg->Fill(minY);
3643 fHMinD->Fill(minD[0], minD[1]);
3656 xmean /= Float_t(nfound); // Middle x
3657 track2.PropagateTo(xmean); // Propagate track to the center
3660 // Choose one of the variants
3665 Double_t sumdy = 0.0;
3666 Double_t sumdy2 = 0.0;
3667 Double_t sumx = 0.0;
3668 Double_t sumxy = 0.0;
3669 Double_t sumx2 = 0.0;
3670 Double_t mpads = 0.0;
3676 Double_t moffset[10]; // Mean offset
3677 Double_t mean[10]; // Mean value
3678 Double_t angle[10]; // Angle
3680 Double_t smoffset[10]; // Sigma of mean offset
3681 Double_t smean[10]; // Sigma of mean value
3682 Double_t sangle[10]; // Sigma of angle
3683 Double_t smeanangle[10]; // Correlation
3685 Double_t sigmas[10];
3686 Double_t tchi2s[10]; // Chi2s for tracklet
3688 for (Int_t it = 0; it < 10; it++) {
3694 moffset[it] = 0.0; // Mean offset
3695 mean[it] = 0.0; // Mean value
3696 angle[it] = 0.0; // Angle
3698 smoffset[it] = 1.0e5; // Sigma of mean offset
3699 smean[it] = 1.0e5; // Sigma of mean value
3700 sangle[it] = 1.0e5; // Sigma of angle
3701 smeanangle[it] = 0.0; // Correlation
3704 tchi2s[it] = 1.0e5; // Chi2s for tracklet
3711 for (Int_t it = 0; it < t1 - t0; it++) {
3715 for (Int_t dt = -3; dt <= 3; dt++) {
3719 if (it+dt > t1-t0) {
3722 if (!cl[0][it+dt]) {
3725 zmean[it] += cl[0][it+dt]->GetZ();
3728 zmean[it] /= nmean[it];
3731 for (Int_t it = 0; it < t1 - t0; it++) {
3735 for (Int_t ih = 0; ih < 10; ih++) {
3736 dz[ih][it] = -100.0;
3737 dy[ih][it] = -100.0;
3741 Double_t xcluster = cl[ih][it]->GetX();
3744 track2.GetProlongation(xcluster,ytrack,ztrack );
3745 dz[ih][it] = cl[ih][it]->GetZ()- ztrack; // Calculate distance from track in z
3746 dy[ih][it] = cl[ih][it]->GetY() + dz[ih][it]*h01 - ytrack; // and in y
3753 if ((TMath::Abs(cl[0][it]->GetZ()-zmean[it]) > padlength * 0.8) &&
3755 if (TMath::Abs(cl[1][it]->GetZ()-zmean[it]) < padlength * 0.5) {
3763 // Iterative choice of "best path"
3765 Int_t label = TMath::Abs(track->GetLabel());
3768 for (Int_t iter = 0; iter < 9; iter++) {
3783 for (Int_t it = 0; it < t1 - t0; it++) {
3785 if (!cl[best[iter][it]][it]) {
3789 // Calculates pad-row changes
3790 Double_t zbefore = cl[best[iter][it]][it]->GetZ();
3791 Double_t zafter = cl[best[iter][it]][it]->GetZ();
3792 for (Int_t itd = it - 1; itd >= 0; itd--) {
3793 if (cl[best[iter][itd]][itd]) {
3794 zbefore = cl[best[iter][itd]][itd]->GetZ();
3798 for (Int_t itd = it + 1; itd < t1 - t0; itd++) {
3799 if (cl[best[iter][itd]][itd]) {
3800 zafter = cl[best[iter][itd]][itd]->GetZ();
3804 if ((TMath::Abs(cl[best[iter][it]][it]->GetZ()-zbefore) > 0.1) &&
3805 (TMath::Abs(cl[best[iter][it]][it]->GetZ()- zafter) > 0.1)) {
3809 Double_t dx = x[it]-xmean; // Distance to reference x
3810 sumz += cl[best[iter][it]][it]->GetZ();
3812 sumdy += dy[best[iter][it]][it];
3813 sumdy2 += dy[best[iter][it]][it]*dy[best[iter][it]][it];
3816 sumxy += dx*dy[best[iter][it]][it];
3817 mpads += cl[best[iter][it]][it]->GetNPads();
3818 if ((cl[best[iter][it]][it]->GetLabel(0) == label) ||
3819 (cl[best[iter][it]][it]->GetLabel(1) == label) ||
3820 (cl[best[iter][it]][it]->GetLabel(2) == label)) {
3830 // calculates line parameters
3832 Double_t det = sum*sumx2 - sumx*sumx;
3833 angle[iter] = (sum*sumxy - sumx*sumdy) / det;
3834 mean[iter] = (sumx2*sumdy - sumx*sumxy) / det;
3835 meanz[iter] = sumz / sum;
3836 moffset[iter] = sumdy / sum;
3837 mpads /= sum; // Mean number of pads
3839 Double_t sigma2 = 0.0; // Normalized residuals - for line fit
3840 Double_t sigma1 = 0.0; // Normalized residuals - constant fit
3842 for (Int_t it = 0; it < t1 - t0; it++) {
3843 if (!cl[best[iter][it]][it]) {
3846 Double_t dx = x[it] - xmean;
3847 Double_t ytr = mean[iter] + angle[iter] * dx;
3848 sigma2 += (dy[best[iter][it]][it] - ytr)
3849 * (dy[best[iter][it]][it] - ytr);
3850 sigma1 += (dy[best[iter][it]][it] - moffset[iter])
3851 * (dy[best[iter][it]][it] - moffset[iter]);
3854 sigma2 /= (sum - 2); // Normalized residuals
3855 sigma1 /= (sum - 1); // Normalized residuals
3856 smean[iter] = sigma2 * (sumx2 / det); // Estimated error2 of mean
3857 sangle[iter] = sigma2 * ( sum / det); // Estimated error2 of angle
3858 smeanangle[iter] = sigma2 * (-sumx / det); // Correlation
3859 sigmas[iter] = TMath::Sqrt(sigma1);
3860 smoffset[iter] = (sigma1 / sum) + 0.01*0.01; // Sigma of mean offset + unisochronity sigma
3863 // Iterative choice of "better path"
3865 for (Int_t it = 0; it < t1 - t0; it++) {
3867 if (!cl[best[iter][it]][it]) {
3871 // Add unisochronity + angular effect contribution
3872 Double_t sigmatr2 = smoffset[iter] + 0.5*tany*tany;
3873 Double_t sweight = 1.0/sigmatr2 + 1.0/track->GetSigmaY2();
3874 Double_t weighty = (moffset[iter] / sigmatr2) / sweight; // Weighted mean
3875 Double_t sigmacl = TMath::Sqrt(sigma1*sigma1 + track->GetSigmaY2());
3876 Double_t mindist = 100000.0;
3879 for (Int_t ih = 0; ih < 10; ih++) {
3883 Double_t dist2 = (dy[ih][it] - weighty) / sigmacl;
3884 dist2 *= dist2; // Chi2 distance
3885 if (dist2 < mindist) {
3890 best[iter+1][it] = ihbest;
3894 // Update best hypothesy if better chi2 according tracklet position and angle
3896 Double_t sy2 = smean[iter] + track->GetSigmaY2();
3897 Double_t sa2 = sangle[iter] + track->GetSigmaSnp2(); // track->fCee;
3898 Double_t say = track->GetSigmaSnpY(); // track->fCey;
3899 //Double_t chi20 = mean[bestiter]*mean[bestiter ] / sy2+angle[bestiter]*angle[bestiter]/sa2;
3900 //Double_t chi21 = mean[iter]*mean[iter] / sy2+angle[iter]*angle[iter]/sa2;
3902 Double_t detchi = sy2*sa2 - say*say;
3903 Double_t invers[3] = {sa2/detchi,sy2/detchi,-say/detchi}; // Inverse value of covariance matrix
3905 Double_t chi20 = mean[bestiter] * mean[bestiter] * invers[0]
3906 + angle[bestiter] * angle[bestiter] * invers[1]
3907 + 2.0 * mean[bestiter] * angle[bestiter] * invers[2];
3908 Double_t chi21 = mean[iter] * mean[iter] * invers[0]
3909 + angle[iter] * angle[iter] * invers[1]
3910 + 2.0 * mean[iter] * angle[iter] * invers[2];
3911 tchi2s[iter] = chi21;
3913 if ((changes[iter] <= changes[bestiter]) &&
3923 Double_t sigma2 = sigmas[0]; // Choose as sigma from 0 iteration
3924 Short_t maxpos = -1;
3925 Float_t maxcharge = 0.0;
3926 Short_t maxpos4 = -1;
3927 Float_t maxcharge4 = 0.0;
3928 Short_t maxpos5 = -1;
3929 Float_t maxcharge5 = 0.0;
3931 //if (tchi2s[bestiter]>25.) sigma2*=tchi2s[bestiter]/25.;
3932 //if (tchi2s[bestiter]>25.) sigma2=1000.; // dont'accept
3934 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(AliTRDcalibDB::Instance()->GetVdrift(0,0,0)
3935 ,-AliTracker::GetBz()*0.1);
3936 Double_t expectederr = sigma2*sigma2 + 0.01*0.01;
3938 expectederr += (mpads - 3.5) * 0.04;
3940 if (changes[bestiter] > 1) {
3941 expectederr += changes[bestiter] * 0.01;
3943 expectederr += (0.03 * (tany-exB)*(tany-exB)) * 15.0;
3944 //if (tchi2s[bestiter]>18.) expectederr*= tchi2s[bestiter]/18.;
3945 //expectederr+=10000;
3947 for (Int_t it = 0; it < t1 - t0; it++) {
3949 if (!cl[best[bestiter][it]][it]) {
3953 cl[best[bestiter][it]][it]->SetSigmaY2(expectederr); // Set cluster error
3954 if (!cl[best[bestiter][it]][it]->IsUsed()) {
3955 cl[best[bestiter][it]][it]->SetY(cl[best[bestiter][it]][it]->GetY());
3956 //cl[best[bestiter][it]][it]->Use();
3959 // Time bins with maximal charge
3960 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge) {
3961 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3962 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3965 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge4) {
3966 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 4) {
3967 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3968 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3972 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge5) {
3973 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 5) {
3974 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3975 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3979 // Time bins with maximal charge
3980 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge) {
3981 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3982 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3985 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge4) {
3986 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 4) {
3987 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3988 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3992 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge5) {
3993 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 5) {
3994 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3995 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3999 clusters[it+t0] = indexes[best[bestiter][it]][it];
4001 // Still needed ????
4002 //if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>4 &&
4003 //cl[best[bestiter][it]][it]->GetLocalTimeBin()<18) clusters[it+t0]
4004 // = indexes[best[bestiter][it]][it]; //Test
4009 // Set tracklet parameters
4011 Double_t trackleterr2 = smoffset[bestiter] + 0.01*0.01;
4013 trackleterr2 += (mpads - 3.5) * 0.04;
4015 trackleterr2 += changes[bestiter] * 0.01;
4016 trackleterr2 *= TMath::Max(14.0 - nfound,1.0);
4017 trackleterr2 += 0.2 * (tany-exB)*(tany-exB);
4019 // Set tracklet parameters
4021 ,track2.GetY() + moffset[bestiter]
4025 tracklet.SetTilt(h01);
4026 tracklet.SetP0(mean[bestiter]);
4027 tracklet.SetP1(angle[bestiter]);
4028 tracklet.SetN(nfound);
4029 tracklet.SetNCross(changes[bestiter]);
4030 tracklet.SetPlane(plane);
4031 tracklet.SetSigma2(expectederr);
4032 tracklet.SetChi2(tchi2s[bestiter]);
4033 tracklet.SetMaxPos(maxpos,maxpos4,maxpos5);
4034 track->SetTracklets(plane,tracklet);
4035 track->SetNWrong(track->GetNWrong() + nbad[0]);
4040 TClonesArray array0("AliTRDcluster");
4041 TClonesArray array1("AliTRDcluster");
4042 array0.ExpandCreateFast(t1 - t0 + 1);
4043 array1.ExpandCreateFast(t1 - t0 + 1);
4044 TTreeSRedirector &cstream = *fDebugStreamer;
4045 AliTRDcluster dummy;
4049 for (Int_t it = 0; it < t1 - t0; it++) {
4050 dy0[it] = dy[0][it];
4051 dyb[it] = dy[best[bestiter][it]][it];
4053 new(array0[it]) AliTRDcluster(*cl[0][it]);
4056 new(array0[it]) AliTRDcluster(dummy);
4058 if(cl[best[bestiter][it]][it]) {
4059 new(array1[it]) AliTRDcluster(*cl[best[bestiter][it]][it]);
4062 new(array1[it]) AliTRDcluster(dummy);
4065 TGraph graph0(t1-t0,x,dy0);
4066 TGraph graph1(t1-t0,x,dyb);
4067 TGraph graphy(t1-t0,x,yt);
4068 TGraph graphz(t1-t0,x,zt);
4070 if (AliTRDReconstructor::StreamLevel() > 0) {
4071 cstream << "tracklet"
4072 << "track.=" << track // Track parameters
4073 << "tany=" << tany // Tangent of the local track angle
4074 << "xmean=" << xmean // Xmean - reference x of tracklet
4075 << "tilt=" << h01 // Tilt angle
4076 << "nall=" << nall // Number of foundable clusters
4077 << "nfound=" << nfound // Number of found clusters
4078 << "clfound=" << clfound // Total number of found clusters in road
4079 << "mpads=" << mpads // Mean number of pads per cluster
4080 << "plane=" << plane // Plane number
4081 << "detector=" << detector // Detector number
4082 << "road=" << road // The width of the used road
4083 << "graph0.=" << &graph0 // x - y = dy for closest cluster
4084 << "graph1.=" << &graph1 // x - y = dy for second closest cluster
4085 << "graphy.=" << &graphy // y position of the track
4086 << "graphz.=" << &graphz // z position of the track
4087 //<< "fCl.=" << &array0 // closest cluster
4088 //<< "fCl2.=" << &array1 // second closest cluster
4089 << "maxpos=" << maxpos // Maximal charge postion
4090 << "maxcharge=" << maxcharge // Maximal charge
4091 << "maxpos4=" << maxpos4 // Maximal charge postion - after bin 4
4092 << "maxcharge4=" << maxcharge4 // Maximal charge - after bin 4
4093 << "maxpos5=" << maxpos5 // Maximal charge postion - after bin 5
4094 << "maxcharge5=" << maxcharge5 // Maximal charge - after bin 5
4095 << "bestiter=" << bestiter // Best iteration number
4096 << "tracklet.=" << &tracklet // Corrspond to the best iteration
4097 << "tchi20=" << tchi2s[0] // Chi2 of cluster in the 0 iteration
4098 << "tchi2b=" << tchi2s[bestiter] // Chi2 of cluster in the best iteration
4099 << "sigmas0=" << sigmas[0] // Residuals sigma
4100 << "sigmasb=" << sigmas[bestiter] // Residulas sigma
4101 << "ngood0=" << ngood[0] // Number of good clusters in 0 iteration
4102 << "nbad0=" << nbad[0] // Number of bad clusters in 0 iteration
4103 << "ngoodb=" << ngood[bestiter] // in best iteration
4104 << "nbadb=" << nbad[bestiter] // in best iteration
4105 << "changes0=" << changes[0] // Changes of pardrows in iteration number 0
4106 << "changesb=" << changes[bestiter] // Changes of pardrows in best iteration
4107 << "moffset0=" << moffset[0] // Offset fixing angle in iter=0
4108 << "smoffset0=" << smoffset[0] // Sigma of offset fixing angle in iter=0
4109 << "moffsetb=" << moffset[bestiter] // Offset fixing angle in iter=best
4110 << "smoffsetb=" << smoffset[bestiter] // Sigma of offset fixing angle in iter=best
4111 << "mean0=" << mean[0] // Mean dy in iter=0;
4112 << "smean0=" << smean[0] // Sigma of mean dy in iter=0
4113 << "meanb=" << mean[bestiter] // Mean dy in iter=best
4114 << "smeanb=" << smean[bestiter] // Sigma of mean dy in iter=best
4115 << "angle0=" << angle[0] // Angle deviation in the iteration number 0
4116 << "sangle0=" << sangle[0] // Sigma of angular deviation in iteration number 0
4117 << "angleb=" << angle[bestiter] // Angle deviation in the best iteration
4118 << "sangleb=" << sangle[bestiter] // Sigma of angle deviation in the best iteration
4119 << "expectederr=" << expectederr // Expected error of cluster position
4127 //_____________________________________________________________________________
4128 Int_t AliTRDtracker::Freq(Int_t n, const Int_t *inlist
4129 , Int_t *outlist, Bool_t down)
4132 // Sort eleements according occurancy
4133 // The size of output array has is 2*n
4136 Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
4137 Int_t *sindexF = new Int_t[2*n];
4138 for (Int_t i = 0; i < n; i++) {
4142 TMath::Sort(n,inlist,sindexS,down);
4144 Int_t last = inlist[sindexS[0]];
4147 sindexF[0+n] = last;
4151 for (Int_t i = 1; i < n; i++) {
4152 val = inlist[sindexS[i]];
4154 sindexF[countPos]++;
4158 sindexF[countPos+n] = val;
4159 sindexF[countPos]++;
4167 // Sort according frequency
4168 TMath::Sort(countPos,sindexF,sindexS,kTRUE);
4170 for (Int_t i = 0; i < countPos; i++) {
4171 outlist[2*i ] = sindexF[sindexS[i]+n];
4172 outlist[2*i+1] = sindexF[sindexS[i]];
4182 //_____________________________________________________________________________
4183 AliTRDtrack *AliTRDtracker::RegisterSeed(AliTRDseed *seeds, Double_t *params)
4189 Double_t alpha = AliTRDgeometry::GetAlpha();
4190 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
4194 c[ 1] = 0.0; c[ 2] = 2.0;
4195 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
4196 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
4197 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
4200 AliTRDcluster *cl = 0;
4202 for (Int_t ilayer = 0; ilayer < 6; ilayer++) {
4203 if (seeds[ilayer].IsOK()) {
4204 for (Int_t itime = 22; itime > 0; itime--) {
4205 if (seeds[ilayer].GetIndexes(itime) > 0) {
4206 index = seeds[ilayer].GetIndexes(itime);
4207 cl = seeds[ilayer].GetClusters(itime);
4220 AliTRDtrack *track = new AliTRDtrack(cl
4225 ,params[6]*alpha+shift);
4226 // SetCluster(cl, 0); // A. Bercuci
4227 track->PropagateTo(params[0]-5.0);
4228 track->ResetCovariance(1);
4230 Int_t rc = FollowBackProlongation(*track);
4237 track->CookdEdxTimBin();
4238 CookLabel(track,0.9);
4244 //////////////////////////////////////////////////////////////////////////////////////////
4246 void AliTRDtracker::InitLogHists() {
4248 fHBackfit = new TH1D("logTRD_backfit", "", 40, -0.5, 39.5);
4249 fHRefit = new TH1D("logTRD_refit", "", 40, -0.5, 39.5);
4250 fHClSearch = new TH1D("logTRD_clSearch", "", 60, -0.5, 59.5);
4252 fHX = new TH1D("logTRD_X", ";x (cm)", 200, 50, 400);
4253 fHNCl = new TH1D("logTRD_ncl", "", 40, -0.5, 39.5);
4254 fHNClTrack = new TH1D("logTRD_nclTrack", "", 180, -0.5, 179.5);
4256 fHMinYPos = new TH1D("logTRD_minYPos", ";#delta Y (cm)", 400, -6, 6);
4257 fHMinYNeg = new TH1D("logTRD_minYNeg", ";#delta Y (cm)", 400, -6, 6);
4258 fHMinZ = new TH1D("logTRD_minZ", ";#delta Z (cm)", 400, -20, 20);
4259 fHMinD = new TH2D("logTRD_minD", ";#delta Y (cm);#delta Z (cm)", 100, -6, 6, 100, -50, 50);
4261 fHDeltaX = new TH1D("logTRD_deltaX", ";#delta X (cm)", 100, -5, 5);
4262 fHXCl = new TH1D("logTRD_xCl", ";cluster x position (cm)", 1000, 280, 380);
4264 const char *nameFindCl[4] = {"logTRD_clY", "logTRD_clZ", "logTRD_clB", "logTRD_clG"};
4266 for(int i=0; i<4; i++) {
4267 fHFindCl[i] = new TH1D(nameFindCl[i], "", 30, -0.5, 29.5);
4271 //////////////////////////////////////////////////////////////////////////////////////////
4273 void AliTRDtracker::SaveLogHists() {
4275 TDirectory *sav = gDirectory;
4278 TSeqCollection *col = gROOT->GetListOfFiles();
4279 int N = col->GetEntries();
4280 for(int i=0; i<N; i++) {
4281 logFile = (TFile*)col->At(i);
4282 if (strstr(logFile->GetName(), "AliESDs.root")) break;
4286 fHBackfit->Write(fHBackfit->GetName(), TObject::kOverwrite);
4287 fHRefit->Write(fHRefit->GetName(), TObject::kOverwrite);
4288 fHClSearch->Write(fHClSearch->GetName(), TObject::kOverwrite);
4289 fHX->Write(fHX->GetName(), TObject::kOverwrite);
4290 fHNCl->Write(fHNCl->GetName(), TObject::kOverwrite);
4291 fHNClTrack->Write(fHNClTrack->GetName(), TObject::kOverwrite);
4293 fHMinYPos->Write(fHMinYPos->GetName(), TObject::kOverwrite);
4294 fHMinYNeg->Write(fHMinYNeg->GetName(), TObject::kOverwrite);
4295 fHMinD->Write(fHMinD->GetName(), TObject::kOverwrite);
4296 fHMinZ->Write(fHMinZ->GetName(), TObject::kOverwrite);
4298 fHDeltaX->Write(fHDeltaX->GetName(), TObject::kOverwrite);
4299 fHXCl->Write(fHXCl->GetName(), TObject::kOverwrite);
4302 for(int i=0; i<4; i++)
4303 fHFindCl[i]->Write(fHFindCl[i]->GetName(), TObject::kOverwrite);
4310 //////////////////////////////////////////////////////////////////////////////////////////