]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TPC/AliTPCtrackerMI.h
Add cell rejection based on energy in sorrounding cells, 4 cells in corners with...
[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 "AliTracker.h"
17#include "AliTPCreco.h"
18#include "AliTPCclusterMI.h"
19#include "AliTPCtrackerSector.h"
20
21
22class TFile;
23class AliTPCParam;
24class AliTPCseed;
25class AliTPCTrackerPoint;
26class AliESDEvent;
27class TTree;
28class AliESDkink;
29class TTreeSRedirector;
30class AliTrackPoint;
31
32
33
34class AliTPCtrackerMI : public AliTracker {
35public:
36 AliTPCtrackerMI();
37 AliTPCtrackerMI(const AliTPCParam *par);
38 virtual ~AliTPCtrackerMI();
39 //
40 void SetIteration(Int_t iteration){fIteration = iteration;}
41 virtual Int_t Clusters2Tracks (AliESDEvent *const esd);
42 virtual Int_t RefitInward (AliESDEvent *esd);
43 virtual Int_t LoadClusters (TTree * const tree);
44 virtual Int_t LoadClusters (const TObjArray * arr); // another input
45 virtual Int_t LoadClusters (const TClonesArray * arr); // another input
46 Int_t LoadClusters();
47 void UnloadClusters();
48 Int_t LoadInnerSectors();
49 Int_t LoadOuterSectors();
50 virtual void FillClusterArray(TObjArray* array) const;
51 void Transform(AliTPCclusterMI * cluster);
52 //
53 void FillESD(const TObjArray* arr);
54 void DeleteSeeds();
55 void SetDebug(Int_t debug){ fDebug = debug;}
56 void FindKinks(TObjArray * array, AliESDEvent * esd);
57 //
58 void FindCurling(const TObjArray * array, AliESDEvent * esd, Int_t iter);
59 void FindSplitted(TObjArray * array, AliESDEvent * esd, Int_t iter);
60 void FindMultiMC(const TObjArray * array, AliESDEvent * esd, Int_t iter);
61 //
62 void UpdateKinkQualityM(AliTPCseed * seed);
63 void UpdateKinkQualityD(AliTPCseed * seed);
64 Int_t CheckKinkPoint(AliTPCseed*seed, AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
65 Int_t RefitKink(AliTPCseed &mother, AliTPCseed &daughter, const AliESDkink &kink);
66 Int_t ReadSeeds(const TFile *in);
67 TObjArray * GetSeeds() const {return fSeeds;}
68 //
69 AliCluster * GetCluster(Int_t index) const {return (AliCluster*)GetClusterMI(index);}
70 AliTPCclusterMI *GetClusterMI(Int_t index) const;
71 Int_t Clusters2Tracks();
72 virtual void CookLabel(AliKalmanTrack *tk,Float_t wrong) const;
73 virtual Int_t CookLabel(AliTPCseed *const t,Float_t wrong, Int_t first,Int_t last ) const;
74
75 void RotateToLocal(AliTPCseed *seed);
76
77 Int_t FollowProlongation(AliTPCseed& t, Int_t rf=0, Int_t step=1);
78 Bool_t GetTrackPoint(Int_t index, AliTrackPoint &p ) const;
79
80 Int_t FollowBackProlongation(AliTPCseed& t, Int_t rf);
81 Int_t FollowToNext(AliTPCseed& t, Int_t nr);
82 Int_t UpdateClusters(AliTPCseed& t, Int_t nr);
83 Int_t FollowToNextCluster( AliTPCseed& t, Int_t nr);
84
85 Int_t PropagateBack(const TObjArray *const arr);
86 Int_t PropagateBack(AliESDEvent * event);
87 Int_t PropagateBack(AliTPCseed *const pt, Int_t row0, Int_t row1);
88 Int_t PropagateForward();
89 Int_t PropagateForward2(const TObjArray *const arr);
90
91 void SortTracks(TObjArray * arr, Int_t mode) const;
92
93 virtual Double_t ErrY2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
94 virtual Double_t ErrZ2(AliTPCseed* seed, const AliTPCclusterMI * cl = 0);
95
96 Double_t F1(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
97 Double_t F1old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
98 Double_t F2(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
99 Double_t F2old(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t x3,Double_t y3) const;
100
101 Double_t F3(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2) const;
102 Double_t F3n(Double_t x1,Double_t y1, Double_t x2,Double_t y2, Double_t z1,Double_t z2,
103 Double_t c) const;
104 Bool_t GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z) const;
105
106 public:
107
108
109 Float_t OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t &sum2);
110 void SignShared(AliTPCseed * s1, AliTPCseed * s2);
111 void SignShared(TObjArray * arr);
112
113 void RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal);
114
115 void StopNotActive(const TObjArray * arr, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const;
116 void StopNotActive(AliTPCseed * seed, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const;
117 Int_t AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster);
118
119private:
120 Bool_t IsFindable(AliTPCseed & t);
121 AliTPCtrackerMI(const AliTPCtrackerMI& r); //dummy copy constructor
122 AliTPCtrackerMI &operator=(const AliTPCtrackerMI& r);//dummy assignment operator
123 void AddCovariance(AliTPCseed * seed); // add covariance
124 void AddCovarianceAdd(AliTPCseed * seed); // add covariance
125
126 inline AliTPCtrackerRow &GetRow(Int_t sec, Int_t row);
127 inline Bool_t IsActive(Int_t sec, Int_t row);
128 inline Double_t GetXrow(Int_t row) const;
129 inline Double_t GetMaxY(Int_t row) const;
130 inline Int_t GetRowNumber(Double_t x) const;
131 Int_t GetRowNumber(Double_t x[3]) const;
132 inline Double_t GetPadPitchLength(Double_t x) const;
133 inline Double_t GetPadPitchLength(Int_t row) const;
134
135 void GetShape(AliTPCseed * seed, Int_t row);
136
137 void ReadSeeds(const AliESDEvent *const event, Int_t direction); //read seeds from the event
138
139 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);
140 void MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4], Float_t deltay = -1);
141
142 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);
143
144
145 AliTPCseed *MakeSeed(AliTPCseed *const track, Float_t r0, Float_t r1, Float_t r2); //reseed
146 AliTPCseed *ReSeed(const AliTPCseed *t, Float_t r0, Float_t r1, Float_t r2); //reseed
147 AliTPCseed *ReSeed(AliTPCseed *t, Int_t r0, Bool_t forward); //reseed
148
149
150
151 AliTPCseed * ReSeed(AliTPCseed *t);
152 //Int_t LoadInnerSectors();
153 //Int_t LoadOuterSectors();
154 void DumpClusters(Int_t iter, TObjArray *trackArray);
155 void UnsignClusters();
156 void SignClusters(const TObjArray * arr, Float_t fnumber=3., Float_t fdensity=2.);
157
158 void ParallelTracking(TObjArray *const arr, Int_t rfirst, Int_t rlast);
159 void Tracking(TObjArray * arr);
160 TObjArray * Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy=-1, Int_t dsec=0);
161 TObjArray * Tracking();
162 TObjArray * TrackingSpecial();
163 void SumTracks(TObjArray *arr1,TObjArray *&arr2) const;
164 void PrepareForBackProlongation(const TObjArray *const arr, Float_t fac) const;
165 void PrepareForProlongation(TObjArray *const arr, Float_t fac) const;
166
167 Int_t UpdateTrack(AliTPCseed *t, Int_t accept); //update trackinfo
168
169 void MakeBitmaps(AliTPCseed *t);
170
171 const Int_t fkNIS; //number of inner sectors
172 AliTPCtrackerSector *fInnerSec; //array of inner sectors;
173 const Int_t fkNOS; //number of outer sectors
174 AliTPCtrackerSector *fOuterSec; //array of outer sectors;
175
176 Int_t fN; //number of loaded sectors
177 AliTPCtrackerSector *fSectors; //pointer to loaded sectors;
178 //
179 TTree * fInput; // input tree with clusters
180 TTree * fOutput; // output tree with tracks
181 TTree * fSeedTree; // output tree with seeds - filled in debug mode 1
182 TTree * fTreeDebug; // output with a debug information about track
183 AliESDEvent * fEvent; // output with esd tracks
184 Int_t fDebug; // debug option
185 Bool_t fNewIO; // indicated if we have data using New IO
186 Int_t fNtracks; //current number of tracks
187 TObjArray *fSeeds; //array of track seeds
188 Int_t fIteration; // indicate iteration - 0 - froward -1 back - 2forward - back->forward
189 // TObjArray * fTrackPointPool; // ! pool with track points
190 // TObjArray * fSeedPool; //! pool with seeds
191 Double_t fXRow[200]; // radius of the pad row
192 Double_t fYMax[200]; // max y for given pad row
193 Double_t fPadLength[200]; // max y for given pad row
194 const AliTPCParam *fkParam; //pointer to the parameters
195 TTreeSRedirector *fDebugStreamer; //!debug streamer
196 ClassDef(AliTPCtrackerMI,2)
197};
198
199
200AliTPCtrackerRow & AliTPCtrackerMI::GetRow(Int_t sec, Int_t row)
201{
202 //
203 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()]:fInnerSec[sec][row];
204}
205
206Bool_t AliTPCtrackerMI::IsActive(Int_t sec, Int_t row)
207{
208 //
209 // check if the given sector row is active
210 //
211 return (row>=fInnerSec->GetNRows()) ? fOuterSec[sec][row-fInnerSec->GetNRows()].GetN()>0:fInnerSec[sec][row].GetN()>0;
212}
213
214
215Double_t AliTPCtrackerMI::GetXrow(Int_t row) const {
216 // return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetX(row-fInnerSec->GetNRows()):fInnerSec->GetX(row);
217 return fXRow[row];
218}
219
220Double_t AliTPCtrackerMI::GetMaxY(Int_t row) const {
221 //return (row>=fInnerSec->GetNRows()) ? fOuterSec->GetMaxY(row-fInnerSec->GetNRows()):fInnerSec->GetMaxY(row);
222 return fYMax[row];
223}
224
225Int_t AliTPCtrackerMI::GetRowNumber(Double_t x) const
226{
227 //
228 return (x>133.) ? fOuterSec->GetRowNumber(x)+fInnerSec->GetNRows():fInnerSec->GetRowNumber(x);
229}
230
231Double_t AliTPCtrackerMI::GetPadPitchLength(Double_t x) const
232{
233 //
234 return (x>133.) ? fOuterSec->GetPadPitchLength(x):fInnerSec->GetPadPitchLength(x);
235 //return fPadLength[row];
236}
237
238Double_t AliTPCtrackerMI::GetPadPitchLength(Int_t row) const
239{
240 //
241 return fPadLength[row];
242}
243
244
245
246#endif
247
248