]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCFCVoltError3D.h
Possibility to vary the cos(PA) cut
[u/mrichter/AliRoot.git] / TPC / AliTPCFCVoltError3D.h
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 ////////////////////////////////////////////////////////////////////////////
8 // AliTPCFCVoltError3D class                                              //
9 // Authors: Jim Thomas, Stefan Rossegger                                  //
10 ////////////////////////////////////////////////////////////////////////////
11
12 #include "AliTPCCorrection.h"
13
14
15 class AliTPCFCVoltError3D : public AliTPCCorrection {
16 public:
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
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;}
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)
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;}
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 ...
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;}
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
58
59
60   void InitFCVoltError3D(); // Fill the lookup tables
61   void ForceInitFCVoltError3D() { fInitLookUp=kFALSE; InitFCVoltError3D(); }
62
63   virtual void Print(const Option_t* option="") const;
64
65
66
67 protected:
68   virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
69
70 private:
71
72   AliTPCFCVoltError3D(const AliTPCFCVoltError3D &);               // not implemented
73   AliTPCFCVoltError3D &operator=(const AliTPCFCVoltError3D &);    // not implemented
74
75   Float_t fC0; // coefficient C0           (compare Jim Thomas's notes for definitions)
76   Float_t fC1; // coefficient C1           (compare Jim Thomas's notes for definitions)
77   Float_t fRodVoltShiftA[36];      // Rod (&strips) shift in Volt (40V~1mm) 
78   Float_t fRodVoltShiftC[36];      // Rod (&strips) shift in Volt (40V~1mm) 
79   Float_t fRotatedClipVoltA[2];    // rotated clips at HV rod
80   Float_t fRotatedClipVoltC[2];    // rotated clips at HV rod
81   Float_t fCopperRodShiftA[36];    // only Rod shift 
82   Float_t fCopperRodShiftC[36];    // only Rod shift 
83
84   Bool_t fInitLookUp;           // flag to check if the Look Up table was created (SUM)
85   Bool_t fInitLookUpBasic[6];   // ! flag if the basic lookup was created (shifted Rod (IFC,OFC) or rotated clip (IFC,OFC))
86
87
88   TMatrixF *fLookUpErOverEz[kNPhi];   // Array to store electric field integral (int Er/Ez)
89   TMatrixF *fLookUpEphiOverEz[kNPhi]; // Array to store electric field integral (int Er/Ez)
90   TMatrixF *fLookUpDeltaEz[kNPhi];    // Array to store electric field integral (int Er/Ez)
91
92   // basic numbers for the poisson relaxation //can be set individually in each class
93   enum {kRows   =257}; // grid size in r direction used in the poisson relaxation // ( 2**n + 1 ) eg. 65, 129, 257 etc.
94   enum {kColumns=129}; // grid size in z direction used in the poisson relaxation // ( 2**m + 1 ) eg. 65, 129, 257 etc.
95   enum {kPhiSlicesPerSector = 10 }; // number of points in phi slices
96   enum {kPhiSlices = 1+kPhiSlicesPerSector*3 };      // number of points in phi for the basic lookup tables
97   enum {kIterations=100}; // Number of iterations within the poisson relaxation 
98
99   // ugly way to store "partial" look up tables
100   // needed for the faster calculation of the final distortion table
101
102   // for Rod and Strip shift
103   TMatrixD *fLookUpBasic1ErOverEz[kPhiSlices];   // ! Array to store electric field integral (int Er/Ez)
104   TMatrixD *fLookUpBasic1EphiOverEz[kPhiSlices]; // ! Array to store electric field integral (int Ephi/Ez)
105   TMatrixD *fLookUpBasic1DeltaEz[kPhiSlices];    // ! Array to store electric field integral (int Ez)
106
107   TMatrixD *fLookUpBasic2ErOverEz[kPhiSlices];   // ! Array to store electric field integral 
108   TMatrixD *fLookUpBasic2EphiOverEz[kPhiSlices]; // ! Array to store electric field integral 
109   TMatrixD *fLookUpBasic2DeltaEz[kPhiSlices];    // ! Array to store electric field integral 
110
111   // for rotated clips
112   TMatrixD *fLookUpBasic3ErOverEz[kPhiSlices];   // ! Array to store electric field integral 
113   TMatrixD *fLookUpBasic3EphiOverEz[kPhiSlices]; // ! Array to store electric field integral 
114   TMatrixD *fLookUpBasic3DeltaEz[kPhiSlices];    // ! Array to store electric field integral 
115
116   TMatrixD *fLookUpBasic4ErOverEz[kPhiSlices];   // ! Array to store electric field integral 
117   TMatrixD *fLookUpBasic4EphiOverEz[kPhiSlices]; // ! Array to store electric field integral 
118   TMatrixD *fLookUpBasic4DeltaEz[kPhiSlices];    // ! Array to store electric field integral 
119
120   // for (only rod) shift (copper rods)
121   TMatrixD *fLookUpBasic5ErOverEz[kPhiSlices];   // ! Array to store electric field integral 
122   TMatrixD *fLookUpBasic5EphiOverEz[kPhiSlices]; // ! Array to store electric field integral 
123   TMatrixD *fLookUpBasic5DeltaEz[kPhiSlices];    // ! Array to store electric field integral 
124
125   TMatrixD *fLookUpBasic6ErOverEz[kPhiSlices];   // ! Array to store electric field integral 
126   TMatrixD *fLookUpBasic6EphiOverEz[kPhiSlices]; // ! Array to store electric field integral 
127   TMatrixD *fLookUpBasic6DeltaEz[kPhiSlices];    // ! Array to store electric field integral 
128
129
130   ClassDef(AliTPCFCVoltError3D,3); //
131 };
132
133 #endif