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