ROOT logo
#ifndef ALIITSTRACKERMI_H
#define ALIITSTRACKERMI_H
/* Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

/* $Id$ */

//-------------------------------------------------------------------------
//                          ITS tracker
//     reads AliITSclusterMI clusters and creates AliITStrackMI tracks
//           Origin: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
//           Current support and development: 
//                     Andrea Dainese, andrea.dainese@lnl.infn.it
//-------------------------------------------------------------------------

class TTree;
class TTreeSRedirector;
class AliESDEvent;

class AliITSPlaneEff;
class AliITSChannelStatus;
class AliITSDetTypeRec;
class AliPlaneEff;

#include <TObjArray.h>

#include "AliITStrackMI.h"
#include "AliITSRecPoint.h"
#include "AliTracker.h"
#include "AliRefArray.h"
#include "AliITSPIDResponse.h"

//-------------------------------------------------------------------------
class AliITStrackerMI : public AliTracker {
public:
  AliITStrackerMI();
  AliITStrackerMI(const Char_t *geom);
  virtual ~AliITStrackerMI();
  AliCluster *GetCluster(Int_t index) const;
  virtual Bool_t GetTrackPoint(Int_t index, AliTrackPoint& p) const;
  virtual Bool_t GetTrackPointTrackingError(Int_t index, 
			AliTrackPoint& p, const AliESDtrack *t);
  AliITSRecPoint *GetClusterLayer(Int_t layn, Int_t ncl) const
                        {return fgLayers[layn].GetCluster(ncl);}
  Int_t GetNumberOfClustersLayer(Int_t layn) const 
                        {return fgLayers[layn].GetNumberOfClusters();}
  Int_t LoadClusters(TTree *cf);
  void UnloadClusters();
  void FillClusterArray(TObjArray* array) const;
  Int_t Clusters2Tracks(AliESDEvent *event);
  Int_t PropagateBack(AliESDEvent *event);
  Int_t RefitInward(AliESDEvent *event);
  Bool_t RefitAt(Double_t x, AliITStrackMI *track, 
		 const AliITStrackMI *clusters, Bool_t extra=kFALSE, Bool_t planeeff=kFALSE);
  Bool_t RefitAt(Double_t x, AliITStrackMI *track, 
		 const Int_t *clusters, Bool_t extra=kFALSE, Bool_t planeeff=kFALSE);
  void SetupFirstPass(const Int_t *flags,const Double_t *cuts=0);
  void SetupSecondPass(const Int_t *flags,const Double_t *cuts=0);

  void SetLastLayerToTrackTo(Int_t l=0) {fLastLayerToTrackTo=l;} 
  void UseClusters(const AliKalmanTrack *t, Int_t from=0) const;

  void  GetDCASigma(const AliITStrackMI* track, Float_t & sigmarfi, Float_t &sigmaz);
  Double_t GetPredictedChi2MI(AliITStrackMI* track, const AliITSRecPoint *cluster,Int_t layer);
  Int_t UpdateMI(AliITStrackMI* track, const AliITSRecPoint* cl,Double_t chi2,Int_t layer) const;
  AliPlaneEff *GetPlaneEff() {return (AliPlaneEff*)fPlaneEff;}   // return the pointer to AliPlaneEff
  void SetDetTypeRec(const AliITSDetTypeRec *detTypeRec) {fkDetTypeRec = detTypeRec; ReadBadFromDetTypeRec(); }
  TObjArray* GetTrackHypothesys()  {return &fTrackHypothesys;}
  TObjArray* GetBestHypothesys()   {return &fBestHypothesys;}
  TObjArray* GetOriginal()         {return &fOriginal;}
  TTreeSRedirector *GetDebugStreamer() const {return fDebugStreamer;}
  static Int_t CorrectForTPCtoITSDeadZoneMaterial(AliITStrackMI *t);
  void  SetForceSkippingOfLayer();
  Int_t ForceSkippingOfLayer(Int_t l) const { return fForceSkippingOfLayer[l]; }
  //
  // methods for debugging (RS) >>
  Int_t FindClusterOfTrack(int label, int lr, int* store) const;
  //  Int_t GetPattern(const AliITStrackMI* track, char* patt);
  // methods for debugging (RS) <<
  //
  class AliITSdetector { 
  public:
    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) {}
    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) {}
    ~AliITSdetector() {if(fChipIsBad) delete [] fChipIsBad;}
    inline void GetGlobalXYZ( const AliITSRecPoint *cl, Double_t xyz[3]) const;
    Double_t GetR()   const {return fR;}
    Double_t GetRmisal()   const {return fRmisal;}
    Double_t GetPhi() const {return fPhi;}
    Double_t GetYmin() const {return fYmin;}
    Double_t GetYmax() const {return fYmax;}
    Double_t GetZmin() const {return fZmin;}
    Double_t GetZmax() const {return fZmax;}
    Bool_t   IsBad() const {return fIsBad;}
    Int_t    GetNChips() const {return fNChips;}
    Bool_t   IsChipBad(Int_t iChip) const {return (fChipIsBad ? fChipIsBad[iChip] : kFALSE);}
    void SetRmisal(Double_t rmisal) {fRmisal = rmisal;}
    void SetYmin(Double_t min) {fYmin = min;}
    void SetYmax(Double_t max) {fYmax = max;}
    void SetZmin(Double_t min) {fZmin = min;}
    void SetZmax(Double_t max) {fZmax = max;}
    void SetBad() {fIsBad = kTRUE;}
    void ReadBadDetectorAndChips(Int_t ilayer,Int_t idet,const AliITSDetTypeRec *detTypeRec);
  private:
    AliITSdetector(const AliITSdetector& det);
    AliITSdetector & operator=(const AliITSdetector& det){
      this->~AliITSdetector();new(this) AliITSdetector(det);
      return *this;}
    Double_t fR;    // polar coordinates: r 
    Double_t fRmisal;    // polar coordinates: r, with misalignment 
    Double_t fPhi;  // polar coordinates: phi
    Double_t fSinPhi; // sin of phi;
    Double_t fCosPhi; // cos of phi 
    Double_t fYmin;   //  local y minimal
    Double_t fYmax;   //  local max y
    Double_t fZmin;   //  local z min
    Double_t fZmax;   //  local z max
    Bool_t fIsBad;    // is detector dead or noisy?
    Int_t fNChips;    // number of chips
    Bool_t *fChipIsBad; //[fNChips] is chip dead or noisy? 
  };

  class AliITSlayer {
  public:
    AliITSlayer();
    AliITSlayer(Double_t r, Double_t p, Double_t z, Int_t nl, Int_t nd);
    ~AliITSlayer();
    Int_t InsertCluster(AliITSRecPoint *c);
    void  SortClusters();
    void ResetClusters();
    void ResetWeights();
    void SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin,Double_t ymax);
    const AliITSRecPoint *GetNextCluster(Int_t &ci,Bool_t test=kFALSE);
    void ResetRoad();
    Double_t GetRoad() const {return fRoad;}
    Double_t GetR() const {return fR;}
    Int_t FindClusterIndex(Float_t z) const;
    AliITSRecPoint *GetCluster(Int_t i) const {return i<fN ? fClusters[i]:0;} 
    Float_t  *GetWeight(Int_t i)  {return i<fN ? &fClusterWeight[i]:0;}
    AliITSdetector &GetDetector(Int_t n) const { return fDetectors[n]; }
    Int_t FindDetectorIndex(Double_t phi, Double_t z) const;
    Double_t GetThickness(Double_t y, Double_t z, Double_t &x0) const;
    Int_t InRoad() const ;
    Int_t GetNumberOfClusters() const {return fN;}
    Int_t GetNladders() const {return fNladders;}
    Int_t GetNdetectors() const {return fNdetectors;}
    Int_t GetSkip() const {return fSkip;}
    void  SetSkip(Int_t skip){fSkip=skip;}
    void IncAccepted(){fAccepted++;}
    Int_t GetAccepted() const {return fAccepted;}    
    Int_t GetClusterTracks(Int_t i, Int_t j) const {return fClusterTracks[i][j];}
    void SetClusterTracks(Int_t i, Int_t j, Int_t c) {fClusterTracks[i][j]=c;}
    Int_t FindClusterForLabel(Int_t label, Int_t *store) const; //RS
  protected:
    AliITSlayer(const AliITSlayer& layer);
    AliITSlayer & operator=(const AliITSlayer& layer){
      this->~AliITSlayer();new(this) AliITSlayer(layer);
      return *this;}
    Double_t fR;                // mean radius of this layer
    Double_t fPhiOffset;        // offset of the first detector in Phi
    Int_t fNladders;            // number of ladders
    Double_t fZOffset;          // offset of the first detector in Z
    Int_t fNdetectors;          // detectors/ladder
    AliITSdetector *fDetectors; // array of detectors
    Int_t fN;                   // number of clusters
    AliITSRecPoint *fClusters[AliITSRecoParam::kMaxClusterPerLayer]; // pointers to clusters
    Int_t        fClusterIndex[AliITSRecoParam::kMaxClusterPerLayer]; // pointers to clusters
    Float_t fY[AliITSRecoParam::kMaxClusterPerLayer];                // y position of the clusters      
    Float_t fZ[AliITSRecoParam::kMaxClusterPerLayer];                // z position of the clusters      
    Float_t fYB[2];                                       // ymin and ymax
    //
    AliITSRecPoint *fClusters5[6][AliITSRecoParam::kMaxClusterPerLayer5]; // pointers to clusters -     slice in y
    Int_t        fClusterIndex5[6][AliITSRecoParam::kMaxClusterPerLayer5]; // pointers to clusters -     slice in y    
    Float_t fY5[6][AliITSRecoParam::kMaxClusterPerLayer5];                // y position of the clusters  slice in y    
    Float_t fZ5[6][AliITSRecoParam::kMaxClusterPerLayer5];                // z position of the clusters  slice in y 
    Int_t fN5[6];                                       // number of cluster in slice
    Float_t fDy5;                                       //delta y
    Float_t fBy5[6][2];                                    //slice borders
    //
    AliITSRecPoint *fClusters10[11][AliITSRecoParam::kMaxClusterPerLayer10]; // pointers to clusters -     slice in y
    Int_t        fClusterIndex10[11][AliITSRecoParam::kMaxClusterPerLayer10]; // pointers to clusters -     slice in y    
    Float_t fY10[11][AliITSRecoParam::kMaxClusterPerLayer10];                // y position of the clusters  slice in y    
    Float_t fZ10[11][AliITSRecoParam::kMaxClusterPerLayer10];                // z position of the clusters  slice in y 
    Int_t fN10[11];                                       // number of cluster in slice
    Float_t fDy10;                                        // delta y
    Float_t fBy10[11][2];                                 // slice borders
    //
    AliITSRecPoint *fClusters20[21][AliITSRecoParam::kMaxClusterPerLayer20]; // pointers to clusters -     slice in y
    Int_t        fClusterIndex20[21][AliITSRecoParam::kMaxClusterPerLayer20]; // pointers to clusters -     slice in y    
    Float_t fY20[21][AliITSRecoParam::kMaxClusterPerLayer20];                // y position of the clusters  slice in y    
    Float_t fZ20[21][AliITSRecoParam::kMaxClusterPerLayer20];                // z position of the clusters  slice in y 
    Int_t fN20[21];                                       // number of cluster in slice
    Float_t fDy20;                                        //delta y 
    Float_t fBy20[21][2];                                 //slice borders
    //
    AliITSRecPoint** fClustersCs;                         //clusters table in current slice
    Int_t   *fClusterIndexCs;                             //cluster index in current slice 
    Float_t *fYcs;                                        //y position in current slice
    Float_t *fZcs;                                        //z position in current slice
    Int_t    fNcs;                                        //number of clusters in current slice    
    Int_t fCurrentSlice;                                  //current slice
    //
    Float_t  fClusterWeight[AliITSRecoParam::kMaxClusterPerLayer]; // probabilistic weight of the cluster
    Int_t    fClusterTracks[4][AliITSRecoParam::kMaxClusterPerLayer]; //tracks registered to given cluster
    Float_t fZmin;      //    the
    Float_t fZmax;      //    edges
    Float_t fYmin;      //   of  the
    Float_t fYmax;      //   "window"
    Int_t fI;            // index of the current cluster within the "window"
    Int_t fImax;            // index of the last cluster within the "window"    
    Int_t fSkip;     // indicates possibility to skip cluster
    Int_t fAccepted;     // accept indicator 
    Double_t fRoad;      // road defined by the cluster density
    Double_t fMaxSigmaClY; // maximum cluster error Y (to enlarge road)
    Double_t fMaxSigmaClZ; // maximum cluster error Z (to enlarge road)
    Double_t fNMaxSigmaCl; // number of sigma for road enlargement
  };
  AliITStrackerMI::AliITSlayer    & GetLayer(Int_t layer) const;
  AliITStrackerMI::AliITSdetector & GetDetector(Int_t layer, Int_t n) const {return GetLayer(layer).GetDetector(n); }
  Int_t GetNearestLayer(const Double_t *xr) const;  //get nearest upper layer close to the point xr
  void SetCurrentEsdTrack(Int_t i) {fCurrentEsdTrack=i;}
  void FollowProlongationTree(AliITStrackMI * otrack, Int_t esdindex, Bool_t constrain);
  //
  void   FlagFakes(const TObjArray &itsTracks);
  //
protected:
  Bool_t ComputeRoad(AliITStrackMI* track,Int_t ilayer,Int_t idet,Double_t &zmin,Double_t &zmax,Double_t &ymin,Double_t &ymax) const;
    
  void CookLabel(AliKalmanTrack *t,Float_t wrong) const;
  void CookLabel(AliITStrackMI *t,Float_t wrong) const;
  Double_t GetEffectiveThickness();
  Int_t    GetEffectiveThicknessLbyL(Double_t* xMS, Double_t* x2x0MS);
  void ResetBestTrack() {
     fBestTrack.~AliITStrackMI();
     new(&fBestTrack) AliITStrackMI(fTrackToFollow);
  }
  void ResetTrackToFollow(const AliITStrackMI &t) {
     fTrackToFollow.~AliITStrackMI();
     new(&fTrackToFollow) AliITStrackMI(t);
  }
  void CookdEdx(AliITStrackMI* track);

  Int_t GetParticleId(const AliESDtrack* track) const{
    ULong_t trStatus=track->GetStatus();  
    Bool_t isSA=kTRUE; if(trStatus&AliESDtrack::kTPCin) isSA=kFALSE;
    return fITSPid->GetParticleIdFromdEdxVsP(track->P(),track->GetITSsignal(),isSA);
  }
  Int_t GetParticleId(const AliITStrackV2* track) const{
    if(track->GetESDtrack()) return GetParticleId(track->GetESDtrack());
    return fITSPid->GetParticleIdFromdEdxVsP(track->P(),track->GetdEdx(),kFALSE);
  }

  Double_t GetNormalizedChi2(AliITStrackMI * track, Int_t mode);
  Double_t GetTruncatedChi2(const AliITStrackMI * track, Float_t fac);
  Double_t NormalizedChi2(AliITStrackMI * track, Int_t layer);
  Double_t GetInterpolatedChi2(const AliITStrackMI * forwardtrack,const AliITStrackMI * backtrack);  
  Double_t GetMatchingChi2(const AliITStrackMI * track1,const AliITStrackMI * track2);
  Double_t GetSPDDeadZoneProbability(Double_t zpos, Double_t zerr) const;

  Float_t    *GetWeight(Int_t index);
  void AddTrackHypothesys(AliITStrackMI * track, Int_t esdindex);
  void SortTrackHypothesys(Int_t esdindex, Int_t maxcut, Int_t mode);
  AliITStrackMI * GetBestHypothesys(Int_t esdindex, AliITStrackMI * original, Int_t checkmax); 
  void  GetBestHypothesysMIP(TObjArray &itsTracks); 
  void RegisterClusterTracks(const AliITStrackMI* track, Int_t id);
  void UnRegisterClusterTracks(const AliITStrackMI* track, Int_t id);
  Float_t GetNumberOfSharedClusters(AliITStrackMI* track,Int_t id, Int_t list[6], AliITSRecPoint *clist[6]);
  Int_t GetOverlapTrack(const AliITStrackMI *track, Int_t trackID, Int_t &shared, Int_t clusterlist[6], Int_t overlist[6]);
  AliITStrackMI * GetBest2Tracks(Int_t trackID1, Int_t treackID2, Float_t th0, Float_t th1,AliITStrackMI* original);
  Float_t  * GetErrY(Int_t trackindex) const {return &fCoefficients[trackindex*48];}
  Float_t  * GetErrZ(Int_t trackindex) const {return &fCoefficients[trackindex*48+12];}
  Float_t  * GetNy(Int_t trackindex) const {return &fCoefficients[trackindex*48+24];}
  Float_t  * GetNz(Int_t trackindex) const {return &fCoefficients[trackindex*48+36];}
  void       SignDeltas(const TObjArray *clusterArray, Float_t zv);
  void MakeCoefficients(Int_t ntracks);
  void BuildMaterialLUT(TString material);
  void MakeTrksMaterialLUT(Int_t ntracks);
  void DeleteTrksMaterialLUT();
  Int_t CorrectForPipeMaterial(AliITStrackMI *t, TString direction="inward");
  Int_t CorrectForShieldMaterial(AliITStrackMI *t, TString shield, TString direction="inward");
  Int_t CorrectForLayerMaterial(AliITStrackMI *t, Int_t layerindex, Double_t oldGlobXYZ[3], TString direction="inward");
  void UpdateESDtrack(AliITStrackMI* track, ULong_t flags) const;
  void ReadBadFromDetTypeRec();
  Int_t CheckSkipLayer(const AliITStrackMI *track,Int_t ilayer,Int_t idet) const;
  Int_t CheckDeadZone(AliITStrackMI *track,Int_t ilayer,Int_t idet,Double_t dz,Double_t dy,Bool_t noClusters=kFALSE) const;
  Bool_t LocalModuleCoord(Int_t ilayer,Int_t idet,const AliITStrackMI *track,
			  Float_t &xloc,Float_t &zloc) const;
// method to be used for Plane Efficiency evaluation
  Bool_t IsOKForPlaneEff(const AliITStrackMI* track, const Int_t *clusters, Int_t ilayer); // Check if a track is usable
                                                                                           // for Plane Eff evaluation
  void UseTrackForPlaneEff(const AliITStrackMI* track, Int_t ilayer);                            // Use this track for Plane Eff
// 
  Int_t fI;                              // index of the current layer
  static AliITSlayer fgLayers[AliITSgeomTGeo::kNLayers];// ITS layers
  AliITStrackMI fTracks[AliITSgeomTGeo::kNLayers];      // track estimations at the ITS layers
  AliITStrackMI fBestTrack;              // "best" track 
  AliITStrackMI fTrackToFollow;          // followed track
  TObjArray     fTrackHypothesys;        // ! array with track hypothesys- ARRAY is the owner of tracks- MI
  TObjArray     fBestHypothesys;         // ! array with track hypothesys- ARRAY is the owner of tracks- MI
  TObjArray     fOriginal;               // ! array with seeds from the TPC
  Int_t         fBestTrackIndex[100000]; // ! index of the best track
  Int_t         fCurrentEsdTrack;        // ! current esd track           - MI
  Int_t fPass;                           // current pass through the data 
  Int_t fConstraint[2];                  // constraint flags
  Bool_t fAfterV0;                       //indicates V0 founded
  Int_t fForceSkippingOfLayer[AliITSgeomTGeo::kNLayers]; // layers to be skipped
  Int_t fLastLayerToTrackTo;             // the innermost layer to track to
  Float_t * fCoefficients;               //! working array with errors and mean cluster shape
  AliESDEvent  * fEsd;                   //! pointer to the ESD event
  Double_t fSPDdetzcentre[4];            // centres of SPD modules in z
  TString fTrackingPhase;                // current tracking phase
  Int_t fUseTGeo;                        // use TGeo to get material budget
  Int_t   fNtracks;                      // number of tracks to prolong
  Bool_t  fFlagFakes;                    // request fakes flagging
  Bool_t  fSelectBestMIP03;              // use Chi2MIP[0]*Chi2MIP[3] in hypothesis analysis instead of Chi2MIP[0]
  Bool_t  fUseImproveKalman;             // use Kalman version of Improve
  Float_t fxOverX0Pipe;                  // material budget
  Float_t fxTimesRhoPipe;                // material budget
  Float_t fxOverX0Shield[2];             // material budget
  Float_t fxTimesRhoShield[2];           // material budget
  Float_t fxOverX0Layer[6];              // material budget
  Float_t fxTimesRhoLayer[6];            // material budget
  Float_t *fxOverX0PipeTrks;             //! material budget
  Float_t *fxTimesRhoPipeTrks;           //! material budget
  Float_t *fxOverX0ShieldTrks;           //! material budget
  Float_t *fxTimesRhoShieldTrks;         //! material budget
  Float_t *fxOverX0LayerTrks;            //! material budget
  Float_t *fxTimesRhoLayerTrks;          //! material budget
  TTreeSRedirector *fDebugStreamer;      //!debug streamer
  AliITSChannelStatus *fITSChannelStatus;//! bitmaps with channel status for SPD and SDD
  const AliITSDetTypeRec *fkDetTypeRec;         //! ITS det type rec, from AliITSReconstructor
  AliITSPlaneEff *fPlaneEff;             //! Pointer to the ITS plane efficicency
  Bool_t* fSPDChipIntPlaneEff;      //! Map of the SPD chips already intersected by a track (for FO studies)
  AliITSPIDResponse *fITSPid;            //! parameters for ITS pid 
  //
private:
  AliITStrackerMI(const AliITStrackerMI &tracker);
  AliITStrackerMI & operator=(const AliITStrackerMI &tracker);
  ClassDef(AliITStrackerMI,11)   //ITS tracker MI
};




/////////////////////////////////////////////////////////
/////////////////////////////////////////////////////////
/////////////////////////////////////////////////////////





inline void AliITStrackerMI::SetupFirstPass(const Int_t *flags,const Double_t *cuts) {
  // This function sets up flags and cuts for the first tracking pass   
  //
  //   flags[0] - vertex constaint flag                                
  //              negative means "skip the pass"                        
  //              0        means "no constraint"                        
  //              positive means "normal constraint"                    

   fConstraint[0]=flags[0];
   if (!cuts) return;
}

inline void AliITStrackerMI::SetupSecondPass(const Int_t *flags,const Double_t *cuts) {
  // This function sets up flags and cuts for the second tracking pass   
  //
  //   flags[0] - vertex constaint flag                                
  //              negative means "skip the pass"                        
  //              0        means "no constraint"                        
  //              positive means "normal constraint"                    

   fConstraint[1]=flags[0];
   if (!cuts) return;
}

inline void AliITStrackerMI::CookLabel(AliKalmanTrack *t,Float_t wrong) const {
  //--------------------------------------------------------------------
  //This function "cooks" a track label. If label<0, this track is fake.
  //--------------------------------------------------------------------
   Int_t tpcLabel=t->GetLabel();
   if (tpcLabel<0) return;
   AliTracker::CookLabel(t,wrong);
   if (tpcLabel!=TMath::Abs(t->GetLabel())){
     t->SetFakeRatio(1.);
   }
   if (tpcLabel !=t->GetLabel()) {
     t->SetLabel(-tpcLabel);      
   }
}

inline Double_t AliITStrackerMI::NormalizedChi2(AliITStrackMI * track, Int_t layer)
{
  //--------------------------------------------------------------------
  //get normalize chi2
  //--------------------------------------------------------------------
  track->SetNormChi2(layer,2.*track->GetNSkipped()+0.25*track->GetNDeadZone()+track->GetdEdxMismatch()+track->GetChi2()/
  //track->fNormChi2[layer] = 2.*track->fNSkipped+0.25*track->fNDeadZone+track->fdEdxMismatch+track->fChi22/
    TMath::Max(double(track->GetNumberOfClusters()-track->GetNSkipped()),
	       1./(1.+track->GetNSkipped())));
  return track->GetNormChi2(layer);
}
inline void  AliITStrackerMI::AliITSdetector::GetGlobalXYZ(const AliITSRecPoint *cl, Double_t xyz[3]) const
{
  //
  // get cluster coordinates in global cooordinate 
  //
  xyz[2] = cl->GetZ();
  xyz[0] = fR*fCosPhi - cl->GetY()*fSinPhi;
  xyz[1] = fR*fSinPhi + cl->GetY()*fCosPhi;
}
#endif

 AliITStrackerMI.h:1
 AliITStrackerMI.h:2
 AliITStrackerMI.h:3
 AliITStrackerMI.h:4
 AliITStrackerMI.h:5
 AliITStrackerMI.h:6
 AliITStrackerMI.h:7
 AliITStrackerMI.h:8
 AliITStrackerMI.h:9
 AliITStrackerMI.h:10
 AliITStrackerMI.h:11
 AliITStrackerMI.h:12
 AliITStrackerMI.h:13
 AliITStrackerMI.h:14
 AliITStrackerMI.h:15
 AliITStrackerMI.h:16
 AliITStrackerMI.h:17
 AliITStrackerMI.h:18
 AliITStrackerMI.h:19
 AliITStrackerMI.h:20
 AliITStrackerMI.h:21
 AliITStrackerMI.h:22
 AliITStrackerMI.h:23
 AliITStrackerMI.h:24
 AliITStrackerMI.h:25
 AliITStrackerMI.h:26
 AliITStrackerMI.h:27
 AliITStrackerMI.h:28
 AliITStrackerMI.h:29
 AliITStrackerMI.h:30
 AliITStrackerMI.h:31
 AliITStrackerMI.h:32
 AliITStrackerMI.h:33
 AliITStrackerMI.h:34
 AliITStrackerMI.h:35
 AliITStrackerMI.h:36
 AliITStrackerMI.h:37
 AliITStrackerMI.h:38
 AliITStrackerMI.h:39
 AliITStrackerMI.h:40
 AliITStrackerMI.h:41
 AliITStrackerMI.h:42
 AliITStrackerMI.h:43
 AliITStrackerMI.h:44
 AliITStrackerMI.h:45
 AliITStrackerMI.h:46
 AliITStrackerMI.h:47
 AliITStrackerMI.h:48
 AliITStrackerMI.h:49
 AliITStrackerMI.h:50
 AliITStrackerMI.h:51
 AliITStrackerMI.h:52
 AliITStrackerMI.h:53
 AliITStrackerMI.h:54
 AliITStrackerMI.h:55
 AliITStrackerMI.h:56
 AliITStrackerMI.h:57
 AliITStrackerMI.h:58
 AliITStrackerMI.h:59
 AliITStrackerMI.h:60
 AliITStrackerMI.h:61
 AliITStrackerMI.h:62
 AliITStrackerMI.h:63
 AliITStrackerMI.h:64
 AliITStrackerMI.h:65
 AliITStrackerMI.h:66
 AliITStrackerMI.h:67
 AliITStrackerMI.h:68
 AliITStrackerMI.h:69
 AliITStrackerMI.h:70
 AliITStrackerMI.h:71
 AliITStrackerMI.h:72
 AliITStrackerMI.h:73
 AliITStrackerMI.h:74
 AliITStrackerMI.h:75
 AliITStrackerMI.h:76
 AliITStrackerMI.h:77
 AliITStrackerMI.h:78
 AliITStrackerMI.h:79
 AliITStrackerMI.h:80
 AliITStrackerMI.h:81
 AliITStrackerMI.h:82
 AliITStrackerMI.h:83
 AliITStrackerMI.h:84
 AliITStrackerMI.h:85
 AliITStrackerMI.h:86
 AliITStrackerMI.h:87
 AliITStrackerMI.h:88
 AliITStrackerMI.h:89
 AliITStrackerMI.h:90
 AliITStrackerMI.h:91
 AliITStrackerMI.h:92
 AliITStrackerMI.h:93
 AliITStrackerMI.h:94
 AliITStrackerMI.h:95
 AliITStrackerMI.h:96
 AliITStrackerMI.h:97
 AliITStrackerMI.h:98
 AliITStrackerMI.h:99
 AliITStrackerMI.h:100
 AliITStrackerMI.h:101
 AliITStrackerMI.h:102
 AliITStrackerMI.h:103
 AliITStrackerMI.h:104
 AliITStrackerMI.h:105
 AliITStrackerMI.h:106
 AliITStrackerMI.h:107
 AliITStrackerMI.h:108
 AliITStrackerMI.h:109
 AliITStrackerMI.h:110
 AliITStrackerMI.h:111
 AliITStrackerMI.h:112
 AliITStrackerMI.h:113
 AliITStrackerMI.h:114
 AliITStrackerMI.h:115
 AliITStrackerMI.h:116
 AliITStrackerMI.h:117
 AliITStrackerMI.h:118
 AliITStrackerMI.h:119
 AliITStrackerMI.h:120
 AliITStrackerMI.h:121
 AliITStrackerMI.h:122
 AliITStrackerMI.h:123
 AliITStrackerMI.h:124
 AliITStrackerMI.h:125
 AliITStrackerMI.h:126
 AliITStrackerMI.h:127
 AliITStrackerMI.h:128
 AliITStrackerMI.h:129
 AliITStrackerMI.h:130
 AliITStrackerMI.h:131
 AliITStrackerMI.h:132
 AliITStrackerMI.h:133
 AliITStrackerMI.h:134
 AliITStrackerMI.h:135
 AliITStrackerMI.h:136
 AliITStrackerMI.h:137
 AliITStrackerMI.h:138
 AliITStrackerMI.h:139
 AliITStrackerMI.h:140
 AliITStrackerMI.h:141
 AliITStrackerMI.h:142
 AliITStrackerMI.h:143
 AliITStrackerMI.h:144
 AliITStrackerMI.h:145
 AliITStrackerMI.h:146
 AliITStrackerMI.h:147
 AliITStrackerMI.h:148
 AliITStrackerMI.h:149
 AliITStrackerMI.h:150
 AliITStrackerMI.h:151
 AliITStrackerMI.h:152
 AliITStrackerMI.h:153
 AliITStrackerMI.h:154
 AliITStrackerMI.h:155
 AliITStrackerMI.h:156
 AliITStrackerMI.h:157
 AliITStrackerMI.h:158
 AliITStrackerMI.h:159
 AliITStrackerMI.h:160
 AliITStrackerMI.h:161
 AliITStrackerMI.h:162
 AliITStrackerMI.h:163
 AliITStrackerMI.h:164
 AliITStrackerMI.h:165
 AliITStrackerMI.h:166
 AliITStrackerMI.h:167
 AliITStrackerMI.h:168
 AliITStrackerMI.h:169
 AliITStrackerMI.h:170
 AliITStrackerMI.h:171
 AliITStrackerMI.h:172
 AliITStrackerMI.h:173
 AliITStrackerMI.h:174
 AliITStrackerMI.h:175
 AliITStrackerMI.h:176
 AliITStrackerMI.h:177
 AliITStrackerMI.h:178
 AliITStrackerMI.h:179
 AliITStrackerMI.h:180
 AliITStrackerMI.h:181
 AliITStrackerMI.h:182
 AliITStrackerMI.h:183
 AliITStrackerMI.h:184
 AliITStrackerMI.h:185
 AliITStrackerMI.h:186
 AliITStrackerMI.h:187
 AliITStrackerMI.h:188
 AliITStrackerMI.h:189
 AliITStrackerMI.h:190
 AliITStrackerMI.h:191
 AliITStrackerMI.h:192
 AliITStrackerMI.h:193
 AliITStrackerMI.h:194
 AliITStrackerMI.h:195
 AliITStrackerMI.h:196
 AliITStrackerMI.h:197
 AliITStrackerMI.h:198
 AliITStrackerMI.h:199
 AliITStrackerMI.h:200
 AliITStrackerMI.h:201
 AliITStrackerMI.h:202
 AliITStrackerMI.h:203
 AliITStrackerMI.h:204
 AliITStrackerMI.h:205
 AliITStrackerMI.h:206
 AliITStrackerMI.h:207
 AliITStrackerMI.h:208
 AliITStrackerMI.h:209
 AliITStrackerMI.h:210
 AliITStrackerMI.h:211
 AliITStrackerMI.h:212
 AliITStrackerMI.h:213
 AliITStrackerMI.h:214
 AliITStrackerMI.h:215
 AliITStrackerMI.h:216
 AliITStrackerMI.h:217
 AliITStrackerMI.h:218
 AliITStrackerMI.h:219
 AliITStrackerMI.h:220
 AliITStrackerMI.h:221
 AliITStrackerMI.h:222
 AliITStrackerMI.h:223
 AliITStrackerMI.h:224
 AliITStrackerMI.h:225
 AliITStrackerMI.h:226
 AliITStrackerMI.h:227
 AliITStrackerMI.h:228
 AliITStrackerMI.h:229
 AliITStrackerMI.h:230
 AliITStrackerMI.h:231
 AliITStrackerMI.h:232
 AliITStrackerMI.h:233
 AliITStrackerMI.h:234
 AliITStrackerMI.h:235
 AliITStrackerMI.h:236
 AliITStrackerMI.h:237
 AliITStrackerMI.h:238
 AliITStrackerMI.h:239
 AliITStrackerMI.h:240
 AliITStrackerMI.h:241
 AliITStrackerMI.h:242
 AliITStrackerMI.h:243
 AliITStrackerMI.h:244
 AliITStrackerMI.h:245
 AliITStrackerMI.h:246
 AliITStrackerMI.h:247
 AliITStrackerMI.h:248
 AliITStrackerMI.h:249
 AliITStrackerMI.h:250
 AliITStrackerMI.h:251
 AliITStrackerMI.h:252
 AliITStrackerMI.h:253
 AliITStrackerMI.h:254
 AliITStrackerMI.h:255
 AliITStrackerMI.h:256
 AliITStrackerMI.h:257
 AliITStrackerMI.h:258
 AliITStrackerMI.h:259
 AliITStrackerMI.h:260
 AliITStrackerMI.h:261
 AliITStrackerMI.h:262
 AliITStrackerMI.h:263
 AliITStrackerMI.h:264
 AliITStrackerMI.h:265
 AliITStrackerMI.h:266
 AliITStrackerMI.h:267
 AliITStrackerMI.h:268
 AliITStrackerMI.h:269
 AliITStrackerMI.h:270
 AliITStrackerMI.h:271
 AliITStrackerMI.h:272
 AliITStrackerMI.h:273
 AliITStrackerMI.h:274
 AliITStrackerMI.h:275
 AliITStrackerMI.h:276
 AliITStrackerMI.h:277
 AliITStrackerMI.h:278
 AliITStrackerMI.h:279
 AliITStrackerMI.h:280
 AliITStrackerMI.h:281
 AliITStrackerMI.h:282
 AliITStrackerMI.h:283
 AliITStrackerMI.h:284
 AliITStrackerMI.h:285
 AliITStrackerMI.h:286
 AliITStrackerMI.h:287
 AliITStrackerMI.h:288
 AliITStrackerMI.h:289
 AliITStrackerMI.h:290
 AliITStrackerMI.h:291
 AliITStrackerMI.h:292
 AliITStrackerMI.h:293
 AliITStrackerMI.h:294
 AliITStrackerMI.h:295
 AliITStrackerMI.h:296
 AliITStrackerMI.h:297
 AliITStrackerMI.h:298
 AliITStrackerMI.h:299
 AliITStrackerMI.h:300
 AliITStrackerMI.h:301
 AliITStrackerMI.h:302
 AliITStrackerMI.h:303
 AliITStrackerMI.h:304
 AliITStrackerMI.h:305
 AliITStrackerMI.h:306
 AliITStrackerMI.h:307
 AliITStrackerMI.h:308
 AliITStrackerMI.h:309
 AliITStrackerMI.h:310
 AliITStrackerMI.h:311
 AliITStrackerMI.h:312
 AliITStrackerMI.h:313
 AliITStrackerMI.h:314
 AliITStrackerMI.h:315
 AliITStrackerMI.h:316
 AliITStrackerMI.h:317
 AliITStrackerMI.h:318
 AliITStrackerMI.h:319
 AliITStrackerMI.h:320
 AliITStrackerMI.h:321
 AliITStrackerMI.h:322
 AliITStrackerMI.h:323
 AliITStrackerMI.h:324
 AliITStrackerMI.h:325
 AliITStrackerMI.h:326
 AliITStrackerMI.h:327
 AliITStrackerMI.h:328
 AliITStrackerMI.h:329
 AliITStrackerMI.h:330
 AliITStrackerMI.h:331
 AliITStrackerMI.h:332
 AliITStrackerMI.h:333
 AliITStrackerMI.h:334
 AliITStrackerMI.h:335
 AliITStrackerMI.h:336
 AliITStrackerMI.h:337
 AliITStrackerMI.h:338
 AliITStrackerMI.h:339
 AliITStrackerMI.h:340
 AliITStrackerMI.h:341
 AliITStrackerMI.h:342
 AliITStrackerMI.h:343
 AliITStrackerMI.h:344
 AliITStrackerMI.h:345
 AliITStrackerMI.h:346
 AliITStrackerMI.h:347
 AliITStrackerMI.h:348
 AliITStrackerMI.h:349
 AliITStrackerMI.h:350
 AliITStrackerMI.h:351
 AliITStrackerMI.h:352
 AliITStrackerMI.h:353
 AliITStrackerMI.h:354
 AliITStrackerMI.h:355
 AliITStrackerMI.h:356
 AliITStrackerMI.h:357
 AliITStrackerMI.h:358
 AliITStrackerMI.h:359
 AliITStrackerMI.h:360
 AliITStrackerMI.h:361
 AliITStrackerMI.h:362
 AliITStrackerMI.h:363
 AliITStrackerMI.h:364
 AliITStrackerMI.h:365
 AliITStrackerMI.h:366
 AliITStrackerMI.h:367
 AliITStrackerMI.h:368
 AliITStrackerMI.h:369
 AliITStrackerMI.h:370
 AliITStrackerMI.h:371
 AliITStrackerMI.h:372
 AliITStrackerMI.h:373
 AliITStrackerMI.h:374
 AliITStrackerMI.h:375
 AliITStrackerMI.h:376
 AliITStrackerMI.h:377
 AliITStrackerMI.h:378
 AliITStrackerMI.h:379
 AliITStrackerMI.h:380
 AliITStrackerMI.h:381
 AliITStrackerMI.h:382
 AliITStrackerMI.h:383
 AliITStrackerMI.h:384
 AliITStrackerMI.h:385
 AliITStrackerMI.h:386
 AliITStrackerMI.h:387
 AliITStrackerMI.h:388
 AliITStrackerMI.h:389
 AliITStrackerMI.h:390
 AliITStrackerMI.h:391
 AliITStrackerMI.h:392
 AliITStrackerMI.h:393
 AliITStrackerMI.h:394
 AliITStrackerMI.h:395
 AliITStrackerMI.h:396
 AliITStrackerMI.h:397
 AliITStrackerMI.h:398
 AliITStrackerMI.h:399
 AliITStrackerMI.h:400
 AliITStrackerMI.h:401
 AliITStrackerMI.h:402
 AliITStrackerMI.h:403
 AliITStrackerMI.h:404
 AliITStrackerMI.h:405
 AliITStrackerMI.h:406
 AliITStrackerMI.h:407
 AliITStrackerMI.h:408
 AliITStrackerMI.h:409
 AliITStrackerMI.h:410
 AliITStrackerMI.h:411
 AliITStrackerMI.h:412
 AliITStrackerMI.h:413