]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCFCVoltError3D.h
Remove obsolete drawing methods and macros
[u/mrichter/AliRoot.git] / TPC / AliTPCFCVoltError3D.h
CommitLineData
c9cbd2f2 1#ifndef ALITPCFCVOLTERROR3D_H
2#define ALITPCFCVOLTERROR3D_H
3
4/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
6
7////////////////////////////////////////////////////////////////////////////
c9cbd2f2 8// AliTPCFCVoltError3D class //
c9cbd2f2 9// Authors: Jim Thomas, Stefan Rossegger //
10////////////////////////////////////////////////////////////////////////////
11
12#include "AliTPCCorrection.h"
13
14
15class AliTPCFCVoltError3D : public AliTPCCorrection {
16public:
17 AliTPCFCVoltError3D();
18 virtual ~AliTPCFCVoltError3D();
19
20 // initialization and update functions
21 virtual void Init();
22 virtual void Update(const TTimeStamp &timeStamp);
23
24 // common setters and getters for tangled ExB effect
25 virtual void SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2) {
26 fT1=t1; fT2=t2;
27 const Double_t wt0=t2*omegaTau; fC0=1./(1.+wt0*wt0);
28 const Double_t wt1=t1*omegaTau; fC1=wt1/(1.+wt1*wt1);
29 };
30 void SetC0C1(Float_t c0,Float_t c1) {fC0=c0;fC1=c1;} // CAUTION: USE WITH CARE
31 Float_t GetC0() const {return fC0;}
32 Float_t GetC1() const {return fC1;}
33
34 // setters and getters
35
36 // Set rod shift in Voltage equivalents (40V ~ 1mm)
37 // rod numbers: 0-17 (IFC), 18-35 (OFC)
38 // note: strips move accordingly
f2efa2c7 39 void SetRodVoltShiftA(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fRodVoltShiftA[rod]=voltOffset; fInitLookUp=doInit;}
40 void SetRodVoltShiftC(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fRodVoltShiftC[rod]=voltOffset; fInitLookUp=doInit;}
c9cbd2f2 41 Float_t GetRodVoltShiftA(Int_t i) const {return fRodVoltShiftA[i];}// 0-17: IFC, 18-35; OFC
42 Float_t GetRodVoltShiftC(Int_t i) const {return fRodVoltShiftC[i];}// 0-17: IFC, 18-35; OFC
43
44 // Set rotated clip (just at High Voltage RODs) in Voltage equivalents (40V ~ 1mm)
45 // rod number: 0 (IFC), 1 (OFC)
f2efa2c7 46 void SetRotatedClipVoltA(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fRotatedClipVoltA[rod]=voltOffset; fInitLookUp=doInit;}
47 void SetRotatedClipVoltC(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fRotatedClipVoltC[rod]=voltOffset; fInitLookUp=doInit;}
c9cbd2f2 48 Float_t GetRotatedClipVoltA(Int_t i) const {return fRotatedClipVoltA[i];}// (0,1):(IFC,OFC)
49 Float_t GetRotatedClipVoltC(Int_t i) const {return fRotatedClipVoltC[i];}// (0,1):(IFC,OFC)
50
51 // Set rod shift in Voltage equivalents (40V ~ 1mm)
52 // rod numbers: 0-17 (OFC)
53 // note: strips DO NOT move, only the copper rods do ...
f2efa2c7 54 void SetCopperRodShiftA(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fCopperRodShiftA[rod]=voltOffset; fInitLookUp=doInit;}
55 void SetCopperRodShiftC(Int_t rod, Float_t voltOffset, Bool_t doInit=kTRUE) {fCopperRodShiftC[rod]=voltOffset; fInitLookUp=doInit;}
25732bff 56 Float_t GetCopperRodShiftA(Int_t i) const {return fCopperRodShiftA[i];}// 0-17: IFC, 18-35; OFC
57 Float_t GetCopperRodShiftC(Int_t i) const {return fCopperRodShiftC[i];}// 0-17: IFC, 18-35; OFC
c9cbd2f2 58
59
60 void InitFCVoltError3D(); // Fill the lookup tables
61
62 virtual void Print(const Option_t* option="") const;
63
64protected:
65 virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
66
67private:
68
69 AliTPCFCVoltError3D(const AliTPCFCVoltError3D &); // not implemented
70 AliTPCFCVoltError3D &operator=(const AliTPCFCVoltError3D &); // not implemented
71
72 Float_t fC0; // coefficient C0 (compare Jim Thomas's notes for definitions)
73 Float_t fC1; // coefficient C1 (compare Jim Thomas's notes for definitions)
74 Float_t fRodVoltShiftA[36]; // Rod (&strips) shift in Volt (40V~1mm)
75 Float_t fRodVoltShiftC[36]; // Rod (&strips) shift in Volt (40V~1mm)
76 Float_t fRotatedClipVoltA[2]; // rotated clips at HV rod
77 Float_t fRotatedClipVoltC[2]; // rotated clips at HV rod
25732bff 78 Float_t fCopperRodShiftA[36]; // only Rod shift
79 Float_t fCopperRodShiftC[36]; // only Rod shift
c9cbd2f2 80
81 Bool_t fInitLookUp; // flag to check it the Look Up table was created (SUM)
25732bff 82 Bool_t fInitLookUpBasic[6]; // ! flag if the basic lookup was created (shifted Rod (IFC,OFC) or rotated clip (IFC,OFC))
c9cbd2f2 83
84
85 TMatrixD *fLookUpErOverEz[kNPhi]; // Array to store electric field integral (int Er/Ez)
86 TMatrixD *fLookUpEphiOverEz[kNPhi]; // Array to store electric field integral (int Er/Ez)
87 TMatrixD *fLookUpDeltaEz[kNPhi]; // Array to store electric field integral (int Er/Ez)
88
89 // basic numbers for the poisson relaxation //can be set individually in each class
90 enum {kRows =257}; // grid size in r direction used in the poisson relaxation // ( 2**n + 1 ) eg. 65, 129, 257 etc.
91 enum {kColumns=129}; // grid size in z direction used in the poisson relaxation // ( 2**m + 1 ) eg. 65, 129, 257 etc.
92 enum {kPhiSlicesPerSector = 10 }; // number of points in phi slices
93 enum {kPhiSlices = 1+kPhiSlicesPerSector*3 }; // number of points in phi for the basic lookup tables
94 enum {kIterations=100}; // Number of iterations within the poisson relaxation
95
96 // ugly way to store "partial" look up tables
97 // needed for the faster calculation of the final distortion table
98
99 // for Rod and Strip shift
25732bff 100 TMatrixD *fLookUpBasic1ErOverEz[kPhiSlices]; // ! Array to store electric field integral (int Er/Ez)
101 TMatrixD *fLookUpBasic1EphiOverEz[kPhiSlices]; // ! Array to store electric field integral (int Ephi/Ez)
102 TMatrixD *fLookUpBasic1DeltaEz[kPhiSlices]; // ! Array to store electric field integral (int Ez)
c9cbd2f2 103
25732bff 104 TMatrixD *fLookUpBasic2ErOverEz[kPhiSlices]; // ! Array to store electric field integral
105 TMatrixD *fLookUpBasic2EphiOverEz[kPhiSlices]; // ! Array to store electric field integral
106 TMatrixD *fLookUpBasic2DeltaEz[kPhiSlices]; // ! Array to store electric field integral
c9cbd2f2 107
108 // for rotated clips
25732bff 109 TMatrixD *fLookUpBasic3ErOverEz[kPhiSlices]; // ! Array to store electric field integral
110 TMatrixD *fLookUpBasic3EphiOverEz[kPhiSlices]; // ! Array to store electric field integral
111 TMatrixD *fLookUpBasic3DeltaEz[kPhiSlices]; // ! Array to store electric field integral
c9cbd2f2 112
25732bff 113 TMatrixD *fLookUpBasic4ErOverEz[kPhiSlices]; // ! Array to store electric field integral
114 TMatrixD *fLookUpBasic4EphiOverEz[kPhiSlices]; // ! Array to store electric field integral
115 TMatrixD *fLookUpBasic4DeltaEz[kPhiSlices]; // ! Array to store electric field integral
c9cbd2f2 116
25732bff 117 // for (only rod) shift (copper rods)
118 TMatrixD *fLookUpBasic5ErOverEz[kPhiSlices]; // ! Array to store electric field integral
119 TMatrixD *fLookUpBasic5EphiOverEz[kPhiSlices]; // ! Array to store electric field integral
120 TMatrixD *fLookUpBasic5DeltaEz[kPhiSlices]; // ! Array to store electric field integral
c9cbd2f2 121
25732bff 122 TMatrixD *fLookUpBasic6ErOverEz[kPhiSlices]; // ! Array to store electric field integral
123 TMatrixD *fLookUpBasic6EphiOverEz[kPhiSlices]; // ! Array to store electric field integral
124 TMatrixD *fLookUpBasic6DeltaEz[kPhiSlices]; // ! Array to store electric field integral
c9cbd2f2 125
25732bff 126
127 ClassDef(AliTPCFCVoltError3D,2); //
c9cbd2f2 128};
129
130#endif