]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCCorrection.h
d1d55dd92f45b6ce3a47c647af488e28a09fc996
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCCorrection.h
1 #ifndef ALI_TPC_CORRECTION_H
2 #define ALI_TPC_CORRECTION_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7
8 ////////////////////////////////////////////////////////////////////////////////
9 // AliTPCCorrection class                                                     //
10 ////////////////////////////////////////////////////////////////////////////////
11
12
13 #include <TNamed.h>
14 #include "TMatrixD.h"
15 #include "TMatrixF.h"
16 class TH2F;
17 class TTimeStamp;
18 class TCollection;
19 class TTreeSRedirector;
20 class AliExternalTrackParam;
21 class TTree;
22 class THnSparse;
23 class AliESDVertex;
24
25
26 class AliTPCCorrection : public TNamed {
27 public:
28   enum CompositionType {kParallel,kQueue};
29
30   AliTPCCorrection();
31   AliTPCCorrection(const char *name,const char *title);
32   virtual ~AliTPCCorrection();
33   
34
35   // functions to correct a space point
36           void CorrectPoint (      Float_t x[],const Short_t roc);
37           void CorrectPointLocal(Float_t x[],const Short_t roc);
38           void CorrectPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
39           virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
40
41   virtual void GetCorrectionDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
42   virtual void GetCorrectionIntegralDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
43
44   // functions to distort a space point
45           void DistortPoint (      Float_t x[],const Short_t roc);
46           void DistortPointLocal(Float_t x[],const Short_t roc);
47           void DistortPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
48   virtual void GetDistortion(const Float_t x[],const Short_t roc,Float_t dx[]);
49
50   // initialization and update functions
51   virtual void Init();
52   virtual void Update(const TTimeStamp &timeStamp);
53
54   // convenience functions
55   virtual void Print(Option_t* option="") const;
56  
57   TH2F* CreateHistoDRinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
58   TH2F* CreateHistoDRPhiinXY(Float_t z=10.,Int_t nx=100,Int_t nphi=100);
59   TH2F* CreateHistoDZinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
60
61   TH2F* CreateHistoDRinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
62   TH2F* CreateHistoDRPhiinZR(Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
63   TH2F* CreateHistoDZinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
64
65
66   TTree* CreateDistortionTree(Double_t step=5);
67   static void  MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run,  Float_t refX, Int_t type, Int_t integ=1);
68   static void  MakeDistortionMapCosmic(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run,  Float_t refX, Int_t type);
69   static void  MakeDistortionMapSector(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run, Int_t type);
70   // normally called directly in the correction classes which inherit from this class
71   virtual void SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2);
72   AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir,TTreeSRedirector *pcstream);
73   void StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment=0);
74   static void MakeTrackDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step=1, Int_t offset=0, Bool_t debug=0);
75   static void MakeSectorDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step=1, Int_t offset=0, Bool_t debug=0);
76   static void MakeLaserDistortionTreeOld(TTree* tree, TObjArray *corrArray, Int_t itype);
77   static void MakeLaserDistortionTree(TTree* tree, TObjArray *corrArray, Int_t itype);
78
79   void FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts);
80
81   static void AddVisualCorrection(AliTPCCorrection* corr, Int_t position);
82   static Double_t GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0);
83   static Double_t GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
84   //
85   static Double_t GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0,Double_t delta=5);
86   static Double_t GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5);
87   
88
89 protected:
90   TH2F* CreateTH2F(const char *name,const char *title,
91                    const char *xlabel,const char *ylabel,const char *zlabel,
92                    Int_t nbinsx,Double_t xlow,Double_t xup,
93                    Int_t nbinsy,Double_t ylow,Double_t yup);
94  
95   static const Double_t fgkTPCZ0;       // nominal gating grid position 
96   static const Double_t fgkIFCRadius;   // Mean Radius of the Inner Field Cage ( 82.43 min,  83.70 max) (cm)
97   static const Double_t fgkOFCRadius;   // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
98   static const Double_t fgkZOffSet;     // Offset from CE: calculate all distortions closer to CE as if at this point
99   static const Double_t fgkCathodeV;    // Cathode Voltage (volts)
100   static const Double_t fgkGG;          // Gating Grid voltage (volts)
101   static const Double_t fgkdvdE;        // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
102   static const Double_t fgkEM;          // charge/mass in [C/kg]
103   static const Double_t fgke0;          // vacuum permittivity [A·s/(V·m)]
104
105
106   enum {kNR=   72  };     // Number of R points in the table for interpolating distortion data
107   enum {kNPhi= 18*10+1};  // Number of Phi points in the table for interpolating distortion data ( plus one extra for 360 == 0 ) 
108   enum {kNZ=   166};      // Number of Z points in the table for interpolating distortion data
109
110   Double_t fgkRList[kNR]; // points in the radial direction (for the lookup table)
111   Double_t fgkPhiList[kNPhi]; // points in the phi direction (for the lookup table)
112   Double_t fgkZList[kNZ]; // points in the z direction (for the lookup table)
113
114   // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
115   Int_t fILow, fJLow, fKLow;          // variable to help in the interpolation 
116   // Double_t versions
117   void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z, 
118                                  const Double_t er[kNZ][kNR], Double_t &erValue );
119   void Interpolate3DEdistortion( const Int_t order, const Double_t r, const Float_t phi, const Double_t z, 
120                                  const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], 
121                                  const Double_t ez[kNZ][kNPhi][kNR],
122                                  Double_t &erValue, Double_t &ephiValue, Double_t &ezValue);
123   // TMatrixD versions (for e.g. Poisson relaxation)
124   Double_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y, 
125                               const Int_t nx,  const Int_t ny, const Double_t xv[], const Double_t yv[], 
126                               const TMatrixD &array );
127   Double_t Interpolate3DTable( const Int_t order, const Double_t x,   const Double_t y,   const Double_t z,
128                               const Int_t  nx,    const Int_t  ny,    const Int_t  nz,
129                               const Double_t xv[], const Double_t yv[], const Double_t zv[],
130                               TMatrixD **arrayofArrays );
131   Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], 
132                         const Int_t order, const Double_t x );
133   void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
134  
135   // TMatrixF versions (smaller size, e.g. for final look up table)
136   Float_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y, 
137                               const Int_t nx,  const Int_t ny, const Double_t xv[], const Double_t yv[], 
138                               const TMatrixF &array );
139   Float_t Interpolate3DTable( const Int_t order, const Double_t x,   const Double_t y,   const Double_t z,
140                               const Int_t  nx,    const Int_t  ny,    const Int_t  nz,
141                               const Double_t xv[], const Double_t yv[], const Double_t zv[],
142                                TMatrixF **arrayofArrays ); 
143   Float_t Interpolate( const Double_t xArray[], const Float_t yArray[], 
144                         const Int_t order, const Double_t x );
145
146   virtual Int_t IsPowerOfTwo ( Int_t i ) const  ;
147
148   
149  
150   // Algorithms to solve the laplace or poisson equation 
151   void PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity, 
152                            TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
153                            const Int_t rows, const Int_t columns, const Int_t iterations,
154                            const Bool_t rocDisplacement = kTRUE);
155
156   void PoissonRelaxation3D( TMatrixD **arrayofArrayV, TMatrixD **arrayofChargeDensities, 
157                             TMatrixD **arrayofEroverEz, TMatrixD **arrayofEPhioverEz, TMatrixD **arrayofEz,
158                             const Int_t rows, const Int_t columns,  const Int_t phislices, 
159                             const Float_t deltaphi, const Int_t iterations, const Int_t summetry,
160                             const Bool_t rocDisplacement = kTRUE); 
161     
162 protected:
163   Double_t fT1;         // tensor term of wt - T1
164   Double_t fT2;         // tensor term of wt - T2
165   static TObjArray *fgVisualCorrection;  // array of orrection for visualization
166 private:
167   AliTPCCorrection(const AliTPCCorrection &);               // not implemented
168   AliTPCCorrection &operator=(const AliTPCCorrection &);    // not implemented
169
170   void InitLookUpfulcrums();   // to initialize the grid of the look up table
171
172   ClassDef(AliTPCCorrection,4);
173 };
174
175 #endif