]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCorrection.h
Split task for ITS tracks: 1. ITS tracks with TPC partner, 2. pure ITS standalone
[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 class TH2F;
36 class TTimeStamp;
37 class TCollection;
38 class TTreeSRedirector;
39 class AliExternalTrackParam;
40 class TTree;
41 class THnSparse;
42
43 class AliTPCCorrection : public TNamed {
44 public:
45   enum CompositionType {kParallel,kQueue};
46
47   AliTPCCorrection();
48   AliTPCCorrection(const char *name,const char *title);
49   virtual ~AliTPCCorrection();
50   
51
52   // functions to correct a space point
53           void CorrectPoint (      Float_t x[],const Short_t roc);
54           void CorrectPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
55   virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
56
57   // functions to distort a space point
58           void DistortPoint (      Float_t x[],const Short_t roc);
59           void DistortPoint (const Float_t x[],const Short_t roc,Float_t xp[]);
60   virtual void GetDistortion(const Float_t x[],const Short_t roc,Float_t dx[]);
61
62   // initialization and update functions
63   virtual void Init();
64   virtual void Update(const TTimeStamp &timeStamp);
65
66   // convenience functions
67   virtual void Print(Option_t* option="") const;
68  
69   TH2F* CreateHistoDRinXY   (Float_t z=10.,Int_t nx=100,Int_t ny=100);
70   TH2F* CreateHistoDRPhiinXY(Float_t z=10.,Int_t nx=100,Int_t nphi=100);
71   TH2F* CreateHistoDRinZR   (Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
72   TH2F* CreateHistoDRPhiinZR(Float_t phi=0.,Int_t nZ=100,Int_t nR=100);
73   TTree* CreateDistortionTree(Double_t step=5);
74   static void  MakeDistortionMap(THnSparse * his0, TTreeSRedirector *pcstream, const char* hname, Int_t run);
75   // normally called directly in the correction classes which inherit from this class
76   virtual void SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2);
77   AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam & trackIn, Double_t refX, Int_t dir,TTreeSRedirector *pcstream);
78   void StoreInOCDB(Int_t startRun, Int_t endRun, const char *comment=0);
79   static void MakeTrackDistortionTree(TTree *tinput, Int_t dtype, Int_t ptype, TObjArray * corrArray, Int_t step=1, Bool_t debug=0);
80 protected:
81   TH2F* CreateTH2F(const char *name,const char *title,
82                    const char *xlabel,const char *ylabel,const char *zlabel,
83                    Int_t nbinsx,Double_t xlow,Double_t xup,
84                    Int_t nbinsy,Double_t ylow,Double_t yup);
85  
86   static const Double_t fgkTPC_Z0;      // nominal gating grid position 
87   static const Double_t fgkIFCRadius;   // Mean Radius of the Inner Field Cage ( 82.43 min,  83.70 max) (cm)
88   static const Double_t fgkOFCRadius;   // Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
89   static const Double_t fgkZOffSet;     // Offset from CE: calculate all distortions closer to CE as if at this point
90   static const Double_t fgkCathodeV;    // Cathode Voltage (volts)
91   static const Double_t fgkGG;          // Gating Grid voltage (volts)
92
93   enum {kNR=   92};              // Number of R points in the table for interpolating distortion data
94   enum {kNZ=  270};              // Number of Z points in the table for interpolating distortion data
95   static const Double_t fgkRList[kNR];
96   static const Double_t fgkZList[kNZ];
97
98   // Simple Interpolation functions: e.g. with tricubic interpolation (not yet in TH3)
99   Int_t fJLow; 
100   Int_t fKLow;
101   void Interpolate2DEdistortion( const Int_t order, const Double_t r, const Double_t z, 
102                                  const Double_t er[kNZ][kNR], Double_t &er_value );
103   Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], 
104                         const Int_t order, const Double_t x );
105   void Search( const Int_t n, const Double_t xArray[], const Double_t x, Int_t &low );
106   
107 protected:
108   Double_t fT1;         // tensor term of wt - T1
109   Double_t fT2;         // tensor term of wt - T2
110
111   ClassDef(AliTPCCorrection,2);
112 };
113
114 #endif