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