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