![]() |
AliRoot Core
3dc7879 (3dc7879)
|
AliTPCCorrection class. More...
#include <AliTPCCorrection.h>
Public Types | |
enum | CompositionType { kParallel, kQueue } |
enum | IntegrationType { kIntegral, kDifferential } |
Public Member Functions | |
AliTPCCorrection () | |
AliTPCCorrection (const char *name, const char *title) | |
virtual | ~AliTPCCorrection () |
virtual Bool_t | AddCorrectionCompact (AliTPCCorrection *corr, Double_t weight) |
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 | GetCorrection (const Float_t x[], Short_t roc, Float_t dx[]) |
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 | Init () |
virtual void | Update (const TTimeStamp &timeStamp) |
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) |
TH2F * | CreateHistoDRinPhiR (Float_t z=10., Int_t nPhi=180, Int_t nR=100, Bool_t useSector=kTRUE, Float_t shift=0.) |
TH2F * | CreateHistoDRPhiinPhiR (Float_t z=10., Int_t nPhi=180, Int_t nR=100, Bool_t useSector=kTRUE, Float_t shift=0.) |
TH2F * | CreateHistoDZinPhiR (Float_t z=10., Int_t nPhi=180, Int_t nR=100, Bool_t useSector=kTRUE, Float_t shift=0.) |
TTree * | CreateDistortionTree (Double_t step=5, Int_t type=0) |
virtual void | SetOmegaTauT1T2 (Float_t omegaTau, Float_t t1, Float_t t2) |
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) |
Static Public Member Functions | |
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 AliTPCCorrection * | GetVisualCorrection (Int_t position) |
static AliTPCCorrection * | GetVisualCorrection (const char *corName) |
static TObjArray * | GetVisualCorrections () |
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) |
static Double_t | GetCorrectionSector (Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0) |
static Double_t | GetDistortionSector (Double_t sector, Double_t r, Double_t kZ, Int_t axisType, Int_t corrType=0) |
Protected Types | |
enum | { kNR = 72 } |
enum | { kNPhi = 18*10+1 } |
enum | { kNZ = 166 } |
Protected Member Functions | |
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, IntegrationType integrationType=kIntegral) |
void | SetIsLocal (Bool_t isLocal) |
Bool_t | IsLocal () const |
Protected Attributes | |
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... | |
IntegrationType | fIntegrationType |
Presentation of the underlying corrections, integrated, or differential. More... | |
Static Protected Attributes | |
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)/hera/alice/wiechula/calib/guiTrees. 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 TObjArray * | fgVisualCorrection =0 |
array of orrection for visualization More... | |
Private Member Functions | |
AliTPCCorrection (const AliTPCCorrection &) | |
AliTPCCorrection & | operator= (const AliTPCCorrection &) |
void | InitLookUpfulcrums () |
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.
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,$$,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).
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,$$)-plane and for the other corrections like $dr$ and $rd$ are available
Note: This class is normally used via the class AliTPCComposedCorrection
Definition at line 26 of file AliTPCCorrection.h.
|
protected |
Enumerator | |
---|---|
kNR |
Definition at line 130 of file AliTPCCorrection.h.
|
protected |
Enumerator | |
---|---|
kNPhi |
Definition at line 131 of file AliTPCCorrection.h.
|
protected |
Enumerator | |
---|---|
kNZ |
Definition at line 132 of file AliTPCCorrection.h.
Enumerator | |
---|---|
kParallel | |
kQueue |
Definition at line 28 of file AliTPCCorrection.h.
Enumerator | |
---|---|
kIntegral | |
kDifferential |
Definition at line 29 of file AliTPCCorrection.h.
AliTPCCorrection::AliTPCCorrection | ( | ) |
default constructor
Definition at line 126 of file AliTPCCorrection.cxx.
AliTPCCorrection::AliTPCCorrection | ( | const char * | name, |
const char * | title | ||
) |
default constructor, that set the name and title
Definition at line 137 of file AliTPCCorrection.cxx.
|
virtual |
virtual destructor
Definition at line 148 of file AliTPCCorrection.cxx.
|
private |
|
virtual |
Add correction and make them compact Assumptions:
Reimplemented in AliTPCComposedCorrection, AliTPCExBTwist, AliTPCExBBShape, AliTPCCalibGlobalMisalignment, AliTPCBoundaryVoltError, AliTPCFCVoltError3D, and AliTPCROCVoltError3D.
Definition at line 153 of file AliTPCCorrection.cxx.
Referenced by AliTPCComposedCorrection::AddCorrectionCompact().
|
static |
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
Definition at line 3219 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), LoadModels(), MakeAlignCorrection(), MakeAlignFunctionGlobal(), MakeAlignFunctionSector(), MakeEfieldCorrection(), MakeLaserDistortionTree(), MakeQA(), MakeQuadrantCorrection(), RegisterAlignFunction(), RegisterAliTPCBoundaryVoltError(), RegisterAliTPCCalibGlobalMisalignment(), RegisterAliTPCExBShape(), RegisterAliTPCExBTwist(), RegisterAliTPCFCVoltError3D(), RegisterAliTPCFCVoltError3DRodFCSideRadiusType(), RegisterAliTPCROCVoltError3D(), RegisterAliTPCROCVoltError3DSector(), TestCorrection_AliTPCComposedCorrectionAddCorrectionCompact(), TestCorrection_AliTPCComposedCorrectionAddCorrectionCompact_TPCCalibCorrection(), and UpdateEffSectorOCDB().
void AliTPCCorrection::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)
Definition at line 168 of file AliTPCCorrection.cxx.
Referenced by CorrectPointLocal(), AliTPCCorrectionLookupTable::CreateResidual(), GetCorrectionSector(), GetCorrXYZ(), MakeLaserDistortionTreeOld(), and AliTPCTransform::Transform().
void AliTPCCorrection::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)
Definition at line 178 of file AliTPCCorrection.cxx.
void AliTPCCorrection::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)
Definition at line 214 of file AliTPCCorrection.cxx.
TTree * AliTPCCorrection::CreateDistortionTree | ( | Double_t | step = 5 , |
Int_t | type = 0 |
||
) |
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 type - 0: Call GetDistortion() 1: Call GetDistortionIntegralDz()
Definition at line 1993 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::CreateHistoDRinPhiR | ( | Float_t | z = 10. , |
Int_t | nPhi = 180 , |
||
Int_t | nR = 100 , |
||
Bool_t | useSector = kTRUE , |
||
Float_t | shift = 0. |
||
) |
Simple plot functionality. Returns a 2d hisogram which represents the corrections in radial direction (dr) in respect to position z within the PhiR plane. The histogramm has nPhi times nR entries.
Definition at line 717 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::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.
Definition at line 538 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::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.
Definition at line 636 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::CreateHistoDRPhiinPhiR | ( | Float_t | z = 10. , |
Int_t | nPhi = 180 , |
||
Int_t | nR = 100 , |
||
Bool_t | useSector = kTRUE , |
||
Float_t | shift = 0. |
||
) |
Simple plot functionality. Returns a 2d hisogram which represents the corrections in rphi direction (drphi) in respect to position z within the PhiR plane. The histogramm has nPhi times nR entries.
Definition at line 750 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::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.
Definition at line 569 of file AliTPCCorrection.cxx.
Referenced by FitRodShift(), GetCorrScaleFactor(), and MakeComposedCorrection().
TH2F * AliTPCCorrection::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.
Definition at line 662 of file AliTPCCorrection.cxx.
Referenced by AliTPCCorrectionDemo(), GetCorrScaleFactor(), and MakeComposedCorrection().
TH2F * AliTPCCorrection::CreateHistoDZinPhiR | ( | Float_t | z = 10. , |
Int_t | nPhi = 180 , |
||
Int_t | nR = 100 , |
||
Bool_t | useSector = kTRUE , |
||
Float_t | shift = 0. |
||
) |
Simple plot functionality. Returns a 2d hisogram which represents the corrections in longitudinal direction (dz) in respect to position z within the PhiR plane. The histogramm has nPhi times nR entries.
Definition at line 789 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::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.
Definition at line 606 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
TH2F * AliTPCCorrection::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.
Definition at line 693 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
|
protected |
Helper function to create a 2d histogramm of given size
Definition at line 822 of file AliTPCCorrection.cxx.
Referenced by CreateHistoDRinPhiR(), CreateHistoDRinXY(), CreateHistoDRinZR(), CreateHistoDRPhiinPhiR(), CreateHistoDRPhiinXY(), CreateHistoDRPhiinZR(), CreateHistoDZinPhiR(), CreateHistoDZinXY(), CreateHistoDZinZR(), AliTPCSpaceCharge3DDriftLine::CreateHistogramCorrDRInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramCorrDRPhiInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramCorrDZInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramDistDRInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramDistDRPhiInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramDistDZInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramSCInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramSCInZR(), AliTPCSpaceCharge2D2D::CreateHistoSCinXY(), AliTPCSpaceCharge3D::CreateHistoSCinXY(), AliTPCSpaceCharge2D2D::CreateHistoSCinZR(), AliTPCSpaceCharge3D::CreateHistoSCinZR(), and GetVisualCorrections().
void AliTPCCorrection::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)
Definition at line 188 of file AliTPCCorrection.cxx.
Referenced by AliTPCCorrectionLookupTable::CreateResidual(), DistortPointLocal(), FitDistortedTrack(), GetCorr(), GetCorrSector(), GetDistortionSector(), GetDistXYZ(), and MakeLaserDistortionTreeOld().
void AliTPCCorrection::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)
Definition at line 231 of file AliTPCCorrection.cxx.
void AliTPCCorrection::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)
Definition at line 198 of file AliTPCCorrection.cxx.
void AliTPCCorrection::FastSimDistortedVertex | ( | Double_t | orgVertex[3], |
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
Definition at line 3131 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and TestVertex().
AliExternalTrackParam * AliTPCCorrection::FitDistortedTrack | ( | AliExternalTrackParam & | trackIn, |
Double_t | refX, | ||
Int_t | dir, | ||
TTreeSRedirector * | pcstream | ||
) |
Fit the track parameters - without and with distortion
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
Definition at line 1820 of file AliTPCCorrection.cxx.
Referenced by FastSimDistortedVertex(), GetCorrScaleFactor(), MakeDistortionMap(), MakeSectorDistortionTree(), and MakeTrackDistortionTree().
|
virtual |
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 ...
Reimplemented in AliTPCFCVoltError3D, AliTPCCalibGlobalMisalignment, AliTPCGGVoltError, AliTPCSpaceCharge3D, AliTPCROCVoltError3D, AliTPCSpaceCharge3DDriftLine, AliTPCExBTwist, AliTPCBoundaryVoltError, AliTPCSpaceCharge2D2D, AliTPCExBBShape, AliTPCComposedCorrection, AliTPCExBConical, AliTPCSpaceCharge, AliTPCExBEffective, AliTPCExBEffectiveSector, AliTPCInverseCorrection, AliTPCCorrectionDrift, and AliTPCCorrectionLookupTable.
Definition at line 241 of file AliTPCCorrection.cxx.
Referenced by CorrectPoint(), CreateDistortionTree(), CreateHistoDRinPhiR(), CreateHistoDRinXY(), CreateHistoDRinZR(), CreateHistoDRPhiinPhiR(), CreateHistoDRPhiinXY(), CreateHistoDRPhiinZR(), CreateHistoDZinPhiR(), CreateHistoDZinXY(), CreateHistoDZinZR(), AliTPCCorrectionLookupTable::CreateLookupTablePhiBin(), AliTPCComposedCorrection::GetCorrection(), GetCorrectionDz(), AliTPCInverseCorrection::GetDistortion(), AliTPCComposedCorrection::GetDistortion(), and GetDistortion().
|
virtual |
author: maria n.iv anov@ 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 }
Definition at line 257 of file AliTPCCorrection.cxx.
Referenced by GetCorrectionIntegralDz(), and GetCorrXYZDz().
|
virtual |
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) }
Definition at line 415 of file AliTPCCorrection.cxx.
Referenced by CreateDistortionTree(), AliTPCCorrectionLookupTable::CreateLookupTablePhiBin(), and GetCorrXYZIntegrateZ().
|
static |
calculate the correction at given position - check the geffCorr
corrType return values 0 - delta R 1 - delta RPhi 2 - delta Z 3 - delta RPHI
Definition at line 3275 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
inlinevirtual |
Reimplemented in AliTPCCorrectionLookupTable.
Definition at line 60 of file AliTPCCorrection.h.
|
static |
calculate the correction at given position - check the geffCorr
corrType return values 0 - delta R 1 - delta RPhi 2 - delta Z 3 - delta RPHI
Definition at line 3241 of file AliTPCCorrection.cxx.
Referenced by DeltaLookup(), and GetVisualCorrections().
|
static |
return correction at given x,y,z
Definition at line 3345 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections(), MakeLaserDistortionTree(), TestCorrection_AliTPCComposedCorrectionAddCorrectionCompact(), and TestCorrection_AliTPCComposedCorrectionAddCorrectionCompact_TPCCalibCorrection().
|
static |
return correction at given x,y,z
Definition at line 3364 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
static |
return correction at given x,y,z
Definition at line 3388 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
virtual |
This function delivers the distortion values dx in respect to the inital coordinates x roc represents the TPC read out chamber (offline numbering convention)
Reimplemented in AliTPCSpaceCharge3DDriftLine, AliTPCComposedCorrection, AliTPCInverseCorrection, and AliTPCCorrectionLookupTable.
Definition at line 249 of file AliTPCCorrection.cxx.
Referenced by CreateDistortionTree(), AliTPCCorrectionLookupTable::CreateLookupTablePhiBin(), DistortPoint(), AliTPCInverseCorrection::GetCorrection(), AliTPCComposedCorrection::GetDistortion(), and GetDistortionDz().
|
virtual |
author: maria n.iv anov@ 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 }
Definition at line 343 of file AliTPCCorrection.cxx.
Referenced by GetDistortionIntegralDz(), and GetDistXYZDz().
|
virtual |
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) }
Definition at line 466 of file AliTPCCorrection.cxx.
Referenced by CreateDistortionTree(), AliTPCCorrectionLookupTable::CreateLookupTablePhiBin(), and GetDistXYZIntegrateZ().
|
static |
calculate the correction at given position - check the geffCorr
corrType return values 0 - delta R 1 - delta RPhi 2 - delta Z 3 - delta RPHI
Definition at line 3310 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
static |
return correction at given x,y,z
Definition at line 3413 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
static |
return correction at given x,y,z
Definition at line 3432 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
static |
return correction at given x,y,z
Definition at line 3456 of file AliTPCCorrection.cxx.
Referenced by GetVisualCorrections().
|
static |
Get visula correction registered at index=position
Definition at line 3233 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and TestParExample().
|
inlinestatic |
Definition at line 94 of file AliTPCCorrection.h.
|
inlinestatic |
Definition at line 95 of file AliTPCCorrection.h.
|
virtual |
Initialization funtion (not used at the moment)
Reimplemented in AliTPCSpaceCharge3DDriftLine, AliTPCComposedCorrection, AliTPCGGVoltError, AliTPCInverseCorrection, AliTPCExBTwist, AliTPCExBBShape, AliTPCBoundaryVoltError, AliTPCExBEffectiveSector, AliTPCSpaceCharge3D, AliTPCExBConical, AliTPCFCVoltError3D, AliTPCSpaceCharge, AliTPCCorrectionDrift, AliTPCROCVoltError3D, AliTPCExBEffective, and AliTPCSpaceCharge2D2D.
Definition at line 510 of file AliTPCCorrection.cxx.
Referenced by AliTPCcalibDB::GetTPCComposedCorrectionDelta(), AliTPCInverseCorrection::Init(), AliTPCComposedCorrection::Init(), AliTPCcalibDB::Update(), and AliTPCcalibDB::UpdateRunInformations().
|
private |
Initialization of interpolation points - for main look up table (course grid in the middle, fine grid on the borders)
Definition at line 1115 of file AliTPCCorrection.cxx.
Referenced by AliTPCCorrection().
|
protected |
Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation.
Definition at line 972 of file AliTPCCorrection.cxx.
Referenced by Interpolate2DEdistortion(), Interpolate2DTable(), Interpolate3DEdistortion(), and Interpolate3DTable().
|
protected |
Interpolate function Y(x) using linear (order=1) or quadratic (order=2) interpolation. Float version (in order to decrease the OCDB size)
Definition at line 1050 of file AliTPCCorrection.cxx.
|
protected |
Interpolate table - 2D interpolation
Definition at line 850 of file AliTPCCorrection.cxx.
Referenced by AliTPCSpaceCharge::GetCorrection(), AliTPCBoundaryVoltError::GetCorrection(), AliTPCGGVoltError::GetCorrection(), and AliTPCGGVoltError::GetIntErOverEz().
|
protected |
Interpolate table (TMatrix format) - 2D interpolation
Definition at line 912 of file AliTPCCorrection.cxx.
Referenced by AliTPCSpaceCharge2D2D::GetSpaceChargeDensity(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), and AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion().
|
protected |
Interpolate table (TMatrix format) - 2D interpolation Float version (in order to decrease the OCDB size)
Definition at line 989 of file AliTPCCorrection.cxx.
|
protected |
Interpolate table - 3D interpolation
Definition at line 870 of file AliTPCCorrection.cxx.
|
protected |
Interpolate table (TMatrix format) - 3D interpolation
Definition at line 937 of file AliTPCCorrection.cxx.
Referenced by AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCROCVoltError3D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), AliTPCFCVoltError3D::GetCorrection(), AliTPCCorrectionLookupTable::GetInterpolation(), AliTPCSpaceCharge2D2D::GetSpaceChargeDensity(), AliTPCSpaceCharge3D::GetSpaceChargeDensity(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), and AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson().
|
protected |
Interpolate table (TMatrix format) - 3D interpolation Float version (in order to decrease the OCDB size)
Definition at line 1015 of file AliTPCCorrection.cxx.
|
inlineprotected |
Definition at line 186 of file AliTPCCorrection.h.
|
protectedvirtual |
Helperfunction: Check if integer is a power of 2
Definition at line 1810 of file AliTPCCorrection.cxx.
Referenced by PoissonRelaxation2D(), and PoissonRelaxation3D().
|
static |
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
Definition at line 2704 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and MakeLookup().
|
static |
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 maria n.iv anov@ cern .ch
Histo axeses Collection name='TObjArray', class='TObjArray', size=16 0. OBJ: TAxis #Delta #Delta
Definition at line 2796 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
|
static |
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 maria n.iv anov@ cern .ch
Definition at line 2907 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
|
static |
Make a laser fit tree for global minimization
Definition at line 3482 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and MakeLaserTree().
|
static |
Make a laser fit tree for global minimization
Definition at line 2509 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
|
static |
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
Definition at line 2337 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor().
|
static |
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
Definition at line 2129 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and MakeFitTree().
|
private |
|
protected |
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)
Definition at line 1163 of file AliTPCCorrection.cxx.
Referenced by AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), and AliTPCSpaceCharge::InitSpaceChargeDistortion().
|
protected |
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
Definition at line 1398 of file AliTPCCorrection.cxx.
Referenced by AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCROCVoltError3D::InitROCVoltError3D(), and AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson().
|
virtual |
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
Reimplemented in AliTPCComposedCorrection, AliTPCExBBShape, AliTPCCalibGlobalMisalignment, and AliTPCInverseCorrection.
Definition at line 520 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), AliTPCInverseCorrection::Print(), and AliTPCComposedCorrection::Print().
|
protected |
Search an ordered table by starting at the most recently used point
Definition at line 1070 of file AliTPCCorrection.cxx.
Referenced by AliTPCCorrectionLookupTable::BuildExactInverse(), AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), Interpolate2DEdistortion(), Interpolate2DTable(), Interpolate3DEdistortion(), and Interpolate3DTable().
|
inlinevirtual |
Reimplemented in AliTPCCorrectionLookupTable.
Definition at line 59 of file AliTPCCorrection.h.
|
inlineprotected |
Definition at line 185 of file AliTPCCorrection.h.
|
virtual |
Virtual funtion to pass the wt values (might become event dependent) to the inherited classes t1 and t2 represent the "effective omegaTau" corrections and were measured in a dedicated calibration run
Reimplemented in AliTPCSpaceCharge3DDriftLine, AliTPCComposedCorrection, AliTPCGGVoltError, AliTPCExBTwist, AliTPCInverseCorrection, AliTPCExBBShape, AliTPCBoundaryVoltError, AliTPCSpaceCharge, AliTPCExBConical, AliTPCSpaceCharge3D, AliTPCExBEffectiveSector, AliTPCFCVoltError3D, AliTPCROCVoltError3D, AliTPCSpaceCharge2D2D, and AliTPCExBEffective.
Definition at line 528 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), MakeComposedCorrection(), RegisterAliTPCBoundaryVoltError(), RegisterAliTPCROCVoltError3D(), RegisterAliTPCROCVoltError3DSector(), AliTPCInverseCorrection::SetOmegaTauT1T2(), AliTPCComposedCorrection::SetOmegaTauT1T2(), TestCorrection_AliTPCComposedCorrectionAddCorrectionCompact_TPCCalibCorrection(), and TestParExample().
void AliTPCCorrection::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
Definition at line 3107 of file AliTPCCorrection.cxx.
Referenced by GetCorrScaleFactor(), and LoadModels().
|
virtual |
Update function
Reimplemented in AliTPCComposedCorrection, AliTPCGGVoltError, AliTPCInverseCorrection, AliTPCExBTwist, AliTPCExBBShape, AliTPCBoundaryVoltError, AliTPCExBEffectiveSector, AliTPCSpaceCharge3D, AliTPCExBConical, AliTPCFCVoltError3D, AliTPCSpaceCharge, AliTPCCorrectionDrift, AliTPCROCVoltError3D, AliTPCExBEffective, and AliTPCSpaceCharge2D2D.
Definition at line 515 of file AliTPCCorrection.cxx.
Referenced by AliTPCInverseCorrection::Update(), and AliTPCComposedCorrection::Update().
|
staticprotected |
Cathode Voltage (volts)
Definition at line 123 of file AliTPCCorrection.h.
Referenced by AliTPCROCVoltError3D::CreateHistoOfZAlignment(), AliTPCROCVoltError3D::GetROCVoltOffset(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), PoissonRelaxation2D(), and PoissonRelaxation3D().
|
staticprotected |
[cm/V] drift velocity dependency on the E field (from Magboltz for NeCO2N2 at standard environment)
Definition at line 125 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), PoissonRelaxation2D(), and PoissonRelaxation3D().
|
staticprotected |
vacuum permittivity [A·s/(V·m)]
Definition at line 127 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge3DDriftLine::CreateHistogramSCInXY(), AliTPCSpaceCharge3DDriftLine::CreateHistogramSCInZR(), AliTPCSpaceCharge2D2D::CreateHistoSCinXY(), AliTPCSpaceCharge3D::CreateHistoSCinXY(), AliTPCSpaceCharge2D2D::CreateHistoSCinZR(), AliTPCSpaceCharge3D::CreateHistoSCinZR(), AliTPCSpaceCharge2D2D::WriteChargeDistributionToFile(), and AliTPCSpaceCharge3D::WriteChargeDistributionToFile().
|
staticprotected |
charge/mass in [C/kg]
Definition at line 126 of file AliTPCCorrection.h.
|
staticprotected |
Gating Grid voltage (volts)
Definition at line 124 of file AliTPCCorrection.h.
Referenced by AliTPCROCVoltError3D::CreateHistoOfZAlignment(), AliTPCROCVoltError3D::GetROCVoltOffset(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), PoissonRelaxation2D(), and PoissonRelaxation3D().
|
staticprotected |
Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)/hera/alice/wiechula/calib/guiTrees.
radius which renders the "18 rod manifold" best -> compare calc. of Jim Thomas
Definition at line 120 of file AliTPCCorrection.h.
Referenced by AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), PoissonRelaxation2D(), PoissonRelaxation3D(), AliTPCSpaceCharge2D2D::WriteChargeDistributionToFile(), and AliTPCSpaceCharge3D::WriteChargeDistributionToFile().
|
staticprotected |
Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
Definition at line 121 of file AliTPCCorrection.h.
Referenced by AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), PoissonRelaxation2D(), PoissonRelaxation3D(), AliTPCSpaceCharge2D2D::WriteChargeDistributionToFile(), and AliTPCSpaceCharge3D::WriteChargeDistributionToFile().
|
protected |
points in the phi direction (for the lookup table)
Definition at line 135 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCROCVoltError3D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), AliTPCFCVoltError3D::GetCorrection(), AliTPCSpaceCharge2D2D::GetSpaceChargeDensity(), AliTPCSpaceCharge3D::GetSpaceChargeDensity(), AliTPCFCVoltError3D::InitFCVoltError3D(), InitLookUpfulcrums(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), Interpolate3DEdistortion(), AliTPCSpaceCharge3D::SetInputSpaceCharge(), AliTPCSpaceCharge2D2D::SetSCDataFileName(), and AliTPCCorrectionLookupTable::SetupDefaultLimits().
|
protected |
points in the radial direction (for the lookup table)
Definition at line 134 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCROCVoltError3D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), AliTPCFCVoltError3D::GetCorrection(), AliTPCSpaceCharge2D2D::GetSpaceChargeDensity(), AliTPCSpaceCharge3D::GetSpaceChargeDensity(), AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), InitLookUpfulcrums(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), Interpolate2DEdistortion(), Interpolate3DEdistortion(), AliTPCSpaceCharge3D::SetInputSpaceCharge(), AliTPCSpaceCharge2D2D::SetSCDataFileName(), and AliTPCCorrectionLookupTable::SetupDefaultLimits().
|
staticprotected |
nominal gating grid position
Definition at line 119 of file AliTPCCorrection.h.
Referenced by AliTPCROCVoltError3D::CreateHistoOfZAlignment(), AliTPCExBBShape::GetBxAndByOverBz(), AliTPCExBBShape::GetCorrection(), AliTPCExBTwist::GetCorrection(), GetCorrectionIntegralDz(), GetDistortionIntegralDz(), AliTPCROCVoltError3D::GetROCVoltOffset(), AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), PoissonRelaxation2D(), and PoissonRelaxation3D().
|
protected |
points in the z direction (for the lookup table)
Definition at line 136 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCROCVoltError3D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), AliTPCFCVoltError3D::GetCorrection(), AliTPCSpaceCharge2D2D::GetSpaceChargeDensity(), AliTPCSpaceCharge3D::GetSpaceChargeDensity(), AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion(), AliTPCFCVoltError3D::InitFCVoltError3D(), AliTPCGGVoltError::InitGGVoltErrorDistortion(), InitLookUpfulcrums(), AliTPCROCVoltError3D::InitROCVoltError3D(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortion(), AliTPCSpaceCharge3D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge2D2D::InitSpaceCharge3DDistortionCourse(), AliTPCSpaceCharge3D::InitSpaceCharge3DPoisson(), AliTPCSpaceCharge::InitSpaceChargeDistortion(), Interpolate2DEdistortion(), Interpolate3DEdistortion(), AliTPCSpaceCharge3D::SetInputSpaceCharge(), AliTPCSpaceCharge2D2D::SetSCDataFileName(), and AliTPCCorrectionLookupTable::SetupDefaultLimits().
|
staticprotected |
Offset from CE: calculate all distortions closer to CE as if at this point.
Definition at line 122 of file AliTPCCorrection.h.
Referenced by AliTPCSpaceCharge::GetCorrection(), AliTPCSpaceCharge2D2D::GetCorrection(), AliTPCBoundaryVoltError::GetCorrection(), AliTPCROCVoltError3D::GetCorrection(), AliTPCSpaceCharge3D::GetCorrection(), AliTPCGGVoltError::GetCorrection(), AliTPCFCVoltError3D::GetCorrection(), AliTPCGGVoltError::GetIntErOverEz(), and AliTPCCorrectionLookupTable::GetInterpolation().
|
staticprotected |
array of orrection for visualization
Definition at line 191 of file AliTPCCorrection.h.
Referenced by AddVisualCorrection(), AliTPCCorrection(), GetCorrectionSector(), GetCorrSector(), GetCorrXYZ(), GetCorrXYZDz(), GetCorrXYZIntegrateZ(), GetDistortionSector(), GetDistXYZ(), GetDistXYZDz(), GetDistXYZIntegrateZ(), GetVisualCorrection(), and GetVisualCorrections().
|
protected |
Definition at line 139 of file AliTPCCorrection.h.
Referenced by Interpolate3DEdistortion().
|
protected |
Presentation of the underlying corrections, integrated, or differential.
Definition at line 192 of file AliTPCCorrection.h.
Referenced by GetCorrectionDz(), GetDistortionDz(), and PoissonRelaxation3D().
|
protected |
switch to indicate that the distortion is a local vector drphi/dz, dr/dz
Definition at line 190 of file AliTPCCorrection.h.
Referenced by IsLocal(), and SetIsLocal().
|
protected |
Definition at line 139 of file AliTPCCorrection.h.
Referenced by Interpolate2DEdistortion(), and Interpolate3DEdistortion().
|
protected |
variable to help in the interpolation
Definition at line 139 of file AliTPCCorrection.h.
Referenced by Interpolate2DEdistortion(), and Interpolate3DEdistortion().
|
protected |
tensor term of wt - T1
Definition at line 188 of file AliTPCCorrection.h.
Referenced by AliTPCExBEffective::Init(), AliTPCSpaceCharge2D2D::Init(), AliTPCROCVoltError3D::Init(), AliTPCSpaceCharge::Init(), AliTPCExBConical::Init(), AliTPCFCVoltError3D::Init(), AliTPCExBEffectiveSector::Init(), AliTPCSpaceCharge3D::Init(), AliTPCBoundaryVoltError::Init(), AliTPCExBBShape::Init(), AliTPCExBTwist::Init(), AliTPCGGVoltError::Init(), AliTPCSpaceCharge3DDriftLine::Init(), AliTPCExBEffective::Print(), AliTPCExBEffectiveSector::Print(), AliTPCSpaceCharge::Print(), AliTPCExBConical::Print(), AliTPCSpaceCharge2D2D::Print(), AliTPCBoundaryVoltError::Print(), AliTPCExBBShape::Print(), AliTPCSpaceCharge3D::Print(), AliTPCExBTwist::Print(), AliTPCROCVoltError3D::Print(), AliTPCGGVoltError::Print(), AliTPCFCVoltError3D::Print(), AliTPCExBEffective::SetOmegaTauT1T2(), AliTPCSpaceCharge2D2D::SetOmegaTauT1T2(), AliTPCROCVoltError3D::SetOmegaTauT1T2(), AliTPCExBEffectiveSector::SetOmegaTauT1T2(), AliTPCFCVoltError3D::SetOmegaTauT1T2(), AliTPCSpaceCharge3D::SetOmegaTauT1T2(), AliTPCExBConical::SetOmegaTauT1T2(), AliTPCBoundaryVoltError::SetOmegaTauT1T2(), AliTPCSpaceCharge::SetOmegaTauT1T2(), AliTPCExBBShape::SetOmegaTauT1T2(), AliTPCExBTwist::SetOmegaTauT1T2(), AliTPCGGVoltError::SetOmegaTauT1T2(), SetOmegaTauT1T2(), AliTPCSpaceCharge3DDriftLine::SetOmegaTauT1T2(), AliTPCExBEffective::Update(), AliTPCSpaceCharge2D2D::Update(), AliTPCROCVoltError3D::Update(), AliTPCExBConical::Update(), AliTPCFCVoltError3D::Update(), AliTPCSpaceCharge::Update(), AliTPCSpaceCharge3D::Update(), AliTPCExBEffectiveSector::Update(), AliTPCBoundaryVoltError::Update(), AliTPCExBBShape::Update(), AliTPCExBTwist::Update(), and AliTPCGGVoltError::Update().
|
protected |
tensor term of wt - T2
Definition at line 189 of file AliTPCCorrection.h.
Referenced by AliTPCExBEffective::Init(), AliTPCSpaceCharge2D2D::Init(), AliTPCROCVoltError3D::Init(), AliTPCSpaceCharge::Init(), AliTPCExBConical::Init(), AliTPCFCVoltError3D::Init(), AliTPCExBEffectiveSector::Init(), AliTPCSpaceCharge3D::Init(), AliTPCBoundaryVoltError::Init(), AliTPCExBBShape::Init(), AliTPCExBTwist::Init(), AliTPCGGVoltError::Init(), AliTPCSpaceCharge3DDriftLine::Init(), AliTPCExBEffective::Print(), AliTPCExBEffectiveSector::Print(), AliTPCSpaceCharge::Print(), AliTPCExBConical::Print(), AliTPCSpaceCharge2D2D::Print(), AliTPCBoundaryVoltError::Print(), AliTPCExBBShape::Print(), AliTPCSpaceCharge3D::Print(), AliTPCExBTwist::Print(), AliTPCROCVoltError3D::Print(), AliTPCGGVoltError::Print(), AliTPCFCVoltError3D::Print(), AliTPCExBEffective::SetOmegaTauT1T2(), AliTPCSpaceCharge2D2D::SetOmegaTauT1T2(), AliTPCROCVoltError3D::SetOmegaTauT1T2(), AliTPCExBEffectiveSector::SetOmegaTauT1T2(), AliTPCFCVoltError3D::SetOmegaTauT1T2(), AliTPCSpaceCharge3D::SetOmegaTauT1T2(), AliTPCExBConical::SetOmegaTauT1T2(), AliTPCBoundaryVoltError::SetOmegaTauT1T2(), AliTPCSpaceCharge::SetOmegaTauT1T2(), AliTPCExBBShape::SetOmegaTauT1T2(), AliTPCExBTwist::SetOmegaTauT1T2(), AliTPCGGVoltError::SetOmegaTauT1T2(), SetOmegaTauT1T2(), AliTPCSpaceCharge3DDriftLine::SetOmegaTauT1T2(), AliTPCExBEffective::Update(), AliTPCSpaceCharge2D2D::Update(), AliTPCROCVoltError3D::Update(), AliTPCExBConical::Update(), AliTPCFCVoltError3D::Update(), AliTPCSpaceCharge::Update(), AliTPCSpaceCharge3D::Update(), AliTPCExBEffectiveSector::Update(), AliTPCBoundaryVoltError::Update(), AliTPCExBBShape::Update(), AliTPCExBTwist::Update(), and AliTPCGGVoltError::Update().