]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerMI.h
Avoiding geometry lock/unlock (R. Grosso)
[u/mrichter/AliRoot.git] / ITS / AliITStrackerMI.h
1 #ifndef ALIITSTRACKERMI_H
2 #define ALIITSTRACKERMI_H
3 /* Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 //-------------------------------------------------------------------------
9 //                          ITS tracker
10 //     reads AliITSclusterMI clusters and creates AliITStrackMI tracks
11 //           Origin: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
12 //           Current support and development: 
13 //                     Andrea Dainese, andrea.dainese@lnl.infn.it
14 //-------------------------------------------------------------------------
15
16 class TTree;
17 class TTreeSRedirector;
18 class AliESDEvent;
19
20 class AliITSPlaneEff;
21 class AliITSChannelStatus;
22 class AliITSDetTypeRec;
23 class AliPlaneEff;
24
25 #include <TObjArray.h>
26
27 #include "AliITStrackMI.h"
28 #include "AliITSRecPoint.h"
29 #include "AliTracker.h"
30
31
32 //-------------------------------------------------------------------------
33 class AliITStrackerMI : public AliTracker {
34 public:
35   AliITStrackerMI();
36   AliITStrackerMI(const Char_t *geom);
37   virtual ~AliITStrackerMI();
38   AliCluster *GetCluster(Int_t index) const;
39   virtual Bool_t GetTrackPoint(Int_t index, AliTrackPoint& p) const;
40   virtual Bool_t GetTrackPointTrackingError(Int_t index, 
41                         AliTrackPoint& p, const AliESDtrack *t);
42   AliITSRecPoint *GetClusterLayer(Int_t layn, Int_t ncl) const
43                         {return fgLayers[layn].GetCluster(ncl);}
44   Int_t GetNumberOfClustersLayer(Int_t layn) const 
45                         {return fgLayers[layn].GetNumberOfClusters();}
46   Int_t LoadClusters(TTree *cf);
47   void UnloadClusters();
48   void FillClusterArray(TObjArray* array) const;
49   Int_t Clusters2Tracks(AliESDEvent *event);
50   Int_t PropagateBack(AliESDEvent *event);
51   Int_t RefitInward(AliESDEvent *event);
52   Bool_t RefitAt(Double_t x, AliITStrackMI *track, 
53                  const AliITStrackMI *clusters, Bool_t extra=kFALSE, Bool_t planeeff=kFALSE);
54   Bool_t RefitAt(Double_t x, AliITStrackMI *track, 
55                  const Int_t *clusters, Bool_t extra=kFALSE, Bool_t planeeff=kFALSE);
56   void SetupFirstPass(const Int_t *flags,const Double_t *cuts=0);
57   void SetupSecondPass(const Int_t *flags,const Double_t *cuts=0);
58
59   void SetLastLayerToTrackTo(Int_t l=0) {fLastLayerToTrackTo=l;} 
60   void SetLayersNotToSkip(const Int_t *l);
61   void UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
62
63   void  GetDCASigma(const AliITStrackMI* track, Float_t & sigmarfi, Float_t &sigmaz);
64   Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSRecPoint *cluster,Int_t layer);
65   Int_t UpdateMI(AliITStrackMI* track, const AliITSRecPoint* cl,Double_t chi2,Int_t layer) const;
66   AliPlaneEff *GetPlaneEff() {return (AliPlaneEff*)fPlaneEff;}   // return the pointer to AliPlaneEff
67   void SetDetTypeRec(const AliITSDetTypeRec *detTypeRec) {fkDetTypeRec = detTypeRec; ReadBadFromDetTypeRec(); }
68   TObjArray* GetTrackHypothesys() {return &fTrackHypothesys;}
69   TObjArray* GetBestHypothesys()  {return &fBestHypothesys;}
70   TObjArray* GetOriginal()        {return &fOriginal;}
71   TTreeSRedirector *GetDebugStreamer() {return fDebugStreamer;}
72   static Int_t CorrectForTPCtoITSDeadZoneMaterial(AliITStrackMI *t);
73
74   class AliITSdetector { 
75   public:
76     AliITSdetector():fR(0),fRmisal(0),fPhi(0),fSinPhi(0),fCosPhi(0),fYmin(0),fYmax(0),fZmin(0),fZmax(0),fIsBad(kFALSE),fNChips(0),fChipIsBad(0) {}
77     AliITSdetector(Double_t r,Double_t phi):fR(r),fRmisal(r),fPhi(phi),fSinPhi(TMath::Sin(phi)),fCosPhi(TMath::Cos(phi)),fYmin(10000),fYmax(-1000),fZmin(10000),fZmax(-1000),fIsBad(kFALSE),fNChips(0),fChipIsBad(0) {}
78     ~AliITSdetector() {if(fChipIsBad) delete [] fChipIsBad;}
79     inline void GetGlobalXYZ( const AliITSRecPoint *cl, Double_t xyz[3]) const;
80     Double_t GetR()   const {return fR;}
81     Double_t GetRmisal()   const {return fRmisal;}
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     Bool_t   IsBad() const {return fIsBad;}
88     Int_t    GetNChips() const {return fNChips;}
89     Bool_t   IsChipBad(Int_t iChip) const {return (fChipIsBad ? fChipIsBad[iChip] : kFALSE);}
90     void SetRmisal(Double_t rmisal) {fRmisal = rmisal;}
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     void SetBad() {fIsBad = kTRUE;}
96     void ReadBadDetectorAndChips(Int_t ilayer,Int_t idet,const AliITSDetTypeRec *detTypeRec);
97   private:
98     AliITSdetector(const AliITSdetector& det);
99     AliITSdetector & operator=(const AliITSdetector& det){
100       this->~AliITSdetector();new(this) AliITSdetector(det);
101       return *this;}
102     Double_t fR;    // polar coordinates: r 
103     Double_t fRmisal;    // polar coordinates: r, with misalignment 
104     Double_t fPhi;  // polar coordinates: phi
105     Double_t fSinPhi; // sin of phi;
106     Double_t fCosPhi; // cos of phi 
107     Double_t fYmin;   //  local y minimal
108     Double_t fYmax;   //  local max y
109     Double_t fZmin;   //  local z min
110     Double_t fZmax;   //  local z max
111     Bool_t fIsBad;    // is detector dead or noisy?
112     Int_t fNChips;    // number of chips
113     Bool_t *fChipIsBad; //[fNChips] is chip dead or noisy? 
114   };
115
116   class AliITSlayer {
117   public:
118     AliITSlayer();
119     AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
120     ~AliITSlayer();
121     Int_t InsertCluster(AliITSRecPoint *c);
122     void  SortClusters();
123     void ResetClusters();
124     void ResetWeights();
125     void SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin,Double_t ymax);
126     const AliITSRecPoint *GetNextCluster(Int_t &ci,Bool_t test=kFALSE);
127     void ResetRoad();
128     Double_t GetRoad() const {return fRoad;}
129     Double_t GetR() const {return fR;}
130     Int_t FindClusterIndex(Float_t z) const;
131     AliITSRecPoint *GetCluster(Int_t i) const {return i<fN? fClusters[i]:0;} 
132     Float_t  *GetWeight(Int_t i)  {return i<fN ?&fClusterWeight[i]:0;}
133     AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
134     Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
135     Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
136     Int_t InRoad() const ;
137     Int_t GetNumberOfClusters() const {return fN;}
138     Int_t GetNladders() const {return fNladders;}
139     Int_t GetNdetectors() const {return fNdetectors;}
140     Int_t GetSkip() const {return fSkip;}
141     void  SetSkip(Int_t skip){fSkip=skip;}
142     void IncAccepted(){fAccepted++;}
143     Int_t GetAccepted() const {return fAccepted;}    
144     Int_t GetClusterTracks(Int_t i, Int_t j) const {return fClusterTracks[i][j];}
145     void SetClusterTracks(Int_t i, Int_t j, Int_t c) {fClusterTracks[i][j]=c;}
146   protected:
147     AliITSlayer(const AliITSlayer& layer);
148     AliITSlayer & operator=(const AliITSlayer& layer){
149       this->~AliITSlayer();new(this) AliITSlayer(layer);
150       return *this;}
151     Double_t fR;                // mean radius of this layer
152     Double_t fPhiOffset;        // offset of the first detector in Phi
153     Int_t fNladders;            // number of ladders
154     Double_t fZOffset;          // offset of the first detector in Z
155     Int_t fNdetectors;          // detectors/ladder
156     AliITSdetector *fDetectors; // array of detectors
157     Int_t fN;                   // number of clusters
158     AliITSRecPoint *fClusters[AliITSRecoParam::fgkMaxClusterPerLayer]; // pointers to clusters
159     Int_t        fClusterIndex[AliITSRecoParam::fgkMaxClusterPerLayer]; // pointers to clusters
160     Float_t fY[AliITSRecoParam::fgkMaxClusterPerLayer];                // y position of the clusters      
161     Float_t fZ[AliITSRecoParam::fgkMaxClusterPerLayer];                // z position of the clusters      
162     Float_t fYB[2];                                       // ymin and ymax
163     //
164     AliITSRecPoint *fClusters5[6][AliITSRecoParam::fgkMaxClusterPerLayer5]; // pointers to clusters -     slice in y
165     Int_t        fClusterIndex5[6][AliITSRecoParam::fgkMaxClusterPerLayer5]; // pointers to clusters -     slice in y    
166     Float_t fY5[6][AliITSRecoParam::fgkMaxClusterPerLayer5];                // y position of the clusters  slice in y    
167     Float_t fZ5[6][AliITSRecoParam::fgkMaxClusterPerLayer5];                // z position of the clusters  slice in y 
168     Int_t fN5[6];                                       // number of cluster in slice
169     Float_t fDy5;                                       //delta y
170     Float_t fBy5[6][2];                                    //slice borders
171     //
172     AliITSRecPoint *fClusters10[11][AliITSRecoParam::fgkMaxClusterPerLayer10]; // pointers to clusters -     slice in y
173     Int_t        fClusterIndex10[11][AliITSRecoParam::fgkMaxClusterPerLayer10]; // pointers to clusters -     slice in y    
174     Float_t fY10[11][AliITSRecoParam::fgkMaxClusterPerLayer10];                // y position of the clusters  slice in y    
175     Float_t fZ10[11][AliITSRecoParam::fgkMaxClusterPerLayer10];                // z position of the clusters  slice in y 
176     Int_t fN10[11];                                       // number of cluster in slice
177     Float_t fDy10;                                        // delta y
178     Float_t fBy10[11][2];                                 // slice borders
179     //
180     AliITSRecPoint *fClusters20[21][AliITSRecoParam::fgkMaxClusterPerLayer20]; // pointers to clusters -     slice in y
181     Int_t        fClusterIndex20[21][AliITSRecoParam::fgkMaxClusterPerLayer20]; // pointers to clusters -     slice in y    
182     Float_t fY20[21][AliITSRecoParam::fgkMaxClusterPerLayer20];                // y position of the clusters  slice in y    
183     Float_t fZ20[21][AliITSRecoParam::fgkMaxClusterPerLayer20];                // z position of the clusters  slice in y 
184     Int_t fN20[21];                                       // number of cluster in slice
185     Float_t fDy20;                                        //delta y 
186     Float_t fBy20[21][2];                                 //slice borders
187     //
188     AliITSRecPoint** fClustersCs;                         //clusters table in current slice
189     Int_t   *fClusterIndexCs;                             //cluster index in current slice 
190     Float_t *fYcs;                                        //y position in current slice
191     Float_t *fZcs;                                        //z position in current slice
192     Int_t    fNcs;                                        //number of clusters in current slice    
193     Int_t fCurrentSlice;                                  //current slice
194     //
195     Float_t  fClusterWeight[AliITSRecoParam::fgkMaxClusterPerLayer]; // probabilistic weight of the cluster
196     Int_t    fClusterTracks[4][AliITSRecoParam::fgkMaxClusterPerLayer]; //tracks registered to given cluster
197     Float_t fZmax;      //    edges
198     Float_t fYmin;      //   of  the
199     Float_t fYmax;      //   "window"
200     Int_t fI;            // index of the current cluster within the "window"
201     Int_t fImax;            // index of the last cluster within the "window"    
202     Int_t fSkip;     // indicates possibility to skip cluster
203     Int_t fAccepted;     // accept indicator 
204     Double_t fRoad;      // road defined by the cluster density
205   };
206   AliITStrackerMI::AliITSlayer    & GetLayer(Int_t layer) const;
207   AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
208   Int_t GetNearestLayer(const Double_t *xr) const;  //get nearest upper layer close to the point xr
209   void SetCurrentEsdTrack(Int_t i) {fCurrentEsdTrack=i;}
210   void FollowProlongationTree(AliITStrackMI * otrack, Int_t esdindex, Bool_t constrain);
211
212 protected:
213   Bool_t ComputeRoad(AliITStrackMI* track,Int_t ilayer,Int_t idet,Double_t &zmin,Double_t &zmax,Double_t &ymin,Double_t &ymax) const;
214   
215   void FindV02(AliESDEvent *event);  //try to find V0
216   void RefitV02(const AliESDEvent *event);  //try to refit  V0's
217   void UpdateTPCV0(const AliESDEvent *event);  //try to update, or reject TPC  V0s
218   
219   void CookLabel(AliKalmanTrack *t,Float_t wrong) const;
220   void CookLabel(AliITStrackMI *t,Float_t wrong) const;
221   Double_t GetEffectiveThickness();
222   void ResetBestTrack() {
223      fBestTrack.~AliITStrackMI();
224      new(&fBestTrack) AliITStrackMI(fTrackToFollow);
225   }
226   void ResetTrackToFollow(const AliITStrackMI &t) {
227      fTrackToFollow.~AliITStrackMI();
228      new(&fTrackToFollow) AliITStrackMI(t);
229   }
230   void CookdEdx(AliITStrackMI* track);
231   Double_t GetNormalizedChi2(AliITStrackMI * track, Int_t mode);
232   Double_t GetTruncatedChi2(const AliITStrackMI * track, Float_t fac);
233   Double_t NormalizedChi2(AliITStrackMI * track, Int_t layer);
234   Double_t GetInterpolatedChi2(const AliITStrackMI * forwardtrack,const AliITStrackMI * backtrack);  
235   Double_t GetMatchingChi2(const AliITStrackMI * track1,const AliITStrackMI * track2);
236   Double_t GetSPDDeadZoneProbability(Double_t zpos, Double_t zerr) const;
237
238   Float_t    *GetWeight(Int_t index);
239   void AddTrackHypothesys(AliITStrackMI * track, Int_t esdindex);
240   void SortTrackHypothesys(Int_t esdindex, Int_t maxcut, Int_t mode);
241   AliITStrackMI * GetBestHypothesys(Int_t esdindex, AliITStrackMI * original, Int_t checkmax); 
242   void  GetBestHypothesysMIP(TObjArray &itsTracks); 
243   void RegisterClusterTracks(const AliITStrackMI* track, Int_t id);
244   void UnRegisterClusterTracks(const AliITStrackMI* track, Int_t id);
245   Float_t GetNumberOfSharedClusters(AliITStrackMI* track,Int_t id, Int_t list[6], AliITSRecPoint *clist[6]);
246   Int_t GetOverlapTrack(const AliITStrackMI *track, Int_t trackID, Int_t &shared, Int_t clusterlist[6], Int_t overlist[6]);
247   AliITStrackMI * GetBest2Tracks(Int_t trackID1, Int_t treackID2, Float_t th0, Float_t th1);
248   Float_t  * GetErrY(Int_t trackindex) const {return &fCoefficients[trackindex*48];}
249   Float_t  * GetErrZ(Int_t trackindex) const {return &fCoefficients[trackindex*48+12];}
250   Float_t  * GetNy(Int_t trackindex) const {return &fCoefficients[trackindex*48+24];}
251   Float_t  * GetNz(Int_t trackindex) const {return &fCoefficients[trackindex*48+36];}
252   void       SignDeltas(const TObjArray *clusterArray, Float_t zv);
253   void MakeCoefficients(Int_t ntracks);
254   void BuildMaterialLUT(TString material);
255   void MakeTrksMaterialLUT(Int_t ntracks);
256   void DeleteTrksMaterialLUT();
257   Int_t CorrectForPipeMaterial(AliITStrackMI *t, TString direction="inward");
258   Int_t CorrectForShieldMaterial(AliITStrackMI *t, TString shield, TString direction="inward");
259   Int_t CorrectForLayerMaterial(AliITStrackMI *t, Int_t layerindex, Double_t oldGlobXYZ[3], TString direction="inward");
260   void UpdateESDtrack(AliITStrackMI* track, ULong_t flags) const;
261   void ReadBadFromDetTypeRec();
262   Int_t CheckSkipLayer(const AliITStrackMI *track,Int_t ilayer,Int_t idet) const;
263   Int_t CheckDeadZone(AliITStrackMI *track,Int_t ilayer,Int_t idet,Double_t dz,Double_t dy,Bool_t noClusters=kFALSE) const;
264   Bool_t LocalModuleCoord(Int_t ilayer,Int_t idet,const AliITStrackMI *track,
265                           Float_t &xloc,Float_t &zloc) const;
266 // method to be used for Plane Efficiency evaluation
267   Bool_t IsOKForPlaneEff(const AliITStrackMI* track, const Int_t *clusters, Int_t ilayer) const; // Check if a track is usable
268                                                                                            // for Plane Eff evaluation
269   void UseTrackForPlaneEff(const AliITStrackMI* track, Int_t ilayer);                            // Use this track for Plane Eff
270 // 
271   Int_t fI;                              // index of the current layer
272   static AliITSlayer fgLayers[AliITSgeomTGeo::kNLayers];// ITS layers
273   AliITStrackMI fTracks[AliITSgeomTGeo::kNLayers];      // track estimations at the ITS layers
274   AliITStrackMI fBestTrack;              // "best" track 
275   AliITStrackMI fTrackToFollow;          // followed track
276   TObjArray     fTrackHypothesys;        // ! array with track hypothesys- ARRAY is the owner of tracks- MI
277   TObjArray     fBestHypothesys;         // ! array with track hypothesys- ARRAY is the owner of tracks- MI
278   TObjArray     fOriginal;               // ! array with seeds from the TPC
279   Int_t         fBestTrackIndex[100000]; // ! index of the best track
280   Int_t         fCurrentEsdTrack;        // ! current esd track           - MI
281   Int_t fPass;                           // current pass through the data 
282   Int_t fConstraint[2];                  // constraint flags
283   Bool_t fAfterV0;                       //indicates V0 founded
284   Int_t fLayersNotToSkip[AliITSgeomTGeo::kNLayers];     // layer masks
285   Int_t fLastLayerToTrackTo;             // the innermost layer to track to
286   Float_t * fCoefficients;               //! working array with errors and mean cluster shape
287   AliESDEvent  * fEsd;                   //! pointer to the ESD event
288   Double_t fSPDdetzcentre[4];            // centres of SPD modules in z
289   TString fTrackingPhase;                // current tracking phase
290   Int_t fUseTGeo;                        // use TGeo to get material budget
291   Int_t   fNtracks;                      // number of tracks to prolong
292   Float_t fxOverX0Pipe;                  // material budget
293   Float_t fxTimesRhoPipe;                // material budget
294   Float_t fxOverX0Shield[2];             // material budget
295   Float_t fxTimesRhoShield[2];           // material budget
296   Float_t fxOverX0Layer[6];              // material budget
297   Float_t fxTimesRhoLayer[6];            // material budget
298   Float_t *fxOverX0PipeTrks;             //! material budget
299   Float_t *fxTimesRhoPipeTrks;           //! material budget
300   Float_t *fxOverX0ShieldTrks;           //! material budget
301   Float_t *fxTimesRhoShieldTrks;         //! material budget
302   Float_t *fxOverX0LayerTrks;            //! material budget
303   Float_t *fxTimesRhoLayerTrks;          //! material budget
304   TTreeSRedirector *fDebugStreamer;      //!debug streamer
305   AliITSChannelStatus *fITSChannelStatus;//! bitmaps with channel status for SPD and SDD
306   const AliITSDetTypeRec *fkDetTypeRec;         //! ITS det type rec, from AliITSReconstructor
307   AliITSPlaneEff *fPlaneEff;             //! Pointer to the ITS plane efficicency
308 private:
309   AliITStrackerMI(const AliITStrackerMI &tracker);
310   AliITStrackerMI & operator=(const AliITStrackerMI &tracker);
311   ClassDef(AliITStrackerMI,7)   //ITS tracker MI
312 };
313
314
315
316
317 /////////////////////////////////////////////////////////
318 /////////////////////////////////////////////////////////
319 /////////////////////////////////////////////////////////
320
321
322
323
324
325 inline void AliITStrackerMI::SetupFirstPass(const Int_t *flags,const Double_t *cuts) {
326   // This function sets up flags and cuts for the first tracking pass   
327   //
328   //   flags[0] - vertex constaint flag                                
329   //              negative means "skip the pass"                        
330   //              0        means "no constraint"                        
331   //              positive means "normal constraint"                    
332
333    fConstraint[0]=flags[0];
334    if (cuts==0) return;
335 }
336
337 inline void AliITStrackerMI::SetupSecondPass(const Int_t *flags,const Double_t *cuts) {
338   // This function sets up flags and cuts for the second tracking pass   
339   //
340   //   flags[0] - vertex constaint flag                                
341   //              negative means "skip the pass"                        
342   //              0        means "no constraint"                        
343   //              positive means "normal constraint"                    
344
345    fConstraint[1]=flags[0];
346    if (cuts==0) return;
347 }
348
349 inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
350   //--------------------------------------------------------------------
351   //This function "cooks" a track label. If label<0, this track is fake.
352   //--------------------------------------------------------------------
353    Int_t tpcLabel=t->GetLabel();
354    if (tpcLabel<0) return;
355    AliTracker::CookLabel(t,wrong);
356    if (tpcLabel!=TMath::Abs(t->GetLabel())){
357      t->SetFakeRatio(1.);
358    }
359    if (tpcLabel !=t->GetLabel()) {
360      t->SetLabel(-tpcLabel);      
361    }
362 }
363
364 inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
365 {
366   //--------------------------------------------------------------------
367   //get normalize chi2
368   //--------------------------------------------------------------------
369   track->SetNormChi2(layer,2.*track->GetNSkipped()+0.25*track->GetNDeadZone()+track->GetdEdxMismatch()+track->GetChi2()/
370   //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
371     TMath::Max(double(track->GetNumberOfClusters()-track->GetNSkipped()),
372                1./(1.+track->GetNSkipped())));
373   return track->GetNormChi2(layer);
374 }
375 inline void  AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSRecPoint *cl, Double_t xyz[3]) const
376 {
377   //
378   // get cluster coordinates in global cooordinate 
379   //
380   xyz[2] = cl->GetZ();
381   xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
382   xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;
383 }
384 #endif