]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCCorrection.h
The speed-up of the Poisson relaxation was lost adding it back
[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   virtual void GetDistortionDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
51   virtual void GetDistortionIntegralDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
52   
53   // initialization and update functions
54   virtual void Init();
55   virtual void Update(const TTimeStamp &timeStamp);
56
57   // convenience functions
58   virtual void Print(Option_t* option="") const;
59  
60   TH2F* CreateHistoDRinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
61   TH2F* CreateHistoDRPhiinXY(Float_t z=10.,Int_t nx=100,Int_t nphi=100);
62   TH2F* CreateHistoDZinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
63
64   TH2F* CreateHistoDRinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
65   TH2F* CreateHistoDRPhiinZR(Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
66   TH2F* CreateHistoDZinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
67
68
69   TTree* CreateDistortionTree(Double_t step=5);
70   static void  MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run,  Float_t refX, Int_t type, Int_t integ=1);
71   static void  MakeDistortionMapCosmic(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run,  Float_t refX, Int_t type);
72   static void  MakeDistortionMapSector(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run, Int_t type);
73   // normally called directly in the correction classes which inherit from this class
74   virtual void SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2);
75   AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir,TTreeSRedirector *pcstream);
76   void StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment=0);
77   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);
78   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);
79   static void MakeLaserDistortionTreeOld(TTree* tree, TObjArray *corrArray, Int_t itype);
80   static void MakeLaserDistortionTree(TTree* tree, TObjArray *corrArray, Int_t itype);
81
82   void FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts);
83
84   static void AddVisualCorrection(AliTPCCorrection* corr, Int_t position);
85   static Double_t GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0);
86   static Double_t GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
87   //
88   static Double_t GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0,Double_t delta=5);
89   static Double_t GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5);
90
91   static Double_t GetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
92   //
93   static Double_t GetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0,Double_t delta=5);
94   static Double_t GetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5);
95   
96
97 protected:
98   TH2F* CreateTH2F(const char *name,const char *title,
99                    const char *xlabel,const char *ylabel,const char *zlabel,
100                    Int_t nbinsx,Double_t xlow,Double_t xup,
101                    Int_t nbinsy,Double_t ylow,Double_t yup);
102  
103   static const Double_t fgkTPCZ0;       // nominal gating grid position 
104   static const Double_t fgkIFCRadius;   // Mean Radius of the Inner Field Cage ( 82.43 min,  83.70 max) (cm)
105   static const Double_t fgkOFCRadius;   // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
106   static const Double_t fgkZOffSet;     // Offset from CE: calculate all distortions closer to CE as if at this point
107   static const Double_t fgkCathodeV;    // Cathode Voltage (volts)
108   static const Double_t fgkGG;          // Gating Grid voltage (volts)
109   static const Double_t fgkdvdE;        // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
110   static const Double_t fgkEM;          // charge/mass in [C/kg]
111   static const Double_t fgke0;          // vacuum permittivity [A·s/(V·m)]
112
113
114   enum {kNR=   72  };     // Number of R points in the table for interpolating distortion data
115   enum {kNPhi= 18*10+1};  // Number of Phi points in the table for interpolating distortion data ( plus one extra for 360 == 0 ) 
116   enum {kNZ=   166};      // Number of Z points in the table for interpolating distortion data
117
118   Double_t fgkRList[kNR]; // points in the radial direction (for the lookup table)
119   Double_t fgkPhiList[kNPhi]; // points in the phi direction (for the lookup table)
120   Double_t fgkZList[kNZ]; // points in the z direction (for the lookup table)
121
122   // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
123   Int_t fILow, fJLow, fKLow;          // variable to help in the interpolation 
124   // Double_t versions
125   void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z, 
126                                  const Double_t er[kNZ][kNR], Double_t &erValue );
127   void Interpolate3DEdistortion( const Int_t order, const Double_t r, const Float_t phi, const Double_t z, 
128                                  const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], 
129                                  const Double_t ez[kNZ][kNPhi][kNR],
130                                  Double_t &erValue, Double_t &ephiValue, Double_t &ezValue);
131   // TMatrixD versions (for e.g. Poisson relaxation)
132   Double_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y, 
133                               const Int_t nx,  const Int_t ny, const Double_t xv[], const Double_t yv[], 
134                               const TMatrixD &array );
135   Double_t Interpolate3DTable( const Int_t order, const Double_t x,   const Double_t y,   const Double_t z,
136                               const Int_t  nx,    const Int_t  ny,    const Int_t  nz,
137                               const Double_t xv[], const Double_t yv[], const Double_t zv[],
138                               TMatrixD **arrayofArrays );
139   Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], 
140                         const Int_t order, const Double_t x );
141   void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
142  
143   // TMatrixF versions (smaller size, e.g. for final look up table)
144   Float_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y, 
145                               const Int_t nx,  const Int_t ny, const Double_t xv[], const Double_t yv[], 
146                               const TMatrixF &array );
147   Float_t Interpolate3DTable( const Int_t order, const Double_t x,   const Double_t y,   const Double_t z,
148                               const Int_t  nx,    const Int_t  ny,    const Int_t  nz,
149                               const Double_t xv[], const Double_t yv[], const Double_t zv[],
150                                TMatrixF **arrayofArrays ); 
151   Float_t Interpolate( const Double_t xArray[], const Float_t yArray[], 
152                         const Int_t order, const Double_t x );
153
154   virtual Int_t IsPowerOfTwo ( Int_t i ) const  ;
155
156   
157  
158   // Algorithms to solve the laplace or poisson equation 
159   void PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity, 
160                            TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
161                            const Int_t rows, const Int_t columns, const Int_t iterations,
162                            const Bool_t rocDisplacement = kTRUE);
163
164   void PoissonRelaxation3D( TMatrixD **arrayofArrayV, TMatrixD **arrayofChargeDensities, 
165                             TMatrixD **arrayofEroverEz, TMatrixD **arrayofEPhioverEz, TMatrixD **arrayofEz,
166                             const Int_t rows, const Int_t columns,  const Int_t phislices, 
167                             const Float_t deltaphi, const Int_t iterations, const Int_t summetry,
168                             const Bool_t rocDisplacement = kTRUE); 
169   void   SetIsLocal(Bool_t isLocal){fIsLocal=isLocal;}
170   Bool_t IsLocal() const  { return fIsLocal;}
171 protected:
172   Double_t fT1;         // tensor term of wt - T1
173   Double_t fT2;         // tensor term of wt - T2
174   Bool_t fIsLocal;      // switch to indicate that the distortion is a local vector drphi/dz, dr/dz
175   static TObjArray *fgVisualCorrection;  // array of orrection for visualization
176 private:
177   AliTPCCorrection(const AliTPCCorrection &);               // not implemented
178   AliTPCCorrection &operator=(const AliTPCCorrection &);    // not implemented
179
180   void InitLookUpfulcrums();   // to initialize the grid of the look up table
181
182   ClassDef(AliTPCCorrection,4);
183 };
184
185 #endif