1 #ifndef ALITPCTRACKER_H
2 #define ALITPCTRACKER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
9 //-------------------------------------------------------
14 //-------------------------------------------------------
17 #include "AliTracker.h"
18 #include "AliTPCreco.h"
19 #include "AliTPCclusterMI.h"
20 #include "AliTPCtrackerSector.h"
26 class AliTPCTrackerPoint;
31 class TTreeSRedirector;
33 class AliDCSSensorArray;
38 class AliTPCtracker : public AliTracker {
41 AliTPCtracker(const AliTPCParam *par);
42 virtual ~AliTPCtracker();
44 void SetIteration(Int_t iteration){fIteration = iteration;}
45 virtual Int_t Clusters2TracksHLT(AliESDEvent *const esd, const AliESDEvent *hltEvent);
46 virtual Int_t Clusters2Tracks (AliESDEvent *const esd);
47 virtual Int_t RefitInward (AliESDEvent *esd);
48 virtual Int_t LoadClusters (TTree * const tree);
49 virtual Int_t LoadClusters (const TObjArray * arr); // another input
50 virtual Int_t LoadClusters (const TClonesArray * arr); // another input
51 virtual Int_t PostProcess(AliESDEvent *esd);
53 void UnloadClusters();
54 Int_t LoadInnerSectors();
55 Int_t LoadOuterSectors();
56 virtual void FillClusterArray(TObjArray* array) const;
57 void Transform(AliTPCclusterMI * cluster);
58 void ApllyTailCancellation();
60 void FillESD(const TObjArray* arr);
62 void SetDebug(Int_t debug){ fDebug = debug;}
63 void FindKinks(TObjArray * array, AliESDEvent * esd);
65 void FindCurling(const TObjArray * array, AliESDEvent * esd, Int_t iter);
66 void FindSplitted(TObjArray * array, AliESDEvent * esd, Int_t iter);
67 void FindMultiMC(const TObjArray * array, AliESDEvent * esd, Int_t iter);
69 void UpdateKinkQualityM(AliTPCseed * seed);
70 void UpdateKinkQualityD(AliTPCseed * seed);
71 Int_t CheckKinkPoint(AliTPCseed*seed, AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
72 Int_t RefitKink(AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
73 Int_t ReadSeeds(const TFile *in);
74 TObjArray * GetSeeds() const {return fSeeds;}
76 AliCluster * GetCluster(Int_t index) const {return (AliCluster*)GetClusterMI(index);}
77 AliTPCclusterMI *GetClusterMI(Int_t index) const;
78 Int_t Clusters2Tracks();
79 virtual void CookLabel(AliKalmanTrack *tk,Float_t wrong) const;
80 virtual Int_t CookLabel(AliTPCseed *const t,Float_t wrong, Int_t first,Int_t last ) const;
82 void RotateToLocal(AliTPCseed *seed);
84 Int_t FollowProlongation(AliTPCseed& t, Int_t rf=0, Int_t step=1, Bool_t fromSeeds=0);
85 Bool_t GetTrackPoint(Int_t index, AliTrackPoint &p ) const;
87 Int_t FollowBackProlongation(AliTPCseed& t, Int_t rf, Bool_t fromSeeds=0);
88 Int_t FollowToNext(AliTPCseed& t, Int_t nr);
89 Int_t UpdateClusters(AliTPCseed& t, Int_t nr);
90 Int_t FollowToNextCluster( AliTPCseed& t, Int_t nr);
92 Int_t PropagateBack(const TObjArray *const arr);
93 Int_t PropagateBack(AliESDEvent * event);
94 Int_t PropagateBack(AliTPCseed *const pt, Int_t row0, Int_t row1);
95 Int_t PropagateForward();
96 Int_t PropagateForward2(const TObjArray *const arr);
98 void SortTracks(TObjArray * arr, Int_t mode) const;
100 virtual Double_t ErrY2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
101 virtual Double_t ErrZ2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
103 Double_t F1(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
104 Double_t F1old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
105 Double_t F2(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
106 Double_t F2old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
108 Double_t F3(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2) const;
109 Double_t F3n(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2,
111 Bool_t GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z) const;
113 void ResetSeedsPool();
114 void MarkSeedFree( TObject* seed );
115 TObject *&NextFreeSeed();
118 void SetUseHLTClusters(Int_t useHLTClusters) {fUseHLTClusters = useHLTClusters;} // set usage from HLT clusters from rec.C options
120 inline void SetTPCtrackerSectors(AliTPCtrackerSector *innerSec, AliTPCtrackerSector *outerSec); // set the AliTPCtrackerSector arrays from outside (toy MC)
122 Float_t OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t &sum2);
123 void SignShared(AliTPCseed * s1, AliTPCseed * s2);
124 void SignShared(TObjArray * arr);
126 void RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal);
128 Int_t AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster);
130 Bool_t IsTPCHVDipEvent(AliESDEvent const *esdEvent);
132 // public for ToyMC usage
133 void MakeSeeds2(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1, Bool_t bconstrain=kTRUE);
134 void MakeSeeds3(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1, Int_t ddsec=0);
135 void SumTracks(TObjArray *arr1,TObjArray *&arr2);
136 void SignClusters(const TObjArray * arr, Float_t fnumber=3., Float_t fdensity=2.);
139 Bool_t IsFindable(AliTPCseed & t);
140 AliTPCtracker(const AliTPCtracker& r); //dummy copy constructor
141 AliTPCtracker &operator=(const AliTPCtracker& r);//dummy assignment operator
142 void AddCovariance(AliTPCseed * seed); // add covariance
143 void AddCovarianceAdd(AliTPCseed * seed); // add covariance
145 inline AliTPCtrackerRow &GetRow(Int_t sec, Int_t row);
146 inline Bool_t IsActive(Int_t sec, Int_t row);
147 inline Double_t GetXrow(Int_t row) const;
148 inline Double_t GetMaxY(Int_t row) const;
149 inline Int_t GetRowNumber(Double_t x) const;
150 Int_t GetRowNumber(Double_t x[3]) const;
151 inline Double_t GetPadPitchLength(Double_t x) const;
152 inline Double_t GetPadPitchLength(Int_t row) const;
154 void GetShape(AliTPCseed * seed, Int_t row);
156 void ReadSeeds(const AliESDEvent *const event, Int_t direction); //read seeds from the event
158 void MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1);
161 AliTPCseed *MakeSeed(AliTPCseed *const track, Float_t r0, Float_t r1, Float_t r2); //reseed
162 AliTPCseed *ReSeed(const AliTPCseed *t, Float_t r0, Float_t r1, Float_t r2); //reseed
163 AliTPCseed *ReSeed(AliTPCseed *t, Int_t r0, Bool_t forward); //reseed
167 AliTPCseed * ReSeed(AliTPCseed *t);
168 //Int_t LoadInnerSectors();
169 //Int_t LoadOuterSectors();
170 void DumpClusters(Int_t iter, TObjArray *trackArray);
171 void UnsignClusters();
173 void ParallelTracking(TObjArray *const arr, Int_t rfirst, Int_t rlast);
174 void Tracking(TObjArray * arr);
175 TObjArray * Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy=-1, Int_t dsec=0);
176 TObjArray * Tracking();
177 TObjArray * TrackingSpecial();
178 void PrepareForBackProlongation(const TObjArray *const arr, Float_t fac) const;
179 void PrepareForProlongation(TObjArray *const arr, Float_t fac) const;
181 Int_t UpdateTrack(AliTPCseed *t, Int_t accept); //update trackinfo
183 void MakeESDBitmaps(AliTPCseed *t, AliESDtrack *esd);
185 Int_t PropagateToRowHLT(AliTPCseed *pt, int nrow);
186 void TrackFollowingHLT(TObjArray *const arr);
187 TObjArray * MakeSeedsHLT(const AliESDEvent *hltEvent);
189 const Int_t fkNIS; //number of inner sectors
190 AliTPCtrackerSector *fInnerSec; //array of inner sectors;
191 const Int_t fkNOS; //number of outer sectors
192 AliTPCtrackerSector *fOuterSec; //array of outer sectors;
194 Int_t fN; //number of loaded sectors
195 AliTPCtrackerSector *fSectors; //pointer to loaded sectors;
197 TTree * fInput; // input tree with clusters
198 TTree * fOutput; // output tree with tracks
199 TTree * fSeedTree; // output tree with seeds - filled in debug mode 1
200 TTree * fTreeDebug; // output with a debug information about track
201 AliESDEvent * fEvent; // output with esd tracks
202 const AliESDEvent * fEventHLT; // input with HLT tracks
203 Int_t fDebug; // debug option
204 Bool_t fNewIO; // indicated if we have data using New IO
205 Int_t fNtracks; //current number of tracks
206 TObjArray *fSeeds; //array of track seeds
207 Int_t fIteration; // indicate iteration - 0 - froward -1 back - 2forward - back->forward
208 // TObjArray * fTrackPointPool; // ! pool with track points
209 Double_t fXRow[200]; // radius of the pad row
210 Double_t fYMax[200]; // max y for given pad row
211 Double_t fPadLength[200]; // max y for given pad row
212 const AliTPCParam *fkParam; //pointer to the parameters
213 TTreeSRedirector *fDebugStreamer; //!debug streamer
214 Int_t fUseHLTClusters; // use HLT clusters instead of offline clusters
216 TClonesArray* fSeedsPool; //! pool of seeds
217 TArrayI fFreeSeedsID; //! array of ID's of freed seeds
218 Int_t fNFreeSeeds; //! number of seeds freed in the pool
219 Int_t fLastSeedID; //! id of the pool seed on which is returned by the NextFreeSeed method
221 ClassDef(AliTPCtracker,4)
225 AliTPCtrackerRow & AliTPCtracker::GetRow(Int_t sec, Int_t row)
228 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()]:fInnerSec[sec][row];
231 Bool_t AliTPCtracker::IsActive(Int_t sec, Int_t row)
234 // check if the given sector row is active
236 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()].GetN()>0:fInnerSec[sec][row].GetN()>0;
240 Double_t AliTPCtracker::GetXrow(Int_t row) const {
241 // return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetX(row-fInnerSec->GetNRows()):fInnerSec->GetX(row);
245 Double_t AliTPCtracker::GetMaxY(Int_t row) const {
246 //return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetMaxY(row-fInnerSec->GetNRows()):fInnerSec->GetMaxY(row);
250 Int_t AliTPCtracker::GetRowNumber(Double_t x) const
253 return (x>133.) ? fOuterSec->GetRowNumber(x)+fInnerSec->GetNRows():fInnerSec->GetRowNumber(x);
256 Double_t AliTPCtracker::GetPadPitchLength(Double_t x) const
259 return (x>133.) ? fOuterSec->GetPadPitchLength(x):fInnerSec->GetPadPitchLength(x);
260 //return fPadLength[row];
263 Double_t AliTPCtracker::GetPadPitchLength(Int_t row) const
266 return fPadLength[row];
269 void AliTPCtracker::SetTPCtrackerSectors(AliTPCtrackerSector *innerSec, AliTPCtrackerSector *outerSec)
272 fInnerSec = innerSec;
273 fOuterSec = outerSec;