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