]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCorrection.h
Warning removal
[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* CreateHistoDRinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
75   TH2F* CreateHistoDRPhiinZR(Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
76   TTree* CreateDistortionTree(Double_t step=5);
77   static void  MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run);
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, Bool_t debug=0);
83   static void MakeLaserDistortionTree(TTree* tree, TObjArray *corrArray, Int_t itype);
84   void  FastSimDistortedVertex(Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector * const pcstream, Double_t etaCuts);
85
86 protected:
87   TH2F* CreateTH2F(const char *name,const char *title,
88                    const char *xlabel,const char *ylabel,const char *zlabel,
89                    Int_t nbinsx,Double_t xlow,Double_t xup,
90                    Int_t nbinsy,Double_t ylow,Double_t yup);
91  
92   static const Double_t fgkTPCZ0;      // nominal gating grid position 
93   static const Double_t fgkIFCRadius;   // Mean Radius of the Inner Field Cage ( 82.43 min,  83.70 max) (cm)
94   static const Double_t fgkOFCRadius;   // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
95   static const Double_t fgkZOffSet;     // Offset from CE: calculate all distortions closer to CE as if at this point
96   static const Double_t fgkCathodeV;    // Cathode Voltage (volts)
97   static const Double_t fgkGG;          // Gating Grid voltage (volts)
98
99   enum {kNR=   92};              // Number of R points in the table for interpolating distortion data
100   enum {kNZ=  270};              // Number of Z points in the table for interpolating distortion data
101   static const Double_t fgkRList[kNR]; // points in the radial direction (for the lookup table)
102   static const Double_t fgkZList[kNZ]; // points in the z direction (for the lookup table)
103
104   // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
105   Int_t fJLow;         // variable to help in the interpolation 
106   Int_t fKLow;         // variable to help in the interpolation 
107   void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z, 
108                                  const Double_t er[kNZ][kNR], Double_t &erValue );
109   Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], 
110                         const Int_t order, const Double_t x );
111   void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
112   virtual Int_t IsPowerOfTwo ( Int_t i ) const  ;
113     
114   // Algorithms to solve the laplace or possion equation 
115   void PoissonRelaxation2D(TMatrixD &arrayV, const TMatrixD &chargeDensity, TMatrixD &arrayErOverEz, const Int_t rows, const Int_t columns, const Int_t iterations );
116
117
118 protected:
119   Double_t fT1;         // tensor term of wt - T1
120   Double_t fT2;         // tensor term of wt - T2
121 private:
122   ClassDef(AliTPCCorrection,2);
123 };
124
125 #endif