ROOT logo
AliRoot » TPC » BASE » AliTPCCorrection

class AliTPCCorrection: public TNamed



AliTPCCorrection class

The AliTPCCorrection class provides a general framework to deal with space point distortions. An correction class which inherits from here is for example AliTPCExBBShape or AliTPCExBTwist.
General virtual functions are (for example) CorrectPoint(x,roc) where x is the vector of initial positions in cartesian coordinates and roc represents the read-out chamber number according to the offline numbering convention. The vector x is overwritten with the corrected coordinates.
An alternative usage would be CorrectPoint(x,roc,dx), which leaves the vector x untouched, but returns the distortions via the vector dx.
This class is normally used via the general class AliTPCComposedCorrection.

Furthermore, the class contains basic geometrical descriptions like field cage radii (fgkIFCRadius, fgkOFCRadius) and length (fgkTPCZ0) plus the voltages. Also, the definitions of size and widths of the fulcrums building the grid of the final look-up table, which is then interpolated, is defined in kNX and fgkXList).

All physics-model classes below are derived from this class in order to not duplicate code and to allow a uniform treatment of all physics models.

Poisson solver

A numerical solver of the Poisson equation (relaxation technique) is implemented for 2-dimensional geometries (r,z) as well as for 3-dimensional problems (r,$\phi$,z). The corresponding function names are PoissonRelaxation?D. The relevant function arguments are the arrays of the boundary and initial conditions (ArrayofArrayV, ArrayofChargeDensities) as well as the grid granularity which is used during the calculation. These inputs can be chosen according to the needs of the physical effect which is supposed to be simulated. In the 3D version, different symmetry conditions can be set in order to reduce the calculation time (used in AliTPCFCVoltError3D).

Unified plotting functionality

Generic plot functions were implemented. They return a histogram pointer in the chosen plane of the TPC drift volume with a selectable grid granularity and the magnitude of the correction vector. For example, the function CreateHistoDZinXY(z,nx,ny) returns a 2-dimensional histogram which contains the longitudinal corrections $dz$ in the (x,y)-plane at the given z position with the granularity of nx and ny. The magnitude of the corrections is defined by the class from which this function is called. In the same manner, standard plots for the (r,$\phi$)-plane and for the other corrections like $dr$ and $rd\phi$ are available

Note: This class is normally used via the class AliTPCComposedCorrection



output of MACRO_AliTPCCorrection_1_cAliTPCCorrection
   {
   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
   TCanvas *c2 = new TCanvas("cAliTPCCorrection","cAliTPCCorrection",700,1050);  c2->Divide(2,3);
   AliTPCROCVoltError3D roc; // EXAMPLE PLOTS - SEE BELOW
   SetROCDataFileName("$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root");
   roc.SetOmegaTauT1T2(0,1,1); 
   Float_t z0 = 1; 
   c2->cd(1); roc.CreateHistoDRinXY(1.,300,300)->Draw("cont4z"); 
   c2->cd(3);roc.CreateHistoDRPhiinXY(1.,300,300)->Draw("cont4z"); 
   c2->cd(5);roc.CreateHistoDZinXY(1.,300,300)->Draw("cont4z"); 
   Float_t phi0=0.5;
   c2->cd(2);roc.CreateHistoDRinZR(phi0)->Draw("surf2"); 
   c2->cd(4);roc.CreateHistoDRPhiinZR(phi0)->Draw("surf2"); 
   c2->cd(6);roc.CreateHistoDZinZR(phi0)->Draw("surf2"); 
   return c2;
   } 
 
   

Date: 27/04/2010
Authors: Magnus Mager, Stefan Rossegger, Jim Thomas


Function Members (Methods)

public:
AliTPCCorrection()
AliTPCCorrection(const char* name, const char* title)
virtual~AliTPCCorrection()
voidTObject::AbstractMethod(const char* method) const
virtual Bool_tAddCorrectionCompact(AliTPCCorrection* corr, Double_t weight)
static voidAddVisualCorrection(AliTPCCorrection* corr, Int_t position)
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
virtual voidTNamed::Copy(TObject& named) const
voidCorrectPoint(Float_t* x, Short_t roc)
voidCorrectPoint(const Float_t* x, Short_t roc, Float_t* xp)
voidCorrectPointLocal(Float_t* x, Short_t roc)
TTree*CreateDistortionTree(Double_t step = 5)
TH2F*CreateHistoDRinXY(Float_t z = 10., Int_t nx = 100, Int_t ny = 100)
TH2F*CreateHistoDRinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
TH2F*CreateHistoDRPhiinXY(Float_t z = 10., Int_t nx = 100, Int_t nphi = 100)
TH2F*CreateHistoDRPhiinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
TH2F*CreateHistoDZinXY(Float_t z = 10., Int_t nx = 100, Int_t ny = 100)
TH2F*CreateHistoDZinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
voidDistortPoint(Float_t* x, Short_t roc)
voidDistortPoint(const Float_t* x, Short_t roc, Float_t* xp)
voidDistortPointLocal(Float_t* x, Short_t roc)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
voidFastSimDistortedVertex(Double_t* orgVertex, Int_t nTracks, AliESDVertex& aV, AliESDVertex& avOrg, AliESDVertex& cV, AliESDVertex& cvOrg, TTreeSRedirector *const pcstream, Double_t etaCuts)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
AliExternalTrackParam*FitDistortedTrack(AliExternalTrackParam& trackIn, Double_t refX, Int_t dir, TTreeSRedirector* pcstream)
virtual voidGetCorrection(const Float_t* x, Short_t roc, Float_t* dx)
virtual voidGetCorrectionDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
virtual voidGetCorrectionIntegralDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
virtual Float_tGetCorrScaleFactor() const
static Double_tGetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType = 0)
static Double_tGetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0)
static Double_tGetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
static Double_tGetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
virtual voidGetDistortion(const Float_t* x, Short_t roc, Float_t* dx)
virtual voidGetDistortionDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
virtual voidGetDistortionIntegralDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
static Double_tGetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0)
static Double_tGetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
static Double_tGetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
virtual const char*TNamed::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
static AliTPCCorrection*GetVisualCorrection(Int_t position)
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidInit()
virtual voidTObject::Inspect() constMENU
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
static voidMakeDistortionMap(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Float_t refX, Int_t type, Int_t integ = 1)
static voidMakeDistortionMapCosmic(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Float_t refX, Int_t type)
static voidMakeDistortionMapSector(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Int_t type)
static voidMakeLaserDistortionTree(TTree* tree, TObjArray* corrArray, Int_t itype)
static voidMakeLaserDistortionTreeOld(TTree* tree, TObjArray* corrArray, Int_t itype)
static voidMakeSectorDistortionTree(TTree* tinput, Int_t dtype, Int_t ptype, const TObjArray* corrArray, Int_t step = 1, Int_t offset = 0, Bool_t debug = 0)
static voidMakeTrackDistortionTree(TTree* tinput, Int_t dtype, Int_t ptype, const TObjArray* corrArray, Int_t step = 1, Int_t offset = 0, Bool_t debug = 0)
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTObject::Notify()
voidTObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidPrint(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(basic_ostream<char,char_traits<char> >& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidSetCorrScaleFactor(Float_t)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
virtual voidTNamed::SetName(const char* name)MENU
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
virtual voidSetOmegaTauT1T2(Float_t omegaTau, Float_t t1, Float_t t2)
virtual voidTNamed::SetTitle(const char* title = "")MENU
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector&)
virtual Int_tTNamed::Sizeof() const
voidStoreInOCDB(Int_t startRun, Int_t endRun, const char* comment = 0)
virtual voidStreamer(TBuffer&)
voidStreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidUpdate(const TTimeStamp& timeStamp)
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
TH2F*CreateTH2F(const char* name, const char* title, const char* xlabel, const char* ylabel, const char* zlabel, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
Double_tInterpolate(const Double_t* xArray, const Double_t* yArray, Int_t order, Double_t x)
Float_tInterpolate(const Double_t* xArray, const Float_t* yArray, Int_t order, Double_t x)
voidInterpolate2DEdistortion(Int_t order, Double_t r, Double_t z, const Double_t er[][kNR], Double_t& erValue)
Double_tInterpolate2DTable(Int_t order, Double_t x, Double_t y, Int_t nx, Int_t ny, const Double_t* xv, const Double_t* yv, const TMatrixD& array)
Float_tInterpolate2DTable(Int_t order, Double_t x, Double_t y, Int_t nx, Int_t ny, const Double_t* xv, const Double_t* yv, const TMatrixF& array)
voidInterpolate3DEdistortion(Int_t order, Double_t r, Float_t phi, Double_t z, const Double_t er[][kNPhi][kNR], const Double_t ephi[][kNPhi][kNR], const Double_t ez[][kNPhi][kNR], Double_t& erValue, Double_t& ephiValue, Double_t& ezValue)
Double_tInterpolate3DTable(Int_t order, Double_t x, Double_t y, Double_t z, Int_t nx, Int_t ny, Int_t nz, const Double_t* xv, const Double_t* yv, const Double_t* zv, TMatrixD** arrayofArrays)
Float_tInterpolate3DTable(Int_t order, Double_t x, Double_t y, Double_t z, Int_t nx, Int_t ny, Int_t nz, const Double_t* xv, const Double_t* yv, const Double_t* zv, TMatrixF** arrayofArrays)
Bool_tIsLocal() const
virtual Int_tIsPowerOfTwo(Int_t i) const
voidTObject::MakeZombie()
voidPoissonRelaxation2D(TMatrixD& arrayV, TMatrixD& chargeDensity, TMatrixD& arrayErOverEz, TMatrixD& arrayDeltaEz, Int_t rows, Int_t columns, Int_t iterations, Bool_t rocDisplacement = kTRUE)
voidPoissonRelaxation3D(TMatrixD** arrayofArrayV, TMatrixD** arrayofChargeDensities, TMatrixD** arrayofEroverEz, TMatrixD** arrayofEPhioverEz, TMatrixD** arrayofEz, Int_t rows, Int_t columns, Int_t phislices, Float_t deltaphi, Int_t iterations, Int_t summetry, Bool_t rocDisplacement = kTRUE)
voidSearch(Int_t n, const Double_t* xArray, Double_t x, Int_t& low)
voidSetIsLocal(Bool_t isLocal)

Data Members

public:
enum CompositionType { kParallel
kQueue
};
enum { kNR
kNPhi
kNZ
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
Int_tfILow
Bool_tfIsLocalswitch to indicate that the distortion is a local vector drphi/dz, dr/dz
Int_tfJLow
Int_tfKLowvariable to help in the interpolation
TStringTNamed::fNameobject identifier
Double_tfT1tensor term of wt - T1
Double_tfT2tensor term of wt - T2
TStringTNamed::fTitleobject title
static TObjArray*fgVisualCorrectionarray of orrection for visualization
static const Double_tfgkCathodeVCathode Voltage (volts)
static const Double_tfgkEMcharge/mass in [C/kg]
static const Double_tfgkGGGating Grid voltage (volts)
static const Double_tfgkIFCRadiusMean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)
static const Double_tfgkOFCRadiusMean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
Double_tfgkPhiList[181]points in the phi direction (for the lookup table)
Double_tfgkRList[72]points in the radial direction (for the lookup table)
static const Double_tfgkTPCZ0nominal gating grid position
Double_tfgkZList[166]points in the z direction (for the lookup table)
static const Double_tfgkZOffSetOffset from CE: calculate all distortions closer to CE as if at this point
static const Double_tfgkdvdE[cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
static const Double_tfgke0vacuum permittivity [A·s/(V·m)]

Class Charts

Inheritance Chart:
TNamed
AliTPCCorrection
AliTPCBoundaryVoltError
AliTPCCalibGlobalMisalignment
AliTPCComposedCorrection
AliTPCCorrectionDrift
AliTPCCorrectionLookupTable
AliTPCExBBShape
AliTPCExBEffective
AliTPCExBEffectiveSector
AliTPCExBTwist
AliTPCFCVoltError3D
 [more...]

Function documentation

AliTPCCorrection()
 default constructor

AliTPCCorrection(const char* name, const char* title)
 default constructor, that set the name and title

~AliTPCCorrection()
 virtual destructor

Bool_t AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight)
 Add correction  and make them compact
 Assumptions:
  - origin of distortion/correction are additive
  - only correction ot the same type supported ()
void CorrectPoint(Float_t* x, Short_t roc)
 Corrects the initial coordinates x (cartesian coordinates)
 according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void CorrectPoint(const Float_t* x, Short_t roc, Float_t* xp)
 Corrects the initial coordinates x (cartesian coordinates) and stores the new
 (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void DistortPoint(Float_t* x, Short_t roc)
 Distorts the initial coordinates x (cartesian coordinates)
 according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void DistortPointLocal(Float_t* x, Short_t roc)
 Distorts the initial coordinates x (cartesian coordinates)
 according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void CorrectPointLocal(Float_t* x, Short_t roc)
 Distorts the initial coordinates x (cartesian coordinates)
 according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void DistortPoint(const Float_t* x, Short_t roc, Float_t* xp)
 Distorts the initial coordinates x (cartesian coordinates) and stores the new
 (distorted) coordinates in xp. The distortion is set according to the given effect (inherited classes)
 roc represents the TPC read out chamber (offline numbering convention)

void GetCorrection(const Float_t* x, Short_t roc, Float_t* dx)
 This function delivers the correction values dx in respect to the inital coordinates x
 roc represents the TPC read out chamber (offline numbering convention)
 Note: The dx is overwritten by the inherited effectice class ...

void GetDistortion(const Float_t* x, Short_t roc, Float_t* dx)
 This function delivers the distortion values dx in respect to the inital coordinates x
 roc represents the TPC read out chamber (offline numbering convention)

void GetCorrectionDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
 author: marian.ivanov@cern.ch

 In this (virtual)function calculates the dx'/dz,  dy'/dz  and dz'/dz at given point (x,y,z)
 Generic implementation. Better precision can be acchieved knowing the internal structure
 of underlying trasnformation. Derived classes can reimplement it.
 To calculate correction is fitted in small neighberhood:
 (x+-delta,y+-delta,z+-delta) where delta is an argument

 Input parameters:
   x[]   - space point corrdinate
   roc   - readout chamber identifier (important e.g to do not miss the side of detector)
   delta - define the size of neighberhood
 Output parameter:
   dx[] - array {dx'/dz,  dy'/dz ,  dz'/dz }
void GetDistortionDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
 author: marian.ivanov@cern.ch

 In this (virtual)function calculates the dx'/dz,  dy'/dz  and dz'/dz at given point (x,y,z)
 Generic implementation. Better precision can be acchieved knowing the internal structure
 of underlying trasnformation. Derived classes can reimplement it.
 To calculate distortion is fitted in small neighberhood:
 (x+-delta,y+-delta,z+-delta) where delta is an argument

 Input parameters:
   x[]   - space point corrdinate
   roc   - readout chamber identifier (important e.g to do not miss the side of detector)
   delta - define the size of neighberhood
 Output parameter:
   dx[] - array {dx'/dz,  dy'/dz ,  dz'/dz }
void GetCorrectionIntegralDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
 Integrate 3D distortion along drift lines starting from the roc plane
   to the expected z position of the point, this assumes that dz is small
   and the error propagating to z' instead of the correct z is negligible
 To define the drift lines virtual function  AliTPCCorrection::GetCorrectionDz is used

 Input parameters:
   x[]   - space point corrdinate
   roc   - readout chamber identifier (important e.g to do not miss the side of detector)
   delta - define the size of neighberhood
 Output parameter:
   dx[] - array { integral(dx'/dz),  integral(dy'/dz) ,  integral(dz'/dz) }
void GetDistortionIntegralDz(const Float_t* x, Short_t roc, Float_t* dx, Float_t delta)
 Integrate 3D distortion along drift lines
 To define the drift lines virtual function  AliTPCCorrection::GetCorrectionDz is used

 Input parameters:
   x[]   - space point corrdinate
   roc   - readout chamber identifier (important e.g to do not miss the side of detector)
   delta - define the size of neighberhood
 Output parameter:
   dx[] - array { integral(dx'/dz),  integral(dy'/dz) ,  integral(dz'/dz) }
void Init()
 Initialization funtion (not used at the moment)

void Update(const TTimeStamp& timeStamp)
 Update function

void Print(Option_t* option = "") const
 Print function to check which correction classes are used
 option=="d" prints details regarding the setted magnitude
 option=="a" prints the C0 and C1 coefficents for calibration purposes

TH2F* CreateHistoDRinXY(Float_t z = 10., Int_t nx = 100, Int_t ny = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in radial direction (dr)
 in respect to position z within the XY plane.
 The histogramm has nx times ny entries.

TH2F* CreateHistoDRPhiinXY(Float_t z = 10., Int_t nx = 100, Int_t nphi = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
 in respect to position z within the XY plane.
 The histogramm has nx times ny entries.

TH2F* CreateHistoDZinXY(Float_t z = 10., Int_t nx = 100, Int_t ny = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
 in respect to position z within the XY plane.
 The histogramm has nx times ny entries.

TH2F* CreateHistoDRinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in r direction (dr)
 in respect to angle phi within the ZR plane.
 The histogramm has nx times ny entries.

TH2F* CreateHistoDRPhiinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in rphi direction (drphi)
 in respect to angle phi within the ZR plane.
 The histogramm has nx times ny entries.

TH2F* CreateHistoDZinZR(Float_t phi = 0., Int_t nZ = 100, Int_t nR = 100)
 Simple plot functionality.
 Returns a 2d hisogram which represents the corrections in longitudinal direction (dz)
 in respect to angle phi within the ZR plane.
 The histogramm has nx times ny entries.

TH2F* CreateTH2F(const char* name, const char* title, const char* xlabel, const char* ylabel, const char* zlabel, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup)
 Helper function to create a 2d histogramm of given size

void Interpolate2DEdistortion(Int_t order, Double_t r, Double_t z, const Double_t er[][kNR], Double_t& erValue)
 Interpolate table - 2D interpolation

{0,0,0,0,0}
void Interpolate3DEdistortion(Int_t order, Double_t r, Float_t phi, Double_t z, const Double_t er[][kNPhi][kNR], const Double_t ephi[][kNPhi][kNR], const Double_t ez[][kNPhi][kNR], Double_t& erValue, Double_t& ephiValue, Double_t& ezValue)
 Interpolate table - 3D interpolation

{0,0,0,0,0}
Double_t Interpolate2DTable(Int_t order, Double_t x, Double_t y, Int_t nx, Int_t ny, const Double_t* xv, const Double_t* yv, const TMatrixD& array)
 Interpolate table (TMatrix format) - 2D interpolation

Double_t Interpolate3DTable(Int_t order, Double_t x, Double_t y, Double_t z, Int_t nx, Int_t ny, Int_t nz, const Double_t* xv, const Double_t* yv, const Double_t* zv, TMatrixD** arrayofArrays)
 Interpolate table (TMatrix format) - 3D interpolation

Double_t Interpolate( const Double_t xArray[], const Double_t yArray[], Int_t order, Double_t x )
 Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.

Float_t Interpolate2DTable(Int_t order, Double_t x, Double_t y, Int_t nx, Int_t ny, const Double_t* xv, const Double_t* yv, const TMatrixF& array)
 Interpolate table (TMatrix format) - 2D interpolation
 Float version (in order to decrease the OCDB size)

Float_t Interpolate3DTable(Int_t order, Double_t x, Double_t y, Double_t z, Int_t nx, Int_t ny, Int_t nz, const Double_t* xv, const Double_t* yv, const Double_t* zv, TMatrixF** arrayofArrays)
 Interpolate table (TMatrix format) - 3D interpolation
 Float version (in order to decrease the OCDB size)

Float_t Interpolate( const Double_t xArray[], const Float_t yArray[], Int_t order, Double_t x )
 Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
 Float version (in order to decrease the OCDB size)

void Search(Int_t n, const Double_t* xArray, Double_t x, Int_t& low)
 Search an ordered table by starting at the most recently used point

void InitLookUpfulcrums()
 Initialization of interpolation points - for main look up table
   (course grid in the middle, fine grid on the borders)

void PoissonRelaxation2D(TMatrixD& arrayV, TMatrixD& chargeDensity, TMatrixD& arrayErOverEz, TMatrixD& arrayDeltaEz, Int_t rows, Int_t columns, Int_t iterations, Bool_t rocDisplacement = kTRUE)
 Solve Poisson's Equation by Relaxation Technique in 2D (assuming cylindrical symmetry)

 Solve Poissons equation in a cylindrical coordinate system. The arrayV matrix must be filled with the
 boundary conditions on the first and last rows, and the first and last columns.  The remainder of the
 array can be blank or contain a preliminary guess at the solution.  The Charge density matrix contains
 the enclosed spacecharge density at each point. The charge density matrix can be full of zero's if
 you wish to solve Laplaces equation however it should not contain random numbers or you will get
 random numbers back as a solution.
 Poisson's equation is solved by iteratively relaxing the matrix to the final solution.  In order to
 speed up the convergence to the best solution, this algorithm does a binary expansion of the solution
 space.  First it solves the problem on a very sparse grid by skipping rows and columns in the original
 matrix.  Then it doubles the number of points and solves the problem again.  Then it doubles the
 number of points and solves the problem again.  This happens several times until the maximum number
 of points has been included in the array.

 NOTE: In order for this algorithmto work, the number of rows and columns must be a power of 2 plus one.
 So rows == 2**M + 1 and columns == 2**N + 1.  The number of rows and columns can be different.

 NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation

 Original code by Jim Thomas (STAR TPC Collaboration)

void PoissonRelaxation3D(TMatrixD** arrayofArrayV, TMatrixD** arrayofChargeDensities, TMatrixD** arrayofEroverEz, TMatrixD** arrayofEPhioverEz, TMatrixD** arrayofEz, Int_t rows, Int_t columns, Int_t phislices, Float_t deltaphi, Int_t iterations, Int_t summetry, Bool_t rocDisplacement = kTRUE)
 3D - Solve Poisson's Equation in 3D by Relaxation Technique

    NOTE: In order for this algorith to work, the number of rows and columns must be a power of 2 plus one.
    The number of rows and COLUMNS can be different.

    ROWS       ==  2**M + 1
    COLUMNS    ==  2**N + 1
    PHISLICES  ==  Arbitrary but greater than 3

    DeltaPhi in Radians

    SYMMETRY = 0 if no phi symmetries, and no phi boundary conditions
             = 1 if we have reflection symmetry at the boundaries (eg. sector symmetry or half sector symmetries).

 NOTE: rocDisplacement is used to include (or ignore) the ROC misalignment in the dz calculation
Int_t IsPowerOfTwo(Int_t i) const
 Helperfunction: Check if integer is a power of 2

AliExternalTrackParam * FitDistortedTrack(AliExternalTrackParam& trackIn, Double_t refX, Int_t dir, TTreeSRedirector* pcstream)
 Fit the track parameters - without and with distortion
 1. Space points in the TPC are simulated along the trajectory
 2. Space points distorted
 3. Fits  the non distorted and distroted track to the reference plane at refX
 4. For visualization and debugging  purposes the space points and tracks can be stored  in the tree - using the TTreeSRedirector functionality

 trackIn   - input track parameters
 refX     - reference X to fit the track
 dir      - direction - out=1 or in=-1
 pcstream -  debug streamer to check the results

 see AliExternalTrackParam.h documentation:
 track1.fP[0] - local y (rphi)
 track1.fP[1] - z
 track1.fP[2] - sinus of local inclination angle
 track1.fP[3] - tangent of deep angle
 track1.fP[4] - 1/pt
TTree* CreateDistortionTree(Double_t step = 5)
 create the distortion tree on a mesh with granularity given by step
 return the tree with distortions at given position
 Map is created on the mesh with given step size

void MakeTrackDistortionTree(TTree* tinput, Int_t dtype, Int_t ptype, const TObjArray* corrArray, Int_t step = 1, Int_t offset = 0, Bool_t debug = 0)
 Make a fit tree:
 For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
 calculates partial distortions
 Partial distortion is stored in the resulting tree
 Output is storred in the file distortion_<dettype>_<partype>.root
 Partial  distortion is stored with the name given by correction name


 Parameters of function:
 input     - input tree
 dtype     - distortion type 0 - ITSTPC,  1 -TPCTRD, 2 - TPCvertex , 3 - TPC-TOF,  4 - TPCTPC track crossing
 ppype     - parameter type
 corrArray - array with partial corrections
 step      - skipe entries  - if 1 all entries processed - it is slow
 debug     0 if debug on also space points dumped - it is slow
void MakeSectorDistortionTree(TTree* tinput, Int_t dtype, Int_t ptype, const TObjArray* corrArray, Int_t step = 1, Int_t offset = 0, Bool_t debug = 0)
 Make a fit tree:
 For each partial correction (specified in array) and given track topology (phi, theta, snp, refX)
 calculates partial distortions
 Partial distortion is stored in the resulting tree
 Output is storred in the file distortion_<dettype>_<partype>.root
 Partial  distortion is stored with the name given by correction name


 Parameters of function:
 input     - input tree
 dtype     - distortion type 10 - IROC-OROC
 ppype     - parameter type
 corrArray - array with partial corrections
 step      - skipe entries  - if 1 all entries processed - it is slow
 debug     0 if debug on also space points dumped - it is slow
void MakeLaserDistortionTreeOld(TTree* tree, TObjArray* corrArray, Int_t itype)
 Make a laser fit tree for global minimization

void MakeDistortionMap(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Float_t refX, Int_t type, Int_t integ = 1)
 make a distortion map out ou fthe residual histogram
 Results are written to the debug streamer - pcstream
 Parameters:
   his0       - input (4D) residual histogram
   pcstream   - file to write the tree
   run        - run number
   refX       - track matching reference X
   type       - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
 THnSparse axes:
 OBJ: TAxis     #Delta  #Delta
 OBJ: TAxis     tanTheta        tan(#Theta)
 OBJ: TAxis     phi     #phi
 OBJ: TAxis     snp     snp
void MakeDistortionMapCosmic(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Float_t refX, Int_t type)
 make a distortion map out ou fthe residual histogram
 Results are written to the debug streamer - pcstream
 Parameters:
   his0       - input (4D) residual histogram
   pcstream   - file to write the tree
   run        - run number
   refX       - track matching reference X
   type       - 0- y 1-z,2 -snp, 3-theta, 4=1/pt
 marian.ivanov@cern.ch

  Histo axeses
   Collection name='TObjArray', class='TObjArray', size=16
  0. OBJ: TAxis     #Delta  #Delta
  1. OBJ: TAxis     N_{cl}  N_{cl}
  2. OBJ: TAxis     dca_{r} (cm)    dca_{r} (cm)
  3. OBJ: TAxis     z (cm)  z (cm)
  4. OBJ: TAxis     sin(#phi)       sin(#phi)
  5. OBJ: TAxis     tan(#theta)     tan(#theta)
  6. OBJ: TAxis     1/pt (1/GeV)    1/pt (1/GeV)
  7. OBJ: TAxis     pt (GeV)        pt (GeV)
  8. OBJ: TAxis     alpha   alpha
void MakeDistortionMapSector(THnSparse* his0, TTreeSRedirector* pcstream, const char* hname, Int_t run, Int_t type)
 make a distortion map out of the residual histogram
 Results are written to the debug streamer - pcstream
 Parameters:
   his0       - input (4D) residual histogram
   pcstream   - file to write the tree
   run        - run number
   type       - 0- y 1-z,2 -snp, 3-theta
 marian.ivanov@cern.ch
void StoreInOCDB(Int_t startRun, Int_t endRun, const char* comment = 0)
 Store object in the OCDB
 By default the object is stored in the current directory
 default comment consit of user name and the date

void FastSimDistortedVertex(Double_t* orgVertex, Int_t nTracks, AliESDVertex& aV, AliESDVertex& avOrg, AliESDVertex& cV, AliESDVertex& cvOrg, TTreeSRedirector *const pcstream, Double_t etaCuts)
 Fast method to simulate the influence of the given distortion on the vertex reconstruction

void AddVisualCorrection(AliTPCCorrection* corr, Int_t position)
 make correction available for visualization using
 TFormula, TFX and TTree::Draw
 important in order to check corrections and also compute dervied variables
 e.g correction partial derivatives

 NOTE - class is not owner of correction

AliTPCCorrection* GetVisualCorrection(Int_t position)
 Get visula correction registered at index=position

Double_t GetCorrSector(Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType = 0)
 calculate the correction at given position - check the geffCorr

 corrType return values
 0 - delta R
 1 - delta RPhi
 2 - delta Z
 3 - delta RPHI

Double_t GetCorrXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0)
 return correction at given x,y,z

Double_t GetCorrXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
 return correction at given x,y,z

Double_t GetCorrXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
 return correction at given x,y,z

Double_t GetDistXYZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0)
 return correction at given x,y,z

Double_t GetDistXYZDz(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
 return correction at given x,y,z

Double_t GetDistXYZIntegrateZ(Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType = 0, Double_t delta = 5)
 return correction at given x,y,z

void MakeLaserDistortionTree(TTree* tree, TObjArray* corrArray, Int_t itype)
 Make a laser fit tree for global minimization

AliTPCCorrection()
void SetCorrScaleFactor(Float_t )
 map scaling
{ ; }
Float_t GetCorrScaleFactor() const
{ return 1.; }
void SetOmegaTauT1T2(Float_t omegaTau, Float_t t1, Float_t t2)
 normally called directly in the correction classes which inherit from this class
void SetIsLocal(Bool_t isLocal)
{fIsLocal=isLocal;}
Bool_t IsLocal() const
{ return fIsLocal;}
AliTPCCorrection & operator=(const AliTPCCorrection& )