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