]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCorrection.h
Missing ; added.
[u/mrichter/AliRoot.git] / TPC / AliTPCCorrection.h
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 //                                                                            //
9 // AliTPCCorrection class                                                     //
10 //                                                                            //
11 // This class provides a general framework to deal with space point           //
12 // distortions. An correction class which inherits from here is for example   //
13 // AliTPCExBBShape or AliTPCExBTwist                                          //
14 //                                                                            //
15 // General functions are (for example):                                       //
16 //   CorrectPoint(x,roc) where x is the vector of inital positions in         //
17 //   cartesian coordinates and roc represents the Read Out chamber number     //
18 //   according to the offline naming convention. The vector x is overwritten  //
19 //   with the corrected coordinates.                                          //
20 //                                                                            //
21 // An alternative usage would be CorrectPoint(x,roc,dx), which leaves the     //
22 //   vector x untouched, put returns the distortions via the vector dx        //
23 //                                                                            //
24 // The class allows "effective Omega Tau" corrections to be shifted to the    //
25 // single distortion classes.                                                 //
26 //                                                                            //
27 // Note: This class is normally used via the class AliTPCComposedCorrection   //
28 //                                                                            //
29 // date: 27/04/2010                                                           //
30 // Authors: Magnus Mager, Stefan Rossegger, Jim Thomas                        //
31 ////////////////////////////////////////////////////////////////////////////////
32
33
34 #include <TNamed.h>
35 #include "TMatrixD.h"
36 class TH2F;
37 class TTimeStamp;
38 class TCollection;
39 class TTreeSRedirector;
40 class AliExternalTrackParam;
41 class TTree;
42 class THnSparse;
43 class AliESDVertex;
44
45
46 class AliTPCCorrection : public TNamed {
47 public:
48   enum CompositionType {kParallel,kQueue};
49
50   AliTPCCorrection();
51   AliTPCCorrection(const char *name,const char *title);
52   virtual ~AliTPCCorrection();
53   
54
55   // functions to correct a space point
56           void CorrectPoint (      Float_t x[],const Short_t roc);
57           void CorrectPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
58   virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
59
60   // functions to distort a space point
61           void DistortPoint (      Float_t x[],const Short_t roc);
62           void DistortPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
63   virtual void GetDistortion(const Float_t x[],const Short_t roc,Float_t dx[]);
64
65   // initialization and update functions
66   virtual void Init();
67   virtual void Update(const TTimeStamp &timeStamp);
68
69   // convenience functions
70   virtual void Print(Option_t* option="") const;
71  
72   TH2F* CreateHistoDRinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
73   TH2F* CreateHistoDRPhiinXY(Float_t z=10.,Int_t nx=100,Int_t nphi=100);
74   TH2F* CreateHistoDZinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
75
76   TH2F* CreateHistoDRinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
77   TH2F* CreateHistoDRPhiinZR(Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
78   TH2F* CreateHistoDZinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
79
80
81   TTree* CreateDistortionTree(Double_t step=5);
82   static void  MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run);
83   // normally called directly in the correction classes which inherit from this class
84   virtual void SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2);
85   AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir,TTreeSRedirector *pcstream);
86   void StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment=0);
87   static void MakeTrackDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, const TObjArray * corrArray, Int_t step=1, Bool_t debug=0);
88   static void MakeLaserDistortionTree(TTree* tree, TObjArray *corrArray, Int_t itype);
89
90   void FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts);
91
92   static void AddVisualCorrection(AliTPCCorrection* corr, Int_t position);
93   static Double_t GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0);
94   static Double_t GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0);
95
96 protected:
97   TH2F* CreateTH2F(const char *name,const char *title,
98                    const char *xlabel,const char *ylabel,const char *zlabel,
99                    Int_t nbinsx,Double_t xlow,Double_t xup,
100                    Int_t nbinsy,Double_t ylow,Double_t yup);
101  
102   static const Double_t fgkTPCZ0;       // nominal gating grid position 
103   static const Double_t fgkIFCRadius;   // Mean Radius of the Inner Field Cage ( 82.43 min,  83.70 max) (cm)
104   static const Double_t fgkOFCRadius;   // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
105   static const Double_t fgkZOffSet;     // Offset from CE: calculate all distortions closer to CE as if at this point
106   static const Double_t fgkCathodeV;    // Cathode Voltage (volts)
107   static const Double_t fgkGG;          // Gating Grid voltage (volts)
108   static const Double_t fgkdvdE;        // [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
109   static const Double_t fgkEM;          // charge/mass in [C/kg]
110   static const Double_t fgke0;          // vacuum permittivity [A·s/(V·m)]
111
112
113   enum {kNR=   72  };     // Number of R points in the table for interpolating distortion data
114   enum {kNPhi= 18*10+1};  // Number of Phi points in the table for interpolating distortion data ( plus one extra for 360 == 0 ) 
115   enum {kNZ=   166};      // Number of Z points in the table for interpolating distortion data
116
117   Double_t fgkRList[kNR]; // points in the radial direction (for the lookup table)
118   Double_t fgkPhiList[kNPhi]; // points in the phi direction (for the lookup table)
119   Double_t fgkZList[kNZ]; // points in the z direction (for the lookup table)
120
121   // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
122   Int_t fILow, fJLow, fKLow;          // variable to help in the interpolation 
123   void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z, 
124                                  const Double_t er[kNZ][kNR], Double_t &erValue );
125   void Interpolate3DEdistortion( const Int_t order, const Double_t r, const Float_t phi, const Double_t z, 
126                                  const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], 
127                                  const Double_t ez[kNZ][kNPhi][kNR],
128                                  Double_t &erValue, Double_t &ephiValue, Double_t &ezValue);
129   Double_t Interpolate2DTable( const Int_t order, const Double_t x, const Double_t y, 
130                               const Int_t nx,  const Int_t ny, const Double_t xv[], const Double_t yv[], 
131                               const TMatrixD &array );
132   Double_t Interpolate3DTable( const Int_t order, const Double_t x,   const Double_t y,   const Double_t z,
133                               const Int_t  nx,    const Int_t  ny,    const Int_t  nz,
134                               const Double_t xv[], const Double_t yv[], const Double_t zv[],
135                               TMatrixD **arrayofArrays );
136   Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], 
137                         const Int_t order, const Double_t x );
138   void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
139   virtual Int_t IsPowerOfTwo ( Int_t i ) const  ;
140
141   
142  
143   // Algorithms to solve the laplace or poisson equation 
144   void PoissonRelaxation2D(TMatrixD &arrayV, TMatrixD &chargeDensity, 
145                            TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz,
146                            const Int_t rows, const Int_t columns, const Int_t iterations,
147                            const Bool_t rocDisplacement = kTRUE);
148
149   void PoissonRelaxation3D( TMatrixD **arrayofArrayV, TMatrixD **arrayofChargeDensities, 
150                             TMatrixD **arrayofEroverEz, TMatrixD **arrayofEPhioverEz, TMatrixD **arrayofEz,
151                             const Int_t rows, const Int_t columns,  const Int_t phislices, 
152                             const Float_t deltaphi, const Int_t iterations, const Int_t summetry,
153                             const Bool_t rocDisplacement = kTRUE); 
154     
155 protected:
156   Double_t fT1;         // tensor term of wt - T1
157   Double_t fT2;         // tensor term of wt - T2
158   static TObjArray *fgVisualCorrection;  // array of orrection for visualization
159 private:
160   AliTPCCorrection(const AliTPCCorrection &);               // not implemented
161   AliTPCCorrection &operator=(const AliTPCCorrection &);    // not implemented
162
163   void InitLookUpfulcrums();   // to initialize the grid of the look up table
164
165   ClassDef(AliTPCCorrection,4);
166 };
167
168 #endif