]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/Base/AliTPCCorrection.h
skip clusters when the roc voltage is funny
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCCorrection.h
CommitLineData
0116859c 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
b4caed64 7
0116859c 8////////////////////////////////////////////////////////////////////////////////
0116859c 9// AliTPCCorrection class //
0116859c 10////////////////////////////////////////////////////////////////////////////////
11
12
13#include <TNamed.h>
1b923461 14#include "TMatrixD.h"
2bf29b72 15#include "TMatrixF.h"
0116859c 16class TH2F;
e527a1b9 17class TTimeStamp;
0116859c 18class TCollection;
cf5b0aa0 19class TTreeSRedirector;
20class AliExternalTrackParam;
ffab0c37 21class TTree;
8b63d99c 22class THnSparse;
ca58ed4e 23class AliESDVertex;
1b923461 24
25
0116859c 26class AliTPCCorrection : public TNamed {
27public:
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);
46e89793 37 void CorrectPointLocal(Float_t x[],const Short_t roc);
0116859c 38 void CorrectPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
fdbbc146 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);
05da1b4e 43
0116859c 44 // functions to distort a space point
45 void DistortPoint ( Float_t x[],const Short_t roc);
46e89793 46 void DistortPointLocal(Float_t x[],const Short_t roc);
0116859c 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[]);
2d4e971f 49
50 virtual void GetDistortionDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
05da1b4e 51 virtual void GetDistortionIntegralDz(const Float_t x[],const Short_t roc,Float_t dx[], Float_t delta);
52
e527a1b9 53 // initialization and update functions
0116859c 54 virtual void Init();
e527a1b9 55 virtual void Update(const TTimeStamp &timeStamp);
0116859c 56
862220e2 57 // map scaling
58 virtual void SetCorrScaleFactor(Float_t /*val*/) { ; }
59 virtual Float_t GetCorrScaleFactor() const { return 1.; }
60
0116859c 61 // convenience functions
62 virtual void Print(Option_t* option="") const;
63
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);
c9cbd2f2 66 TH2F* CreateHistoDZinXY (Float_t z=10.,Int_t nx=100,Int_t ny=100);
67
0116859c 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);
c9cbd2f2 70 TH2F* CreateHistoDZinZR (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
71
72
ffab0c37 73 TTree* CreateDistortionTree(Double_t step=5);
97d17739 74 static void MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run, Float_t refX, Int_t type, Int_t integ=1);
cfe2c39a 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);
0116859c 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);
be67055b 79 AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir,TTreeSRedirector *pcstream);
ffab0c37 80 void StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment=0);
46e89793 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);
7f4cb119 84 static void MakeLaserDistortionTree(TTree* tree, TObjArray *corrArray, Int_t itype);
f1817479 85
86 void FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts);
87
88 static void AddVisualCorrection(AliTPCCorrection* corr, Int_t position);
287fbdfa 89 static Double_t GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0);
f1817479 90 static Double_t GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
fdbbc146 91 //
92 static Double_t GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0,Double_t delta=5);
93 static Double_t GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5);
2d4e971f 94
95 static Double_t GetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
96 //
97 static Double_t GetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0,Double_t delta=5);
98 static Double_t GetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5);
fdbbc146 99
ca58ed4e 100
0116859c 101protected:
102 TH2F* CreateTH2F(const char *name,const char *title,
103 const char *xlabel,const char *ylabel,const char *zlabel,
104 Int_t nbinsx,Double_t xlow,Double_t xup,
105 Int_t nbinsy,Double_t ylow,Double_t yup);
106
c9cbd2f2 107 static const Double_t fgkTPCZ0; // nominal gating grid position
0116859c 108 static const Double_t fgkIFCRadius; // Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)
109 static const Double_t fgkOFCRadius; // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
110 static const Double_t fgkZOffSet; // Offset from CE: calculate all distortions closer to CE as if at this point
111 static const Double_t fgkCathodeV; // Cathode Voltage (volts)
112 static const Double_t fgkGG; // Gating Grid voltage (volts)
c9cbd2f2 113 static const Double_t fgkdvdE; // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
114 static const Double_t fgkEM; // charge/mass in [C/kg]
115 static const Double_t fgke0; // vacuum permittivity [A·s/(V·m)]
0116859c 116
c9cbd2f2 117
35ae345f 118 enum {kNR= 72 }; // Number of R points in the table for interpolating distortion data
119 enum {kNPhi= 18*10+1}; // Number of Phi points in the table for interpolating distortion data ( plus one extra for 360 == 0 )
120 enum {kNZ= 166}; // Number of Z points in the table for interpolating distortion data
c9cbd2f2 121
35ae345f 122 Double_t fgkRList[kNR]; // points in the radial direction (for the lookup table)
c9cbd2f2 123 Double_t fgkPhiList[kNPhi]; // points in the phi direction (for the lookup table)
35ae345f 124 Double_t fgkZList[kNZ]; // points in the z direction (for the lookup table)
0116859c 125
126 // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
c9cbd2f2 127 Int_t fILow, fJLow, fKLow; // variable to help in the interpolation
2bf29b72 128 // Double_t versions
0116859c 129 void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z,
b1f0a2a5 130 const Double_t er[kNZ][kNR], Double_t &erValue );
c9cbd2f2 131 void Interpolate3DEdistortion( const Int_t order, const Double_t r, const Float_t phi, const Double_t z,
35ae345f 132 const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR],
133 const Double_t ez[kNZ][kNPhi][kNR],
c9cbd2f2 134 Double_t &erValue, Double_t &ephiValue, Double_t &ezValue);
2bf29b72 135 // TMatrixD versions (for e.g. Poisson relaxation)
c9cbd2f2 136 Double_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 TMatrixD &array );
139 Double_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 TMatrixD **arrayofArrays );
0116859c 143 Double_t Interpolate( const Double_t xArray[], const Double_t yArray[],
144 const Int_t order, const Double_t x );
145 void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
2bf29b72 146
147 // TMatrixF versions (smaller size, e.g. for final look up table)
148 Float_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y,
149 const Int_t nx, const Int_t ny, const Double_t xv[], const Double_t yv[],
150 const TMatrixF &array );
151 Float_t Interpolate3DTable( const Int_t order, const Double_t x, const Double_t y, const Double_t z,
152 const Int_t nx, const Int_t ny, const Int_t nz,
153 const Double_t xv[], const Double_t yv[], const Double_t zv[],
154 TMatrixF **arrayofArrays );
155 Float_t Interpolate( const Double_t xArray[], const Float_t yArray[],
156 const Int_t order, const Double_t x );
157
710bda39 158 virtual Int_t IsPowerOfTwo ( Int_t i ) const ;
1b923461 159
c9cbd2f2 160
161
162 // Algorithms to solve the laplace or poisson equation
163 void PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity,
164 TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
165 const Int_t rows, const Int_t columns, const Int_t iterations,
166 const Bool_t rocDisplacement = kTRUE);
167
168 void PoissonRelaxation3D( TMatrixD **arrayofArrayV, TMatrixD **arrayofChargeDensities,
169 TMatrixD **arrayofEroverEz, TMatrixD **arrayofEPhioverEz, TMatrixD **arrayofEz,
170 const Int_t rows, const Int_t columns, const Int_t phislices,
171 const Float_t deltaphi, const Int_t iterations, const Int_t summetry,
172 const Bool_t rocDisplacement = kTRUE);
d9ef0909 173 void SetIsLocal(Bool_t isLocal){fIsLocal=isLocal;}
174 Bool_t IsLocal() const { return fIsLocal;}
534fd34a 175protected:
176 Double_t fT1; // tensor term of wt - T1
177 Double_t fT2; // tensor term of wt - T2
d9ef0909 178 Bool_t fIsLocal; // switch to indicate that the distortion is a local vector drphi/dz, dr/dz
f1817479 179 static TObjArray *fgVisualCorrection; // array of orrection for visualization
7f4cb119 180private:
c9cbd2f2 181 AliTPCCorrection(const AliTPCCorrection &); // not implemented
182 AliTPCCorrection &operator=(const AliTPCCorrection &); // not implemented
183
35ae345f 184 void InitLookUpfulcrums(); // to initialize the grid of the look up table
185
c3fea858 186 ClassDef(AliTPCCorrection,5);
0116859c 187};
188
189#endif