AliRoot Core  v5-06-30 (35d6c57)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliTPCROCVoltError3D Class Reference

The class calculates the space point distortions due to z offsets of the. More...

#include <AliTPCROCVoltError3D.h>

Inheritance diagram for AliTPCROCVoltError3D:

Public Member Functions

 AliTPCROCVoltError3D ()
 
virtual ~AliTPCROCVoltError3D ()
 
virtual Bool_t AddCorrectionCompact (AliTPCCorrection *corr, Double_t weight)
 
virtual void Init ()
 
virtual void Update (const TTimeStamp &timeStamp)
 
virtual void SetOmegaTauT1T2 (Float_t omegaTau, Float_t t1, Float_t t2)
 
void SetC0C1 (Float_t c0, Float_t c1)
 
Float_t GetC0 () const
 
Float_t GetC1 () const
 
void SetROCData (TMatrixD *matrix)
 
void SetROCDataFileName (const char *fname)
 
const Char_t * GetROCDataFileName () const
 
void SetROCDisplacement (Bool_t flag)
 
Bool_t GetROCDisplacement () const
 
void SetElectronArrivalCorrection (Bool_t flag)
 
Bool_t GetElectronArrivalCorrection () const
 
void InitROCVoltError3D ()
 
void ForceInitROCVoltError3D ()
 
Float_t GetROCVoltOffset (Int_t side, Float_t r0, Float_t phi0)
 
TH2F * CreateHistoOfZAlignment (Int_t side, Int_t nx=250, Int_t ny=250)
 
virtual void Print (const Option_t *option="") const
 
TMatrixD * GetMatrix () const
 
- Public Member Functions inherited from AliTPCCorrection
 AliTPCCorrection ()
 
 AliTPCCorrection (const char *name, const char *title)
 
virtual ~AliTPCCorrection ()
 
void CorrectPoint (Float_t x[], Short_t roc)
 
void CorrectPointLocal (Float_t x[], Short_t roc)
 
void CorrectPoint (const Float_t x[], Short_t roc, Float_t xp[])
 
virtual void GetCorrectionDz (const Float_t x[], Short_t roc, Float_t dx[], Float_t delta)
 
virtual void GetCorrectionIntegralDz (const Float_t x[], Short_t roc, Float_t dx[], Float_t delta)
 
void DistortPoint (Float_t x[], Short_t roc)
 
void DistortPointLocal (Float_t x[], Short_t roc)
 
void DistortPoint (const Float_t x[], Short_t roc, Float_t xp[])
 
virtual void GetDistortion (const Float_t x[], Short_t roc, Float_t dx[])
 
virtual void GetDistortionDz (const Float_t x[], Short_t roc, Float_t dx[], Float_t delta)
 
virtual void GetDistortionIntegralDz (const Float_t x[], Short_t roc, Float_t dx[], Float_t delta)
 
virtual void SetCorrScaleFactor (Float_t)
 
virtual Float_t GetCorrScaleFactor () const
 
virtual void Print (Option_t *option="") const
 
TH2F * CreateHistoDRinXY (Float_t z=10., Int_t nx=100, Int_t ny=100)
 
TH2F * CreateHistoDRPhiinXY (Float_t z=10., Int_t nx=100, Int_t nphi=100)
 
TH2F * CreateHistoDZinXY (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 * CreateHistoDRPhiinZR (Float_t phi=0., Int_t nZ=100, Int_t nR=100)
 
TH2F * CreateHistoDZinZR (Float_t phi=0., Int_t nZ=100, Int_t nR=100)
 
TTree * CreateDistortionTree (Double_t step=5)
 
AliExternalTrackParam * FitDistortedTrack (AliExternalTrackParam &trackIn, Double_t refX, Int_t dir, TTreeSRedirector *pcstream)
 
void StoreInOCDB (Int_t startRun, Int_t endRun, const char *comment=0)
 
void FastSimDistortedVertex (Double_t orgVertex[3], Int_t nTracks, AliESDVertex &aV, AliESDVertex &avOrg, AliESDVertex &cV, AliESDVertex &cvOrg, TTreeSRedirector *const pcstream, Double_t etaCuts)
 

Protected Member Functions

virtual void GetCorrection (const Float_t x[], const Short_t roc, Float_t dx[])
 
- Protected Member Functions inherited from AliTPCCorrection
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)
 
void Interpolate2DEdistortion (Int_t order, Double_t r, Double_t z, const Double_t er[kNZ][kNR], Double_t &erValue)
 
void Interpolate3DEdistortion (Int_t order, Double_t r, Float_t phi, Double_t z, const Double_t er[kNZ][kNPhi][kNR], const Double_t ephi[kNZ][kNPhi][kNR], const Double_t ez[kNZ][kNPhi][kNR], Double_t &erValue, Double_t &ephiValue, Double_t &ezValue)
 
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)
 
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)
 
Double_t Interpolate (const Double_t xArray[], const Double_t yArray[], Int_t order, Double_t x)
 
void Search (Int_t n, const Double_t xArray[], Double_t x, Int_t &low)
 
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)
 
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)
 
Float_t Interpolate (const Double_t xArray[], const Float_t yArray[], Int_t order, Double_t x)
 
virtual Int_t IsPowerOfTwo (Int_t i) const
 
void PoissonRelaxation2D (TMatrixD &arrayV, TMatrixD &chargeDensity, TMatrixD &arrayErOverEz, TMatrixD &arrayDeltaEz, Int_t rows, Int_t columns, Int_t iterations, Bool_t rocDisplacement=kTRUE)
 
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)
 
void SetIsLocal (Bool_t isLocal)
 
Bool_t IsLocal () const
 

Private Types

enum  { kRows =257 }
 
enum  { kColumns =129 }
 
enum  { kPhiSlicesPerSector =10 }
 
enum  { kPhiSlices = 18*kPhiSlicesPerSector }
 
enum  { kIterations =100 }
 

Private Member Functions

 AliTPCROCVoltError3D (const AliTPCROCVoltError3D &)
 
AliTPCROCVoltError3Doperator= (const AliTPCROCVoltError3D &)
 

Private Attributes

Float_t fC0
 coefficient C0 (compare Jim Thomas's notes for definitions) More...
 
Float_t fC1
 coefficient C1 (compare Jim Thomas's notes for definitions) More...
 
Bool_t fROCdisplacement
 flag on wheter to consider the ROC displacement More...
 
Bool_t fElectronArrivalCorrection
 flag on wheter to consider the difference More...
 
Bool_t fInitLookUp
 flag to check it the Look Up table was created (SUM) More...
 
TMatrixF * fLookUpErOverEz [kNPhi]
 Array to store electric field integral (int Er/Ez) More...
 
TMatrixF * fLookUpEphiOverEz [kNPhi]
 Array to store electric field integral (int Er/Ez) More...
 
TMatrixF * fLookUpDeltaEz [kNPhi]
 Array to store electric field integral (int Er/Ez) More...
 
TString fROCDataFileName
 filename of the survey data containing the lin Fit values More...
 
TMatrixD * fdzDataLinFit
 Linear fits of dz survey points (each sector=72) (z0,slopeX,slopeY) More...
 

Additional Inherited Members

- Public Types inherited from AliTPCCorrection
enum  CompositionType { kParallel, kQueue }
 
- Static Public Member Functions inherited from AliTPCCorrection
static void MakeDistortionMap (THnSparse *his0, TTreeSRedirector *pcstream, const char *hname, Int_t run, Float_t refX, Int_t type, Int_t integ=1)
 
static void MakeDistortionMapCosmic (THnSparse *his0, TTreeSRedirector *pcstream, const char *hname, Int_t run, Float_t refX, Int_t type)
 
static void MakeDistortionMapSector (THnSparse *his0, TTreeSRedirector *pcstream, const char *hname, Int_t run, Int_t type)
 
static 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)
 
static 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)
 
static void MakeLaserDistortionTreeOld (TTree *tree, TObjArray *corrArray, Int_t itype)
 
static void MakeLaserDistortionTree (TTree *tree, TObjArray *corrArray, Int_t itype)
 
static void AddVisualCorrection (AliTPCCorrection *corr, Int_t position)
 
static AliTPCCorrectionGetVisualCorrection (Int_t position)
 
static AliTPCCorrectionGetVisualCorrection (const char *corName)
 
static TObjArrayGetVisualCorrections ()
 
static Double_t GetCorrSector (Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0)
 
static Double_t GetCorrXYZ (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0)
 
static Double_t GetCorrXYZDz (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5)
 
static Double_t GetCorrXYZIntegrateZ (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5)
 
static Double_t GetDistXYZ (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0)
 
static Double_t GetDistXYZDz (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5)
 
static Double_t GetDistXYZIntegrateZ (Double_t gx, Double_t gy, Double_t gz, Int_t axisType, Int_t corrType=0, Double_t delta=5)
 
- Protected Types inherited from AliTPCCorrection
enum  { kNR = 72 }
 
enum  { kNPhi = 18*10+1 }
 
enum  { kNZ = 166 }
 
- Protected Attributes inherited from AliTPCCorrection
Double_t fgkRList [kNR]
 points in the radial direction (for the lookup table) More...
 
Double_t fgkPhiList [kNPhi]
 points in the phi direction (for the lookup table) More...
 
Double_t fgkZList [kNZ]
 points in the z direction (for the lookup table) More...
 
Int_t fILow
 
Int_t fJLow
 
Int_t fKLow
 variable to help in the interpolation More...
 
Double_t fT1
 tensor term of wt - T1 More...
 
Double_t fT2
 tensor term of wt - T2 More...
 
Bool_t fIsLocal
 switch to indicate that the distortion is a local vector drphi/dz, dr/dz More...
 
- Static Protected Attributes inherited from AliTPCCorrection
static const Double_t fgkTPCZ0 = 249.7
 nominal gating grid position More...
 
static const Double_t fgkIFCRadius = 83.5
 Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm) More...
 
static const Double_t fgkOFCRadius = 254.5
 Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm) More...
 
static const Double_t fgkZOffSet = 0.2
 Offset from CE: calculate all distortions closer to CE as if at this point. More...
 
static const Double_t fgkCathodeV = -100000.0
 Cathode Voltage (volts) More...
 
static const Double_t fgkGG = -70.0
 Gating Grid voltage (volts) More...
 
static const Double_t fgkdvdE = 0.0024
 [cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment) More...
 
static const Double_t fgkEM = -1.602176487e-19/9.10938215e-31
 charge/mass in [C/kg] More...
 
static const Double_t fgke0 = 8.854187817e-12
 vacuum permittivity [A·s/(V·m)] More...
 
static TObjArrayfgVisualCorrection =0
 array of orrection for visualization More...
 

Detailed Description

The class calculates the space point distortions due to z offsets of the.

ROCs via the residual voltage technique (Poisson relaxation) in 3D. Since the GG (part of the ROCs) represents the closure of the FC in z direction, every misalignment in z produces not only dz distortions but also electrical field inhomogeneities throughout the volume, which produces additional dr and rd$$ distortions.

Each ROC can be misaligned (in z direction) in three ways. A general z0 offset, an inclination along the x and an inclination along the y axis. The z-misalignment's can be set via the function SetROCData(TMatrixD *mat) for each single chamber (ROC). The array size has to be (72,3) representing the 72 chambers according to the offline numbering scheme (IROC: roc$<$36; OROC: roc$$36) and the three misalignment's (see the source code for further details).

Internally, these z offsets (unit is cm) are recalculated into residual voltage equivalents in order to make use of the relaxation technique.

One has two possibilities when calculating the $dz$ distortions. The resulting distortions are purely due to the change of the drift velocity (along with the change of the drift field) when the SetROCDisplacement is FALSE. For this class, this is a rather unphysical setting and should be avoided. When the flag is set to TRUE, the corresponding offset in z is added to the dz calculation of the outer ROCs. For the Alice TPC gas, both effects are of similar magnitude. This means, if the drift length is sufficiently large, a z-offset of a chamber appears to have (approx.) twice the magnitude when one looks at the actual dz distortions.

In addition, this class allows a correction regarding the different arrival times of the electrons due to the geometrical difference of the inner and outer chambers. The impact was simulated via Garfield. If the flag is set via the function SetElectronArrivalCorrection, the electron-arrival correction is added to the dz calculation.

AliTPCROCVoltError3D_cxx_99b6733.png
Picture from ROOT macro
Author
Jim Thomas, Stefan Rossegger
Date
08/08/2010
Author
Jim Thomas, Stefan Rossegger

Definition at line 14 of file AliTPCROCVoltError3D.h.

Member Enumeration Documentation

anonymous enum
private
Enumerator
kRows 

Definition at line 88 of file AliTPCROCVoltError3D.h.

anonymous enum
private
Enumerator
kColumns 

Definition at line 89 of file AliTPCROCVoltError3D.h.

anonymous enum
private
Enumerator
kPhiSlicesPerSector 

Definition at line 90 of file AliTPCROCVoltError3D.h.

anonymous enum
private
Enumerator
kPhiSlices 

Definition at line 91 of file AliTPCROCVoltError3D.h.

anonymous enum
private
Enumerator
kIterations 

Definition at line 92 of file AliTPCROCVoltError3D.h.

Constructor & Destructor Documentation

AliTPCROCVoltError3D::AliTPCROCVoltError3D ( )

Definition at line 70 of file AliTPCROCVoltError3D.cxx.

AliTPCROCVoltError3D::~AliTPCROCVoltError3D ( )
virtual

destructor

Definition at line 96 of file AliTPCROCVoltError3D.cxx.

AliTPCROCVoltError3D::AliTPCROCVoltError3D ( const AliTPCROCVoltError3D )
private

Member Function Documentation

Bool_t AliTPCROCVoltError3D::AddCorrectionCompact ( AliTPCCorrection corr,
Double_t  weight 
)
virtual

Add correction and make them compact Assumptions:

  • origin of distortion/correction are additive
  • only correction ot the same type supported ()

Reimplemented from AliTPCCorrection.

Definition at line 110 of file AliTPCROCVoltError3D.cxx.

TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment ( Int_t  side,
Int_t  nx = 250,
Int_t  ny = 250 
)

return a simple histogramm containing the input to the poisson solver (z positions of the Read-out chambers, linearly interpolated)

Definition at line 475 of file AliTPCROCVoltError3D.cxx.

void AliTPCROCVoltError3D::ForceInitROCVoltError3D ( )
inline

Definition at line 53 of file AliTPCROCVoltError3D.h.

Referenced by GetCorrection().

Float_t AliTPCROCVoltError3D::GetC0 ( ) const
inline

Definition at line 30 of file AliTPCROCVoltError3D.h.

Float_t AliTPCROCVoltError3D::GetC1 ( ) const
inline

Definition at line 31 of file AliTPCROCVoltError3D.h.

void AliTPCROCVoltError3D::GetCorrection ( const Float_t  x[],
const Short_t  roc,
Float_t  dx[] 
)
protectedvirtual

Calculates the correction due e.g. residual voltage errors on the TPC boundaries

Reimplemented from AliTPCCorrection.

Definition at line 212 of file AliTPCROCVoltError3D.cxx.

Bool_t AliTPCROCVoltError3D::GetElectronArrivalCorrection ( ) const
inline

Definition at line 49 of file AliTPCROCVoltError3D.h.

TMatrixD* AliTPCROCVoltError3D::GetMatrix ( ) const
inline
const Char_t* AliTPCROCVoltError3D::GetROCDataFileName ( ) const
inline

Definition at line 35 of file AliTPCROCVoltError3D.h.

Bool_t AliTPCROCVoltError3D::GetROCDisplacement ( ) const
inline

Definition at line 42 of file AliTPCROCVoltError3D.h.

Float_t AliTPCROCVoltError3D::GetROCVoltOffset ( Int_t  side,
Float_t  r0,
Float_t  phi0 
)

Returns the dz alignment data (in voltage equivalents) at the given position

Definition at line 438 of file AliTPCROCVoltError3D.cxx.

Referenced by CreateHistoOfZAlignment(), and InitROCVoltError3D().

void AliTPCROCVoltError3D::Init ( )
virtual

Initialization funtion

Reimplemented from AliTPCCorrection.

Definition at line 152 of file AliTPCROCVoltError3D.cxx.

void AliTPCROCVoltError3D::InitROCVoltError3D ( )

Initialization of the Lookup table which contains the solutions of the Dirichlet boundary problem Calculation of the single 3D-Poisson solver is done just if needed (see basic lookup tables in header file)

Definition at line 297 of file AliTPCROCVoltError3D.cxx.

Referenced by ForceInitROCVoltError3D(), GetCorrection(), and Init().

AliTPCROCVoltError3D& AliTPCROCVoltError3D::operator= ( const AliTPCROCVoltError3D )
private
void AliTPCROCVoltError3D::Print ( const Option_t *  option = "") const
virtual

Print function to check the settings of the Rod shifts and the rotated clips option=="a" prints the C0 and C1 coefficents for calibration purposes

Definition at line 515 of file AliTPCROCVoltError3D.cxx.

void AliTPCROCVoltError3D::SetC0C1 ( Float_t  c0,
Float_t  c1 
)
inline

Definition at line 29 of file AliTPCROCVoltError3D.h.

void AliTPCROCVoltError3D::SetElectronArrivalCorrection ( Bool_t  flag)
inline
virtual void AliTPCROCVoltError3D::SetOmegaTauT1T2 ( Float_t  omegaTau,
Float_t  t1,
Float_t  t2 
)
inlinevirtual

Reimplemented from AliTPCCorrection.

Definition at line 24 of file AliTPCROCVoltError3D.h.

Referenced by Init(), and Update().

void AliTPCROCVoltError3D::SetROCData ( TMatrixD *  matrix)

Set a z alignment map of the chambers not via a file, but directly via a TMatrix(72,3), where dz = p0 + p1*(lx-133.4) + p2*ly (all in cm)

Definition at line 143 of file AliTPCROCVoltError3D.cxx.

Referenced by RegisterAliTPCROCVoltError3D(), RegisterAliTPCROCVoltError3DSector(), and TestCorrection_AliTPCRocVoltError3DAddCorrectionCompact().

void AliTPCROCVoltError3D::SetROCDataFileName ( const char *  fname)

Set / load the ROC data (linear fit of ROC misalignments)

Definition at line 185 of file AliTPCROCVoltError3D.cxx.

void AliTPCROCVoltError3D::SetROCDisplacement ( Bool_t  flag)
inline

Definition at line 39 of file AliTPCROCVoltError3D.h.

void AliTPCROCVoltError3D::Update ( const TTimeStamp &  timeStamp)
virtual

Update function

Reimplemented from AliTPCCorrection.

Definition at line 169 of file AliTPCROCVoltError3D.cxx.

Member Data Documentation

Float_t AliTPCROCVoltError3D::fC0
private

coefficient C0 (compare Jim Thomas's notes for definitions)

Definition at line 69 of file AliTPCROCVoltError3D.h.

Referenced by GetC0(), GetCorrection(), Print(), SetC0C1(), and SetOmegaTauT1T2().

Float_t AliTPCROCVoltError3D::fC1
private

coefficient C1 (compare Jim Thomas's notes for definitions)

Definition at line 70 of file AliTPCROCVoltError3D.h.

Referenced by GetC1(), GetCorrection(), Print(), SetC0C1(), and SetOmegaTauT1T2().

TMatrixD* AliTPCROCVoltError3D::fdzDataLinFit
private

Linear fits of dz survey points (each sector=72) (z0,slopeX,slopeY)

Definition at line 85 of file AliTPCROCVoltError3D.h.

Referenced by AddCorrectionCompact(), GetMatrix(), GetROCVoltOffset(), Print(), SetROCData(), SetROCDataFileName(), and ~AliTPCROCVoltError3D().

Bool_t AliTPCROCVoltError3D::fElectronArrivalCorrection
private

flag on wheter to consider the difference

Definition at line 74 of file AliTPCROCVoltError3D.h.

Referenced by AddCorrectionCompact(), GetCorrection(), GetElectronArrivalCorrection(), Print(), and SetElectronArrivalCorrection().

Bool_t AliTPCROCVoltError3D::fInitLookUp
private

flag to check it the Look Up table was created (SUM)

Definition at line 78 of file AliTPCROCVoltError3D.h.

Referenced by ForceInitROCVoltError3D(), GetCorrection(), Init(), InitROCVoltError3D(), Print(), SetElectronArrivalCorrection(), SetROCDataFileName(), and SetROCDisplacement().

TMatrixF* AliTPCROCVoltError3D::fLookUpDeltaEz[kNPhi]
private

Array to store electric field integral (int Er/Ez)

Definition at line 82 of file AliTPCROCVoltError3D.h.

Referenced by GetCorrection(), InitROCVoltError3D(), and ~AliTPCROCVoltError3D().

TMatrixF* AliTPCROCVoltError3D::fLookUpEphiOverEz[kNPhi]
private

Array to store electric field integral (int Er/Ez)

Definition at line 81 of file AliTPCROCVoltError3D.h.

Referenced by GetCorrection(), InitROCVoltError3D(), and ~AliTPCROCVoltError3D().

TMatrixF* AliTPCROCVoltError3D::fLookUpErOverEz[kNPhi]
private

Array to store electric field integral (int Er/Ez)

Definition at line 80 of file AliTPCROCVoltError3D.h.

Referenced by GetCorrection(), InitROCVoltError3D(), and ~AliTPCROCVoltError3D().

TString AliTPCROCVoltError3D::fROCDataFileName
private

filename of the survey data containing the lin Fit values

Definition at line 84 of file AliTPCROCVoltError3D.h.

Referenced by GetROCDataFileName(), Print(), and SetROCDataFileName().

Bool_t AliTPCROCVoltError3D::fROCdisplacement
private

flag on wheter to consider the ROC displacement

Definition at line 72 of file AliTPCROCVoltError3D.h.

Referenced by AddCorrectionCompact(), GetROCDisplacement(), InitROCVoltError3D(), Print(), and SetROCDisplacement().


The documentation for this class was generated from the following files: