2 #ifndef ALITPCTRACKER_H
3 #define ALITPCTRACKER_H
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
10 //-------------------------------------------------------
15 //-------------------------------------------------------
19 #include "AliTracker.h"
20 #include "AliTPCreco.h"
21 #include "AliTPCclusterMI.h"
22 #include "AliTPCtrackerSector.h"
23 #include "AliESDfriend.h"
30 class AliTPCTrackerPoint;
35 class TTreeSRedirector;
37 class AliDCSSensorArray;
42 class AliTPCtracker : public AliTracker {
44 enum EStreamFlags{ // flags to store addition data/code debugging infomation which is not stored in ESD but in specaial TPCdebug.root file
45 kStreamTransform =0x00001, // flag:stream cluster transformation
46 kStreamErrParam =0x00002, // flag:stream in debug mode cluster and track extrapolation at given row together with error nad shape estimate
47 kStreamFilterClusterInfo =0x00004, // flag:stream TPC data ouliers filtering information
48 kStreamClDump =0x00008, // flag:stream clusters at the end of process (signed with useage flags)
49 kStreamRemoveUsed =0x00010, // flag:stream information about TPC tracks which were descarded (double track removal)
50 kStreamRemoveShort =0x00020, // flag:stream information about TPC tracks which were discarded (short track removal)
51 kStreamSplitted2 =0x00040, // flag:stream information about discarded TPC tracks pair algorithm
52 kStreamFillESD =0x00080, // flag: stream track information in FillESD function (after track Iteration 0)
53 kStreamCPropagateBack =0x00100, // flag: stream track information in PropagateBack function (after tracking Iteration 1)
54 kStreamRecoverBack =0x00200, // flag: stream track information for track failing PropagateBack function and recovered back
55 kStreamRefitInward =0x00400, // flag: stream track information in RefitInward function (after tracking Iteration 2)
56 kStreamRecoverIn =0x00800, // flag: stream track information for track failing in RefitInward function and recovered back
57 kStreamUpdateTrack =0x01000, // flag: stream track/cluster infroamtion in track update method
59 kStreamCrosstalkMatrix =0x02000, // flag: stream crosstalk matrix as used in the reconstruction at given region of TPC
60 kStreamXtalk =0x04000, // flag: stream crosstalk correction as applied to cluster
61 kStreamIonTail =0x08000, // flag: stream ion tail correction as applied to cluster
62 kStreamFindMultiMC =0x10000, // flag: stream MC infomation about the multiple find track (ONLY for MC data)
63 kStreamFindCurling =0x20000, // flag: stream track infroamtion in the FindCurling tracks method
64 kStreamFindKinks =0x40000, // flag: stream track infroamtion in the FindKinks method
67 AliTPCtracker(const AliTPCParam *par);
68 virtual ~AliTPCtracker();
70 void SetIteration(Int_t iteration){fIteration = iteration;}
71 virtual Int_t Clusters2TracksHLT(AliESDEvent *const esd, const AliESDEvent *hltEvent);
72 virtual Int_t Clusters2Tracks (AliESDEvent *const esd);
73 virtual Int_t RefitInward (AliESDEvent *esd);
74 virtual Int_t LoadClusters (TTree * const tree);
75 virtual Int_t LoadClusters (const TObjArray * arr); // another input
76 virtual Int_t LoadClusters (const TClonesArray * arr); // another input
77 void FilterOutlierClusters(); // filter outlier clusters
78 virtual Int_t PostProcess(AliESDEvent *esd);
80 void UnloadClusters();
81 Int_t LoadInnerSectors();
82 Int_t LoadOuterSectors();
83 virtual void FillClusterArray(TObjArray* array) const;
84 void Transform(AliTPCclusterMI * cluster);
85 void ApplyTailCancellation();
86 void ApplyXtalkCorrection();
87 void CalculateXtalkCorrection();
88 void GetTailValue(Float_t ampfactor,Double_t &ionTailMax,Double_t &ionTailTotal,TGraphErrors **graphRes,Float_t *indexAmpGraphs,AliTPCclusterMI *cl0,AliTPCclusterMI *cl1);
90 void FillESD(const TObjArray* arr);
92 void SetDebug(Int_t debug){ fDebug = debug;}
93 void FindKinks(TObjArray * array, AliESDEvent * esd);
95 void FindCurling(const TObjArray * array, AliESDEvent * esd, Int_t iter);
96 void FindSplitted(TObjArray * array, AliESDEvent * esd, Int_t iter);
97 void FindMultiMC(const TObjArray * array, AliESDEvent * esd, Int_t iter);
99 void UpdateKinkQualityM(AliTPCseed * seed);
100 void UpdateKinkQualityD(AliTPCseed * seed);
101 Int_t CheckKinkPoint(AliTPCseed*seed, AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
102 Int_t RefitKink(AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
103 Int_t ReadSeeds(const TFile *in);
104 TObjArray * GetSeeds() const {return fSeeds;}
105 void SetSeeds(TObjArray * seeds) { fSeeds = seeds;}
107 AliCluster * GetCluster(Int_t index) const {return (AliCluster*)GetClusterMI(index);}
108 AliTPCclusterMI *GetClusterMI(Int_t index) const;
109 Int_t Clusters2Tracks();
110 virtual void CookLabel(AliKalmanTrack *tk,Float_t wrong) const;
111 virtual Int_t CookLabel(AliTPCseed *const t,Float_t wrong, Int_t first,Int_t last ) const;
113 void RotateToLocal(AliTPCseed *seed);
115 Int_t FollowProlongation(AliTPCseed& t, Int_t rf=0, Int_t step=1, Bool_t fromSeeds=0);
116 Bool_t GetTrackPoint(Int_t index, AliTrackPoint &p ) const;
118 Int_t FollowBackProlongation(AliTPCseed& t, Int_t rf, Bool_t fromSeeds=0);
119 Int_t FollowToNext(AliTPCseed& t, Int_t nr);
120 Int_t UpdateClusters(AliTPCseed& t, Int_t nr);
121 Int_t FollowToNextCluster( AliTPCseed& t, Int_t nr);
123 Int_t PropagateBack(const TObjArray *const arr);
124 Int_t PropagateBack(AliESDEvent * event);
125 Int_t PropagateBack(AliTPCseed *const pt, Int_t row0, Int_t row1);
126 Int_t PropagateForward();
127 Int_t PropagateForward2(const TObjArray *const arr);
129 void SortTracks(TObjArray * arr, Int_t mode) const;
131 virtual Double_t ErrY2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
132 virtual Double_t ErrZ2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
134 Double_t F1(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
135 Double_t F1old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
136 Double_t F2(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
137 Double_t F2old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
139 Double_t F3(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2) const;
140 Double_t F3n(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2,
142 Bool_t GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z) const;
144 void ResetSeedsPool();
145 void MarkSeedFree( TObject* seed );
146 TObject *&NextFreeSeed();
149 void SetUseHLTClusters(Int_t useHLTClusters) {fUseHLTClusters = useHLTClusters;} // set usage from HLT clusters from rec.C options
151 inline void SetTPCtrackerSectors(AliTPCtrackerSector *innerSec, AliTPCtrackerSector *outerSec); // set the AliTPCtrackerSector arrays from outside (toy MC)
153 Float_t OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t &sum2);
154 void SignShared(AliTPCseed * s1, AliTPCseed * s2);
155 void SignShared(TObjArray * arr);
157 void RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal);
159 Int_t AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster);
161 Bool_t IsTPCHVDipEvent(AliESDEvent const *esdEvent);
163 // public for ToyMC usage
164 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);
165 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);
166 void SumTracks(TObjArray *arr1,TObjArray *&arr2);
167 void SignClusters(const TObjArray * arr, Float_t fnumber=3., Float_t fdensity=2.);
170 Bool_t IsFindable(AliTPCseed & t);
171 AliTPCtracker(const AliTPCtracker& r); //dummy copy constructor
172 AliTPCtracker &operator=(const AliTPCtracker& r);//dummy assignment operator
173 void AddCovariance(AliTPCseed * seed); // add covariance
174 void AddCovarianceAdd(AliTPCseed * seed); // add covariance
176 inline AliTPCtrackerRow &GetRow(Int_t sec, Int_t row);
177 inline Bool_t IsActive(Int_t sec, Int_t row);
178 inline Double_t GetXrow(Int_t row) const;
179 inline Double_t GetMaxY(Int_t row) const;
180 inline Int_t GetRowNumber(Double_t x) const;
181 Int_t GetRowNumber(Double_t x[3]) const;
182 inline Double_t GetPadPitchLength(Double_t x) const;
183 inline Double_t GetPadPitchLength(Int_t row) const;
185 void GetShape(AliTPCseed * seed, Int_t row);
187 void ReadSeeds(const AliESDEvent *const event, Int_t direction); //read seeds from the event
189 void MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1);
192 AliTPCseed *MakeSeed(AliTPCseed *const track, Float_t r0, Float_t r1, Float_t r2); //reseed
193 AliTPCseed *ReSeed(const AliTPCseed *t, Float_t r0, Float_t r1, Float_t r2); //reseed
194 AliTPCseed *ReSeed(AliTPCseed *t, Int_t r0, Bool_t forward); //reseed
198 AliTPCseed * ReSeed(AliTPCseed *t);
199 //Int_t LoadInnerSectors();
200 //Int_t LoadOuterSectors();
201 void DumpClusters(Int_t iter, TObjArray *trackArray);
202 void UnsignClusters();
204 void FillClusterOccupancyInfo();
206 void ParallelTracking(TObjArray *const arr, Int_t rfirst, Int_t rlast);
207 void Tracking(TObjArray * arr);
208 TObjArray * Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy=-1, Int_t dsec=0);
209 TObjArray * Tracking();
210 TObjArray * TrackingSpecial();
211 void PrepareForBackProlongation(const TObjArray *const arr, Float_t fac) const;
212 void PrepareForProlongation(TObjArray *const arr, Float_t fac) const;
214 Int_t UpdateTrack(AliTPCseed *t, Int_t accept); //update trackinfo
216 void MakeESDBitmaps(AliTPCseed *t, AliESDtrack *esd);
218 Int_t PropagateToRowHLT(AliTPCseed *pt, int nrow);
219 void TrackFollowingHLT(TObjArray *const arr);
220 TObjArray * MakeSeedsHLT(const AliESDEvent *hltEvent);
222 const Int_t fkNIS; //number of inner sectors
223 AliTPCtrackerSector *fInnerSec; //array of inner sectors;
224 const Int_t fkNOS; //number of outer sectors
225 AliTPCtrackerSector *fOuterSec; //array of outer sectors;
227 Int_t fN; //number of loaded sectors
228 AliTPCtrackerSector *fSectors; //pointer to loaded sectors;
230 TTree * fInput; // input tree with clusters
231 TTree * fOutput; // output tree with tracks
232 TTree * fSeedTree; // output tree with seeds - filled in debug mode 1
233 TTree * fTreeDebug; // output with a debug information about track
234 AliESDEvent * fEvent; // output with esd tracks
235 const AliESDEvent * fEventHLT; // input with HLT tracks
236 Int_t fDebug; // debug option
237 Bool_t fNewIO; // indicated if we have data using New IO
238 Int_t fNtracks; //current number of tracks
239 TObjArray *fSeeds; //array of track seeds
240 Int_t fIteration; // indicate iteration - 0 - froward -1 back - 2forward - back->forward
241 // TObjArray * fTrackPointPool; // ! pool with track points
242 Double_t fXRow[200]; // radius of the pad row
243 Double_t fYMax[200]; // max y for given pad row
244 Double_t fPadLength[200]; // max y for given pad row
245 const AliTPCParam *fkParam; //pointer to the parameters
246 TTreeSRedirector *fDebugStreamer; //!debug streamer
247 Int_t fUseHLTClusters; // use HLT clusters instead of offline clusters
250 TObjArray * fCrossTalkSignalArray; // for 36 sectors
251 TClonesArray* fSeedsPool; //! pool of seeds
252 TArrayI fFreeSeedsID; //! array of ID's of freed seeds
253 Int_t fNFreeSeeds; //! number of seeds freed in the pool
254 Int_t fLastSeedID; //! id of the pool seed on which is returned by the NextFreeSeed method
256 ClassDef(AliTPCtracker,4)
260 AliTPCtrackerRow & AliTPCtracker::GetRow(Int_t sec, Int_t row)
263 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()]:fInnerSec[sec][row];
266 Bool_t AliTPCtracker::IsActive(Int_t sec, Int_t row)
269 // check if the given sector row is active
271 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()].GetN()>0:fInnerSec[sec][row].GetN()>0;
275 Double_t AliTPCtracker::GetXrow(Int_t row) const {
276 // return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetX(row-fInnerSec->GetNRows()):fInnerSec->GetX(row);
280 Double_t AliTPCtracker::GetMaxY(Int_t row) const {
281 //return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetMaxY(row-fInnerSec->GetNRows()):fInnerSec->GetMaxY(row);
285 Int_t AliTPCtracker::GetRowNumber(Double_t x) const
288 return (x>133.) ? fOuterSec->GetRowNumber(x)+fInnerSec->GetNRows():fInnerSec->GetRowNumber(x);
291 Double_t AliTPCtracker::GetPadPitchLength(Double_t x) const
294 return (x>133.) ? fOuterSec->GetPadPitchLength(x):fInnerSec->GetPadPitchLength(x);
295 //return fPadLength[row];
298 Double_t AliTPCtracker::GetPadPitchLength(Int_t row) const
301 return fPadLength[row];
304 void AliTPCtracker::SetTPCtrackerSectors(AliTPCtrackerSector *innerSec, AliTPCtrackerSector *outerSec)
307 fInnerSec = innerSec;
308 fOuterSec = outerSec;