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 */
6 //-------------------------------------------------------------------------
8 // reads AliITSclusterMI clusters and creates AliITStrackMI tracks
9 // Origin: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
10 //-------------------------------------------------------------------------
12 #include <TObjArray.h>
14 #include "AliTracker.h"
15 #include "AliITStrackMI.h"
16 #include "AliITSclusterV2.h"
17 #include "AliV0vertex.h"
29 //-------------------------------------------------------------------------
30 class AliITStrackerMI : public AliTracker {
32 AliITStrackerMI():AliTracker(){}
33 virtual ~AliITStrackerMI() {
34 delete[] fCoeficients;
37 AliITStrackerMI(const AliITSgeom *geom);
38 AliCluster *GetCluster(Int_t index) const;
39 AliITSclusterV2 *GetClusterLayer(Int_t layn, Int_t ncl) const
40 {return fgLayers[layn].GetCluster(ncl);}
41 Int_t GetNumberOfClustersLayer(Int_t layn) const
42 {return fgLayers[layn].GetNumberOfClusters();}
43 Int_t LoadClusters(TTree *cf);
44 void UnloadClusters();
45 Int_t Clusters2Tracks(AliESD *event);
46 Int_t PropagateBack(AliESD *event);
47 Int_t RefitInward(AliESD *event);
48 Bool_t RefitAt(Double_t x, AliITStrackMI *seed, const AliITStrackMI *t);
49 void SetupFirstPass(Int_t *flags, Double_t *cuts=0);
50 void SetupSecondPass(Int_t *flags, Double_t *cuts=0);
52 void SetLastLayerToTrackTo(Int_t l=0) {fLastLayerToTrackTo=l;}
53 void SetLayersNotToSkip(Int_t *l);
54 void UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
55 void GetNTeor(Int_t layer, const AliITSclusterV2* cl, Float_t theta, Float_t phi, Float_t &ny, Float_t &nz);
56 Int_t GetError(Int_t layer, const AliITSclusterV2*cl, Float_t theta, Float_t phi, Float_t expQ, Float_t &erry, Float_t &errz);
57 Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSclusterV2 *cluster,Int_t layer);
58 Int_t UpdateMI(AliITStrackMI* track, const AliITSclusterV2* cl,Double_t chi2,Int_t layer) const;
59 class AliITSdetector {
62 AliITSdetector(Double_t r,Double_t phi) {fR=r; fPhi=phi; fSinPhi = TMath::Sin(phi); fCosPhi = TMath::Cos(phi);
63 fYmin=10000;fYmax=-1000; fZmin=10000;fZmax=-1000;}
64 inline void GetGlobalXYZ( const AliITSclusterV2 *cl, Double_t xyz[3]) const;
65 Double_t GetR() const {return fR;}
66 Double_t GetPhi() const {return fPhi;}
67 Double_t GetYmin() const {return fYmin;}
68 Double_t GetYmax() const {return fYmax;}
69 Double_t GetZmin() const {return fZmin;}
70 Double_t GetZmax() const {return fZmax;}
71 void SetYmin(Double_t min) {fYmin = min;}
72 void SetYmax(Double_t max) {fYmax = max;}
73 void SetZmin(Double_t min) {fZmin = min;}
74 void SetZmax(Double_t max) {fZmax = max;}
76 Double_t fR; // polar coordinates
77 Double_t fPhi; // of this detector
78 Double_t fSinPhi; // sin of phi;
79 Double_t fCosPhi; // cos of phi
80 Double_t fYmin; // local y minimal
81 Double_t fYmax; // local max y
82 Double_t fZmin; // local z min
83 Double_t fZmax; // local z max
87 friend class AliITStrackerMI;
90 AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
92 Int_t InsertCluster(AliITSclusterV2 *c);
96 void SelectClusters(Double_t zmi,Double_t zma,Double_t ymi,Double_t yma);
97 const AliITSclusterV2 *GetNextCluster(Int_t &ci);
99 Double_t GetRoad() const {return fRoad;}
100 Double_t GetR() const {return fR;}
101 Int_t FindClusterIndex(Float_t z) const;
102 AliITSclusterV2 *GetCluster(Int_t i) const {return i<fN? fClusters[i]:0;}
103 Float_t *GetWeight(Int_t i) {return i<fN ?&fClusterWeight[i]:0;}
104 AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
105 Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
106 Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
107 Int_t InRoad() const ;
108 Int_t GetNumberOfClusters() const {return fN;}
109 Int_t GetNladders() const {return fNladders;}
110 Int_t GetNdetectors() const {return fNdetectors;}
111 Int_t GetSkip() const {return fSkip;}
112 void SetSkip(Int_t skip){fSkip=skip;}
113 void IncAccepted(){fAccepted++;}
114 Int_t GetAccepted() const {return fAccepted;}
116 AliITSlayer(const AliITSlayer& layer){;}
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 AliITSclusterV2 *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
130 AliITSclusterV2 *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
138 AliITSclusterV2 *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
146 AliITSclusterV2 *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
154 AliITSclusterV2** 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
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
172 AliITStrackerMI::AliITSlayer & GetLayer(Int_t layer) const;
173 AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
176 void FindV0(AliESD *event); //try to find V0
177 Double_t TestV0(AliHelix *h1, AliHelix *h2, AliESDV0MI *vertex); //try to find V0 - return DCA
178 Double_t FindBestPair(Int_t esdtrack0, Int_t esdtrack1,AliESDV0MI *vertex); // try to find best pair from the tree of track hyp.
179 void CookLabel(AliKalmanTrack *t,Float_t wrong) const;
180 void CookLabel(AliITStrackMI *t,Float_t wrong) const;
181 Double_t GetEffectiveThickness(Double_t y, Double_t z) const;
182 void FollowProlongationTree(AliITStrackMI * otrack, Int_t esdindex);
183 void ResetBestTrack() {
184 fBestTrack.~AliITStrackMI();
185 new(&fBestTrack) AliITStrackMI(fTrackToFollow);
187 void ResetTrackToFollow(const AliITStrackMI &t) {
188 fTrackToFollow.~AliITStrackMI();
189 new(&fTrackToFollow) AliITStrackMI(t);
191 void CookdEdx(AliITStrackMI* track);
192 Double_t GetNormalizedChi2(AliITStrackMI * track, Int_t mode);
193 Double_t GetTruncatedChi2(AliITStrackMI * track, Float_t fac);
194 Double_t NormalizedChi2(AliITStrackMI * track, Int_t layer);
195 Double_t GetInterpolatedChi2(AliITStrackMI * forwardtrack, AliITStrackMI * backtrack);
196 Double_t GetMatchingChi2(AliITStrackMI * track1, AliITStrackMI * track2);
197 Double_t GetDeadZoneProbability(Double_t zpos, Double_t zerr);
199 Float_t *GetWeight(Int_t index);
200 void AddTrackHypothesys(AliITStrackMI * track, Int_t esdindex);
201 void SortTrackHypothesys(Int_t esdindex, Int_t maxcut, Int_t mode);
202 AliITStrackMI * GetBestHypothesys(Int_t esdindex, AliITStrackMI * original, Int_t checkmax);
203 void GetBestHypothesysMIP(TObjArray &itsTracks);
204 void RegisterClusterTracks(AliITStrackMI* track, Int_t id);
205 void UnRegisterClusterTracks(AliITStrackMI* track, Int_t id);
206 Float_t GetNumberOfSharedClusters(AliITStrackMI* track,Int_t id, Int_t list[6], AliITSclusterV2 *clist[6]);
207 Int_t GetOverlapTrack(AliITStrackMI *track, Int_t trackID, Int_t &shared, Int_t clusterlist[6], Int_t overlist[6]);
208 AliITStrackMI * GetBest2Tracks(Int_t trackID1, Int_t treackID2, Float_t th0, Float_t th1);
209 Float_t * GetErrY(Int_t trackindex) const {return &fCoeficients[trackindex*48];}
210 Float_t * GetErrZ(Int_t trackindex) const {return &fCoeficients[trackindex*48+12];}
211 Float_t * GetNy(Int_t trackindex) const {return &fCoeficients[trackindex*48+24];}
212 Float_t * GetNz(Int_t trackindex) const {return &fCoeficients[trackindex*48+36];}
213 void SignDeltas( TObjArray *ClusterArray, Float_t zv);
214 void MakeCoeficients(Int_t ntracks);
215 void UpdateESDtrack(AliITStrackMI* track, ULong_t flags) const;
216 Int_t fI; // index of the current layer
217 static AliITSlayer fgLayers[kMaxLayer];// ITS layers
218 AliITStrackMI fTracks[kMaxLayer]; // track estimations at the ITS layers
219 AliITStrackMI fBestTrack; // "best" track
220 AliITStrackMI fTrackToFollow; // followed track
221 TObjArray fTrackHypothesys; // ! array with track hypothesys- ARRAY is the owner of tracks- MI
222 TObjArray fV0Array; // ! array of V0
223 Int_t fBestTrackIndex[100000]; // ! index of the best track
224 Int_t fCurrentEsdTrack; // ! current esd track - MI
225 Int_t fPass; // current pass through the data
226 Int_t fConstraint[2]; // constraint flags
228 Int_t fLayersNotToSkip[kMaxLayer]; // layer masks
229 Int_t fLastLayerToTrackTo; // the innermost layer to track to
230 Float_t * fCoeficients; //! working array with errors and mean cluser shape
232 AliITStrackerMI(const AliITStrackerMI * tracker){;}
233 ClassDef(AliITStrackerMI,2) //ITS tracker MI
239 /////////////////////////////////////////////////////////
240 /////////////////////////////////////////////////////////
241 /////////////////////////////////////////////////////////
247 inline void AliITStrackerMI::SetupFirstPass(Int_t *flags, Double_t *cuts) {
248 // This function sets up flags and cuts for the first tracking pass
250 // flags[0] - vertex constaint flag
251 // negative means "skip the pass"
252 // 0 means "no constraint"
253 // positive means "normal constraint"
255 fConstraint[0]=flags[0];
259 inline void AliITStrackerMI::SetupSecondPass(Int_t *flags, Double_t *cuts) {
260 // This function sets up flags and cuts for the second tracking pass
262 // flags[0] - vertex constaint flag
263 // negative means "skip the pass"
264 // 0 means "no constraint"
265 // positive means "normal constraint"
267 fConstraint[1]=flags[0];
271 inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
272 //--------------------------------------------------------------------
273 //This function "cooks" a track label. If label<0, this track is fake.
274 //--------------------------------------------------------------------
275 Int_t tpcLabel=t->GetLabel();
276 if (tpcLabel<0) return;
277 AliTracker::CookLabel(t,wrong);
278 if (tpcLabel!=TMath::Abs(t->GetLabel())){
281 if (tpcLabel !=t->GetLabel()) {
282 t->SetLabel(-tpcLabel);
286 inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
288 //--------------------------------------------------------------------
290 //--------------------------------------------------------------------
291 track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->GetChi2()/
292 //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
293 TMath::Max(double(track->GetNumberOfClusters()-track->fNSkipped),
294 1./(1.+track->fNSkipped));
295 return track->fNormChi2[layer];
297 inline void AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSclusterV2 *cl, Double_t xyz[3]) const
300 // get cluster coordinates in global cooordinate
303 xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
304 xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;