79369a33faff27d4b0897082f94d4a635907cdd1
[u/mrichter/AliRoot.git] / ITS / AliITStrackerMI.h
1 #ifndef ALIITSTRACKERMI_H
2 #define ALIITSTRACKERMI_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 //-------------------------------------------------------------------------
7 //                          ITS tracker
8 //     reads AliITSclusterMI clusters and creates AliITStrackMI tracks
9 //           Origin: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
10 //-------------------------------------------------------------------------
11
12 #include <TObjArray.h>
13
14 #include "AliTracker.h"
15 #include "AliITStrackMI.h"
16 #include "AliITSclusterV2.h"
17 #include "AliV0vertex.h"
18
19 class AliESD;
20 class AliITSgeom;
21 class TTree;
22 class AliHelix;
23 class AliV0vertex;
24 class AliESDV0MI;
25
26
27
28
29 //-------------------------------------------------------------------------
30 class AliITStrackerMI : public AliTracker {
31 public:
32   AliITStrackerMI():AliTracker(){}
33   AliITStrackerMI(const AliITSgeom *geom);
34   ~AliITStrackerMI();
35   AliCluster *GetCluster(Int_t index) const;
36   AliITSclusterV2 *GetClusterLayer(Int_t layn, Int_t ncl) const
37                         {return fgLayers[layn].GetCluster(ncl);}
38   Int_t GetNumberOfClustersLayer(Int_t layn) const 
39                         {return fgLayers[layn].GetNumberOfClusters();}
40   Int_t LoadClusters(TTree *cf);
41   void UnloadClusters();
42   Int_t Clusters2Tracks(AliESD *event);
43   Int_t PropagateBack(AliESD *event);
44   Int_t RefitInward(AliESD *event);
45   Bool_t RefitAt(Double_t x, AliITStrackMI *seed, const AliITStrackMI *t);
46   void SetupFirstPass(Int_t *flags, Double_t *cuts=0);
47   void SetupSecondPass(Int_t *flags, Double_t *cuts=0);
48
49   void SetLastLayerToTrackTo(Int_t l=0) {fLastLayerToTrackTo=l;} 
50   void SetLayersNotToSkip(Int_t *l);
51   void UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
52   void GetNTeor(Int_t layer, const AliITSclusterV2* cl, Float_t theta, Float_t phi, Float_t &ny, Float_t &nz);
53   Int_t  GetError(Int_t layer, const AliITSclusterV2*cl, Float_t theta, Float_t phi, Float_t expQ, Float_t &erry, Float_t &errz);
54
55   void  GetDCASigma(AliITStrackMI* track, Float_t & sigmarfi, Float_t &sigmaz);
56   Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSclusterV2 *cluster,Int_t layer);
57   Int_t UpdateMI(AliITStrackMI* track, const AliITSclusterV2* cl,Double_t chi2,Int_t layer) const;
58   class AliITSdetector { 
59   public:
60     AliITSdetector(){}
61     AliITSdetector(Double_t r,Double_t phi) {fR=r; fPhi=phi; fSinPhi = TMath::Sin(phi); fCosPhi = TMath::Cos(phi);
62     fYmin=10000;fYmax=-1000; fZmin=10000;fZmax=-1000;}
63     inline void GetGlobalXYZ( const AliITSclusterV2 *cl, Double_t xyz[3]) const;
64     Double_t GetR()   const {return fR;}
65     Double_t GetPhi() const {return fPhi;}
66     Double_t GetYmin() const {return fYmin;}
67     Double_t GetYmax() const {return fYmax;}
68     Double_t GetZmin() const {return fZmin;}
69     Double_t GetZmax() const {return fZmax;}
70     void SetYmin(Double_t min) {fYmin = min;}
71     void SetYmax(Double_t max) {fYmax = max;}
72     void SetZmin(Double_t min) {fZmin = min;}
73     void SetZmax(Double_t max) {fZmax = max;}
74   private:
75     Double_t fR;    // polar coordinates 
76     Double_t fPhi;  // of this detector
77     Double_t fSinPhi; // sin of phi;
78     Double_t fCosPhi; // cos of phi 
79     Double_t fYmin;   //  local y minimal
80     Double_t fYmax;   //  local max y
81     Double_t fZmin;   //  local z min
82     Double_t fZmax;   //  local z max
83   };
84
85   class AliITSlayer {
86     friend class AliITStrackerMI;
87   public:
88     AliITSlayer();
89     AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
90    ~AliITSlayer();
91     Int_t InsertCluster(AliITSclusterV2 *c);
92     void  SortClusters();
93     void ResetClusters();
94     void ResetWeights();
95     void SelectClusters(Double_t zmi,Double_t zma,Double_t ymi,Double_t yma);
96     const AliITSclusterV2 *GetNextCluster(Int_t &ci);
97     void ResetRoad();
98     Double_t GetRoad() const {return fRoad;}
99     Double_t GetR() const {return fR;}
100     Int_t FindClusterIndex(Float_t z) const;
101     AliITSclusterV2 *GetCluster(Int_t i) const {return i<fN? fClusters[i]:0;} 
102     Float_t         *GetWeight(Int_t i)   {return i<fN ?&fClusterWeight[i]:0;}
103     AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
104     Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
105     Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
106     Int_t InRoad() const ;
107     Int_t GetNumberOfClusters() const {return fN;}
108     Int_t GetNladders() const {return fNladders;}
109     Int_t GetNdetectors() const {return fNdetectors;}
110     Int_t GetSkip() const {return fSkip;}
111     void  SetSkip(Int_t skip){fSkip=skip;}
112     void IncAccepted(){fAccepted++;}
113     Int_t GetAccepted() const {return fAccepted;}
114   protected:
115     AliITSlayer(const AliITSlayer& /*layer*/){;}
116     Double_t fR;                // mean radius of this layer
117     Double_t fPhiOffset;        // offset of the first detector in Phi
118     Int_t fNladders;            // number of ladders
119     Double_t fZOffset;          // offset of the first detector in Z
120     Int_t fNdetectors;          // detectors/ladder
121     AliITSdetector *fDetectors; // array of detectors
122     Int_t fN;                   // number of clusters
123     AliITSclusterV2 *fClusters[kMaxClusterPerLayer]; // pointers to clusters
124     Int_t        fClusterIndex[kMaxClusterPerLayer]; // pointers to clusters
125     Float_t fY[kMaxClusterPerLayer];                // y position of the clusters      
126     Float_t fZ[kMaxClusterPerLayer];                // z position of the clusters      
127     Float_t fYB[2];                                       // ymin and ymax
128     //
129     AliITSclusterV2 *fClusters5[6][kMaxClusterPerLayer5]; // pointers to clusters -     slice in y
130     Int_t        fClusterIndex5[6][kMaxClusterPerLayer5]; // pointers to clusters -     slice in y    
131     Float_t fY5[6][kMaxClusterPerLayer5];                // y position of the clusters  slice in y    
132     Float_t fZ5[6][kMaxClusterPerLayer5];                // z position of the clusters  slice in y 
133     Int_t fN5[6];                                       // number of cluster in slice
134     Float_t fDy5;                                       //delta y
135     Float_t fBy5[6][2];                                    //slice borders
136     //
137     AliITSclusterV2 *fClusters10[11][kMaxClusterPerLayer10]; // pointers to clusters -     slice in y
138     Int_t        fClusterIndex10[11][kMaxClusterPerLayer10]; // pointers to clusters -     slice in y    
139     Float_t fY10[11][kMaxClusterPerLayer10];                // y position of the clusters  slice in y    
140     Float_t fZ10[11][kMaxClusterPerLayer10];                // z position of the clusters  slice in y 
141     Int_t fN10[11];                                       // number of cluster in slice
142     Float_t fDy10;                                        // delta y
143     Float_t fBy10[11][2];                                 // slice borders
144     //
145     AliITSclusterV2 *fClusters20[21][kMaxClusterPerLayer20]; // pointers to clusters -     slice in y
146     Int_t        fClusterIndex20[21][kMaxClusterPerLayer20]; // pointers to clusters -     slice in y    
147     Float_t fY20[21][kMaxClusterPerLayer20];                // y position of the clusters  slice in y    
148     Float_t fZ20[21][kMaxClusterPerLayer20];                // z position of the clusters  slice in y 
149     Int_t fN20[21];                                       // number of cluster in slice
150     Float_t fDy20;                                        //delta y 
151     Float_t fBy20[21][2];                                 //slice borders
152     //
153     AliITSclusterV2** fClustersCs;                         //clusters table in current slice
154     Int_t   *fClusterIndexCs;                             //cluster index in current slice 
155     Float_t *fYcs;                                        //y position in current slice
156     Float_t *fZcs;                                        //z position in current slice
157     Int_t    fNcs;                                        //number of clusters in current slice    
158     Int_t fCurrentSlice;                                  //current slice
159     //
160     Float_t  fClusterWeight[kMaxClusterPerLayer]; // probabilistic weight of the cluster
161     Int_t    fClusterTracks[4][kMaxClusterPerLayer]; //tracks registered to given cluster
162     Float_t fZmax;      //    edges
163     Float_t fYmin;      //   of  the
164     Float_t fYmax;      //   "window"
165     Int_t fI;            // index of the current cluster within the "window"
166     Int_t fImax;            // index of the last cluster within the "window"    
167     Int_t fSkip;     // indicates possibility to skip cluster
168     Int_t fAccepted;     // accept indicator 
169     Double_t fRoad;      // road defined by the cluster density
170   };
171   AliITStrackerMI::AliITSlayer    & GetLayer(Int_t layer) const;
172   AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
173
174 protected:
175   void FindV0(AliESD *event);  //try to find V0
176   Double_t  TestV0(AliHelix *h1, AliHelix *h2, AliESDV0MI *vertex, Double_t &rmin);  //try to find V0 - return DCA
177   Bool_t  FindBestPair(Int_t esdtrack0, Int_t esdtrack1,AliESDV0MI *vertex, Int_t &ibest0, Int_t &ibest1);  // try to find best pair from the tree of track hyp.
178   void CookLabel(AliKalmanTrack *t,Float_t wrong) const;
179   void CookLabel(AliITStrackMI *t,Float_t wrong) const;
180   Double_t GetEffectiveThickness(Double_t y, Double_t z) const;
181   void FollowProlongationTree(AliITStrackMI * otrack, Int_t esdindex);
182   void ResetBestTrack() {
183      fBestTrack.~AliITStrackMI();
184      new(&fBestTrack) AliITStrackMI(fTrackToFollow);
185   }
186   void ResetTrackToFollow(const AliITStrackMI &t) {
187      fTrackToFollow.~AliITStrackMI();
188      new(&fTrackToFollow) AliITStrackMI(t);
189   }
190   void CookdEdx(AliITStrackMI* track);
191   Double_t GetNormalizedChi2(AliITStrackMI * track, Int_t mode);
192   Double_t GetTruncatedChi2(AliITStrackMI * track, Float_t fac);
193   Double_t NormalizedChi2(AliITStrackMI * track, Int_t layer);
194   Double_t GetInterpolatedChi2(AliITStrackMI * forwardtrack, AliITStrackMI * backtrack);  
195   Double_t GetMatchingChi2(AliITStrackMI * track1, AliITStrackMI * track2);
196   Double_t GetDeadZoneProbability(Double_t zpos, Double_t zerr);
197
198   Float_t    *GetWeight(Int_t index);
199   void AddTrackHypothesys(AliITStrackMI * track, Int_t esdindex);
200   void SortTrackHypothesys(Int_t esdindex, Int_t maxcut, Int_t mode);
201   AliITStrackMI * GetBestHypothesys(Int_t esdindex, AliITStrackMI * original, Int_t checkmax); 
202   void  GetBestHypothesysMIP(TObjArray &itsTracks); 
203   void RegisterClusterTracks(AliITStrackMI* track, Int_t id);
204   void UnRegisterClusterTracks(AliITStrackMI* track, Int_t id);
205   Float_t GetNumberOfSharedClusters(AliITStrackMI* track,Int_t id, Int_t list[6], AliITSclusterV2 *clist[6]);
206   Int_t GetOverlapTrack(AliITStrackMI *track, Int_t trackID, Int_t &shared, Int_t clusterlist[6], Int_t overlist[6]);
207   AliITStrackMI * GetBest2Tracks(Int_t trackID1, Int_t treackID2, Float_t th0, Float_t th1);
208   Float_t  * GetErrY(Int_t trackindex) const {return &fCoeficients[trackindex*48];}
209   Float_t  * GetErrZ(Int_t trackindex) const {return &fCoeficients[trackindex*48+12];}
210   Float_t  * GetNy(Int_t trackindex) const {return &fCoeficients[trackindex*48+24];}
211   Float_t  * GetNz(Int_t trackindex) const {return &fCoeficients[trackindex*48+36];}
212   void       SignDeltas( TObjArray *ClusterArray, Float_t zv);
213   void MakeCoeficients(Int_t ntracks);
214   void UpdateESDtrack(AliITStrackMI* track, ULong_t flags) const;
215   Int_t fI;                              // index of the current layer
216   static AliITSlayer fgLayers[kMaxLayer];// ITS layers
217   AliITStrackMI fTracks[kMaxLayer];      // track estimations at the ITS layers
218   AliITStrackMI fBestTrack;              // "best" track 
219   AliITStrackMI fTrackToFollow;          // followed track
220   TObjArray     fTrackHypothesys;        // ! array with track hypothesys- ARRAY is the owner of tracks- MI
221   TObjArray     fOriginal;               // ! array with seeds from the TPC
222   Int_t         fBestTrackIndex[100000]; // ! index of the best track
223   Int_t         fCurrentEsdTrack;        // ! current esd track           - MI
224   Int_t fPass;                           // current pass through the data 
225   Int_t fConstraint[2];                  // constraint flags
226   Bool_t fAfterV0;                       //indicates V0 founded
227   Int_t fLayersNotToSkip[kMaxLayer];     // layer masks
228   Int_t fLastLayerToTrackTo;             // the innermost layer to track to
229   Float_t * fCoeficients;                //! working array with errors and mean cluser shape
230  private:
231   AliITStrackerMI(const AliITStrackerMI * /*tracker*/){;}
232   ClassDef(AliITStrackerMI,2)   //ITS tracker MI
233 };
234
235
236
237
238 /////////////////////////////////////////////////////////
239 /////////////////////////////////////////////////////////
240 /////////////////////////////////////////////////////////
241
242
243
244
245
246 inline void AliITStrackerMI::SetupFirstPass(Int_t *flags, Double_t *cuts) {
247   // This function sets up flags and cuts for the first tracking pass   
248   //
249   //   flags[0] - vertex constaint flag                                
250   //              negative means "skip the pass"                        
251   //              0        means "no constraint"                        
252   //              positive means "normal constraint"                    
253
254    fConstraint[0]=flags[0];
255    if (cuts==0) return;
256 }
257
258 inline void AliITStrackerMI::SetupSecondPass(Int_t *flags, Double_t *cuts) {
259   // This function sets up flags and cuts for the second tracking pass   
260   //
261   //   flags[0] - vertex constaint flag                                
262   //              negative means "skip the pass"                        
263   //              0        means "no constraint"                        
264   //              positive means "normal constraint"                    
265
266    fConstraint[1]=flags[0];
267    if (cuts==0) return;
268 }
269
270 inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
271   //--------------------------------------------------------------------
272   //This function "cooks" a track label. If label<0, this track is fake.
273   //--------------------------------------------------------------------
274    Int_t tpcLabel=t->GetLabel();
275    if (tpcLabel<0) return;
276    AliTracker::CookLabel(t,wrong);
277    if (tpcLabel!=TMath::Abs(t->GetLabel())){
278      t->SetFakeRatio(1.);
279    }
280    if (tpcLabel !=t->GetLabel()) {
281      t->SetLabel(-tpcLabel);      
282    }
283 }
284
285 inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
286 {
287   //--------------------------------------------------------------------
288   //get normalize chi2
289   //--------------------------------------------------------------------
290   track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->GetChi2()/
291   //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
292     TMath::Max(double(track->GetNumberOfClusters()-track->fNSkipped),
293                1./(1.+track->fNSkipped));
294   return track->fNormChi2[layer];
295 }
296 inline void  AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSclusterV2 *cl, Double_t xyz[3]) const
297 {
298   //
299   // get cluster coordinates in global cooordinate 
300   //
301   xyz[2] = cl->GetZ();
302   xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
303   xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;
304 }
305 #endif