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