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