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