]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerMI.h
Method GetPattern removed
[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 #include "AliRefArray.h"
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 UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
61
62   void  GetDCASigma(const AliITStrackMI* track, Float_t & sigmarfi, Float_t &sigmaz);
63   Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSRecPoint *cluster,Int_t layer);
64   Int_t UpdateMI(AliITStrackMI* track, const AliITSRecPoint* cl,Double_t chi2,Int_t layer) const;
65   AliPlaneEff *GetPlaneEff() {return (AliPlaneEff*)fPlaneEff;}   // return the pointer to AliPlaneEff
66   void SetDetTypeRec(const AliITSDetTypeRec *detTypeRec) {fkDetTypeRec = detTypeRec; ReadBadFromDetTypeRec(); }
67   TObjArray* GetTrackHypothesys()  {return &fTrackHypothesys;}
68   TObjArray* GetBestHypothesys()   {return &fBestHypothesys;}
69   TObjArray* GetOriginal()         {return &fOriginal;}
70   TTreeSRedirector *GetDebugStreamer() const {return fDebugStreamer;}
71   static Int_t CorrectForTPCtoITSDeadZoneMaterial(AliITStrackMI *t);
72   void  SetForceSkippingOfLayer();
73   Int_t ForceSkippingOfLayer(Int_t l) const { return fForceSkippingOfLayer[l]; }
74   //
75   // methods for debugging (RS) >>
76   Int_t FindClusterOfTrack(int label, int lr, int* store) const;
77   //  Int_t GetPattern(const AliITStrackMI* track, char* patt);
78   // methods for debugging (RS) <<
79   //
80   class AliITSdetector { 
81   public:
82     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) {}
83     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) {}
84     ~AliITSdetector() {if(fChipIsBad) delete [] fChipIsBad;}
85     inline void GetGlobalXYZ( const AliITSRecPoint *cl, Double_t xyz[3]) const;
86     Double_t GetR()   const {return fR;}
87     Double_t GetRmisal()   const {return fRmisal;}
88     Double_t GetPhi() const {return fPhi;}
89     Double_t GetYmin() const {return fYmin;}
90     Double_t GetYmax() const {return fYmax;}
91     Double_t GetZmin() const {return fZmin;}
92     Double_t GetZmax() const {return fZmax;}
93     Bool_t   IsBad() const {return fIsBad;}
94     Int_t    GetNChips() const {return fNChips;}
95     Bool_t   IsChipBad(Int_t iChip) const {return (fChipIsBad ? fChipIsBad[iChip] : kFALSE);}
96     void SetRmisal(Double_t rmisal) {fRmisal = rmisal;}
97     void SetYmin(Double_t min) {fYmin = min;}
98     void SetYmax(Double_t max) {fYmax = max;}
99     void SetZmin(Double_t min) {fZmin = min;}
100     void SetZmax(Double_t max) {fZmax = max;}
101     void SetBad() {fIsBad = kTRUE;}
102     void ReadBadDetectorAndChips(Int_t ilayer,Int_t idet,const AliITSDetTypeRec *detTypeRec);
103   private:
104     AliITSdetector(const AliITSdetector& det);
105     AliITSdetector & operator=(const AliITSdetector& det){
106       this->~AliITSdetector();new(this) AliITSdetector(det);
107       return *this;}
108     Double_t fR;    // polar coordinates: r 
109     Double_t fRmisal;    // polar coordinates: r, with misalignment 
110     Double_t fPhi;  // polar coordinates: phi
111     Double_t fSinPhi; // sin of phi;
112     Double_t fCosPhi; // cos of phi 
113     Double_t fYmin;   //  local y minimal
114     Double_t fYmax;   //  local max y
115     Double_t fZmin;   //  local z min
116     Double_t fZmax;   //  local z max
117     Bool_t fIsBad;    // is detector dead or noisy?
118     Int_t fNChips;    // number of chips
119     Bool_t *fChipIsBad; //[fNChips] is chip dead or noisy? 
120   };
121
122   class AliITSlayer {
123   public:
124     AliITSlayer();
125     AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
126     ~AliITSlayer();
127     Int_t InsertCluster(AliITSRecPoint *c);
128     void  SortClusters();
129     void ResetClusters();
130     void ResetWeights();
131     void SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin,Double_t ymax);
132     const AliITSRecPoint *GetNextCluster(Int_t &ci,Bool_t test=kFALSE);
133     void ResetRoad();
134     Double_t GetRoad() const {return fRoad;}
135     Double_t GetR() const {return fR;}
136     Int_t FindClusterIndex(Float_t z) const;
137     AliITSRecPoint *GetCluster(Int_t i) const {return i<fN ? fClusters[i]:0;} 
138     Float_t  *GetWeight(Int_t i)  {return i<fN ? &fClusterWeight[i]:0;}
139     AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
140     Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
141     Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
142     Int_t InRoad() const ;
143     Int_t GetNumberOfClusters() const {return fN;}
144     Int_t GetNladders() const {return fNladders;}
145     Int_t GetNdetectors() const {return fNdetectors;}
146     Int_t GetSkip() const {return fSkip;}
147     void  SetSkip(Int_t skip){fSkip=skip;}
148     void IncAccepted(){fAccepted++;}
149     Int_t GetAccepted() const {return fAccepted;}    
150     Int_t GetClusterTracks(Int_t i, Int_t j) const {return fClusterTracks[i][j];}
151     void SetClusterTracks(Int_t i, Int_t j, Int_t c) {fClusterTracks[i][j]=c;}
152     Int_t FindClusterForLabel(Int_t label, Int_t *store); //RS
153   protected:
154     AliITSlayer(const AliITSlayer& layer);
155     AliITSlayer & operator=(const AliITSlayer& layer){
156       this->~AliITSlayer();new(this) AliITSlayer(layer);
157       return *this;}
158     Double_t fR;                // mean radius of this layer
159     Double_t fPhiOffset;        // offset of the first detector in Phi
160     Int_t fNladders;            // number of ladders
161     Double_t fZOffset;          // offset of the first detector in Z
162     Int_t fNdetectors;          // detectors/ladder
163     AliITSdetector *fDetectors; // array of detectors
164     Int_t fN;                   // number of clusters
165     AliITSRecPoint *fClusters[AliITSRecoParam::fgkMaxClusterPerLayer]; // pointers to clusters
166     Int_t        fClusterIndex[AliITSRecoParam::fgkMaxClusterPerLayer]; // pointers to clusters
167     Float_t fY[AliITSRecoParam::fgkMaxClusterPerLayer];                // y position of the clusters      
168     Float_t fZ[AliITSRecoParam::fgkMaxClusterPerLayer];                // z position of the clusters      
169     Float_t fYB[2];                                       // ymin and ymax
170     //
171     AliITSRecPoint *fClusters5[6][AliITSRecoParam::fgkMaxClusterPerLayer5]; // pointers to clusters -     slice in y
172     Int_t        fClusterIndex5[6][AliITSRecoParam::fgkMaxClusterPerLayer5]; // pointers to clusters -     slice in y    
173     Float_t fY5[6][AliITSRecoParam::fgkMaxClusterPerLayer5];                // y position of the clusters  slice in y    
174     Float_t fZ5[6][AliITSRecoParam::fgkMaxClusterPerLayer5];                // z position of the clusters  slice in y 
175     Int_t fN5[6];                                       // number of cluster in slice
176     Float_t fDy5;                                       //delta y
177     Float_t fBy5[6][2];                                    //slice borders
178     //
179     AliITSRecPoint *fClusters10[11][AliITSRecoParam::fgkMaxClusterPerLayer10]; // pointers to clusters -     slice in y
180     Int_t        fClusterIndex10[11][AliITSRecoParam::fgkMaxClusterPerLayer10]; // pointers to clusters -     slice in y    
181     Float_t fY10[11][AliITSRecoParam::fgkMaxClusterPerLayer10];                // y position of the clusters  slice in y    
182     Float_t fZ10[11][AliITSRecoParam::fgkMaxClusterPerLayer10];                // z position of the clusters  slice in y 
183     Int_t fN10[11];                                       // number of cluster in slice
184     Float_t fDy10;                                        // delta y
185     Float_t fBy10[11][2];                                 // slice borders
186     //
187     AliITSRecPoint *fClusters20[21][AliITSRecoParam::fgkMaxClusterPerLayer20]; // pointers to clusters -     slice in y
188     Int_t        fClusterIndex20[21][AliITSRecoParam::fgkMaxClusterPerLayer20]; // pointers to clusters -     slice in y    
189     Float_t fY20[21][AliITSRecoParam::fgkMaxClusterPerLayer20];                // y position of the clusters  slice in y    
190     Float_t fZ20[21][AliITSRecoParam::fgkMaxClusterPerLayer20];                // z position of the clusters  slice in y 
191     Int_t fN20[21];                                       // number of cluster in slice
192     Float_t fDy20;                                        //delta y 
193     Float_t fBy20[21][2];                                 //slice borders
194     //
195     AliITSRecPoint** fClustersCs;                         //clusters table in current slice
196     Int_t   *fClusterIndexCs;                             //cluster index in current slice 
197     Float_t *fYcs;                                        //y position in current slice
198     Float_t *fZcs;                                        //z position in current slice
199     Int_t    fNcs;                                        //number of clusters in current slice    
200     Int_t fCurrentSlice;                                  //current slice
201     //
202     Float_t  fClusterWeight[AliITSRecoParam::fgkMaxClusterPerLayer]; // probabilistic weight of the cluster
203     Int_t    fClusterTracks[4][AliITSRecoParam::fgkMaxClusterPerLayer]; //tracks registered to given cluster
204     Float_t fZmin;      //    the
205     Float_t fZmax;      //    edges
206     Float_t fYmin;      //   of  the
207     Float_t fYmax;      //   "window"
208     Int_t fI;            // index of the current cluster within the "window"
209     Int_t fImax;            // index of the last cluster within the "window"    
210     Int_t fSkip;     // indicates possibility to skip cluster
211     Int_t fAccepted;     // accept indicator 
212     Double_t fRoad;      // road defined by the cluster density
213     Double_t fMaxSigmaClY; // maximum cluster error Y (to enlarge road)
214     Double_t fMaxSigmaClZ; // maximum cluster error Z (to enlarge road)
215     Double_t fNMaxSigmaCl; // number of sigma for road enlargement
216   };
217   AliITStrackerMI::AliITSlayer    & GetLayer(Int_t layer) const;
218   AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
219   Int_t GetNearestLayer(const Double_t *xr) const;  //get nearest upper layer close to the point xr
220   void SetCurrentEsdTrack(Int_t i) {fCurrentEsdTrack=i;}
221   void FollowProlongationTree(AliITStrackMI * otrack, Int_t esdindex, Bool_t constrain);
222   //
223   void   FlagFakes(TObjArray &itsTracks);
224   //
225 protected:
226   Bool_t ComputeRoad(AliITStrackMI* track,Int_t ilayer,Int_t idet,Double_t &zmin,Double_t &zmax,Double_t &ymin,Double_t &ymax) const;
227     
228   void CookLabel(AliKalmanTrack *t,Float_t wrong) const;
229   void CookLabel(AliITStrackMI *t,Float_t wrong) const;
230   Double_t GetEffectiveThickness();
231   Int_t    GetEffectiveThicknessLbyL(Double_t* xMS, Double_t* x2x0MS);
232   void ResetBestTrack() {
233      fBestTrack.~AliITStrackMI();
234      new(&fBestTrack) AliITStrackMI(fTrackToFollow);
235   }
236   void ResetTrackToFollow(const AliITStrackMI &t) {
237      fTrackToFollow.~AliITStrackMI();
238      new(&fTrackToFollow) AliITStrackMI(t);
239   }
240   void CookdEdx(AliITStrackMI* track);
241   Double_t GetNormalizedChi2(AliITStrackMI * track, Int_t mode);
242   Double_t GetTruncatedChi2(const AliITStrackMI * track, Float_t fac);
243   Double_t NormalizedChi2(AliITStrackMI * track, Int_t layer);
244   Double_t GetInterpolatedChi2(const AliITStrackMI * forwardtrack,const AliITStrackMI * backtrack);  
245   Double_t GetMatchingChi2(const AliITStrackMI * track1,const AliITStrackMI * track2);
246   Double_t GetSPDDeadZoneProbability(Double_t zpos, Double_t zerr) const;
247
248   Float_t    *GetWeight(Int_t index);
249   void AddTrackHypothesys(AliITStrackMI * track, Int_t esdindex);
250   void SortTrackHypothesys(Int_t esdindex, Int_t maxcut, Int_t mode);
251   AliITStrackMI * GetBestHypothesys(Int_t esdindex, AliITStrackMI * original, Int_t checkmax); 
252   void  GetBestHypothesysMIP(TObjArray &itsTracks); 
253   void RegisterClusterTracks(const AliITStrackMI* track, Int_t id);
254   void UnRegisterClusterTracks(const AliITStrackMI* track, Int_t id);
255   Float_t GetNumberOfSharedClusters(AliITStrackMI* track,Int_t id, Int_t list[6], AliITSRecPoint *clist[6]);
256   Int_t GetOverlapTrack(const AliITStrackMI *track, Int_t trackID, Int_t &shared, Int_t clusterlist[6], Int_t overlist[6]);
257   AliITStrackMI * GetBest2Tracks(Int_t trackID1, Int_t treackID2, Float_t th0, Float_t th1);
258   Float_t  * GetErrY(Int_t trackindex) const {return &fCoefficients[trackindex*48];}
259   Float_t  * GetErrZ(Int_t trackindex) const {return &fCoefficients[trackindex*48+12];}
260   Float_t  * GetNy(Int_t trackindex) const {return &fCoefficients[trackindex*48+24];}
261   Float_t  * GetNz(Int_t trackindex) const {return &fCoefficients[trackindex*48+36];}
262   void       SignDeltas(const TObjArray *clusterArray, Float_t zv);
263   void MakeCoefficients(Int_t ntracks);
264   void BuildMaterialLUT(TString material);
265   void MakeTrksMaterialLUT(Int_t ntracks);
266   void DeleteTrksMaterialLUT();
267   Int_t CorrectForPipeMaterial(AliITStrackMI *t, TString direction="inward");
268   Int_t CorrectForShieldMaterial(AliITStrackMI *t, TString shield, TString direction="inward");
269   Int_t CorrectForLayerMaterial(AliITStrackMI *t, Int_t layerindex, Double_t oldGlobXYZ[3], TString direction="inward");
270   void UpdateESDtrack(AliITStrackMI* track, ULong_t flags) const;
271   void ReadBadFromDetTypeRec();
272   Int_t CheckSkipLayer(const AliITStrackMI *track,Int_t ilayer,Int_t idet) const;
273   Int_t CheckDeadZone(AliITStrackMI *track,Int_t ilayer,Int_t idet,Double_t dz,Double_t dy,Bool_t noClusters=kFALSE) const;
274   Bool_t LocalModuleCoord(Int_t ilayer,Int_t idet,const AliITStrackMI *track,
275                           Float_t &xloc,Float_t &zloc) const;
276 // method to be used for Plane Efficiency evaluation
277   Bool_t IsOKForPlaneEff(const AliITStrackMI* track, const Int_t *clusters, Int_t ilayer) const; // Check if a track is usable
278                                                                                            // for Plane Eff evaluation
279   void UseTrackForPlaneEff(const AliITStrackMI* track, Int_t ilayer);                            // Use this track for Plane Eff
280 // 
281   Int_t fI;                              // index of the current layer
282   static AliITSlayer fgLayers[AliITSgeomTGeo::kNLayers];// ITS layers
283   AliITStrackMI fTracks[AliITSgeomTGeo::kNLayers];      // track estimations at the ITS layers
284   AliITStrackMI fBestTrack;              // "best" track 
285   AliITStrackMI fTrackToFollow;          // followed track
286   TObjArray     fTrackHypothesys;        // ! array with track hypothesys- ARRAY is the owner of tracks- MI
287   TObjArray     fBestHypothesys;         // ! array with track hypothesys- ARRAY is the owner of tracks- MI
288   TObjArray     fOriginal;               // ! array with seeds from the TPC
289   Int_t         fBestTrackIndex[100000]; // ! index of the best track
290   Int_t         fCurrentEsdTrack;        // ! current esd track           - MI
291   Int_t fPass;                           // current pass through the data 
292   Int_t fConstraint[2];                  // constraint flags
293   Bool_t fAfterV0;                       //indicates V0 founded
294   Int_t fForceSkippingOfLayer[AliITSgeomTGeo::kNLayers]; // layers to be skipped
295   Int_t fLastLayerToTrackTo;             // the innermost layer to track to
296   Float_t * fCoefficients;               //! working array with errors and mean cluster shape
297   AliESDEvent  * fEsd;                   //! pointer to the ESD event
298   Double_t fSPDdetzcentre[4];            // centres of SPD modules in z
299   TString fTrackingPhase;                // current tracking phase
300   Int_t fUseTGeo;                        // use TGeo to get material budget
301   Int_t   fNtracks;                      // number of tracks to prolong
302   Bool_t  fFlagFakes;                    // request fakes flagging
303   Bool_t  fSelectBestMIP03;              // use Chi2MIP[0]*Chi2MIP[3] in hypothesis analysis instead of Chi2MIP[0]
304   Bool_t  fUseImproveKalman;             // use Kalman version of Improve
305   Float_t fxOverX0Pipe;                  // material budget
306   Float_t fxTimesRhoPipe;                // material budget
307   Float_t fxOverX0Shield[2];             // material budget
308   Float_t fxTimesRhoShield[2];           // material budget
309   Float_t fxOverX0Layer[6];              // material budget
310   Float_t fxTimesRhoLayer[6];            // material budget
311   Float_t *fxOverX0PipeTrks;             //! material budget
312   Float_t *fxTimesRhoPipeTrks;           //! material budget
313   Float_t *fxOverX0ShieldTrks;           //! material budget
314   Float_t *fxTimesRhoShieldTrks;         //! material budget
315   Float_t *fxOverX0LayerTrks;            //! material budget
316   Float_t *fxTimesRhoLayerTrks;          //! material budget
317   TTreeSRedirector *fDebugStreamer;      //!debug streamer
318   AliITSChannelStatus *fITSChannelStatus;//! bitmaps with channel status for SPD and SDD
319   const AliITSDetTypeRec *fkDetTypeRec;         //! ITS det type rec, from AliITSReconstructor
320   AliITSPlaneEff *fPlaneEff;             //! Pointer to the ITS plane efficicency
321   //
322 private:
323   AliITStrackerMI(const AliITStrackerMI &tracker);
324   AliITStrackerMI & operator=(const AliITStrackerMI &tracker);
325   ClassDef(AliITStrackerMI,9)   //ITS tracker MI
326 };
327
328
329
330
331 /////////////////////////////////////////////////////////
332 /////////////////////////////////////////////////////////
333 /////////////////////////////////////////////////////////
334
335
336
337
338
339 inline void AliITStrackerMI::SetupFirstPass(const Int_t *flags,const Double_t *cuts) {
340   // This function sets up flags and cuts for the first tracking pass   
341   //
342   //   flags[0] - vertex constaint flag                                
343   //              negative means "skip the pass"                        
344   //              0        means "no constraint"                        
345   //              positive means "normal constraint"                    
346
347    fConstraint[0]=flags[0];
348    if (!cuts) return;
349 }
350
351 inline void AliITStrackerMI::SetupSecondPass(const Int_t *flags,const Double_t *cuts) {
352   // This function sets up flags and cuts for the second tracking pass   
353   //
354   //   flags[0] - vertex constaint flag                                
355   //              negative means "skip the pass"                        
356   //              0        means "no constraint"                        
357   //              positive means "normal constraint"                    
358
359    fConstraint[1]=flags[0];
360    if (!cuts) return;
361 }
362
363 inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
364   //--------------------------------------------------------------------
365   //This function "cooks" a track label. If label<0, this track is fake.
366   //--------------------------------------------------------------------
367    Int_t tpcLabel=t->GetLabel();
368    if (tpcLabel<0) return;
369    AliTracker::CookLabel(t,wrong);
370    if (tpcLabel!=TMath::Abs(t->GetLabel())){
371      t->SetFakeRatio(1.);
372    }
373    if (tpcLabel !=t->GetLabel()) {
374      t->SetLabel(-tpcLabel);      
375    }
376 }
377
378 inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
379 {
380   //--------------------------------------------------------------------
381   //get normalize chi2
382   //--------------------------------------------------------------------
383   track->SetNormChi2(layer,2.*track->GetNSkipped()+0.25*track->GetNDeadZone()+track->GetdEdxMismatch()+track->GetChi2()/
384   //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
385     TMath::Max(double(track->GetNumberOfClusters()-track->GetNSkipped()),
386                1./(1.+track->GetNSkipped())));
387   return track->GetNormChi2(layer);
388 }
389 inline void  AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSRecPoint *cl, Double_t xyz[3]) const
390 {
391   //
392   // get cluster coordinates in global cooordinate 
393   //
394   xyz[2] = cl->GetZ();
395   xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
396   xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;
397 }
398 #endif