]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerMI.h
Adding AliComparisonDraw (Marian)
[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 //-------------------------------------------------------------------------
8 //                          ITS tracker
9 //     reads AliITSclusterMI clusters and creates AliITStrackMI tracks
10 //           Origin: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
11 //-------------------------------------------------------------------------
12
13 class TTree;
14 class TTreeSRedirector;
15 class AliESDEvent;
16
17 #include <TObjArray.h>
18 #include "AliITSRecPoint.h"
19 #include "AliITStrackMI.h"
20 #include "AliTracker.h"
21
22 //-------------------------------------------------------------------------
23 class AliITStrackerMI : public AliTracker {
24 public:
25   AliITStrackerMI();
26   AliITStrackerMI(const Char_t *geom);
27   virtual ~AliITStrackerMI();
28   AliCluster *GetCluster(Int_t index) const;
29   virtual Bool_t GetTrackPoint(Int_t index, AliTrackPoint& p) const;
30   virtual Bool_t GetTrackPointTrackingError(Int_t index, 
31                         AliTrackPoint& p, const AliESDtrack *t);
32   AliITSRecPoint *GetClusterLayer(Int_t layn, Int_t ncl) const
33                         {return fgLayers[layn].GetCluster(ncl);}
34   Int_t GetNumberOfClustersLayer(Int_t layn) const 
35                         {return fgLayers[layn].GetNumberOfClusters();}
36   Int_t LoadClusters(TTree *cf);
37   void UnloadClusters();
38   Int_t Clusters2Tracks(AliESDEvent *event);
39   Int_t PropagateBack(AliESDEvent *event);
40   Int_t RefitInward(AliESDEvent *event);
41   Bool_t RefitAt(Double_t x, AliITStrackMI *seed, 
42                  const AliITStrackMI *t, Bool_t extra=kFALSE);
43   Bool_t RefitAt(Double_t x, AliITStrackMI *seed, const Int_t *clindex);
44   void SetupFirstPass(Int_t *flags, Double_t *cuts=0);
45   void SetupSecondPass(Int_t *flags, Double_t *cuts=0);
46
47   void SetLastLayerToTrackTo(Int_t l=0) {fLastLayerToTrackTo=l;} 
48   void SetLayersNotToSkip(Int_t *l);
49   void UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
50   void GetNTeor(Int_t layer, const AliITSRecPoint* cl, Float_t theta, Float_t phi, Float_t &ny, Float_t &nz);
51   Int_t  GetError(Int_t layer, const AliITSRecPoint*cl, Float_t theta, Float_t phi, Float_t expQ, Float_t &erry, Float_t &errz);
52
53   void  GetDCASigma(AliITStrackMI* track, Float_t & sigmarfi, Float_t &sigmaz);
54   Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSRecPoint *cluster,Int_t layer);
55   Int_t UpdateMI(AliITStrackMI* track, const AliITSRecPoint* cl,Double_t chi2,Int_t layer) const;
56   class AliITSdetector { 
57   public:
58     AliITSdetector():fR(0),fPhi(0),fSinPhi(0),fCosPhi(0),fYmin(0),fYmax(0),fZmin(0),fZmax(0){}
59     AliITSdetector(Double_t r,Double_t phi):fR(r),fPhi(phi),fSinPhi(TMath::Sin(phi)),fCosPhi(TMath::Cos(phi)),fYmin(10000),fYmax(-1000),fZmin(10000),fZmax(-1000) {}
60     inline void GetGlobalXYZ( const AliITSRecPoint *cl, Double_t xyz[3]) const;
61     Double_t GetR()   const {return fR;}
62     Double_t GetPhi() const {return fPhi;}
63     Double_t GetYmin() const {return fYmin;}
64     Double_t GetYmax() const {return fYmax;}
65     Double_t GetZmin() const {return fZmin;}
66     Double_t GetZmax() const {return fZmax;}
67     void SetYmin(Double_t min) {fYmin = min;}
68     void SetYmax(Double_t max) {fYmax = max;}
69     void SetZmin(Double_t min) {fZmin = min;}
70     void SetZmax(Double_t max) {fZmax = max;}
71   private:
72     Double_t fR;    // polar coordinates 
73     Double_t fPhi;  // of this detector
74     Double_t fSinPhi; // sin of phi;
75     Double_t fCosPhi; // cos of phi 
76     Double_t fYmin;   //  local y minimal
77     Double_t fYmax;   //  local max y
78     Double_t fZmin;   //  local z min
79     Double_t fZmax;   //  local z max
80   };
81
82   class AliITSlayer {
83   public:
84     AliITSlayer();
85     AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
86    ~AliITSlayer();
87     Int_t InsertCluster(AliITSRecPoint *c);
88     void  SortClusters();
89     void ResetClusters();
90     void ResetWeights();
91     void SelectClusters(Double_t zmi,Double_t zma,Double_t ymi,Double_t yma);
92     const AliITSRecPoint *GetNextCluster(Int_t &ci);
93     void ResetRoad();
94     Double_t GetRoad() const {return fRoad;}
95     Double_t GetR() const {return fR;}
96     Int_t FindClusterIndex(Float_t z) const;
97     AliITSRecPoint *GetCluster(Int_t i) const {return i<fN? fClusters[i]:0;} 
98     Float_t         *GetWeight(Int_t i) {return i<fN ?&fClusterWeight[i]:0;}
99     AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
100     Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
101     Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
102     Int_t InRoad() const ;
103     Int_t GetNumberOfClusters() const {return fN;}
104     Int_t GetNladders() const {return fNladders;}
105     Int_t GetNdetectors() const {return fNdetectors;}
106     Int_t GetSkip() const {return fSkip;}
107     void  SetSkip(Int_t skip){fSkip=skip;}
108     void IncAccepted(){fAccepted++;}
109     Int_t GetAccepted() const {return fAccepted;}    
110     Int_t GetClusterTracks(Int_t i, Int_t j) const {return fClusterTracks[i][j];}
111     void SetClusterTracks(Int_t i, Int_t j, Int_t c) {fClusterTracks[i][j]=c;}
112   protected:
113     AliITSlayer(const AliITSlayer& layer);
114     AliITSlayer & operator=(const AliITSlayer& layer){
115       this->~AliITSlayer();new(this) AliITSlayer(layer);
116       return *this;}
117     Double_t fR;                // mean radius of this layer
118     Double_t fPhiOffset;        // offset of the first detector in Phi
119     Int_t fNladders;            // number of ladders
120     Double_t fZOffset;          // offset of the first detector in Z
121     Int_t fNdetectors;          // detectors/ladder
122     AliITSdetector *fDetectors; // array of detectors
123     Int_t fN;                   // number of clusters
124     AliITSRecPoint *fClusters[kMaxClusterPerLayer]; // pointers to clusters
125     Int_t        fClusterIndex[kMaxClusterPerLayer]; // pointers to clusters
126     Float_t fY[kMaxClusterPerLayer];                // y position of the clusters      
127     Float_t fZ[kMaxClusterPerLayer];                // z position of the clusters      
128     Float_t fYB[2];                                       // ymin and ymax
129     //
130     AliITSRecPoint *fClusters5[6][kMaxClusterPerLayer5]; // pointers to clusters -     slice in y
131     Int_t        fClusterIndex5[6][kMaxClusterPerLayer5]; // pointers to clusters -     slice in y    
132     Float_t fY5[6][kMaxClusterPerLayer5];                // y position of the clusters  slice in y    
133     Float_t fZ5[6][kMaxClusterPerLayer5];                // z position of the clusters  slice in y 
134     Int_t fN5[6];                                       // number of cluster in slice
135     Float_t fDy5;                                       //delta y
136     Float_t fBy5[6][2];                                    //slice borders
137     //
138     AliITSRecPoint *fClusters10[11][kMaxClusterPerLayer10]; // pointers to clusters -     slice in y
139     Int_t        fClusterIndex10[11][kMaxClusterPerLayer10]; // pointers to clusters -     slice in y    
140     Float_t fY10[11][kMaxClusterPerLayer10];                // y position of the clusters  slice in y    
141     Float_t fZ10[11][kMaxClusterPerLayer10];                // z position of the clusters  slice in y 
142     Int_t fN10[11];                                       // number of cluster in slice
143     Float_t fDy10;                                        // delta y
144     Float_t fBy10[11][2];                                 // slice borders
145     //
146     AliITSRecPoint *fClusters20[21][kMaxClusterPerLayer20]; // pointers to clusters -     slice in y
147     Int_t        fClusterIndex20[21][kMaxClusterPerLayer20]; // pointers to clusters -     slice in y    
148     Float_t fY20[21][kMaxClusterPerLayer20];                // y position of the clusters  slice in y    
149     Float_t fZ20[21][kMaxClusterPerLayer20];                // z position of the clusters  slice in y 
150     Int_t fN20[21];                                       // number of cluster in slice
151     Float_t fDy20;                                        //delta y 
152     Float_t fBy20[21][2];                                 //slice borders
153     //
154     AliITSRecPoint** fClustersCs;                         //clusters table in current slice
155     Int_t   *fClusterIndexCs;                             //cluster index in current slice 
156     Float_t *fYcs;                                        //y position in current slice
157     Float_t *fZcs;                                        //z position in current slice
158     Int_t    fNcs;                                        //number of clusters in current slice    
159     Int_t fCurrentSlice;                                  //current slice
160     //
161     Float_t  fClusterWeight[kMaxClusterPerLayer]; // probabilistic weight of the cluster
162     Int_t    fClusterTracks[4][kMaxClusterPerLayer]; //tracks registered to given cluster
163     Float_t fZmax;      //    edges
164     Float_t fYmin;      //   of  the
165     Float_t fYmax;      //   "window"
166     Int_t fI;            // index of the current cluster within the "window"
167     Int_t fImax;            // index of the last cluster within the "window"    
168     Int_t fSkip;     // indicates possibility to skip cluster
169     Int_t fAccepted;     // accept indicator 
170     Double_t fRoad;      // road defined by the cluster density
171   };
172   AliITStrackerMI::AliITSlayer    & GetLayer(Int_t layer) const;
173   AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
174
175 protected:
176   Int_t GetNearestLayer(const Double_t *xr) const;  //get nearest upper layer close to the point xr
177   void FindV02(AliESDEvent *event);  //try to find V0
178   void RefitV02(AliESDEvent *event);  //try to refit  V0's
179   void UpdateTPCV0(AliESDEvent *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], AliITSRecPoint *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 &fCoefficients[trackindex*48];}
211   Float_t  * GetErrZ(Int_t trackindex) const {return &fCoefficients[trackindex*48+12];}
212   Float_t  * GetNy(Int_t trackindex) const {return &fCoefficients[trackindex*48+24];}
213   Float_t  * GetNz(Int_t trackindex) const {return &fCoefficients[trackindex*48+36];}
214   void       SignDeltas( TObjArray *ClusterArray, Float_t zv);
215   void MakeCoefficients(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 * fCoefficients;                //! working array with errors and mean cluster shape
233   AliESDEvent  * fEsd;                   //! pointer to the ESD event
234   Double_t fSPDdetzcentre[4];            // centres of SPD modules in z
235   Bool_t fUseTGeo;                       // use TGeo to get material budget
236   TTreeSRedirector *fDebugStreamer;     //!debug streamer
237 private:
238   AliITStrackerMI(const AliITStrackerMI &tracker);
239   AliITStrackerMI & operator=(const AliITStrackerMI &tracker);
240   ClassDef(AliITStrackerMI,3)   //ITS tracker MI
241 };
242
243
244
245
246 /////////////////////////////////////////////////////////
247 /////////////////////////////////////////////////////////
248 /////////////////////////////////////////////////////////
249
250
251
252
253
254 inline void AliITStrackerMI::SetupFirstPass(Int_t *flags, Double_t *cuts) {
255   // This function sets up flags and cuts for the first tracking pass   
256   //
257   //   flags[0] - vertex constaint flag                                
258   //              negative means "skip the pass"                        
259   //              0        means "no constraint"                        
260   //              positive means "normal constraint"                    
261
262    fConstraint[0]=flags[0];
263    if (cuts==0) return;
264 }
265
266 inline void AliITStrackerMI::SetupSecondPass(Int_t *flags, Double_t *cuts) {
267   // This function sets up flags and cuts for the second tracking pass   
268   //
269   //   flags[0] - vertex constaint flag                                
270   //              negative means "skip the pass"                        
271   //              0        means "no constraint"                        
272   //              positive means "normal constraint"                    
273
274    fConstraint[1]=flags[0];
275    if (cuts==0) return;
276 }
277
278 inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
279   //--------------------------------------------------------------------
280   //This function "cooks" a track label. If label<0, this track is fake.
281   //--------------------------------------------------------------------
282    Int_t tpcLabel=t->GetLabel();
283    if (tpcLabel<0) return;
284    AliTracker::CookLabel(t,wrong);
285    if (tpcLabel!=TMath::Abs(t->GetLabel())){
286      t->SetFakeRatio(1.);
287    }
288    if (tpcLabel !=t->GetLabel()) {
289      t->SetLabel(-tpcLabel);      
290    }
291 }
292
293 inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
294 {
295   //--------------------------------------------------------------------
296   //get normalize chi2
297   //--------------------------------------------------------------------
298   track->SetNormChi2(layer,2.*track->GetNSkipped()+0.25*track->GetNDeadZone()+track->GetdEdxMismatch()+track->GetChi2()/
299   //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
300     TMath::Max(double(track->GetNumberOfClusters()-track->GetNSkipped()),
301                1./(1.+track->GetNSkipped())));
302   return track->GetNormChi2(layer);
303 }
304 inline void  AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSRecPoint *cl, Double_t xyz[3]) const
305 {
306   //
307   // get cluster coordinates in global cooordinate 
308   //
309   xyz[2] = cl->GetZ();
310   xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
311   xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;
312 }
313 #endif