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