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