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