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 ///////////////////////////////////////////////////////////////////////////////
32 #include <TLinearFitter.h>
33 #include <TObjArray.h>
36 #include <TTreeStream.h>
38 #include "AliESDEvent.h"
39 #include "AliESDtrack.h"
40 #include "AliAlignObj.h"
41 #include "AliRieman.h"
42 #include "AliTrackPointArray.h"
44 #include "AliTRDgeometry.h"
45 #include "AliTRDpadPlane.h"
46 #include "AliTRDgeometry.h"
47 #include "AliTRDcluster.h"
48 #include "AliTRDtrack.h"
49 #include "AliTRDseed.h"
50 #include "AliTRDcalibDB.h"
51 #include "AliTRDCommonParam.h"
52 #include "AliTRDtracker.h"
53 #include "AliTRDReconstructor.h"
54 #include "AliTRDCalibraFillHisto.h"
56 ClassImp(AliTRDtracker)
58 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5; //
59 const Float_t AliTRDtracker::fgkLabelFraction = 0.8; //
60 const Double_t AliTRDtracker::fgkMaxChi2 = 12.0; //
61 const Double_t AliTRDtracker::fgkMaxSnp = 0.95; // Maximum local sine of the azimuthal angle
62 const Double_t AliTRDtracker::fgkMaxStep = 2.0; // Maximal step size in propagation
64 //_____________________________________________________________________________
65 AliTRDtracker::AliTRDtracker()
92 // Default constructor
95 for (Int_t i = 0; i < kTrackingSectors; i++) {
103 //_____________________________________________________________________________
104 AliTRDtracker::AliTRDtracker(const AliTRDtracker &t)
113 ,fTimeBinsPerPlane(0)
114 ,fAddTRDseeds(kFALSE)
136 //_____________________________________________________________________________
137 AliTRDtracker::AliTRDtracker(const TFile */*geomfile*/)
141 ,fClusters(new TObjArray(2000))
143 ,fSeeds(new TObjArray(2000))
145 ,fTracks(new TObjArray(1000))
146 ,fTimeBinsPerPlane(0)
147 ,fAddTRDseeds(kFALSE)
167 TDirectory *savedir = gDirectory;
169 fGeom = new AliTRDgeometry();
171 for (Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
173 fTrSec[trS] = new AliTRDtrackingSector(fGeom,geomS);
176 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(0,0);
177 Float_t tiltAngle = TMath::Abs(padPlane->GetTiltingAngle());
178 if (tiltAngle < 0.1) {
182 if (!AliTRDcalibDB::Instance()) {
183 AliFatal("Could not get calibration object");
185 fTimeBinsPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
187 fDebugStreamer = new TTreeSRedirector("TRDdebug.root");
195 //_____________________________________________________________________________
196 AliTRDtracker::~AliTRDtracker()
199 // Destructor of AliTRDtracker
221 for (Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
222 delete fTrSec[geomS];
225 if (fDebugStreamer) {
226 delete fDebugStreamer;
231 //_____________________________________________________________________________
232 Int_t AliTRDtracker::LocalToGlobalID(Int_t lid)
235 // Transform internal TRD ID to global detector ID
238 Int_t isector = fGeom->GetSector(lid);
239 Int_t ichamber = fGeom->GetChamber(lid);
240 Int_t iplan = fGeom->GetPlane(lid);
242 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
245 iLayer = AliGeomManager::kTRD1;
248 iLayer = AliGeomManager::kTRD2;
251 iLayer = AliGeomManager::kTRD3;
254 iLayer = AliGeomManager::kTRD4;
257 iLayer = AliGeomManager::kTRD5;
260 iLayer = AliGeomManager::kTRD6;
264 Int_t modId = isector * fGeom->Ncham() + ichamber;
265 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer,modId);
271 //_____________________________________________________________________________
272 Int_t AliTRDtracker::GlobalToLocalID(Int_t gid)
275 // Transform global detector ID to local detector ID
279 AliGeomManager::ELayerID layerId = AliGeomManager::VolUIDToLayer(gid,modId);
281 Int_t isector = modId / fGeom->Ncham();
282 Int_t ichamber = modId % fGeom->Ncham();
286 case AliGeomManager::kTRD1:
289 case AliGeomManager::kTRD2:
292 case AliGeomManager::kTRD3:
295 case AliGeomManager::kTRD4:
298 case AliGeomManager::kTRD5:
301 case AliGeomManager::kTRD6:
312 Int_t lid = fGeom->GetDetector(iLayer,ichamber,isector);
318 //_____________________________________________________________________________
319 Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track)
322 // Rotates the track when necessary
325 Double_t alpha = AliTRDgeometry::GetAlpha();
326 Double_t y = track->GetY();
327 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
330 if (!track->Rotate( alpha)) {
334 else if (y < -ymax) {
335 if (!track->Rotate(-alpha)) {
344 //_____________________________________________________________________________
345 AliTRDcluster *AliTRDtracker::GetCluster(AliTRDtrack *track, Int_t plane
346 , Int_t timebin, UInt_t &index)
349 // Try to find cluster in the backup list
352 AliTRDcluster *cl =0;
353 Int_t *indexes = track->GetBackupIndexes();
355 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
356 if (indexes[i] == 0) {
359 AliTRDcluster *cli = (AliTRDcluster *) fClusters->UncheckedAt(indexes[i]);
363 if (cli->GetLocalTimeBin() != timebin) {
366 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
367 if (iplane == plane) {
378 //_____________________________________________________________________________
379 Int_t AliTRDtracker::GetLastPlane(AliTRDtrack *track)
382 // Return last updated plane
386 Int_t *indexes = track->GetBackupIndexes();
388 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
389 AliTRDcluster *cli = (AliTRDcluster *) fClusters->UncheckedAt(indexes[i]);
393 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
394 if (iplane > lastplane) {
403 //_____________________________________________________________________________
404 Int_t AliTRDtracker::PropagateBack(AliESDEvent *event)
407 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
408 // backpropagated by the TPC tracker. Each seed is first propagated
409 // to the TRD, and then its prolongation is searched in the TRD.
410 // If sufficiently long continuation of the track is found in the TRD
411 // the track is updated, otherwise it's stored as originaly defined
412 // by the TPC tracker.
415 Int_t found = 0; // number of tracks found
416 Float_t foundMin = 20.0;
418 Int_t nSeed = event->GetNumberOfTracks();
420 // run stand alone tracking
421 if (AliTRDReconstructor::SeedingOn()) Clusters2Tracks(event);
425 Float_t *quality = new Float_t[nSeed];
426 Int_t *index = new Int_t[nSeed];
427 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
428 AliESDtrack *seed = event->GetTrack(iSeed);
429 Double_t covariance[15];
430 seed->GetExternalCovariance(covariance);
431 quality[iSeed] = covariance[0] + covariance[2];
433 // Sort tracks according to covariance of local Y and Z
434 TMath::Sort(nSeed,quality,index,kFALSE);
436 // Backpropagate all seeds
437 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
439 // Get the seeds in sorted sequence
440 AliESDtrack *seed = event->GetTrack(index[iSeed]);
441 fHBackfit->Fill(0); // All seeds
443 // Check the seed status
444 ULong_t status = seed->GetStatus();
445 if ((status & AliESDtrack::kTPCout) == 0) {
446 fHBackfit->Fill(1); // TPC outer edge reached
449 if ((status & AliESDtrack::kTRDout) != 0) {
450 fHBackfit->Fill(2); // TRD outer edge reached (does this happen ?)
454 // Do the back prolongation
455 Int_t lbl = seed->GetLabel();
456 AliTRDtrack *track = new AliTRDtrack(*seed);
457 track->SetSeedLabel(lbl);
458 seed->UpdateTrackParams(track,AliESDtrack::kTRDbackup); // Make backup
460 Float_t p4 = track->GetC();
461 Int_t expectedClr = FollowBackProlongation(*track);
462 fHBackfit->Fill(3); // Back prolongation done
463 fHX->Fill(track->GetX());
465 if ((TMath::Abs(track->GetC() - p4) / TMath::Abs(p4) < 0.2) ||
466 (track->Pt() > 0.8)) {
471 // Make backup for back propagation
474 Int_t foundClr = track->GetNumberOfClusters();
475 if (foundClr >= foundMin) {
477 track->CookdEdxTimBin(seed->GetID()); // A.Bercuci 25.07.07
478 CookLabel(track,1 - fgkLabelFraction);
479 if (track->GetBackupTrack()) UseClusters(track->GetBackupTrack());
482 // Sign only gold tracks
483 if (track->GetChi2() / track->GetNumberOfClusters() < 4) {
484 if ((seed->GetKinkIndex(0) == 0) &&
485 (track->Pt() < 1.5)) {
489 Bool_t isGold = kFALSE;
492 if (track->GetChi2() / track->GetNumberOfClusters() < 5) {
493 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
494 if (track->GetBackupTrack()) {
495 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
502 if ((!isGold) && (track->GetNCross() == 0) &&
503 (track->GetChi2() / track->GetNumberOfClusters() < 7)) {
504 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
505 if (track->GetBackupTrack()) {
506 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
511 if ((!isGold) && (track->GetBackupTrack())) {
512 if ((track->GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
513 seed->UpdateTrackParams(track->GetBackupTrack(),AliESDtrack::kTRDbackup);
518 if ((track->StatusForTOF() > 0) && (track->GetNCross() == 0) && (Float_t(track->GetNumberOfClusters()) / Float_t(track->GetNExpected()) > 0.4)) {
519 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
526 // Debug part of tracking
527 TTreeSRedirector &cstream = *fDebugStreamer;
528 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.
529 if (AliTRDReconstructor::StreamLevel() > 0) {
530 if (track->GetBackupTrack()) {
532 << "EventNrInFile=" << eventNrInFile
535 << "trdback.=" << track->GetBackupTrack()
540 << "EventNrInFile=" << eventNrInFile
543 << "trdback.=" << track
549 // Propagation to the TOF (I.Belikov)
550 if (track->GetStop() == kFALSE) {
553 Double_t xtof = 371.0;
554 Double_t xTOF0 = 370.0;
556 Double_t c2 = track->GetSnp() + track->GetC() * (xtof - track->GetX());
557 if (TMath::Abs(c2) >= 0.99) {
563 PropagateToX(*track,xTOF0,fgkMaxStep);
565 // Energy losses taken to the account - check one more time
566 c2 = track->GetSnp() + track->GetC() * (xtof - track->GetX());
567 if (TMath::Abs(c2) >= 0.99) {
573 //if (!PropagateToX(*track,xTOF0,fgkMaxStep)) {
574 // fHBackfit->Fill(7);
579 Double_t ymax = xtof * TMath::Tan(0.5 * AliTRDgeometry::GetAlpha());
581 track->GetYAt(xtof,GetBz(),y);
583 if (!track->Rotate( AliTRDgeometry::GetAlpha())) {
589 else if (y < -ymax) {
590 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
597 if (track->PropagateTo(xtof)) {
598 seed->UpdateTrackParams(track,AliESDtrack::kTRDout);
601 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
602 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
603 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
605 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
607 //seed->SetTRDtrack(new AliTRDtrack(*track));
608 if (track->GetNumberOfClusters() > foundMin) {
617 if ((track->GetNumberOfClusters() > 15) &&
618 (track->GetNumberOfClusters() > 0.5*expectedClr)) {
619 seed->UpdateTrackParams(track,AliESDtrack::kTRDout);
622 //seed->SetStatus(AliESDtrack::kTRDStop);
623 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
624 for (Int_t j = 0; j <AliESDtrack::kNSlice; j++) {
625 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
627 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
629 //seed->SetTRDtrack(new AliTRDtrack(*track));
634 seed->SetTRDQuality(track->StatusForTOF());
635 seed->SetTRDBudget(track->GetBudget(0));
641 AliInfo(Form("Number of seeds: %d",fNseeds));
642 AliInfo(Form("Number of back propagated TRD tracks: %d",found));
655 //_____________________________________________________________________________
656 Int_t AliTRDtracker::RefitInward(AliESDEvent *event)
659 // Refits tracks within the TRD. The ESD event is expected to contain seeds
660 // at the outer part of the TRD.
661 // The tracks are propagated to the innermost time bin
662 // of the TRD and the ESD event is updated
663 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
666 //Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
667 //Float_t foundMin = fgkMinClustersInTrack * timeBins;
671 //Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
674 // Calibration fill 2D
675 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
677 AliInfo("Could not get Calibra instance\n");
680 Int_t n = event->GetNumberOfTracks();
681 for (Int_t i = 0; i < n; i++) {
683 AliESDtrack *seed = event->GetTrack(i);
684 new (&seed2) AliTRDtrack(*seed);
687 if (seed2.GetX() < 270.0) {
688 seed->UpdateTrackParams(&seed2,AliESDtrack::kTRDbackup); // Backup TPC track - only update
693 ULong_t status = seed->GetStatus();
694 if ((status & AliESDtrack::kTRDout) == 0) {
698 if ((status & AliESDtrack::kTRDin) != 0) {
706 seed2.ResetCovariance(50.0);
708 AliTRDtrack *pt = new AliTRDtrack(seed2,seed2.GetAlpha());
709 Int_t *indexes2 = seed2.GetIndexes();
710 for (Int_t i = 0; i < AliESDtrack::kNPlane;i++) {
711 for (Int_t j = 0; j < AliESDtrack::kNSlice;j++) {
712 pt->SetPIDsignals(seed2.GetPIDsignals(i,j),i,j);
714 pt->SetPIDTimBin(seed2.GetPIDTimBin(i),i);
717 Int_t *indexes3 = pt->GetBackupIndexes();
718 for (Int_t i = 0; i < 200;i++) {
719 if (indexes2[i] == 0) {
722 indexes3[i] = indexes2[i];
725 FollowProlongation(*pt);
727 pt->CookdEdxTimBin(seed->GetID());
728 pt->SetPIDMethod(AliTRDtrack::kLQ); //switch between TRD PID methods
730 if(calibra->GetHisto2d()) calibra->UpdateHistograms(pt);
733 Double_t xTPC = 250.0;
734 if (PropagateToX(*pt,xTPC,fgkMaxStep)) {
736 seed->UpdateTrackParams(pt,AliESDtrack::kTRDrefit);
739 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
740 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
741 seed->SetTRDsignals(pt->GetPIDsignals(i,j),i,j);
743 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
749 // If not prolongation to TPC - propagate without update
751 AliTRDtrack *seed2 = new AliTRDtrack(*seed);
752 seed2->ResetCovariance(5.0);
753 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
756 if (PropagateToX(*pt2,xTPC,fgkMaxStep)) {
759 pt2->CookdEdxTimBin(seed->GetID());
760 seed->UpdateTrackParams(pt2,AliESDtrack::kTRDrefit);
763 for (Int_t i = 0; i < AliESDtrack::kNPlane; i++) {
764 for (Int_t j = 0; j < AliESDtrack::kNSlice; j++) {
765 seed->SetTRDsignals(pt2->GetPIDsignals(i,j),i,j);
767 seed->SetTRDTimBin(pt2->GetPIDTimBin(i),i);
772 // Add TRD track to ESDfriendTrack - maybe this tracks are
773 // not useful for post-processing - TODO make decision
774 if (AliTRDReconstructor::StreamLevel() > 0) {
775 seed->AddCalibObject(new AliTRDtrack(*pt2/*, kTRUE*/));
781 // Add TRD track to ESDfriendTrack
782 if (AliTRDReconstructor::StreamLevel() > 0) {
783 seed->AddCalibObject(new AliTRDtrack(*pt/*, kTRUE*/));
789 AliInfo(Form("Number of loaded seeds: %d",nseed));
790 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
798 //_____________________________________________________________________________
799 Int_t AliTRDtracker::FollowProlongation(AliTRDtrack &t)
802 // Starting from current position on track=t this function tries
803 // to extrapolate the track up to timeBin=0 and to confirm prolongation
804 // if a close cluster is found. Returns the number of clusters
805 // expected to be found in sensitive layers
806 // GeoManager used to estimate mean density
810 Int_t lastplane = GetLastPlane(&t);
813 Int_t expectedNumberOfClusters = 0;
815 for (Int_t iplane = lastplane; iplane >= 0; iplane--) {
817 Int_t row0 = GetGlobalTimeBin(0,iplane,GetTimeBinsPerPlane()-1);
818 Int_t rowlast = GetGlobalTimeBin(0,iplane,0);
820 // Propagate track close to the plane if neccessary
821 Double_t currentx = fTrSec[0]->GetLayer(rowlast)->GetX();
822 if (currentx < (-fgkMaxStep + t.GetX())) {
823 // Propagate closer to chamber - safety space fgkMaxStep
824 if (!PropagateToX(t,currentx+fgkMaxStep,fgkMaxStep)) {
829 if (!AdjustSector(&t)) {
833 // Get material budget
841 // Starting global position
843 // End global position
844 x = fTrSec[0]->GetLayer(row0)->GetX();
845 if (!t.GetProlongation(x,y,z)) {
848 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
849 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
851 AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
852 xrho= param[0]*param[4];
853 xx0 = param[1]; // Get mean propagation parameters
855 // Flags for marking the track momentum and direction. The track is
856 // marked only if it has at least 1 cluster picked up in the current
858 Bool_t kUPDATED = kFALSE;
859 Bool_t kMARKED = kFALSE;
861 // Propagate and update
862 sector = t.GetSector();
863 //for (Int_t itime=GetTimeBinsPerPlane()-1;itime>=0;itime--) {
864 for (Int_t itime = 0 ; itime < GetTimeBinsPerPlane(); itime++) {
866 // Mark track kinematics
867 if (itime > 10 && kUPDATED && !kMARKED) {
868 t.SetTrackSegmentDirMom(iplane);
872 Int_t ilayer = GetGlobalTimeBin(0,iplane,itime);
873 expectedNumberOfClusters++;
874 t.SetNExpected(t.GetNExpected() + 1);
875 if (t.GetX() > 345.0) {
876 t.SetNExpectedLast(t.GetNExpectedLast() + 1);
878 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(ilayer));
879 AliTRDcluster *cl = 0;
881 Double_t maxChi2 = fgkMaxChi2;
886 AliTRDcluster *cl0 = timeBin[0];
888 // No clusters in given time bin
892 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
893 if (plane > lastplane) {
897 Int_t timebin = cl0->GetLocalTimeBin();
898 AliTRDcluster *cl2 = GetCluster(&t,plane,timebin,index);
902 //Double_t h01 = GetTiltFactor(cl); //I.B's fix
903 //maxChi2=t.GetPredictedChi2(cl,h01);
908 //if (cl->GetNPads()<5)
909 Double_t dxsample = timeBin.GetdX();
910 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
911 Double_t h01 = GetTiltFactor(cl);
912 Int_t det = cl->GetDetector();
913 Int_t plane = fGeom->GetPlane(det);
915 if (t.GetX() > 345.0) {
916 t.SetNLast(t.GetNLast() + 1);
917 t.SetChi2Last(t.GetChi2Last() + maxChi2);
920 Double_t xcluster = cl->GetX();
921 t.PropagateTo(xcluster,xx0,xrho);
922 if (!AdjustSector(&t)) {
926 maxChi2 = t.GetPredictedChi2(cl,h01);
927 if (maxChi2 < 1e+10) {
928 if (!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
932 //SetCluster(cl, GetNumberOfClusters()-1);
945 return expectedNumberOfClusters;
949 //_____________________________________________________________________________
950 Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack &t)
953 // Starting from current radial position of track <t> this function
954 // extrapolates the track up to the outer timebin and in the sensitive
955 // layers confirms prolongation if a close cluster is found.
956 // Returns the number of clusters expected to be found in sensitive layers
957 // Uses the geomanager for material description
959 // return number of assigned clusters ?
967 Float_t ratio0 = 0.0;
969 Int_t expectedNumberOfClusters = 0;
971 AliTRDtracklet tracklet;
973 const Int_t kNclusters = 1000;
974 Int_t clusters[kNclusters];
975 for (Int_t i = 0; i < kNclusters; i++) {
979 // // Calibration fill 2D
980 // AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
982 // AliInfo("Could not get Calibra instance\n");
984 // if (calibra->GetMITracking()) {
985 // calibra->ResetTrack();
988 // Loop through the TRD planes
989 for (Int_t iplane = 0; iplane < AliESDtrack::kNPlane; iplane++) {
991 Int_t hb = iplane * 10;
992 fHClSearch->Fill(hb);
994 // Get the global time bin numbers for the first an last
995 // local time bin of the given plane
996 Int_t row0 = GetGlobalTimeBin(0,iplane,GetTimeBinsPerPlane()-1);
997 Int_t rowlast = GetGlobalTimeBin(0,iplane,0);
999 // Get the X coordinates of the propagation layer for the first time bin
1000 Double_t currentx = fTrSec[0]->GetLayer(row0)->GetX();
1001 if (currentx < t.GetX()) {
1002 fHClSearch->Fill(hb+1);
1006 // Propagate closer to the current chamber if neccessary
1007 if (currentx > (fgkMaxStep + t.GetX())) {
1008 if (!PropagateToX(t,currentx-fgkMaxStep,fgkMaxStep)) {
1009 fHClSearch->Fill(hb+2);
1014 // Rotate track to adjacent sector if neccessary
1015 if (!AdjustSector(&t)) {
1016 fHClSearch->Fill(hb+3);
1020 // Check whether azimuthal angle is getting too large
1021 if (TMath::Abs(t.GetSnp()) > fgkMaxSnp) {
1022 fHClSearch->Fill(hb+4);
1032 // Global start position (beginning of chamber)
1034 // X-position of the end of the chamber
1035 x = fTrSec[0]->GetLayer(rowlast)->GetX();
1036 // Get local Y and Z at the X-position of the end of the chamber
1037 if (!t.GetProlongation(x,y,z)) {
1038 fHClSearch->Fill(hb+5);
1041 // Global end position (end of chamber)
1042 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1043 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1046 // Calculate the mean material budget along the path inside the chamber
1047 AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
1048 // The mean propagation parameters (density*length and radiation length)
1049 xrho = param[0]*param[4];
1052 // Find the clusters and tracklet along the path inside the chamber
1053 sector = t.GetSector();
1054 Float_t ncl = FindClusters(sector,row0,rowlast,&t,clusters,tracklet);
1055 fHNCl->Fill(tracklet.GetN());
1057 // Discard if in less than 1/3 of the available timebins
1058 // clusters are found
1059 if (tracklet.GetN() < GetTimeBinsPerPlane()/3) {
1060 fHClSearch->Fill(hb+6);
1065 // Propagate and update track
1067 for (Int_t itime = GetTimeBinsPerPlane()-1; itime >= 0; itime--) {
1068 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
1069 expectedNumberOfClusters++;
1070 t.SetNExpected(t.GetNExpected() + 1);
1071 if (t.GetX() > 345.0) {
1072 t.SetNExpectedLast(t.GetNExpectedLast() + 1);
1074 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(ilayer));
1075 AliTRDcluster *cl = 0;
1077 Double_t maxChi2 = fgkMaxChi2;
1081 if (clusters[ilayer] > 0) {
1082 index = clusters[ilayer];
1083 cl = (AliTRDcluster *)GetCluster(index);
1084 //Double_t h01 = GetTiltFactor(cl); // I.B's fix
1085 //maxChi2=t.GetPredictedChi2(cl,h01); //
1089 //if (cl->GetNPads() < 5)
1090 Double_t dxsample = timeBin.GetdX();
1091 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
1092 Double_t h01 = GetTiltFactor(cl);
1093 Int_t det = cl->GetDetector();
1094 Int_t plane = fGeom->GetPlane(det);
1095 if (t.GetX() > 345.0) {
1096 t.SetNLast(t.GetNLast() + 1);
1097 t.SetChi2Last(t.GetChi2Last() + maxChi2);
1099 Double_t xcluster = cl->GetX();
1100 t.PropagateTo(xcluster,xx0,xrho);
1101 maxChi2 = t.GetPredictedChi2(cl,h01);
1104 if (!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1105 if (!t.Update(cl,maxChi2,index,h01)) {
1108 } // else SetCluster(cl, GetNumberOfClusters()-1); // A.Bercuci 25.07.07
1111 // if (calibra->GetMITracking()) {
1112 // calibra->UpdateHistograms(cl,&t);
1115 // Reset material budget if 2 consecutive gold
1117 if ((t.GetTracklets(plane).GetN() + t.GetTracklets(plane-1).GetN()) > 20) {
1128 ratio0 = ncl / Float_t(fTimeBinsPerPlane);
1129 Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
1130 if ((tracklet.GetChi2() < 18.0) &&
1133 (ratio0+ratio1 > 1.5) &&
1134 (t.GetNCross() == 0) &&
1135 (TMath::Abs(t.GetSnp()) < 0.85) &&
1136 (t.GetNumberOfClusters() > 20)){
1137 //if (ratio0 > 0.8) {
1138 t.MakeBackupTrack(); // Make backup of the track until is gold
1143 return expectedNumberOfClusters;
1147 //_____________________________________________________________________________
1148 Int_t AliTRDtracker::PropagateToX(AliTRDtrack &t, Double_t xToGo, Double_t maxStep)
1151 // Starting from current X-position of track <t> this function
1152 // extrapolates the track up to radial position <xToGo>.
1153 // Returns 1 if track reaches the plane, and 0 otherwise
1156 const Double_t kEpsilon = 0.00001;
1158 // Current track X-position
1159 Double_t xpos = t.GetX();
1161 // Direction: inward or outward
1162 Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
1164 while (((xToGo - xpos) * dir) > kEpsilon) {
1173 // The next step size
1174 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1176 // Get the global position of the starting point
1179 // X-position after next step
1182 // Get local Y and Z at the X-position of the next step
1183 if (!t.GetProlongation(x,y,z)) {
1184 return 0; // No prolongation possible
1187 // The global position of the end point of this prolongation step
1188 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1189 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1192 // Calculate the mean material budget between start and
1193 // end point of this prolongation step
1194 AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
1196 // Propagate the track to the X-position after the next step
1197 if (!t.PropagateTo(x,param[1],param[0]*param[4])) {
1201 // Rotate the track if necessary
1204 // New track X-position
1213 //_____________________________________________________________________________
1214 Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1217 // Fills clusters into TRD tracking_sectors
1218 // Note that the numbering scheme for the TRD tracking_sectors
1219 // differs from that of TRD sectors
1223 if (ReadClusters(fClusters, cTree)) {
1224 AliError("Problem with reading the clusters !");
1227 Int_t ncl = fClusters->GetEntriesFast();
1229 AliInfo(Form("Sorting %d clusters",ncl));
1234 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(ncl);
1235 Int_t detector = c->GetDetector();
1236 Int_t localTimeBin = c->GetLocalTimeBin();
1237 Int_t sector = fGeom->GetSector(detector);
1238 Int_t plane = fGeom->GetPlane(detector);
1239 Int_t trackingSector = sector;
1241 //if (c->GetQ() > 10) {
1242 // Int_t chamber = fGeom->GetChamber(detector);
1245 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1249 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1253 fHXCl->Fill(c->GetX());
1255 fTrSec[trackingSector]->GetLayer(layer)->SetX(c->GetX());
1256 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
1264 //_____________________________________________________________________________
1265 void AliTRDtracker::UnloadClusters()
1268 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1274 nentr = fClusters->GetEntriesFast();
1275 for (i = 0; i < nentr; i++) {
1276 delete fClusters->RemoveAt(i);
1280 nentr = fSeeds->GetEntriesFast();
1281 for (i = 0; i < nentr; i++) {
1282 delete fSeeds->RemoveAt(i);
1285 nentr = fTracks->GetEntriesFast();
1286 for (i = 0; i < nentr; i++) {
1287 delete fTracks->RemoveAt(i);
1290 Int_t nsec = AliTRDgeometry::kNsect;
1291 for (i = 0; i < nsec; i++) {
1292 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1293 fTrSec[i]->GetLayer(pl)->Clear();
1299 //_____________________________________________________________________________
1300 Int_t AliTRDtracker::Clusters2Tracks(AliESDEvent *esd)
1303 // Creates seeds using clusters between position inner plane and outer plane
1306 const Double_t kMaxTheta = 1.0;
1307 const Double_t kMaxPhi = 2.0;
1309 const Double_t kRoad0y = 6.0; // Road for middle cluster
1310 const Double_t kRoad0z = 8.5; // Road for middle cluster
1312 const Double_t kRoad1y = 2.0; // Road in y for seeded cluster
1313 const Double_t kRoad1z = 20.0; // Road in z for seeded cluster
1315 const Double_t kRoad2y = 3.0; // Road in y for extrapolated cluster
1316 const Double_t kRoad2z = 20.0; // Road in z for extrapolated cluster
1317 const Int_t kMaxSeed = 3000;
1319 Int_t maxSec = AliTRDgeometry::kNsect;
1321 // Linear fitters in planes
1322 TLinearFitter fitterTC(2,"hyp2"); // Fitting with tilting pads - kz fixed - kz= Z/x, + vertex const
1323 TLinearFitter fitterT2(4,"hyp4"); // Fitting with tilting pads - kz not fixed
1324 fitterTC.StoreData(kTRUE);
1325 fitterT2.StoreData(kTRUE);
1326 AliRieman rieman(1000); // Rieman fitter
1327 AliRieman rieman2(1000); // Rieman fitter
1329 // Find the maximal and minimal layer for the planes
1331 AliTRDpropagationLayer *reflayers[6];
1332 for (Int_t i = 0; i < 6; i++) {
1333 layers[i][0] = 10000;
1336 for (Int_t ns = 0; ns < maxSec; ns++) {
1337 for (Int_t ilayer = 0; ilayer < fTrSec[ns]->GetNumberOfLayers(); ilayer++) {
1338 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(ilayer));
1342 Int_t det = layer[0]->GetDetector();
1343 Int_t plane = fGeom->GetPlane(det);
1344 if (ilayer < layers[plane][0]) {
1345 layers[plane][0] = ilayer;
1347 if (ilayer > layers[plane][1]) {
1348 layers[plane][1] = ilayer;
1353 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(0,0);
1354 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
1355 Double_t hL[6]; // Tilting angle
1356 Double_t xcl[6]; // X - position of reference cluster
1357 Double_t ycl[6]; // Y - position of reference cluster
1358 Double_t zcl[6]; // Z - position of reference cluster
1360 AliTRDcluster *cl[6] = { 0, 0, 0, 0, 0, 0 }; // Seeding clusters
1361 Float_t padlength[6] = { 10.0, 10.0, 10.0, 10.0, 10.0, 10.0 }; // Current pad-length
1363 Double_t chi2R = 0.0;
1364 Double_t chi2Z = 0.0;
1365 Double_t chi2RF = 0.0;
1366 Double_t chi2ZF = 0.0;
1368 Int_t nclusters; // Total number of clusters
1369 for (Int_t i = 0; i < 6; i++) {
1377 AliTRDseed *pseed = new AliTRDseed[kMaxSeed*6];
1378 AliTRDseed *seed[kMaxSeed];
1379 for (Int_t iseed = 0; iseed < kMaxSeed; iseed++) {
1380 seed[iseed]= &pseed[iseed*6];
1382 AliTRDseed *cseed = seed[0];
1384 Double_t seedquality[kMaxSeed];
1385 Double_t seedquality2[kMaxSeed];
1386 Double_t seedparams[kMaxSeed][7];
1387 Int_t seedlayer[kMaxSeed];
1388 Int_t registered = 0;
1389 Int_t sort[kMaxSeed];
1394 for (Int_t ns = 0; ns < maxSec; ns++) { // Loop over sectors
1395 //for (Int_t ns = 0; ns < 5; ns++) { // Loop over sectors
1397 registered = 0; // Reset registerd seed counter
1398 cseed = seed[registered];
1401 for (Int_t sLayer = 2; sLayer >= 0; sLayer--) {
1402 //for (Int_t dseed = 5; dseed < 15; dseed += 3) {
1405 Int_t dseed = 5 + Int_t(iter) * 3;
1407 // Initialize seeding layers
1408 for (Int_t ilayer = 0; ilayer < 6; ilayer++) {
1409 reflayers[ilayer] = fTrSec[ns]->GetLayer(layers[ilayer][1]-dseed);
1410 xcl[ilayer] = reflayers[ilayer]->GetX();
1413 Double_t xref = (xcl[sLayer+1] + xcl[sLayer+2]) * 0.5;
1414 AliTRDpropagationLayer &layer0 = *reflayers[sLayer+0];
1415 AliTRDpropagationLayer &layer1 = *reflayers[sLayer+1];
1416 AliTRDpropagationLayer &layer2 = *reflayers[sLayer+2];
1417 AliTRDpropagationLayer &layer3 = *reflayers[sLayer+3];
1419 Int_t maxn3 = layer3;
1420 for (Int_t icl3 = 0; icl3 < maxn3; icl3++) {
1422 AliTRDcluster *cl3 = layer3[icl3];
1426 padlength[sLayer+3] = TMath::Sqrt(cl3->GetSigmaZ2() * 12.0);
1427 ycl[sLayer+3] = cl3->GetY();
1428 zcl[sLayer+3] = cl3->GetZ();
1429 Float_t yymin0 = ycl[sLayer+3] - 1.0 - kMaxPhi * (xcl[sLayer+3]-xcl[sLayer+0]);
1430 Float_t yymax0 = ycl[sLayer+3] + 1.0 + kMaxPhi * (xcl[sLayer+3]-xcl[sLayer+0]);
1431 Int_t maxn0 = layer0;
1433 for (Int_t icl0 = layer0.Find(yymin0); icl0 < maxn0; icl0++) {
1435 AliTRDcluster *cl0 = layer0[icl0];
1439 if (cl3->IsUsed() && cl0->IsUsed()) {
1442 ycl[sLayer+0] = cl0->GetY();
1443 zcl[sLayer+0] = cl0->GetZ();
1444 if (ycl[sLayer+0] > yymax0) {
1447 Double_t tanphi = (ycl[sLayer+3]-ycl[sLayer+0]) / (xcl[sLayer+3]-xcl[sLayer+0]);
1448 if (TMath::Abs(tanphi) > kMaxPhi) {
1451 Double_t tantheta = (zcl[sLayer+3]-zcl[sLayer+0]) / (xcl[sLayer+3]-xcl[sLayer+0]);
1452 if (TMath::Abs(tantheta) > kMaxTheta) {
1455 padlength[sLayer+0] = TMath::Sqrt(cl0->GetSigmaZ2() * 12.0);
1457 // Expected position in 1 layer
1458 Double_t y1exp = ycl[sLayer+0] + (tanphi) * (xcl[sLayer+1]-xcl[sLayer+0]);
1459 Double_t z1exp = zcl[sLayer+0] + (tantheta) * (xcl[sLayer+1]-xcl[sLayer+0]);
1460 Float_t yymin1 = y1exp - kRoad0y - tanphi;
1461 Float_t yymax1 = y1exp + kRoad0y + tanphi;
1462 Int_t maxn1 = layer1;
1464 for (Int_t icl1 = layer1.Find(yymin1); icl1 < maxn1; icl1++) {
1466 AliTRDcluster *cl1 = layer1[icl1];
1471 if (cl3->IsUsed()) nusedCl++;
1472 if (cl0->IsUsed()) nusedCl++;
1473 if (cl1->IsUsed()) nusedCl++;
1477 ycl[sLayer+1] = cl1->GetY();
1478 zcl[sLayer+1] = cl1->GetZ();
1479 if (ycl[sLayer+1] > yymax1) {
1482 if (TMath::Abs(ycl[sLayer+1]-y1exp) > kRoad0y+tanphi) {
1485 if (TMath::Abs(zcl[sLayer+1]-z1exp) > kRoad0z) {
1488 padlength[sLayer+1] = TMath::Sqrt(cl1->GetSigmaZ2() * 12.0);
1490 Double_t y2exp = ycl[sLayer+0]+(tanphi) * (xcl[sLayer+2]-xcl[sLayer+0]) + (ycl[sLayer+1]-y1exp);
1491 Double_t z2exp = zcl[sLayer+0]+(tantheta) * (xcl[sLayer+2]-xcl[sLayer+0]);
1492 Int_t index2 = layer2.FindNearestCluster(y2exp,z2exp,kRoad1y,kRoad1z);
1496 AliTRDcluster *cl2 = (AliTRDcluster *) GetCluster(index2);
1497 padlength[sLayer+2] = TMath::Sqrt(cl2->GetSigmaZ2() * 12.0);
1498 ycl[sLayer+2] = cl2->GetY();
1499 zcl[sLayer+2] = cl2->GetZ();
1500 if (TMath::Abs(cl2->GetZ()-z2exp) > kRoad0z) {
1505 rieman.AddPoint(xcl[sLayer+0],ycl[sLayer+0],zcl[sLayer+0],1,10);
1506 rieman.AddPoint(xcl[sLayer+1],ycl[sLayer+1],zcl[sLayer+1],1,10);
1507 rieman.AddPoint(xcl[sLayer+3],ycl[sLayer+3],zcl[sLayer+3],1,10);
1508 rieman.AddPoint(xcl[sLayer+2],ycl[sLayer+2],zcl[sLayer+2],1,10);
1512 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1513 cseed[iLayer].Reset();
1518 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1519 cseed[sLayer+iLayer].SetZref(0,rieman.GetZat(xcl[sLayer+iLayer]));
1520 chi2Z += (cseed[sLayer+iLayer].GetZref(0)- zcl[sLayer+iLayer])
1521 * (cseed[sLayer+iLayer].GetZref(0)- zcl[sLayer+iLayer]);
1522 cseed[sLayer+iLayer].SetZref(1,rieman.GetDZat(xcl[sLayer+iLayer]));
1523 cseed[sLayer+iLayer].SetYref(0,rieman.GetYat(xcl[sLayer+iLayer]));
1524 chi2R += (cseed[sLayer+iLayer].GetYref(0)- ycl[sLayer+iLayer])
1525 * (cseed[sLayer+iLayer].GetYref(0)- ycl[sLayer+iLayer]);
1526 cseed[sLayer+iLayer].SetYref(1,rieman.GetDYat(xcl[sLayer+iLayer]));
1528 if (TMath::Sqrt(chi2R) > 1.0/iter) {
1531 if (TMath::Sqrt(chi2Z) > 7.0/iter) {
1535 Float_t minmax[2] = { -100.0, 100.0 };
1536 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1537 Float_t max = zcl[sLayer+iLayer]+padlength[sLayer+iLayer] * 0.5
1538 + 1.0 - cseed[sLayer+iLayer].GetZref(0);
1539 if (max < minmax[1]) {
1542 Float_t min = zcl[sLayer+iLayer]-padlength[sLayer+iLayer] * 0.5
1543 - 1.0 - cseed[sLayer+iLayer].GetZref(0);
1544 if (min > minmax[0]) {
1549 Bool_t isFake = kFALSE;
1550 if (cl0->GetLabel(0) != cl3->GetLabel(0)) {
1553 if (cl1->GetLabel(0) != cl3->GetLabel(0)) {
1556 if (cl2->GetLabel(0) != cl3->GetLabel(0)) {
1560 if (AliTRDReconstructor::StreamLevel() > 0) {
1561 if ((!isFake) || ((icl3%10) == 0)) { // Debugging print
1562 TTreeSRedirector &cstream = *fDebugStreamer;
1564 << "isFake=" << isFake
1570 << "X0=" << xcl[sLayer+0]
1571 << "X1=" << xcl[sLayer+1]
1572 << "X2=" << xcl[sLayer+2]
1573 << "X3=" << xcl[sLayer+3]
1574 << "Y2exp=" << y2exp
1575 << "Z2exp=" << z2exp
1576 << "Chi2R=" << chi2R
1577 << "Chi2Z=" << chi2Z
1578 << "Seed0.=" << &cseed[sLayer+0]
1579 << "Seed1.=" << &cseed[sLayer+1]
1580 << "Seed2.=" << &cseed[sLayer+2]
1581 << "Seed3.=" << &cseed[sLayer+3]
1582 << "Zmin=" << minmax[0]
1583 << "Zmax=" << minmax[1]
1588 ////////////////////////////////////////////////////////////////////////////////////
1592 ////////////////////////////////////////////////////////////////////////////////////
1598 Bool_t isOK = kTRUE;
1600 for (Int_t jLayer = 0; jLayer < 4; jLayer++) {
1602 cseed[sLayer+jLayer].SetTilt(hL[sLayer+jLayer]);
1603 cseed[sLayer+jLayer].SetPadLength(padlength[sLayer+jLayer]);
1604 cseed[sLayer+jLayer].SetX0(xcl[sLayer+jLayer]);
1606 for (Int_t iter = 0; iter < 2; iter++) {
1609 // In iteration 0 we try only one pad-row
1610 // If quality not sufficient we try 2 pad-rows - about 5% of tracks cross 2 pad-rows
1612 AliTRDseed tseed = cseed[sLayer+jLayer];
1613 Float_t roadz = padlength[sLayer+jLayer] * 0.5;
1615 roadz = padlength[sLayer+jLayer];
1618 Float_t quality = 10000.0;
1620 for (Int_t iTime = 2; iTime < 20; iTime++) {
1622 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[sLayer+jLayer][1]-iTime));
1623 Double_t dxlayer = layer.GetX() - xcl[sLayer+jLayer];
1624 Double_t zexp = cl[sLayer+jLayer]->GetZ();
1627 // Try 2 pad-rows in second iteration
1628 zexp = tseed.GetZref(0) + tseed.GetZref(1) * dxlayer;
1629 if (zexp > cl[sLayer+jLayer]->GetZ()) {
1630 zexp = cl[sLayer+jLayer]->GetZ() + padlength[sLayer+jLayer]*0.5;
1632 if (zexp < cl[sLayer+jLayer]->GetZ()) {
1633 zexp = cl[sLayer+jLayer]->GetZ() - padlength[sLayer+jLayer]*0.5;
1637 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer;
1638 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y,roadz);
1642 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
1644 tseed.SetIndexes(iTime,index);
1645 tseed.SetClusters(iTime,cl); // Register cluster
1646 tseed.SetX(iTime,dxlayer); // Register cluster
1647 tseed.SetY(iTime,cl->GetY()); // Register cluster
1648 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
1654 // Count the number of clusters and distortions into quality
1655 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
1656 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0 + TMath::Abs(dangle) / 0.1
1657 + TMath::Abs(tseed.GetYfit(0) - tseed.GetYref(0)) / 0.2
1658 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
1659 if ((iter == 0) && tseed.IsOK()) {
1660 cseed[sLayer+jLayer] = tseed;
1666 if (tseed.IsOK() && (tquality < quality)) {
1667 cseed[sLayer+jLayer] = tseed;
1672 if (!cseed[sLayer+jLayer].IsOK()) {
1677 cseed[sLayer+jLayer].CookLabels();
1678 cseed[sLayer+jLayer].UpdateUsed();
1679 nusedCl += cseed[sLayer+jLayer].GetNUsed();
1691 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1692 if (cseed[sLayer+iLayer].IsOK()) {
1693 nclusters += cseed[sLayer+iLayer].GetN2();
1699 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1700 rieman.AddPoint(xcl[sLayer+iLayer]
1701 ,cseed[sLayer+iLayer].GetYfitR(0)
1702 ,cseed[sLayer+iLayer].GetZProb()
1711 for (Int_t iLayer = 0; iLayer < 4; iLayer++) {
1712 cseed[sLayer+iLayer].SetYref(0,rieman.GetYat(xcl[sLayer+iLayer]));
1713 chi2R += (cseed[sLayer+iLayer].GetYref(0) - cseed[sLayer+iLayer].GetYfitR(0))
1714 * (cseed[sLayer+iLayer].GetYref(0) - cseed[sLayer+iLayer].GetYfitR(0));
1715 cseed[sLayer+iLayer].SetYref(1,rieman.GetDYat(xcl[sLayer+iLayer]));
1716 cseed[sLayer+iLayer].SetZref(0,rieman.GetZat(xcl[sLayer+iLayer]));
1717 chi2Z += (cseed[sLayer+iLayer].GetZref(0) - cseed[sLayer+iLayer].GetMeanz())
1718 * (cseed[sLayer+iLayer].GetZref(0) - cseed[sLayer+iLayer].GetMeanz());
1719 cseed[sLayer+iLayer].SetZref(1,rieman.GetDZat(xcl[sLayer+iLayer]));
1721 Double_t curv = rieman.GetC();
1726 Double_t sumda = TMath::Abs(cseed[sLayer+0].GetYfitR(1) - cseed[sLayer+0].GetYref(1))
1727 + TMath::Abs(cseed[sLayer+1].GetYfitR(1) - cseed[sLayer+1].GetYref(1))
1728 + TMath::Abs(cseed[sLayer+2].GetYfitR(1) - cseed[sLayer+2].GetYref(1))
1729 + TMath::Abs(cseed[sLayer+3].GetYfitR(1) - cseed[sLayer+3].GetYref(1));
1730 Double_t likea = TMath::Exp(-sumda*10.6);
1731 Double_t likechi2 = 0.0000000001;
1733 likechi2 += TMath::Exp(-TMath::Sqrt(chi2R) * 7.73);
1735 Double_t likechi2z = TMath::Exp(-chi2Z * 0.088) / TMath::Exp(-chi2Z * 0.019);
1736 Double_t likeN = TMath::Exp(-(72 - nclusters) * 0.19);
1737 Double_t like = likea * likechi2 * likechi2z * likeN;
1738 Double_t likePrimY = TMath::Exp(-TMath::Abs(cseed[sLayer+0].GetYref(1) - 130.0*curv) * 1.9);
1739 Double_t likePrimZ = TMath::Exp(-TMath::Abs(cseed[sLayer+0].GetZref(1)
1740 - cseed[sLayer+0].GetZref(0) / xcl[sLayer+0]) * 5.9);
1741 Double_t likePrim = TMath::Max(likePrimY*likePrimZ,0.0005);
1743 seedquality[registered] = like;
1744 seedlayer[registered] = sLayer;
1745 if (TMath::Log(0.000000000000001 + like) < -15) {
1748 AliTRDseed seedb[6];
1749 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1750 seedb[iLayer] = cseed[iLayer];
1753 ////////////////////////////////////////////////////////////////////////////////////
1755 // Full track fit part
1757 ////////////////////////////////////////////////////////////////////////////////////
1764 // Add new layers - avoid long extrapolation
1766 Int_t tLayer[2] = { 0, 0 };
1780 for (Int_t iLayer = 0; iLayer < 2; iLayer++) {
1781 Int_t jLayer = tLayer[iLayer]; // Set tracking layer
1782 cseed[jLayer].Reset();
1783 cseed[jLayer].SetTilt(hL[jLayer]);
1784 cseed[jLayer].SetPadLength(padlength[jLayer]);
1785 cseed[jLayer].SetX0(xcl[jLayer]);
1786 // Get pad length and rough cluster
1787 Int_t indexdummy = reflayers[jLayer]->FindNearestCluster(cseed[jLayer].GetYref(0)
1788 ,cseed[jLayer].GetZref(0)
1791 if (indexdummy <= 0) {
1794 AliTRDcluster *cldummy = (AliTRDcluster *) GetCluster(indexdummy);
1795 padlength[jLayer] = TMath::Sqrt(cldummy->GetSigmaZ2() * 12.0);
1797 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
1799 for (Int_t iLayer = 0; iLayer < 2; iLayer++) {
1801 Int_t jLayer = tLayer[iLayer]; // set tracking layer
1802 if ((jLayer == 0) && !(cseed[1].IsOK())) {
1803 continue; // break not allowed
1805 if ((jLayer == 5) && !(cseed[4].IsOK())) {
1806 continue; // break not allowed
1808 Float_t zexp = cseed[jLayer].GetZref(0);
1809 Double_t zroad = padlength[jLayer] * 0.5 + 1.0;
1811 for (Int_t iter = 0; iter < 2; iter++) {
1813 AliTRDseed tseed = cseed[jLayer];
1814 Float_t quality = 10000.0;
1816 for (Int_t iTime = 2; iTime < 20; iTime++) {
1817 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[jLayer][1]-iTime));
1818 Double_t dxlayer = layer.GetX()-xcl[jLayer];
1819 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer;
1820 Float_t yroad = kRoad1y;
1821 Int_t index = layer.FindNearestCluster(yexp,zexp,yroad,zroad);
1825 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
1826 tseed.SetIndexes(iTime,index);
1827 tseed.SetClusters(iTime,cl); // Register cluster
1828 tseed.SetX(iTime,dxlayer); // Register cluster
1829 tseed.SetY(iTime,cl->GetY()); // Register cluster
1830 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
1835 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
1836 Float_t tquality = (18.0 - tseed.GetN2())/2.0 + TMath::Abs(dangle) / 0.1
1837 + TMath::Abs(tseed.GetYfit(0) - tseed.GetYref(0)) / 0.2
1838 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
1839 if (tquality < quality) {
1840 cseed[jLayer] = tseed;
1849 if ( cseed[jLayer].IsOK()) {
1850 cseed[jLayer].CookLabels();
1851 cseed[jLayer].UpdateUsed();
1852 nusedf += cseed[jLayer].GetNUsed();
1853 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
1859 AliTRDseed bseed[6];
1860 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
1861 bseed[jLayer] = cseed[jLayer];
1863 Float_t lastquality = 10000.0;
1864 Float_t lastchi2 = 10000.0;
1865 Float_t chi2 = 1000.0;
1867 for (Int_t iter = 0; iter < 4; iter++) {
1869 // Sort tracklets according "quality", try to "improve" 4 worst
1870 Float_t sumquality = 0.0;
1871 Float_t squality[6];
1872 Int_t sortindexes[6];
1874 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
1875 if (bseed[jLayer].IsOK()) {
1876 AliTRDseed &tseed = bseed[jLayer];
1877 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
1878 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
1879 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0 + TMath::Abs(dangle) / 0.1
1880 + TMath::Abs(tseed.GetYfit(0) - (tseed.GetYref(0) - zcor)) / 0.2
1881 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
1882 squality[jLayer] = tquality;
1885 squality[jLayer] = -1.0;
1887 sumquality +=squality[jLayer];
1890 if ((sumquality >= lastquality) ||
1891 (chi2 > lastchi2)) {
1894 lastquality = sumquality;
1897 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
1898 cseed[jLayer] = bseed[jLayer];
1901 TMath::Sort(6,squality,sortindexes,kFALSE);
1903 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
1905 Int_t bLayer = sortindexes[jLayer];
1906 AliTRDseed tseed = bseed[bLayer];
1908 for (Int_t iTime = 2; iTime < 20; iTime++) {
1910 AliTRDpropagationLayer &layer = *(fTrSec[ns]->GetLayer(layers[bLayer][1]-iTime));
1911 Double_t dxlayer = layer.GetX() - xcl[bLayer];
1912 Double_t zexp = tseed.GetZref(0);
1913 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
1914 Float_t roadz = padlength[bLayer] + 1;
1915 if (TMath::Abs(tseed.GetZProb() - zexp) > 0.5*padlength[bLayer]) {
1916 roadz = padlength[bLayer] * 0.5;
1918 if (tseed.GetZfit(1)*tseed.GetZref(1) < 0.0) {
1919 roadz = padlength[bLayer] * 0.5;
1921 if (TMath::Abs(tseed.GetZProb() - zexp) < 0.1*padlength[bLayer]) {
1922 zexp = tseed.GetZProb();
1923 roadz = padlength[bLayer] * 0.5;
1926 Double_t yexp = tseed.GetYref(0) + tseed.GetYref(1) * dxlayer - zcor;
1927 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y,roadz);
1931 AliTRDcluster *cl = (AliTRDcluster *) GetCluster(index);
1933 tseed.SetIndexes(iTime,index);
1934 tseed.SetClusters(iTime,cl); // Register cluster
1935 tseed.SetX(iTime,dxlayer); // Register cluster
1936 tseed.SetY(iTime,cl->GetY()); // Register cluster
1937 tseed.SetZ(iTime,cl->GetZ()); // Register cluster
1943 Float_t dangle = tseed.GetYfit(1) - tseed.GetYref(1);
1944 Double_t zcor = tseed.GetTilt() * (tseed.GetZProb() - tseed.GetZref(0));
1945 Float_t tquality = (18.0 - tseed.GetN2()) / 2.0
1946 + TMath::Abs(dangle) / 0.1
1947 + TMath::Abs(tseed.GetYfit(0) - (tseed.GetYref(0) - zcor)) / 0.2
1948 + 2.0 * TMath::Abs(tseed.GetMeanz() - tseed.GetZref(0)) / padlength[jLayer];
1949 if (tquality<squality[bLayer]) {
1950 bseed[bLayer] = tseed;
1956 chi2 = AliTRDseed::FitRiemanTilt(bseed,kTRUE);
1963 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1964 if (TMath::Abs(cseed[iLayer].GetYref(0) / cseed[iLayer].GetX0()) < 0.15) {
1967 if (cseed[iLayer].IsOK()) {
1968 nclusters += cseed[iLayer].GetN2();
1976 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1977 if (cseed[iLayer].IsOK()) {
1978 rieman.AddPoint(xcl[iLayer]
1979 ,cseed[iLayer].GetYfitR(0)
1980 ,cseed[iLayer].GetZProb()
1989 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1990 if (cseed[iLayer].IsOK()) {
1991 cseed[iLayer].SetYref(0,rieman.GetYat(xcl[iLayer]));
1992 chi2RF += (cseed[iLayer].GetYref(0) - cseed[iLayer].GetYfitR(0))
1993 * (cseed[iLayer].GetYref(0) - cseed[iLayer].GetYfitR(0));
1994 cseed[iLayer].SetYref(1,rieman.GetDYat(xcl[iLayer]));
1995 cseed[iLayer].SetZref(0,rieman.GetZat(xcl[iLayer]));
1996 chi2ZF += (cseed[iLayer].GetZref(0) - cseed[iLayer].GetMeanz())
1997 * (cseed[iLayer].GetZref(0) - cseed[iLayer].GetMeanz());
1998 cseed[iLayer].SetZref(1,rieman.GetDZat(xcl[iLayer]));
2001 chi2RF /= TMath::Max((nlayers - 3.0),1.0);
2002 chi2ZF /= TMath::Max((nlayers - 3.0),1.0);
2003 curv = rieman.GetC();
2005 Double_t xref2 = (xcl[2] + xcl[3]) * 0.5; // Middle of the chamber
2006 Double_t dzmf = rieman.GetDZat(xref2);
2007 Double_t zmf = rieman.GetZat(xref2);
2013 fitterTC.ClearPoints();
2014 fitterT2.ClearPoints();
2017 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2019 if (!cseed[iLayer].IsOK()) {
2023 for (Int_t itime = 0; itime < 25; itime++) {
2025 if (!cseed[iLayer].IsUsable(itime)) {
2028 // X relative to the middle chamber
2029 Double_t x = cseed[iLayer].GetX(itime) + cseed[iLayer].GetX0() - xref2;
2030 Double_t y = cseed[iLayer].GetY(itime);
2031 Double_t z = cseed[iLayer].GetZ(itime);
2032 // ExB correction to the correction
2036 Double_t x2 = cseed[iLayer].GetX(itime) + cseed[iLayer].GetX0();
2037 Double_t t = 1.0 / (x2*x2 + y*y);
2039 uvt[0] = 2.0 * x2 * uvt[1]; // u
2040 uvt[2] = 2.0 * hL[iLayer] * uvt[1];
2041 uvt[3] = 2.0 * hL[iLayer] * x * uvt[1];
2042 uvt[4] = 2.0 * (y + hL[iLayer]*z) * uvt[1];
2043 Double_t error = 2.0 * 0.2 * uvt[1];
2044 fitterT2.AddPoint(uvt,uvt[4],error);
2047 // Constrained rieman
2049 z = cseed[iLayer].GetZ(itime);
2050 uvt[0] = 2.0 * x2 * t; // u
2051 uvt[1] = 2.0 * hL[iLayer] * x2 * uvt[1];
2052 uvt[2] = 2.0 * (y + hL[iLayer] * (z - GetZ())) * t;
2053 fitterTC.AddPoint(uvt,uvt[2],error);
2054 rieman2.AddPoint(x2,y,z,1,10);
2064 Double_t rpolz0 = fitterT2.GetParameter(3);
2065 Double_t rpolz1 = fitterT2.GetParameter(4);
2068 // Linear fitter - not possible to make boundaries
2069 // Do not accept non possible z and dzdx combinations
2071 Bool_t acceptablez = kTRUE;
2072 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2073 if (cseed[iLayer].IsOK()) {
2074 Double_t zT2 = rpolz0 + rpolz1 * (xcl[iLayer] - xref2);
2075 if (TMath::Abs(cseed[iLayer].GetZProb() - zT2) > padlength[iLayer] * 0.5 + 1.0) {
2076 acceptablez = kFALSE;
2081 fitterT2.FixParameter(3,zmf);
2082 fitterT2.FixParameter(4,dzmf);
2084 fitterT2.ReleaseParameter(3);
2085 fitterT2.ReleaseParameter(4);
2086 rpolz0 = fitterT2.GetParameter(3);
2087 rpolz1 = fitterT2.GetParameter(4);
2090 Double_t chi2TR = fitterT2.GetChisquare() / Float_t(npointsT);
2091 Double_t chi2TC = fitterTC.GetChisquare() / Float_t(npointsT);
2092 Double_t polz1c = fitterTC.GetParameter(2);
2093 Double_t polz0c = polz1c * xref2;
2094 Double_t aC = fitterTC.GetParameter(0);
2095 Double_t bC = fitterTC.GetParameter(1);
2096 Double_t cC = aC / TMath::Sqrt(bC * bC + 1.0); // Curvature
2097 Double_t aR = fitterT2.GetParameter(0);
2098 Double_t bR = fitterT2.GetParameter(1);
2099 Double_t dR = fitterT2.GetParameter(2);
2100 Double_t cR = 1.0 + bR*bR - dR*aR;
2103 dca = -dR / (TMath::Sqrt(1.0 + bR*bR - dR*aR) + TMath::Sqrt(1.0 + bR*bR));
2104 cR = aR / TMath::Sqrt(cR);
2107 Double_t chi2ZT2 = 0.0;
2108 Double_t chi2ZTC = 0.0;
2109 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2110 if (cseed[iLayer].IsOK()) {
2111 Double_t zT2 = rpolz0 + rpolz1 * (xcl[iLayer] - xref2);
2112 Double_t zTC = polz0c + polz1c * (xcl[iLayer] - xref2);
2113 chi2ZT2 += TMath::Abs(cseed[iLayer].GetMeanz() - zT2);
2114 chi2ZTC += TMath::Abs(cseed[iLayer].GetMeanz() - zTC);
2117 chi2ZT2 /= TMath::Max((nlayers - 3.0),1.0);
2118 chi2ZTC /= TMath::Max((nlayers - 3.0),1.0);
2120 AliTRDseed::FitRiemanTilt(cseed,kTRUE);
2121 Float_t sumdaf = 0.0;
2122 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2123 if (cseed[iLayer].IsOK()) {
2124 sumdaf += TMath::Abs((cseed[iLayer].GetYfit(1) - cseed[iLayer].GetYref(1))
2125 / cseed[iLayer].GetSigmaY2());
2128 sumdaf /= Float_t (nlayers - 2.0);
2131 // Likelihoods for full track
2133 Double_t likezf = TMath::Exp(-chi2ZF * 0.14);
2134 Double_t likechi2C = TMath::Exp(-chi2TC * 0.677);
2135 Double_t likechi2TR = TMath::Exp(-chi2TR * 0.78);
2136 Double_t likeaf = TMath::Exp(-sumdaf * 3.23);
2137 seedquality2[registered] = likezf * likechi2TR * likeaf;
2139 // Still needed ????
2140 // Bool_t isGold = kFALSE;
2142 // if (nlayers == 6 && TMath::Log(0.000000001+seedquality2[index])<-5.) isGold =kTRUE; // gold
2143 // if (nlayers == findable && TMath::Log(0.000000001+seedquality2[index])<-4.) isGold =kTRUE; // gold
2144 // if (isGold &&nusedf<10){
2145 // for (Int_t jLayer=0;jLayer<6;jLayer++){
2146 // if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.1)
2147 // seed[index][jLayer].UseClusters(); //sign gold
2152 if (!cseed[0].IsOK()) {
2154 if (!cseed[1].IsOK()) {
2158 seedparams[registered][0] = cseed[index0].GetX0();
2159 seedparams[registered][1] = cseed[index0].GetYref(0);
2160 seedparams[registered][2] = cseed[index0].GetZref(0);
2161 seedparams[registered][5] = cR;
2162 seedparams[registered][3] = cseed[index0].GetX0() * cR - TMath::Sin(TMath::ATan(cseed[0].GetYref(1)));
2163 seedparams[registered][4] = cseed[index0].GetZref(1)
2164 / TMath::Sqrt(1.0 + cseed[index0].GetYref(1) * cseed[index0].GetYref(1));
2165 seedparams[registered][6] = ns;
2170 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2171 if (!cseed[iLayer].IsOK()) {
2174 if (cseed[iLayer].GetLabels(0) >= 0) {
2175 labels[nlab] = cseed[iLayer].GetLabels(0);
2178 if (cseed[iLayer].GetLabels(1) >= 0) {
2179 labels[nlab] = cseed[iLayer].GetLabels(1);
2183 Freq(nlab,labels,outlab,kFALSE);
2184 Int_t label = outlab[0];
2185 Int_t frequency = outlab[1];
2186 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2187 cseed[iLayer].SetFreq(frequency);
2188 cseed[iLayer].SetC(cR);
2189 cseed[iLayer].SetCC(cC);
2190 cseed[iLayer].SetChi2(chi2TR);
2191 cseed[iLayer].SetChi2Z(chi2ZF);
2195 if (1 || (!isFake)) {
2196 Float_t zvertex = GetZ();
2197 TTreeSRedirector &cstream = *fDebugStreamer;
2198 if (AliTRDReconstructor::StreamLevel() > 0) {
2200 << "isFake=" << isFake
2201 << "Vertex=" << zvertex
2202 << "Rieman2.=" << &rieman2
2203 << "Rieman.=" << &rieman
2211 << "Chi2R=" << chi2R
2212 << "Chi2Z=" << chi2Z
2213 << "Chi2RF=" << chi2RF // Chi2 of trackletes on full track
2214 << "Chi2ZF=" << chi2ZF // Chi2 z on tracklets on full track
2215 << "Chi2ZT2=" << chi2ZT2 // Chi2 z on tracklets on full track - rieman tilt
2216 << "Chi2ZTC=" << chi2ZTC // Chi2 z on tracklets on full track - rieman tilt const
2217 << "Chi2TR=" << chi2TR // Chi2 without vertex constrain
2218 << "Chi2TC=" << chi2TC // Chi2 with vertex constrain
2219 << "C=" << curv // Non constrained - no tilt correction
2220 << "DR=" << dR // DR parameter - tilt correction
2221 << "DCA=" << dca // DCA - tilt correction
2222 << "CR=" << cR // Non constrained curvature - tilt correction
2223 << "CC=" << cC // Constrained curvature
2224 << "Polz0=" << polz0c
2225 << "Polz1=" << polz1c
2226 << "RPolz0=" << rpolz0
2227 << "RPolz1=" << rpolz1
2228 << "Ncl=" << nclusters
2229 << "Nlayers=" << nlayers
2230 << "NUsedS=" << nusedCl
2231 << "NUsed=" << nusedf
2232 << "Findable=" << findable
2234 << "LikePrim=" << likePrim
2235 << "Likechi2C=" << likechi2C
2236 << "Likechi2TR=" << likechi2TR
2237 << "Likezf=" << likezf
2238 << "LikeF=" << seedquality2[registered]
2239 << "S0.=" << &cseed[0]
2240 << "S1.=" << &cseed[1]
2241 << "S2.=" << &cseed[2]
2242 << "S3.=" << &cseed[3]
2243 << "S4.=" << &cseed[4]
2244 << "S5.=" << &cseed[5]
2245 << "SB0.=" << &seedb[0]
2246 << "SB1.=" << &seedb[1]
2247 << "SB2.=" << &seedb[2]
2248 << "SB3.=" << &seedb[3]
2249 << "SB4.=" << &seedb[4]
2250 << "SB5.=" << &seedb[5]
2251 << "Label=" << label
2252 << "Freq=" << frequency
2253 << "sLayer=" << sLayer
2258 if (registered<kMaxSeed - 1) {
2260 cseed = seed[registered];
2263 } // End of loop over layer 1
2265 } // End of loop over layer 0
2267 } // End of loop over layer 3
2269 } // End of loop over seeding time bins
2275 TMath::Sort(registered,seedquality2,sort,kTRUE);
2276 Bool_t signedseed[kMaxSeed];
2277 for (Int_t i = 0; i < registered; i++) {
2278 signedseed[i] = kFALSE;
2281 for (Int_t iter = 0; iter < 5; iter++) {
2283 for (Int_t iseed = 0; iseed < registered; iseed++) {
2285 Int_t index = sort[iseed];
2286 if (signedseed[index]) {
2289 Int_t labelsall[1000];
2290 Int_t nlabelsall = 0;
2291 Int_t naccepted = 0;;
2292 Int_t sLayer = seedlayer[index];
2298 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2300 if (TMath::Abs(seed[index][jLayer].GetYref(0) / xcl[jLayer]) < 0.15) {
2303 if (seed[index][jLayer].IsOK()) {
2304 seed[index][jLayer].UpdateUsed();
2305 ncl +=seed[index][jLayer].GetN2();
2306 nused +=seed[index][jLayer].GetNUsed();
2309 for (Int_t itime = 0; itime < 25; itime++) {
2310 if (seed[index][jLayer].IsUsable(itime)) {
2312 for (Int_t ilab = 0; ilab < 3; ilab++) {
2313 Int_t tindex = seed[index][jLayer].GetClusters(itime)->GetLabel(ilab);
2315 labelsall[nlabelsall] = tindex;
2333 if (TMath::Log(0.000000001+seedquality2[index]) < -5.0) {
2339 if (nlayers < findable) {
2342 if (TMath::Log(0.000000001+seedquality2[index]) < -4.0) {
2348 if ((nlayers == findable) ||
2352 if (TMath::Log(0.000000001+seedquality2[index]) < -6.0) {
2358 if (TMath::Log(0.000000001+seedquality2[index]) < -5.0) {
2364 if (TMath::Log(0.000000001+seedquality2[index]) - nused/(nlayers-3.0) < -15.0) {
2369 signedseed[index] = kTRUE;
2374 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2375 if (seed[index][iLayer].IsOK()) {
2376 if (seed[index][iLayer].GetLabels(0) >= 0) {
2377 labels[nlab] = seed[index][iLayer].GetLabels(0);
2380 if (seed[index][iLayer].GetLabels(1) >= 0) {
2381 labels[nlab] = seed[index][iLayer].GetLabels(1);
2386 Freq(nlab,labels,outlab,kFALSE);
2387 Int_t label = outlab[0];
2388 Int_t frequency = outlab[1];
2389 Freq(nlabelsall,labelsall,outlab,kFALSE);
2390 Int_t label1 = outlab[0];
2391 Int_t label2 = outlab[2];
2392 Float_t fakeratio = (naccepted - outlab[1]) / Float_t(naccepted);
2393 Float_t ratio = Float_t(nused) / Float_t(ncl);
2395 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2396 if ((seed[index][jLayer].IsOK()) &&
2397 (TMath::Abs(seed[index][jLayer].GetYfit(1) - seed[index][jLayer].GetYfit(1)) < 0.2)) {
2398 seed[index][jLayer].UseClusters(); // Sign gold
2403 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.
2404 TTreeSRedirector &cstream = *fDebugStreamer;
2409 AliTRDtrack *track = RegisterSeed(seed[index],seedparams[index]);
2415 AliESDtrack esdtrack;
2416 esdtrack.UpdateTrackParams(track,AliESDtrack::kTRDout);
2417 esdtrack.SetLabel(label);
2418 esd->AddTrack(&esdtrack);
2419 TTreeSRedirector &cstream = *fDebugStreamer;
2420 if (AliTRDReconstructor::StreamLevel() > 0) {
2422 << "EventNrInFile=" << eventNrInFile
2423 << "ESD.=" << &esdtrack
2425 << "trdback.=" << track
2430 if (AliTRDReconstructor::StreamLevel() > 0) {
2433 << "Track.=" << track
2434 << "Like=" << seedquality[index]
2435 << "LikeF=" << seedquality2[index]
2436 << "S0.=" << &seed[index][0]
2437 << "S1.=" << &seed[index][1]
2438 << "S2.=" << &seed[index][2]
2439 << "S3.=" << &seed[index][3]
2440 << "S4.=" << &seed[index][4]
2441 << "S5.=" << &seed[index][5]
2442 << "Label=" << label
2443 << "Label1=" << label1
2444 << "Label2=" << label2
2445 << "FakeRatio=" << fakeratio
2446 << "Freq=" << frequency
2448 << "Nlayers=" << nlayers
2449 << "Findable=" << findable
2450 << "NUsed=" << nused
2451 << "sLayer=" << sLayer
2452 << "EventNrInFile=" << eventNrInFile
2460 } // End of loop over sectors
2467 //_____________________________________________________________________________
2468 Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *clusterTree) const
2471 // Reads AliTRDclusters from the file.
2472 // The names of the cluster tree and branches
2473 // should match the ones used in AliTRDclusterizer::WriteClusters()
2476 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
2477 TObjArray *clusterArray = new TObjArray(nsize+1000);
2479 TBranch *branch = clusterTree->GetBranch("TRDcluster");
2481 AliError("Can't get the branch !");
2484 branch->SetAddress(&clusterArray);
2486 // Loop through all entries in the tree
2487 Int_t nEntries = (Int_t) clusterTree->GetEntries();
2489 AliTRDcluster *c = 0x0;
2490 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
2493 nbytes += clusterTree->GetEvent(iEntry);
2495 // Get the number of points in the detector
2496 Int_t nCluster = clusterArray->GetEntriesFast();
2497 // Loop through all TRD digits
2498 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
2499 if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
2501 //printf("Add cluster 0x%x.\n", c);
2502 clusterArray->RemoveAt(iCluster);
2506 delete clusterArray;
2512 //_____________________________________________________________________________
2513 Bool_t AliTRDtracker::GetTrackPoint(Int_t index, AliTrackPoint &p) const
2516 // Get track space point with index i
2517 // Origin: C.Cheshkov
2520 AliTRDcluster *cl = (AliTRDcluster *) fClusters->UncheckedAt(index);
2521 Int_t idet = cl->GetDetector();
2522 Int_t isector = fGeom->GetSector(idet);
2523 Int_t ichamber = fGeom->GetChamber(idet);
2524 Int_t iplan = fGeom->GetPlane(idet);
2526 local[0] = GetX(isector,iplan,cl->GetLocalTimeBin());
2527 local[1] = cl->GetY();
2528 local[2] = cl->GetZ();
2530 fGeom->RotateBack(idet,local,global);
2531 p.SetXYZ(global[0],global[1],global[2]);
2532 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
2535 iLayer = AliGeomManager::kTRD1;
2538 iLayer = AliGeomManager::kTRD2;
2541 iLayer = AliGeomManager::kTRD3;
2544 iLayer = AliGeomManager::kTRD4;
2547 iLayer = AliGeomManager::kTRD5;
2550 iLayer = AliGeomManager::kTRD6;
2553 Int_t modId = isector * fGeom->Ncham() + ichamber;
2554 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer,modId);
2555 p.SetVolumeID(volid);
2561 //_____________________________________________________________________________
2562 void AliTRDtracker::CookLabel(AliKalmanTrack *pt, Float_t wrong) const
2565 // This cooks a label. Mmmmh, smells good...
2568 Int_t label = 123456789;
2572 Int_t ncl = pt->GetNumberOfClusters();
2574 const Int_t kRange = fTrSec[0]->GetOuterTimeBin() + 1;
2578 Int_t **s = new Int_t* [kRange];
2579 for (i = 0; i < kRange; i++) {
2580 s[i] = new Int_t[2];
2582 for (i = 0; i < kRange; i++) {
2591 for (i = 0; i < ncl; i++) {
2592 index = pt->GetClusterIndex(i);
2593 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2599 for (i = 0; i < ncl; i++) {
2600 index = pt->GetClusterIndex(i);
2601 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2602 for (Int_t k = 0; k < 3; k++) {
2603 label = c->GetLabel(k);
2604 labelAdded = kFALSE;
2607 while ((!labelAdded) && (j < kRange)) {
2608 if ((s[j][0] == label) ||
2611 s[j][1] = s[j][1] + 1;
2623 for (i = 0; i < kRange; i++) {
2624 if (s[i][1] > max) {
2630 for (i = 0; i < kRange; i++) {
2636 if ((1.0 - Float_t(max)/ncl) > wrong) {
2640 pt->SetLabel(label);
2644 //_____________________________________________________________________________
2645 void AliTRDtracker::UseClusters(const AliKalmanTrack *t, Int_t from) const
2648 // Use clusters, but don't abuse them!
2651 const Float_t kmaxchi2 = 18;
2652 const Float_t kmincl = 10;
2654 AliTRDtrack *track = (AliTRDtrack *) t;
2656 Int_t ncl = t->GetNumberOfClusters();
2657 for (Int_t i = from; i < ncl; i++) {
2658 Int_t index = t->GetClusterIndex(i);
2659 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(index);
2660 Int_t iplane = fGeom->GetPlane(c->GetDetector());
2661 if (track->GetTracklets(iplane).GetChi2() > kmaxchi2) {
2664 if (track->GetTracklets(iplane).GetN() < kmincl) {
2667 if (!(c->IsUsed())) {
2674 //_____________________________________________________________________________
2675 Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
2678 // Parametrised "expected" error of the cluster reconstruction in Y
2681 Double_t s = 0.08 * 0.08;
2686 //_____________________________________________________________________________
2687 Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
2690 // Parametrised "expected" error of the cluster reconstruction in Z
2693 Double_t s = 9.0 * 9.0 / 12.0;
2698 //_____________________________________________________________________________
2699 Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
2702 // Returns radial position which corresponds to time bin <localTB>
2703 // in tracking sector <sector> and plane <plane>
2706 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
2707 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
2709 return fTrSec[sector]->GetLayer(pl)->GetX();
2714 //_____________________________________________________________________________
2715 AliTRDtracker::AliTRDtrackingSector
2716 ::AliTRDtrackingSector(AliTRDgeometry *geo, Int_t gs)
2722 // AliTRDtrackingSector Constructor
2725 AliTRDpadPlane *padPlane = 0;
2726 AliTRDpropagationLayer *ppl = 0;
2728 // Get holes description from geometry
2729 Bool_t holes[AliTRDgeometry::kNcham];
2730 for (Int_t icham = 0; icham < AliTRDgeometry::kNcham; icham++) {
2731 holes[icham] = fGeom->IsHole(0,icham,gs);
2734 for (UInt_t i = 0; i < kMaxTimeBinIndex; i++) {
2735 fTimeBinIndex[i] = -1;
2743 // Add layers for each of the planes
2744 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2745 //Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2747 const Int_t kNchambers = AliTRDgeometry::Ncham();
2750 Double_t ymaxsensitive = 0;
2751 Double_t *zc = new Double_t[kNchambers];
2752 Double_t *zmax = new Double_t[kNchambers];
2753 Double_t *zmaxsensitive = new Double_t[kNchambers];
2755 for (Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
2757 ymax = fGeom->GetChamberWidth(plane) / 2.0;
2758 padPlane = fGeom->GetPadPlane(plane,0);
2759 ymaxsensitive = (padPlane->GetColSize(1) * padPlane->GetNcols() - 4.0) / 2.0;
2761 for (Int_t ch = 0; ch < kNchambers; ch++) {
2762 zmax[ch] = fGeom->GetChamberLength(plane,ch) / 2.0;
2763 Float_t pad = padPlane->GetRowSize(1);
2764 Float_t row0 = fGeom->GetRow0(plane,ch,0);
2765 Int_t nPads = fGeom->GetRowMax(plane,ch,0);
2766 zmaxsensitive[ch] = Float_t(nPads) * pad / 2.0;
2767 zc[ch] = -(pad * nPads) / 2.0 + row0;
2770 AliTRDcalibDB *fCalibration = AliTRDcalibDB::Instance();
2771 dx = fCalibration->GetVdrift(0,0,0)
2772 / AliTRDCommonParam::Instance()->GetSamplingFrequency();
2773 rho = 0.00295 * 0.85; //????
2776 Double_t x0 = (Double_t) AliTRDgeometry::GetTime0(plane);
2777 //Double_t xbottom = x0 - dxDrift;
2778 //Double_t xtop = x0 + dxAmp;
2780 //temporary !! (A.Bercuci)
2781 Int_t T0 = (Int_t)fCalibration->GetT0Average(AliTRDgeometry::GetDetector(plane, 2, gs));
2783 Int_t nTimeBins = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
2784 for (Int_t iTime = 0; iTime < nTimeBins; iTime++) {
2786 Double_t xlayer = iTime * dx - dxAmp;
2787 //if (xlayer<0) xlayer = dxAmp / 2.0;
2790 tbIndex = CookTimeBinIndex(plane,iTime);
2791 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex,plane);
2792 ppl->SetYmax(ymax,ymaxsensitive);
2793 ppl->SetZ(zc,zmax,zmaxsensitive);
2794 ppl->SetHoles(holes);
2795 if(iTime == T0) ppl->SetT0();
2807 delete [] zmaxsensitive;
2811 //_____________________________________________________________________________
2812 AliTRDtracker::AliTRDtrackingSector
2813 ::AliTRDtrackingSector(const AliTRDtrackingSector &/*t*/)
2824 //_____________________________________________________________________________
2825 AliTRDtracker::AliTRDtrackingSector
2826 ::~AliTRDtrackingSector()
2832 for (Int_t i = 0; i < fN; i++) {
2838 //_____________________________________________________________________________
2839 Int_t AliTRDtracker::AliTRDtrackingSector
2840 ::CookTimeBinIndex(Int_t plane, Int_t localTB) const
2843 // Depending on the digitization parameters calculates global
2844 // (i.e. for a whole TRD stack of 6 planes) time bin index for
2845 // timebin <localTB> in plane <plane>
2848 Int_t tbPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
2849 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1;
2861 //_____________________________________________________________________________
2862 void AliTRDtracker::AliTRDtrackingSector
2863 ::MapTimeBinLayers()
2866 // For all sensitive time bins sets corresponding layer index
2867 // in the array fTimeBins
2872 for (Int_t i = 0; i < fN; i++) {
2874 index = fLayers[i]->GetTimeBinIndex();
2879 if (index >= (Int_t) kMaxTimeBinIndex) {
2880 //AliWarning(Form("Index %d exceeds allowed maximum of %d!\n"
2881 // ,index,kMaxTimeBinIndex-1));
2885 fTimeBinIndex[index] = i;
2891 //_____________________________________________________________________________
2892 Int_t AliTRDtracker::AliTRDtrackingSector
2893 ::GetLayerNumber(Double_t x) const
2896 // Returns the number of time bin which in radial position is closest to <x>
2899 if (x >= fLayers[fN-1]->GetX()) {
2902 if (x <= fLayers[ 0]->GetX()) {
2908 Int_t m = (b + e) / 2;
2910 for ( ; b < e; m = (b + e) / 2) {
2911 if (x > fLayers[m]->GetX()) {
2919 if (TMath::Abs(x - fLayers[m]->GetX()) > TMath::Abs(x - fLayers[m+1]->GetX())) {
2928 //_____________________________________________________________________________
2929 Int_t AliTRDtracker::AliTRDtrackingSector
2930 ::GetInnerTimeBin() const
2933 // Returns number of the innermost SENSITIVE propagation layer
2936 return GetLayerNumber(0);
2940 //_____________________________________________________________________________
2941 Int_t AliTRDtracker::AliTRDtrackingSector
2942 ::GetOuterTimeBin() const
2945 // Returns number of the outermost SENSITIVE time bin
2948 return GetLayerNumber(GetNumberOfTimeBins() - 1);
2952 //_____________________________________________________________________________
2953 Int_t AliTRDtracker::AliTRDtrackingSector
2954 ::GetNumberOfTimeBins() const
2957 // Returns number of SENSITIVE time bins
2963 for (tb = kMaxTimeBinIndex - 1; tb >= 0; tb--) {
2964 layer = GetLayerNumber(tb);
2974 //_____________________________________________________________________________
2975 void AliTRDtracker::AliTRDtrackingSector
2976 ::InsertLayer(AliTRDpropagationLayer *pl)
2979 // Insert layer <pl> in fLayers array.
2980 // Layers are sorted according to X coordinate.
2983 if (fN == ((Int_t) kMaxLayersPerSector)) {
2984 //AliWarning("Too many layers !\n");
2993 Int_t i = Find(pl->GetX());
2995 memmove(fLayers+i+1,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
3002 //_____________________________________________________________________________
3003 Int_t AliTRDtracker::AliTRDtrackingSector
3004 ::Find(Double_t x) const
3007 // Returns index of the propagation layer nearest to X
3010 if (x <= fLayers[0]->GetX()) {
3014 if (x > fLayers[fN-1]->GetX()) {
3020 Int_t m = (b + e) / 2;
3022 for (; b < e; m = (b + e) / 2) {
3023 if (x > fLayers[m]->GetX()) {
3035 //_____________________________________________________________________________
3036 Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster *c)
3039 // Returns correction factor for tilted pads geometry
3042 Int_t det = c->GetDetector();
3043 Int_t plane = fGeom->GetPlane(det);
3044 AliTRDpadPlane *padPlane = fGeom->GetPadPlane(plane,0);
3045 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
3055 //_____________________________________________________________________________
3056 Int_t AliTRDtracker::FindClusters(Int_t sector, Int_t t0, Int_t t1
3057 , AliTRDtrack *track
3059 , AliTRDtracklet &tracklet)
3062 // Try to find the nearest clusters to the track in the time bins
3063 // between <t0> and <t1>.
3064 // Also the corresponding tracklet is calculated
3065 // Correction coeficients - depend on TRD parameters - to be changed accordingly
3068 const Int_t kN1 = 100;
3069 const Int_t kN2 = 10;
3077 Double_t dz[kN2][kN1];
3078 Double_t dy[kN2][kN1];
3079 Int_t indexes[kN2][kN1]; // Indexes of the clusters in the road
3080 Int_t best[kN2][kN1]; // Index of best matching cluster
3081 AliTRDcluster *cl[kN2][kN1]; // Pointers to the clusters in the road
3083 Double_t xmean = 0.0; // Reference x
3086 // Initialize the arrays
3087 for (Int_t it = 0; it < kN1; it++) {
3096 for (Int_t ih = 0; ih < kN2; ih++) {
3097 indexes[ih][it] = -2;
3099 dz[ih][it] = -100.0;
3100 dy[ih][it] = -100.0;
3106 Double_t x0 = track->GetX();
3107 Double_t sigmaz = TMath::Sqrt(TMath::Abs(track->GetSigmaZ2()));
3112 Int_t detector = -1;
3113 Float_t padlength = 0.0;
3115 AliTRDtrack track2(* track);
3116 Float_t snpy = track->GetSnp();
3117 Float_t tany = TMath::Sqrt(snpy*snpy / (1.0 - snpy*snpy));
3122 Double_t sy2 = ExpectedSigmaY2(x0,track->GetTgl(),track->GetSignedPt());
3123 Double_t sz2 = ExpectedSigmaZ2(x0,track->GetTgl());
3124 Double_t road = 15.0 * TMath::Sqrt(track->GetSigmaY2() + sy2);
3129 for (Int_t it = 0; it < t1-t0; it++) {
3131 Double_t maxChi2[2] = { fgkMaxChi2, fgkMaxChi2 };
3132 AliTRDpropagationLayer &timeBin = *(fTrSec[sector]->GetLayer(it+t0));
3134 continue; // No indexes1
3137 Int_t maxn = timeBin;
3138 x[it] = timeBin.GetX();
3139 track2.PropagateTo(x[it]);
3140 yt[it] = track2.GetY();
3141 zt[it] = track2.GetZ();
3143 Double_t y = yt[it];
3144 Double_t z = zt[it];
3145 Double_t chi2 = 1000000.0;
3149 // Find 2 nearest cluster at given time bin
3151 Int_t checkPoint[4] = { 0, 0, 0, 0 };
3152 Double_t minY = 123456789.0;
3153 Double_t minD[2] = { 1.0, 1.0 };
3155 for (Int_t i = timeBin.Find(y - road); i < maxn; i++) {
3157 AliTRDcluster *c = (AliTRDcluster *) (timeBin[i]);
3158 h01 = GetTiltFactor(c);
3160 Int_t det = c->GetDetector();
3161 plane = fGeom->GetPlane(det);
3162 padlength = TMath::Sqrt(c->GetSigmaZ2() * 12.0);
3165 if (c->GetY() > (y + road)) {
3169 fHDeltaX->Fill(c->GetX() - x[it]);
3171 if (TMath::Abs(c->GetY()-y) < TMath::Abs(minY)) {
3172 minY = c->GetY() - y;
3173 minD[0] = c->GetY() - y;
3174 minD[1] = c->GetZ() - z;
3179 fHMinZ->Fill(c->GetZ() - z);
3180 if ((c->GetZ() - z) * (c->GetZ() - z) > 2 * (12.0 * sz2)) {
3185 Double_t dist = TMath::Abs(c->GetZ() - z);
3186 if (dist > (0.5 * padlength + 6.0 * sigmaz)) {
3187 continue; // 6 sigma boundary cut
3191 // Sigma boundary cost function
3192 Double_t cost = 0.0;
3193 if (dist> (0.5 * padlength - sigmaz)){
3194 cost = (dist - 0.5*padlength) / (2.0 * sigmaz);
3196 cost = (cost + 1.0) * (cost + 1.0);
3202 chi2 = track2.GetPredictedChi2(c,h01) + cost;
3205 if (chi2 > maxChi2[1]) {
3210 // Store the clusters in the road
3211 detector = c->GetDetector();
3212 for (Int_t ih = 2; ih < 9; ih++) {
3213 if (cl[ih][it] == 0) {
3215 indexes[ih][it] = timeBin.GetIndex(i); // Index - 9 - reserved for outliers
3220 if (chi2 < maxChi2[0]) {
3221 maxChi2[1] = maxChi2[0];
3223 indexes[1][it] = indexes[0][it];
3224 cl[1][it] = cl[0][it];
3225 indexes[0][it] = timeBin.GetIndex(i);
3231 indexes[1][it] = timeBin.GetIndex(i);
3235 for(int iCheckPoint = 0; iCheckPoint<4; iCheckPoint++) {
3236 fHFindCl[iCheckPoint]->Fill(checkPoint[iCheckPoint]);
3239 if (checkPoint[3]) {
3240 if (track->GetSignedPt() > 0) {
3241 fHMinYPos->Fill(minY);
3244 fHMinYNeg->Fill(minY);
3246 fHMinD->Fill(minD[0],minD[1]);
3259 xmean /= Float_t(nfound); // Middle x
3260 track2.PropagateTo(xmean); // Propagate track to the center
3263 // Choose one of the variants
3267 Double_t sumdy = 0.0;
3268 Double_t sumdy2 = 0.0;
3269 Double_t sumx = 0.0;
3270 Double_t sumxy = 0.0;
3271 Double_t sumx2 = 0.0;
3272 Double_t mpads = 0.0;
3279 Double_t meanz[kN2];
3280 Double_t moffset[kN2]; // Mean offset
3281 Double_t mean[kN2]; // Mean value
3282 Double_t angle[kN2]; // Angle
3284 Double_t smoffset[kN2]; // Sigma of mean offset
3285 Double_t smean[kN2]; // Sigma of mean value
3286 Double_t sangle[kN2]; // Sigma of angle
3287 Double_t smeanangle[kN2]; // Correlation
3289 Double_t sigmas[kN2];
3290 Double_t tchi2s[kN2]; // Chi2s for tracklet
3292 for (Int_t it = 0; it < kN2; it++) {
3298 moffset[it] = 0.0; // Mean offset
3299 mean[it] = 0.0; // Mean value
3300 angle[it] = 0.0; // Angle
3302 smoffset[it] = 1.0e5; // Sigma of mean offset
3303 smean[it] = 1.0e5; // Sigma of mean value
3304 sangle[it] = 1.0e5; // Sigma of angle
3305 smeanangle[it] = 0.0; // Correlation
3308 tchi2s[it] = 1.0e5; // Chi2s for tracklet
3315 for (Int_t it = 0; it < t1 - t0; it++) {
3319 for (Int_t dt = -3; dt <= 3; dt++) {
3323 if (it+dt > t1-t0) {
3326 if (!cl[0][it+dt]) {
3329 zmean[it] += cl[0][it+dt]->GetZ();
3332 zmean[it] /= nmean[it];
3335 for (Int_t it = 0; it < t1 - t0; it++) {
3339 for (Int_t ih = 0; ih < 10; ih++) {
3340 dz[ih][it] = -100.0;
3341 dy[ih][it] = -100.0;
3345 Double_t xcluster = cl[ih][it]->GetX();
3348 track2.GetProlongation(xcluster,ytrack,ztrack );
3349 dz[ih][it] = cl[ih][it]->GetZ()- ztrack; // Calculate distance from track in z
3350 dy[ih][it] = cl[ih][it]->GetY() + dz[ih][it]*h01 - ytrack; // and in y
3357 if ((TMath::Abs(cl[0][it]->GetZ()-zmean[it]) > padlength * 0.8) &&
3359 if (TMath::Abs(cl[1][it]->GetZ()-zmean[it]) < padlength * 0.5) {
3367 // Iterative choice of "best path"
3369 Int_t label = TMath::Abs(track->GetLabel());
3372 for (Int_t iter = 0; iter < 9; iter++) {
3387 for (Int_t it = 0; it < t1 - t0; it++) {
3389 if (!cl[best[iter][it]][it]) {
3393 // Calculates pad-row changes
3394 Double_t zbefore = cl[best[iter][it]][it]->GetZ();
3395 Double_t zafter = cl[best[iter][it]][it]->GetZ();
3396 for (Int_t itd = it - 1; itd >= 0; itd--) {
3397 if (cl[best[iter][itd]][itd]) {
3398 zbefore = cl[best[iter][itd]][itd]->GetZ();
3402 for (Int_t itd = it + 1; itd < t1 - t0; itd++) {
3403 if (cl[best[iter][itd]][itd]) {
3404 zafter = cl[best[iter][itd]][itd]->GetZ();
3408 if ((TMath::Abs(cl[best[iter][it]][it]->GetZ()-zbefore) > 0.1) &&
3409 (TMath::Abs(cl[best[iter][it]][it]->GetZ()- zafter) > 0.1)) {
3413 // Distance to reference x
3414 Double_t dx = x[it]-xmean;
3415 sumz += cl[best[iter][it]][it]->GetZ();
3417 sumdy += dy[best[iter][it]][it];
3418 sumdy2 += dy[best[iter][it]][it]*dy[best[iter][it]][it];
3421 sumxy += dx*dy[best[iter][it]][it];
3422 mpads += cl[best[iter][it]][it]->GetNPads();
3423 if ((cl[best[iter][it]][it]->GetLabel(0) == label) ||
3424 (cl[best[iter][it]][it]->GetLabel(1) == label) ||
3425 (cl[best[iter][it]][it]->GetLabel(2) == label)) {
3435 // Calculates line parameters
3437 Double_t det = sum*sumx2 - sumx*sumx;
3438 angle[iter] = (sum*sumxy - sumx*sumdy) / det;
3439 mean[iter] = (sumx2*sumdy - sumx*sumxy) / det;
3440 meanz[iter] = sumz / sum;
3441 moffset[iter] = sumdy / sum;
3442 mpads /= sum; // Mean number of pads
3444 Double_t sigma2 = 0.0; // Normalized residuals - for line fit
3445 Double_t sigma1 = 0.0; // Normalized residuals - constant fit
3447 for (Int_t it = 0; it < t1 - t0; it++) {
3448 if (!cl[best[iter][it]][it]) {
3451 Double_t dx = x[it] - xmean;
3452 Double_t ytr = mean[iter] + angle[iter] * dx;
3453 sigma2 += (dy[best[iter][it]][it] - ytr)
3454 * (dy[best[iter][it]][it] - ytr);
3455 sigma1 += (dy[best[iter][it]][it] - moffset[iter])
3456 * (dy[best[iter][it]][it] - moffset[iter]);
3459 sigma2 /= (sum - 2); // Normalized residuals
3460 sigma1 /= (sum - 1); // Normalized residuals
3461 smean[iter] = sigma2 * (sumx2 / det); // Estimated error2 of mean
3462 sangle[iter] = sigma2 * ( sum / det); // Estimated error2 of angle
3463 smeanangle[iter] = sigma2 * (-sumx / det); // Correlation
3464 sigmas[iter] = TMath::Sqrt(sigma1);
3465 smoffset[iter] = (sigma1 / sum) + 0.01*0.01; // Sigma of mean offset + unisochronity sigma
3468 // Iterative choice of "better path"
3470 for (Int_t it = 0; it < t1 - t0; it++) {
3472 if (!cl[best[iter][it]][it]) {
3476 // Add unisochronity + angular effect contribution
3477 Double_t sigmatr2 = smoffset[iter] + 0.5*tany*tany;
3478 Double_t sweight = 1.0/sigmatr2 + 1.0/track->GetSigmaY2();
3479 Double_t weighty = (moffset[iter] / sigmatr2) / sweight; // Weighted mean
3480 Double_t sigmacl = TMath::Sqrt(sigma1*sigma1 + track->GetSigmaY2());
3481 Double_t mindist = 100000.0;
3484 for (Int_t ih = 0; ih < kN2; ih++) {
3488 Double_t dist2 = (dy[ih][it] - weighty) / sigmacl;
3489 dist2 *= dist2; // Chi2 distance
3490 if (dist2 < mindist) {
3496 best[iter+1][it] = ihbest;
3501 // Update best hypothesy if better chi2 according tracklet position and angle
3503 Double_t sy2 = smean[iter] + track->GetSigmaY2();
3504 Double_t sa2 = sangle[iter] + track->GetSigmaSnp2(); // track->fCee;
3505 Double_t say = track->GetSigmaSnpY(); // track->fCey;
3507 Double_t detchi = sy2*sa2 - say*say;
3508 Double_t invers[3] = {sa2/detchi,sy2/detchi,-say/detchi}; // Inverse value of covariance matrix
3510 Double_t chi20 = mean[bestiter] * mean[bestiter] * invers[0]
3511 + angle[bestiter] * angle[bestiter] * invers[1]
3512 + 2.0 * mean[bestiter] * angle[bestiter] * invers[2];
3513 Double_t chi21 = mean[iter] * mean[iter] * invers[0]
3514 + angle[iter] * angle[iter] * invers[1]
3515 + 2.0 * mean[iter] * angle[iter] * invers[2];
3516 tchi2s[iter] = chi21;
3518 if ((changes[iter] <= changes[bestiter]) &&
3528 Double_t sigma2 = sigmas[0]; // Choose as sigma from 0 iteration
3529 Short_t maxpos = -1;
3530 Float_t maxcharge = 0.0;
3531 Short_t maxpos4 = -1;
3532 Float_t maxcharge4 = 0.0;
3533 Short_t maxpos5 = -1;
3534 Float_t maxcharge5 = 0.0;
3536 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(AliTRDcalibDB::Instance()->GetVdrift(0,0,0)
3537 ,-AliTracker::GetBz()*0.1);
3538 Double_t expectederr = sigma2*sigma2 + 0.01*0.01;
3540 expectederr += (mpads - 3.5) * 0.04;
3542 if (changes[bestiter] > 1) {
3543 expectederr += changes[bestiter] * 0.01;
3545 expectederr += (0.03 * (tany-exB)*(tany-exB)) * 15.0;
3547 for (Int_t it = 0; it < t1 - t0; it++) {
3549 if (!cl[best[bestiter][it]][it]) {
3553 // Set cluster error
3554 cl[best[bestiter][it]][it]->SetSigmaY2(expectederr);
3555 if (!cl[best[bestiter][it]][it]->IsUsed()) {
3556 cl[best[bestiter][it]][it]->SetY(cl[best[bestiter][it]][it]->GetY());
3559 // Time bins with maximal charge
3560 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge) {
3561 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3562 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3565 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge4) {
3566 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 4) {
3567 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3568 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3572 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge5) {
3573 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 5) {
3574 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3575 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3579 // Time bins with maximal charge
3580 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge) {
3581 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3582 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3585 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge4) {
3586 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 4) {
3587 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3588 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3592 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ()) > maxcharge5) {
3593 if (cl[best[bestiter][it]][it]->GetLocalTimeBin() >= 5) {
3594 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3595 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3599 clusters[it+t0] = indexes[best[bestiter][it]][it];
3603 // Set tracklet parameters
3604 Double_t trackleterr2 = smoffset[bestiter] + 0.01*0.01;
3606 trackleterr2 += (mpads - 3.5) * 0.04;
3608 trackleterr2 += changes[bestiter] * 0.01;
3609 trackleterr2 *= TMath::Max(14.0 - nfound,1.0);
3610 trackleterr2 += 0.2 * (tany-exB)*(tany-exB);
3613 ,track2.GetY() + moffset[bestiter]
3617 tracklet.SetTilt(h01);
3618 tracklet.SetP0(mean[bestiter]);
3619 tracklet.SetP1(angle[bestiter]);
3620 tracklet.SetN(nfound);
3621 tracklet.SetNCross(changes[bestiter]);
3622 tracklet.SetPlane(plane);
3623 tracklet.SetSigma2(expectederr);
3624 tracklet.SetChi2(tchi2s[bestiter]);
3625 tracklet.SetMaxPos(maxpos,maxpos4,maxpos5);
3626 track->SetTracklets(plane,tracklet);
3627 track->SetNWrong(track->GetNWrong() + nbad[0]);
3632 TClonesArray array0("AliTRDcluster");
3633 TClonesArray array1("AliTRDcluster");
3634 array0.ExpandCreateFast(t1 - t0 + 1);
3635 array1.ExpandCreateFast(t1 - t0 + 1);
3636 TTreeSRedirector &cstream = *fDebugStreamer;
3637 AliTRDcluster dummy;
3641 for (Int_t it = 0; it < t1 - t0; it++) {
3642 dy0[it] = dy[0][it];
3643 dyb[it] = dy[best[bestiter][it]][it];
3645 new(array0[it]) AliTRDcluster(*cl[0][it]);
3648 new(array0[it]) AliTRDcluster(dummy);
3650 if(cl[best[bestiter][it]][it]) {
3651 new(array1[it]) AliTRDcluster(*cl[best[bestiter][it]][it]);
3654 new(array1[it]) AliTRDcluster(dummy);
3658 TGraph graph0(t1-t0,x,dy0);
3659 TGraph graph1(t1-t0,x,dyb);
3660 TGraph graphy(t1-t0,x,yt);
3661 TGraph graphz(t1-t0,x,zt);
3663 if (AliTRDReconstructor::StreamLevel() > 0) {
3664 cstream << "tracklet"
3665 << "track.=" << track // Track parameters
3666 << "tany=" << tany // Tangent of the local track angle
3667 << "xmean=" << xmean // Xmean - reference x of tracklet
3668 << "tilt=" << h01 // Tilt angle
3669 << "nall=" << nall // Number of foundable clusters
3670 << "nfound=" << nfound // Number of found clusters
3671 << "clfound=" << clfound // Total number of found clusters in road
3672 << "mpads=" << mpads // Mean number of pads per cluster
3673 << "plane=" << plane // Plane number
3674 << "detector=" << detector // Detector number
3675 << "road=" << road // The width of the used road
3676 << "graph0.=" << &graph0 // x - y = dy for closest cluster
3677 << "graph1.=" << &graph1 // x - y = dy for second closest cluster
3678 << "graphy.=" << &graphy // y position of the track
3679 << "graphz.=" << &graphz // z position of the track
3680 //<< "fCl.=" << &array0 // closest cluster
3681 //<< "fCl2.=" << &array1 // second closest cluster
3682 << "maxpos=" << maxpos // Maximal charge postion
3683 << "maxcharge=" << maxcharge // Maximal charge
3684 << "maxpos4=" << maxpos4 // Maximal charge postion - after bin 4
3685 << "maxcharge4=" << maxcharge4 // Maximal charge - after bin 4
3686 << "maxpos5=" << maxpos5 // Maximal charge postion - after bin 5
3687 << "maxcharge5=" << maxcharge5 // Maximal charge - after bin 5
3688 << "bestiter=" << bestiter // Best iteration number
3689 << "tracklet.=" << &tracklet // Corrspond to the best iteration
3690 << "tchi20=" << tchi2s[0] // Chi2 of cluster in the 0 iteration
3691 << "tchi2b=" << tchi2s[bestiter] // Chi2 of cluster in the best iteration
3692 << "sigmas0=" << sigmas[0] // Residuals sigma
3693 << "sigmasb=" << sigmas[bestiter] // Residulas sigma
3694 << "ngood0=" << ngood[0] // Number of good clusters in 0 iteration
3695 << "nbad0=" << nbad[0] // Number of bad clusters in 0 iteration
3696 << "ngoodb=" << ngood[bestiter] // in best iteration
3697 << "nbadb=" << nbad[bestiter] // in best iteration
3698 << "changes0=" << changes[0] // Changes of pardrows in iteration number 0
3699 << "changesb=" << changes[bestiter] // Changes of pardrows in best iteration
3700 << "moffset0=" << moffset[0] // Offset fixing angle in iter=0
3701 << "smoffset0=" << smoffset[0] // Sigma of offset fixing angle in iter=0
3702 << "moffsetb=" << moffset[bestiter] // Offset fixing angle in iter=best
3703 << "smoffsetb=" << smoffset[bestiter] // Sigma of offset fixing angle in iter=best
3704 << "mean0=" << mean[0] // Mean dy in iter=0;
3705 << "smean0=" << smean[0] // Sigma of mean dy in iter=0
3706 << "meanb=" << mean[bestiter] // Mean dy in iter=best
3707 << "smeanb=" << smean[bestiter] // Sigma of mean dy in iter=best
3708 << "angle0=" << angle[0] // Angle deviation in the iteration number 0
3709 << "sangle0=" << sangle[0] // Sigma of angular deviation in iteration number 0
3710 << "angleb=" << angle[bestiter] // Angle deviation in the best iteration
3711 << "sangleb=" << sangle[bestiter] // Sigma of angle deviation in the best iteration
3712 << "expectederr=" << expectederr // Expected error of cluster position
3720 //_____________________________________________________________________________
3721 Int_t AliTRDtracker::Freq(Int_t n, const Int_t *inlist
3722 , Int_t *outlist, Bool_t down)
3725 // Sort eleements according occurancy
3726 // The size of output array has is 2*n
3733 Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
3734 Int_t *sindexF = new Int_t[2*n];
3735 for (Int_t i = 0; i < n; i++) {
3739 TMath::Sort(n,inlist,sindexS,down);
3741 Int_t last = inlist[sindexS[0]];
3744 sindexF[0+n] = last;
3748 for (Int_t i = 1; i < n; i++) {
3749 val = inlist[sindexS[i]];
3751 sindexF[countPos]++;
3755 sindexF[countPos+n] = val;
3756 sindexF[countPos]++;
3764 // Sort according frequency
3765 TMath::Sort(countPos,sindexF,sindexS,kTRUE);
3767 for (Int_t i = 0; i < countPos; i++) {
3768 outlist[2*i ] = sindexF[sindexS[i]+n];
3769 outlist[2*i+1] = sindexF[sindexS[i]];
3779 //_____________________________________________________________________________
3780 AliTRDtrack *AliTRDtracker::RegisterSeed(AliTRDseed *seeds, Double_t *params)
3783 // Build a TRD track out of tracklet candidates
3786 // seeds : array of tracklets
3787 // params : track parameters (see MakeSeeds() function body for a detailed description)
3792 // Detailed description
3794 // To be discussed with Marian !!
3797 AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
3798 Int_t nTimeBins = cal->GetNumberOfTimeBins();
3801 Double_t alpha = AliTRDgeometry::GetAlpha();
3802 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
3806 c[ 1] = 0.0; c[ 2] = 2.0;
3807 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
3808 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
3809 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
3812 AliTRDcluster *cl = 0;
3814 for (Int_t ilayer = 0; ilayer < 6; ilayer++) {
3816 if (seeds[ilayer].IsOK()) {
3817 for (Int_t itime = nTimeBins-1; itime > 0; itime--) {
3818 if (seeds[ilayer].GetIndexes(itime) > 0) {
3819 index = seeds[ilayer].GetIndexes(itime);
3820 cl = seeds[ilayer].GetClusters(itime);
3821 //printf("l[%d] index[%d] tb[%d] cptr[%p]\n", ilayer, index, itime, cl);
3834 AliTRDtrack *track = new AliTRDtrack(cl
3839 ,params[6]*alpha+shift);
3840 // SetCluster(cl, 0); // A. Bercuci
3841 track->PropagateTo(params[0]-5.0);
3842 track->ResetCovariance(1);
3844 Int_t rc = FollowBackProlongation(*track);
3851 track->CookdEdxTimBin(-1);
3852 CookLabel(track,0.9);
3859 //_____________________________________________________________________________
3860 void AliTRDtracker::InitLogHists()
3863 // Create the log histograms
3866 fHBackfit = new TH1D("logTRD_backfit" ,""
3868 fHRefit = new TH1D("logTRD_refit" ,""
3870 fHClSearch = new TH1D("logTRD_clSearch",""
3873 fHX = new TH1D("logTRD_X" ,";x (cm)"
3875 fHNCl = new TH1D("logTRD_ncl" ,""
3877 fHNClTrack = new TH1D("logTRD_nclTrack",""
3880 fHMinYPos = new TH1D("logTRD_minYPos" ,";#delta Y (cm)"
3882 fHMinYNeg = new TH1D("logTRD_minYNeg" ,";#delta Y (cm)"
3884 fHMinZ = new TH1D("logTRD_minZ" ,";#delta Z (cm)"
3886 fHMinD = new TH2D("logTRD_minD" ,";#delta Y (cm);#delta Z (cm)"
3890 fHDeltaX = new TH1D("logTRD_deltaX" ,";#delta X (cm)"
3892 fHXCl = new TH1D("logTRD_xCl" ,";cluster x position (cm)"
3895 const Char_t *nameFindCl[4] = { "logTRD_clY"
3900 for (Int_t i = 0; i < 4; i++) {
3901 fHFindCl[i] = new TH1D(nameFindCl[i],"",30,-0.5,29.5);
3906 //_____________________________________________________________________________
3907 void AliTRDtracker::SaveLogHists()
3910 // Save the log histograms in AliESDs.root
3913 TDirectory *sav = gDirectory;
3916 TSeqCollection *col = gROOT->GetListOfFiles();
3917 Int_t nn = col->GetEntries();
3918 for (Int_t i = 0; i < nn; i++) {
3919 logFile = (TFile *) col->At(i);
3920 if (strstr(logFile->GetName(),"AliESDs.root")) {
3927 fHBackfit->Write(fHBackfit->GetName(),TObject::kOverwrite);
3928 fHRefit->Write(fHRefit->GetName(),TObject::kOverwrite);
3929 fHClSearch->Write(fHClSearch->GetName(),TObject::kOverwrite);
3930 fHX->Write(fHX->GetName(),TObject::kOverwrite);
3931 fHNCl->Write(fHNCl->GetName(),TObject::kOverwrite);
3932 fHNClTrack->Write(fHNClTrack->GetName(),TObject::kOverwrite);
3934 fHMinYPos->Write(fHMinYPos->GetName(),TObject::kOverwrite);
3935 fHMinYNeg->Write(fHMinYNeg->GetName(),TObject::kOverwrite);
3936 fHMinD->Write(fHMinD->GetName(),TObject::kOverwrite);
3937 fHMinZ->Write(fHMinZ->GetName(),TObject::kOverwrite);
3939 fHDeltaX->Write(fHDeltaX->GetName(),TObject::kOverwrite);
3940 fHXCl->Write(fHXCl->GetName(),TObject::kOverwrite);
3942 for (Int_t i = 0; i < 4; i++) {
3943 fHFindCl[i]->Write(fHFindCl[i]->GetName(),TObject::kOverwrite);