]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDtracker.h
BlockFilter component added; minor corrections
[u/mrichter/AliRoot.git] / TRD / AliTRDtracker.h
1 #ifndef ALITRDTRACKER_H
2 #define ALITRDTRACKER_H   
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */ 
6
7 /* $Id$ */
8
9 ////////////////////////////////////////////////////////////////////////////
10 //                                                                        //
11 //  The TRD tracker                                                       //
12 //                                                                        //
13 ////////////////////////////////////////////////////////////////////////////
14
15 #include "TObjArray.h" 
16
17 #include "AliTracker.h" 
18
19 class TFile;
20 class TTree;
21 class TH1D;
22 class TH2D;
23 class TParticle;
24 class TParticlePDG;
25 class TTreeSRedirector;
26
27 class AliTRDgeometry;
28 class AliTRDtrack;
29 class AliTRDtracklet;
30 class AliTRDcluster;
31 class AliTRDseed;
32 class AliESDEvent;
33
34 ///////////////////////////////////////////////////////////////////////////////
35 //                                                                           //
36 //  The standard TRD tracker                                                 //  
37 //  Based on Kalman filltering approach                                      //
38 //                                                                           //
39 ///////////////////////////////////////////////////////////////////////////////
40
41 class AliTRDtracker : public AliTracker { 
42
43  public:
44
45   void InitLogHists();
46   void SaveLogHists();
47
48   enum { kMaxLayersPerSector   = 1000
49        , kMaxTimeBinIndex      = 216
50        , kMaxClusterPerTimeBin = 2300
51        , kZones                = 5
52        , kTrackingSectors      = 18   };
53   
54   AliTRDtracker();
55   AliTRDtracker(const AliTRDtracker &t);
56   AliTRDtracker(const TFile *in);
57   virtual         ~AliTRDtracker(); 
58   AliTRDtracker   &operator=(const AliTRDtracker &/*t*/) { return *this;          } 
59   
60   void             SetAddTRDseeds()               { fAddTRDseeds = kTRUE;         }
61   void             SetNoTilt()                    { fNoTilt      = kTRUE;         }
62   
63   Int_t            GetTimeBinsPerPlane() const    { return fTimeBinsPerPlane;     }   
64   Double_t         GetMaxChi2() const             { return fgkMaxChi2;            }
65   Float_t          GetLabelFraction() const       { return fgkLabelFraction;      }
66   Float_t          GetMinClustersInTrack() const  { return fgkMinClustersInTrack; }
67   Int_t            GetLastPlane(AliTRDtrack *track);
68   Double_t         GetTiltFactor(const AliTRDcluster *c);
69   Bool_t           GetTrackPoint(Int_t index, AliTrackPoint& p) const;
70   Double_t         GetX(Int_t sec, Int_t plane, Int_t localTB) const;
71   Double_t         GetX(Int_t sec, Int_t pl) const
72                                                   { return fTrSec[sec]->GetLayer(pl)->GetX();            }
73   Int_t            GetGlobalTimeBin(Int_t sec, Int_t plane, Int_t localTB) const 
74                                                   { return fTrSec[sec]->CookTimeBinIndex(plane,localTB); }
75   Double_t         GetLayerNumber(Int_t sec, Double_t x) const 
76                                                   { return fTrSec[sec]->GetLayerNumber(x);               }
77   AliCluster      *GetCluster(Int_t index) const  { if (index >= fNclusters) return NULL; 
78                                                     return (AliCluster *) fClusters->UncheckedAt(index); }
79   
80   static Int_t     Freq(Int_t n, const Int_t *inlist, Int_t *outlist, Bool_t down);    
81   Int_t            Clusters2Tracks(AliESDEvent *event);
82   Int_t            PropagateBack(AliESDEvent *event);
83   Int_t            RefitInward(AliESDEvent *event);
84   
85   virtual void     CookLabel(AliKalmanTrack *t, Float_t wrong) const;
86   
87   Int_t            LocalToGlobalID(Int_t lid);
88   Int_t            GlobalToLocalID(Int_t gid);
89   
90   Int_t            LoadClusters(TTree *cTree);
91   void             UnloadClusters();
92   virtual void     UseClusters(const AliKalmanTrack *t, Int_t from = 0) const;  
93   Int_t            ReadClusters(TObjArray *array, TTree *in) const;
94   AliTRDcluster   *GetCluster(AliTRDtrack *track, Int_t plane, Int_t timebin, UInt_t &index);
95   Int_t            FindClusters(Int_t sector, Int_t t0, Int_t t1, AliTRDtrack *track
96                               , Int_t *clusters, AliTRDtracklet &tracklet);
97   
98  protected:
99   
100   //__________________________________________________________________________________________________
101   class AliTRDpropagationLayer {
102     
103    public: 
104     
105     AliTRDpropagationLayer(Double_t x, Double_t dx, Double_t rho
106                          , Double_t x0, Int_t tbIndex, Int_t plane); 
107     AliTRDpropagationLayer(const AliTRDpropagationLayer &/*p*/);
108     ~AliTRDpropagationLayer() {  if (fTimeBinIndex >= 0) { 
109                                    delete[] fClusters; 
110                                    delete[] fIndex; 
111                                  } 
112                                }
113     AliTRDpropagationLayer &operator=(const AliTRDpropagationLayer &/*p*/) 
114                                                               { return *this; }
115
116     operator Int_t() const                                    { return fN;    }
117     AliTRDcluster  *operator[](Int_t i)                       { return fClusters[i];          }
118     
119     void     SetZmax(Int_t cham, Double_t center, Double_t w) { fZc[cham]      = center; 
120                                                                 fZmax[cham]    = w;           }
121     void     SetYmax(Double_t w, Double_t wsensitive)         { fYmax          = w;
122                                                                 fYmaxSensitive = wsensitive;  }
123
124     void     SetZ(Double_t* center, Double_t *w, Double_t *wsensitive);
125     void     SetHoles(Bool_t* holes);
126     void     SetHole(Double_t Zmax, Double_t Ymax
127                    , Double_t rho = 1.29e-3, Double_t x0 = 36.66
128                    , Double_t Yc = 0.0, Double_t Zc = 0.0);
129     
130     Double_t GetYmax() const           { return fYmax;                  }
131     Double_t GetZmax(Int_t c) const    { return fZmax[c];               }
132     Double_t GetZc(Int_t c) const      { return fZc[c];                 }
133     UInt_t   GetIndex(Int_t i) const   { return fIndex[i];              }  
134     Double_t GetX() const              { return fX;                     }
135     Double_t GetdX() const             { return fdX;                    }
136     Int_t    GetTimeBinIndex() const   { return fTimeBinIndex;          }     
137     Int_t    GetPlane() const          { return fPlane;                 }
138     Bool_t   IsHole(Int_t zone) const  { return fIsHole[zone];          }              
139     Bool_t   IsSensitive() const       { return (fTimeBinIndex >= 0) ? kTRUE : kFALSE;} 
140     
141     void     Clear() { 
142       for (Int_t i = 0; i < fN; i++) fClusters[i] = NULL;
143       fN = 0;
144     }
145     
146     void     InsertCluster(AliTRDcluster *c, UInt_t index);
147     Int_t    Find(Float_t y) const; 
148     Int_t    FindNearestCluster(Float_t y, Float_t z, Float_t maxroad, Float_t maxroadz) const;
149  
150     void     SetX(Double_t x)          { fX = x; }
151     
152    private:     
153     
154     Int_t                     fN;                            // This is fN
155     Int_t                     fSec;                          // Sector mumber
156     AliTRDcluster           **fClusters;                     // Array of pointers to clusters
157     UInt_t                   *fIndex;                        // Array of cluster indexes
158     Double_t                  fX;                            // X coordinate of the middle plane
159     Double_t                  fdX;                           // Radial thickness of the time bin
160     Double_t                  fRho;                          // Default density of the material 
161     Double_t                  fX0;                           // Default radiation length 
162     Int_t                     fTimeBinIndex;                 // Plane * F(local_tb)  
163     Int_t                     fPlane;                        // Plane number
164     Double_t                  fZc[kZones];                   // Z position of the center for 5 active areas
165     Double_t                  fZmax[kZones];                 // Half of active area length in Z
166     Double_t                  fZmaxSensitive[kZones];        // Sensitive area for detection Z     
167     Bool_t                    fIsHole[kZones];               // Is hole in given sector       
168     Double_t                  fYmax;                         // Half of active area length in Y
169     Double_t                  fYmaxSensitive;                // Half of active area length in Y
170     
171     Bool_t                    fHole;                         // kTRUE if there is a hole in the layer
172     Double_t                  fHoleZc;                       // Z of the center of the hole 
173     Double_t                  fHoleZmax;                     // Half of the hole length in Z
174     Double_t                  fHoleYc;                       // Y of the center of the hole 
175     Double_t                  fHoleYmax;                     // Half of the hole length in Y 
176     Double_t                  fHoleRho;                      // Density of the gas in the hole 
177     Double_t                  fHoleX0;                       // Radiation length of the gas in the hole 
178     
179   };
180   
181   //__________________________________________________________________________________________________
182   class AliTRDtrackingSector {
183     
184    public:
185     
186     AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs);
187     AliTRDtrackingSector(const AliTRDtrackingSector &/*t*/);
188     ~AliTRDtrackingSector();
189     
190     AliTRDtrackingSector &operator=(const AliTRDtrackingSector &/*t*/) { return *this; }
191     
192     Int_t    GetNumberOfLayers() const             { return fN; }
193     Int_t    GetNumberOfTimeBins() const;
194     Int_t    GetLayerNumber(Double_t x) const;
195     Int_t    GetInnerTimeBin() const;
196     Int_t    GetOuterTimeBin() const;
197     Int_t    GetLayerNumber(Int_t tb) const        { return fTimeBinIndex[tb];   }
198     Double_t GetX(Int_t pl) const                  { return fLayers[pl]->GetX(); }
199     AliTRDpropagationLayer* GetLayer(Int_t i)      { return fLayers[i];          }
200     
201     void     MapTimeBinLayers();
202     Int_t    Find(Double_t x) const; 
203     void     InsertLayer(AliTRDpropagationLayer *pl);
204     Int_t    CookTimeBinIndex(Int_t plane, Int_t localTB) const;     
205     
206    private:
207     
208     Int_t                   fN;                              // Total number of layers
209     AliTRDgeometry         *fGeom;                           // Geometry
210     AliTRDpropagationLayer *fLayers[kMaxLayersPerSector];    // Layers   
211     Int_t                   fTimeBinIndex[kMaxTimeBinIndex]; // Time bin index
212     Int_t                   fGeomSector;                     // Sector# in AliTRDgeometry
213     
214   };
215   
216   AliTRDgeometry          *fGeom;                          // Pointer to TRD geometry
217   AliTRDtrackingSector    *fTrSec[kTrackingSectors];       // Array of tracking sectors;    
218   Int_t                    fNclusters;                     // Number of clusters in TRD 
219   TObjArray               *fClusters;                      // List of clusters for all sectors
220   Int_t                    fNseeds;                        // Number of track seeds  
221   TObjArray               *fSeeds;                         // List of track seeds
222   Int_t                    fNtracks;                       // Number of reconstructed tracks 
223   TObjArray               *fTracks;                        // List of reconstructed tracks   
224   Int_t                    fTimeBinsPerPlane;              // Timebins per plane in track prolongation 
225   
226   static const Double_t    fgkMaxChi2;                     // Max increment in track chi2 
227   static const Float_t     fgkMinClustersInTrack;          // Min number of clusters in track
228   static const Float_t     fgkLabelFraction;               // Min fraction of same label
229   static const Double_t    fgkMaxSnp;                      // Maximal snp for tracking
230   static const Double_t    fgkMaxStep;                     // Maximal step for tracking  
231   
232   Bool_t                   fAddTRDseeds;                   // Something else
233   Bool_t                   fNoTilt;                        // No tilt, or what?
234
235   // Histograms
236   TH1D                    *fHBackfit;                      // QA histogram
237   TH1D                    *fHClSearch;                     // QA histogram
238   TH1D                    *fHRefit;                        // QA histogram
239   TH1D                    *fHX;                            // QA histogram
240   TH1D                    *fHNCl;                          // QA histogram
241   TH1D                    *fHNClTrack;                     // QA histogram
242   TH1D                    *fHFindCl[4];                    // QA histogram
243   TH1D                    *fHMinYPos;                      // QA histogram
244   TH1D                    *fHMinYNeg;                      // QA histogram
245   TH1D                    *fHMinZ;                         // QA histogram
246   TH2D                    *fHMinD;                         // QA histogram
247   TH1D                    *fHDeltaX;                       // QA histogram
248   TH1D                    *fHXCl;                          // QA histogram
249   
250   Bool_t   AdjustSector(AliTRDtrack *track); 
251   
252  private:
253   
254   AliTRDtrack *RegisterSeed(AliTRDseed *seeds, Double_t *params);
255   void     MakeSeedsMI(Int_t inner, Int_t outer, AliESDEvent *esd = 0);
256   
257   Int_t    FollowBackProlongation(AliTRDtrack &t);
258   Int_t    FollowProlongation(AliTRDtrack &t);
259   Int_t    PropagateToX(AliTRDtrack &t, Double_t xToGo, Double_t maxStep);
260   Double_t ExpectedSigmaY2(Double_t r, Double_t tgl, Double_t pt) const;
261   Double_t ExpectedSigmaZ2(Double_t r, Double_t tgl) const;
262   
263   TTreeSRedirector        *fDebugStreamer;                 //!Debug streamer
264   
265   ClassDef(AliTRDtracker,3)                                // TRD tracker
266     
267 };
268 #endif