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