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