]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TPC/AliTPCtrackerMI.h
Added switch in recoparam for HLTPreSeeding
[u/mrichter/AliRoot.git] / TPC / AliTPCtrackerMI.h
... / ...
CommitLineData
1#ifndef ALITPCTRACKERMI_H
2#define ALITPCTRACKERMI_H
3/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
5
6
7/* $Id$ */
8
9//-------------------------------------------------------
10// TPC tracker
11// Parallel tracker
12//
13// Origin:
14//-------------------------------------------------------
15
16#include <TArrayI.h>
17#include "AliTracker.h"
18#include "AliTPCreco.h"
19#include "AliTPCclusterMI.h"
20#include "AliTPCtrackerSector.h"
21
22
23class TFile;
24class AliTPCParam;
25class AliTPCseed;
26class AliTPCTrackerPoint;
27class AliESDEvent;
28class AliESDtrack;
29class TTree;
30class AliESDkink;
31class TTreeSRedirector;
32class AliTrackPoint;
33class AliDCSSensorArray;
34class AliDCSSensor;
35
36
37
38class AliTPCtrackerMI : public AliTracker {
39public:
40 AliTPCtrackerMI();
41 AliTPCtrackerMI(const AliTPCParam *par);
42 virtual ~AliTPCtrackerMI();
43 //
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);
52 Int_t LoadClusters();
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();
59 //
60 void FillESD(const TObjArray* arr);
61 void DeleteSeeds();
62 void SetDebug(Int_t debug){ fDebug = debug;}
63 void FindKinks(TObjArray * array, AliESDEvent * esd);
64 //
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);
68 //
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;}
75 //
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;
81
82 void RotateToLocal(AliTPCseed *seed);
83
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;
86
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);
91
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);
97
98 void SortTracks(TObjArray * arr, Int_t mode) const;
99
100 virtual Double_t ErrY2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
101 virtual Double_t ErrZ2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
102
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;
107
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,
110 Double_t c) const;
111 Bool_t GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z) const;
112 //
113 void ResetSeedsPool();
114 void MarkSeedFree( TObject* seed );
115 TObject *&NextFreeSeed();
116 //
117 public:
118 void SetUseHLTClusters(Int_t useHLTClusters) {fUseHLTClusters = useHLTClusters;} // set usage from HLT clusters from rec.C options
119
120 Float_t OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t &sum2);
121 void SignShared(AliTPCseed * s1, AliTPCseed * s2);
122 void SignShared(TObjArray * arr);
123
124 void RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal);
125
126 Int_t AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster);
127
128 Bool_t IsTPCHVDipEvent(AliESDEvent const *esdEvent);
129
130private:
131 Bool_t IsFindable(AliTPCseed & t);
132 AliTPCtrackerMI(const AliTPCtrackerMI& r); //dummy copy constructor
133 AliTPCtrackerMI &operator=(const AliTPCtrackerMI& r);//dummy assignment operator
134 void AddCovariance(AliTPCseed * seed); // add covariance
135 void AddCovarianceAdd(AliTPCseed * seed); // add covariance
136
137 inline AliTPCtrackerRow &GetRow(Int_t sec, Int_t row);
138 inline Bool_t IsActive(Int_t sec, Int_t row);
139 inline Double_t GetXrow(Int_t row) const;
140 inline Double_t GetMaxY(Int_t row) const;
141 inline Int_t GetRowNumber(Double_t x) const;
142 Int_t GetRowNumber(Double_t x[3]) const;
143 inline Double_t GetPadPitchLength(Double_t x) const;
144 inline Double_t GetPadPitchLength(Int_t row) const;
145
146 void GetShape(AliTPCseed * seed, Int_t row);
147
148 void ReadSeeds(const AliESDEvent *const event, Int_t direction); //read seeds from the event
149
150 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);
151 void MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1);
152
153 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);
154
155
156 AliTPCseed *MakeSeed(AliTPCseed *const track, Float_t r0, Float_t r1, Float_t r2); //reseed
157 AliTPCseed *ReSeed(const AliTPCseed *t, Float_t r0, Float_t r1, Float_t r2); //reseed
158 AliTPCseed *ReSeed(AliTPCseed *t, Int_t r0, Bool_t forward); //reseed
159
160
161
162 AliTPCseed * ReSeed(AliTPCseed *t);
163 //Int_t LoadInnerSectors();
164 //Int_t LoadOuterSectors();
165 void DumpClusters(Int_t iter, TObjArray *trackArray);
166 void UnsignClusters();
167 void SignClusters(const TObjArray * arr, Float_t fnumber=3., Float_t fdensity=2.);
168
169 void ParallelTracking(TObjArray *const arr, Int_t rfirst, Int_t rlast);
170 void Tracking(TObjArray * arr);
171 TObjArray * Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy=-1, Int_t dsec=0);
172 TObjArray * Tracking();
173 TObjArray * TrackingSpecial();
174 void SumTracks(TObjArray *arr1,TObjArray *&arr2);
175 void PrepareForBackProlongation(const TObjArray *const arr, Float_t fac) const;
176 void PrepareForProlongation(TObjArray *const arr, Float_t fac) const;
177
178 Int_t UpdateTrack(AliTPCseed *t, Int_t accept); //update trackinfo
179
180 void MakeESDBitmaps(AliTPCseed *t, AliESDtrack *esd);
181
182 Int_t PropagateToRowHLT(AliTPCseed *pt, int nrow);
183 void TrackFollowingHLT(TObjArray *const arr);
184 TObjArray * MakeSeedsHLT(const AliESDEvent *hltEvent);
185
186 const Int_t fkNIS; //number of inner sectors
187 AliTPCtrackerSector *fInnerSec; //array of inner sectors;
188 const Int_t fkNOS; //number of outer sectors
189 AliTPCtrackerSector *fOuterSec; //array of outer sectors;
190
191 Int_t fN; //number of loaded sectors
192 AliTPCtrackerSector *fSectors; //pointer to loaded sectors;
193 //
194 TTree * fInput; // input tree with clusters
195 TTree * fOutput; // output tree with tracks
196 TTree * fSeedTree; // output tree with seeds - filled in debug mode 1
197 TTree * fTreeDebug; // output with a debug information about track
198 AliESDEvent * fEvent; // output with esd tracks
199 const AliESDEvent * fEventHLT; // input with HLT tracks
200 Int_t fDebug; // debug option
201 Bool_t fNewIO; // indicated if we have data using New IO
202 Int_t fNtracks; //current number of tracks
203 TObjArray *fSeeds; //array of track seeds
204 Int_t fIteration; // indicate iteration - 0 - froward -1 back - 2forward - back->forward
205 // TObjArray * fTrackPointPool; // ! pool with track points
206 Double_t fXRow[200]; // radius of the pad row
207 Double_t fYMax[200]; // max y for given pad row
208 Double_t fPadLength[200]; // max y for given pad row
209 const AliTPCParam *fkParam; //pointer to the parameters
210 TTreeSRedirector *fDebugStreamer; //!debug streamer
211 Int_t fUseHLTClusters; // use HLT clusters instead of offline clusters
212 //
213 TClonesArray* fSeedsPool; //! pool of seeds
214 TArrayI fFreeSeedsID; //! array of ID's of freed seeds
215 Int_t fNFreeSeeds; //! number of seeds freed in the pool
216 Int_t fLastSeedID; //! id of the pool seed on which is returned by the NextFreeSeed method
217 //
218 ClassDef(AliTPCtrackerMI,3)
219};
220
221
222AliTPCtrackerRow & AliTPCtrackerMI::GetRow(Int_t sec, Int_t row)
223{
224 //
225 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()]:fInnerSec[sec][row];
226}
227
228Bool_t AliTPCtrackerMI::IsActive(Int_t sec, Int_t row)
229{
230 //
231 // check if the given sector row is active
232 //
233 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()].GetN()>0:fInnerSec[sec][row].GetN()>0;
234}
235
236
237Double_t AliTPCtrackerMI::GetXrow(Int_t row) const {
238 // return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetX(row-fInnerSec->GetNRows()):fInnerSec->GetX(row);
239 return fXRow[row];
240}
241
242Double_t AliTPCtrackerMI::GetMaxY(Int_t row) const {
243 //return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetMaxY(row-fInnerSec->GetNRows()):fInnerSec->GetMaxY(row);
244 return fYMax[row];
245}
246
247Int_t AliTPCtrackerMI::GetRowNumber(Double_t x) const
248{
249 //
250 return (x>133.) ? fOuterSec->GetRowNumber(x)+fInnerSec->GetNRows():fInnerSec->GetRowNumber(x);
251}
252
253Double_t AliTPCtrackerMI::GetPadPitchLength(Double_t x) const
254{
255 //
256 return (x>133.) ? fOuterSec->GetPadPitchLength(x):fInnerSec->GetPadPitchLength(x);
257 //return fPadLength[row];
258}
259
260Double_t AliTPCtrackerMI::GetPadPitchLength(Int_t row) const
261{
262 //
263 return fPadLength[row];
264}
265
266
267
268#endif
269
270