![]() |
Project
|
#include <SpaceCharge.h>
Public Types | |
enum class | IntegrationStrategy { Trapezoidal = 0 , Simpson = 1 , Root = 2 , SimpsonIterative = 3 } |
numerical integration strategys More... | |
enum class | Type { Distortions = 0 , Corrections = 1 } |
enum class | GlobalDistType { Standard = 0 , Fast = 1 , None = 2 } |
enum class | GlobalDistCorrMethod { LocalDistCorr , ElectricalField } |
enum class | MisalignmentType { ShiftedClip , RotatedClip , RodShift } |
enum class | FCType { IFC , OFC } |
Public Member Functions | |
SpaceCharge (const int bfield, const unsigned short nZVertices=129, const unsigned short nRVertices=129, const unsigned short nPhiVertices=360, const bool initBuffers=false) | |
SpaceCharge () | |
constructor for empty object. Can be used when buffers are loaded from file | |
SpaceCharge (SpaceCharge &&)=default | |
default move constructor | |
SpaceCharge & | operator= (SpaceCharge &&)=default |
move assignment | |
void | fillChargeDensityFromHisto (const TH3 &hisSCDensity3D) |
void | fillChargeDensityFromHisto (const TH3 &hisSCDensity3D_A, const TH3 &hisSCDensity3D_C) |
void | fillChargeDensityFromHisto (const char *file, const char *nameA, const char *nameC) |
void | fillChargeDensityFromCalDet (const std::vector< CalDet< float > > &calSCDensity3D) |
void | fillChargeFromIDCs (std::vector< CalDet< float > > &idcZero, const CalDet< float > &mapIBF, const float ionDriftTimeMS=200, const bool normToPadArea=true) |
void | fillChargeFromCalDet (const std::vector< CalDet< float > > &calCharge3D) |
void | fillChargeDensityFromFile (TFile &fInp, const char *name) |
void | calculateDistortionsCorrections (const o2::tpc::Side side, const bool calcVectors=false) |
void | setChargeDensityFromFormula (const AnalyticalFields< DataT > &formulaStruct) |
void | setPotentialBoundaryFromFormula (const AnalyticalFields< DataT > &formulaStruct) |
void | addBoundaryPotential (const SpaceCharge< DataT > &other, const Side side, const float scaling=1) |
void | resetBoundaryPotentialToZeroInRangeZ (float zMin, float zMax, const Side side) |
setting the boundary potential to 0 for z<zMin and z>zMax | |
void | setPotentialFromFormula (const AnalyticalFields< DataT > &formulaStruct) |
void | mirrorPotential (const Side sideRef, const Side sideMirrored) |
void | setSimOneSector () |
void | setDefaultStaticDistortionsGEMFrameChargeUp (const Side side, const DataT deltaPotential=1000) |
void | setPotentialBoundaryGEMFrameAlongR (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameIROCBottomAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameIROCTopAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC1BottomAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC1TopAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC2BottomAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC2TopAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC3BottomAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC3TopAlongPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameOROC3ToOFCPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryGEMFrameIROCToIFCPhi (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
setting the potential from the IROC to the inner field cage | |
void | setPotentialBoundaryInnerRadius (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | setPotentialBoundaryOuterRadius (const std::function< DataT(DataT)> &potentialFunc, const Side side) |
void | poissonSolver (const Side side, const DataT stoppingConvergence=1e-6, const int symmetry=0) |
void | poissonSolver (const DataT stoppingConvergence=1e-6, const int symmetry=0) |
void | calcEField (const Side side) |
void | setEFieldFromFormula (const AnalyticalFields< DataT > &formulaStruct) |
void | scaleChargeDensity (const DataT scalingFactor, const Side side) |
void | scaleChargeDensitySector (const float scalingFactor, const Sector sector) |
void | scaleChargeDensityStack (const float scalingFactor, const Sector sector, const GEMstack stack) |
scaling the space-charge density for given stack | |
void | addChargeDensity (const SpaceCharge< DataT > &otherSC) |
template<typename ElectricFields = AnalyticalFields<DataT>> | |
void | calcLocalDistortionsCorrections (const Type type, const ElectricFields &formulaStruct) |
template<typename ElectricFields = AnalyticalFields<DataT>> | |
void | calcLocalDistortionCorrectionVector (const ElectricFields &formulaStruct) |
template<typename ElectricFields = AnalyticalFields<DataT>> | |
void | calcLocalDistortionsCorrectionsRK4 (const Type type, const Side side) |
template<typename Fields = AnalyticalFields<DataT>> | |
void | calcGlobalCorrections (const Fields &formulaStruct, const int type=3) |
void | calcGlobalCorrectionsEField (const Side side, const int type=3) |
calculate the global corrections using the electric fields (interface for python) | |
template<typename Fields = AnalyticalFields<DataT>> | |
void | calcGlobalDistortions (const Fields &formulaStruct, const int maxIterations=3 *sSteps *129) |
void | init () |
void | calcGlobalDistWithGlobalCorrIterative (const DistCorrInterpolator< DataT > &globCorr, const int maxIter=100, const DataT approachZ=1, const DataT approachR=1, const DataT approachPhi=1, const DataT diffCorr=50e-6, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0) |
void | calcGlobalDistWithGlobalCorrIterative (const Side side, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachZ=1, const DataT approachR=1, const DataT approachPhi=1, const DataT diffCorr=50e-6) |
void | calcGlobalDistWithGlobalCorrIterative (const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachZ=1, const DataT approachR=1, const DataT approachPhi=1, const DataT diffCorr=50e-6) |
void | calcGlobalCorrWithGlobalDistIterative (const Side side, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachZ=1, const DataT approachR=1, const DataT approachPhi=1, const DataT diffCorr=50e-6) |
void | calcGlobalCorrWithGlobalDistIterative (const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachZ=1, const DataT approachR=1, const DataT approachPhi=1, const DataT diffCorr=50e-6) |
void | calcGlobalDistWithGlobalCorrIterativeLinearCartesian (const Side side, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachX=1, const DataT approachY=1, const DataT approachZ=1, const DataT diffCorr=50e-6) |
void | calcGlobalDistWithGlobalCorrIterativeLinearCartesian (const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachX=1, const DataT approachY=1, const DataT approachZ=1, const DataT diffCorr=50e-6) |
void | calcGlobalCorrWithGlobalDistIterativeLinearCartesian (const Side side, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachX=1, const DataT approachY=1, const DataT approachZ=1, const DataT diffCorr=50e-6) |
void | calcGlobalCorrWithGlobalDistIterativeLinearCartesian (const SpaceCharge< DataT > *scSCale=nullptr, float scale=0, const int maxIter=100, const DataT approachX=1, const DataT approachY=1, const DataT approachZ=1, const DataT diffCorr=50e-6) |
unsigned short | getNZVertices () const |
unsigned short | getNRVertices () const |
unsigned short | getNPhiVertices () const |
auto | getC0 () const |
auto | getC1 () const |
auto | getC2 () const |
int | getBField () const |
const auto & | getPotential (const Side side) const & |
void | setPotential (int iz, int ir, int iphi, Side side, float val) |
setting the potential directly for given vertex | |
DataT | getDensityCyl (const DataT z, const DataT r, const DataT phi, const Side side) const |
DataT | getPotentialCyl (const DataT z, const DataT r, const DataT phi, const Side side) const |
std::vector< float > | getPotentialCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side) const |
get the potential for list of given coordinate | |
void | getElectricFieldsCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &eZ, DataT &eR, DataT &ePhi) const |
void | getLocalCorrectionsCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &lcorrZ, DataT &lcorrR, DataT &lcorrRPhi) const |
void | getLocalCorrectionsCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &lcorrZ, std::vector< DataT > &lcorrR, std::vector< DataT > &lcorrRPhi) const |
void | getCorrectionsCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &corrZ, DataT &corrR, DataT &corrRPhi) const |
void | getCorrectionsCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &corrZ, std::vector< DataT > &corrR, std::vector< DataT > &corrRPhi) const |
void | getCorrections (const DataT x, const DataT y, const DataT z, const Side side, DataT &corrX, DataT &corrY, DataT &corrZ) const |
void | getCorrectionsAnalytical (const DataT x, const DataT y, const DataT z, const Side side, DataT &corrX, DataT &corrY, DataT &corrZ) const |
void | getLocalDistortionsCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &ldistZ, DataT &ldistR, DataT &ldistRPhi) const |
void | getLocalDistortionsCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &ldistZ, std::vector< DataT > &ldistR, std::vector< DataT > &ldistRPhi) const |
void | getLocalDistortionVectorCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &lvecdistZ, DataT &lvecdistR, DataT &lvecdistRPhi) const |
void | getLocalDistortionVectorCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &lvecdistZ, std::vector< DataT > &lvecdistR, std::vector< DataT > &lvecdistRPhi) const |
void | getLocalCorrectionVectorCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &lveccorrZ, DataT &lveccorrR, DataT &lveccorrRPhi) const |
void | getLocalCorrectionVectorCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &lveccorrZ, std::vector< DataT > &lveccorrR, std::vector< DataT > &lveccorrRPhi) const |
void | getDistortionsCyl (const DataT z, const DataT r, const DataT phi, const Side side, DataT &distZ, DataT &distR, DataT &distRPhi) const |
void | getDistortionsCyl (const std::vector< DataT > &z, const std::vector< DataT > &r, const std::vector< DataT > &phi, const Side side, std::vector< DataT > &distZ, std::vector< DataT > &distR, std::vector< DataT > &distRPhi) const |
void | getDistortions (const DataT x, const DataT y, const DataT z, const Side side, DataT &distX, DataT &distY, DataT &distZ) const |
void | getDistortionsAnalytical (const DataT x, const DataT y, const DataT z, const Side side, DataT &distX, DataT &distY, DataT &distZ) const |
void | setDistortionsCorrectionsAnalytical (const AnalyticalDistCorr< DataT > &formula) |
set distortions and corrections by an analytical formula | |
const auto & | getDistortionsCorrectionsAnalytical () const & |
void | setUseAnalyticalDistCorr (const bool useAna) |
setting usage of the analytical formula for the distortions and corrections | |
bool | getUseAnalyticalDistCorr () const |
void | correctElectron (GlobalPosition3D &point) |
void | distortElectron (GlobalPosition3D &point, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0) const |
void | setDistortionLookupTables (const DataContainer &distdZ, const DataContainer &distdR, const DataContainer &distdRPhi, const Side side) |
void | setFromFile (std::string_view file, const Side side) |
void | setFromFile (std::string_view file) |
DataT | getGridSpacingR (const Side side) const |
Get grid spacing in r direction. | |
DataT | getGridSpacingZ (const Side side) const |
Get grid spacing in z direction. | |
DataT | getGridSpacingPhi (const Side side) const |
Get grid spacing in phi direction. | |
DataT | getRMin (const Side side) const |
Get inner radius of tpc. | |
DataT | getZMin (const Side side) const |
Get min z position which is used during the calaculations. | |
DataT | getPhiMin (const Side side) const |
Get min phi. | |
DataT | getRMax (const Side side) const |
Get max r. | |
DataT | getZMax (const Side side) const |
Get max z. | |
DataT | getPhiMax (const Side side) const |
Get max phi. | |
const RegularGrid & | getGrid3D (const Side side) const |
Get the grid object. | |
NumericalFields< DataT > | getElectricFieldsInterpolator (const Side side) const |
DistCorrInterpolator< DataT > | getLocalDistInterpolator (const Side side) const |
DistCorrInterpolator< DataT > | getLocalCorrInterpolator (const Side side) const |
DistCorrInterpolator< DataT > | getGlobalDistInterpolator (const Side side) const |
DistCorrInterpolator< DataT > | getGlobalCorrInterpolator (const Side side) const |
DataT | getLocalDistR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalDistZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalDistRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecDistR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecDistZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecDistRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalCorrR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalCorrZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalCorrRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecCorrR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecCorrZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getLocalVecCorrRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalDistR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalDistZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalDistRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalCorrR (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalCorrZ (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getGlobalCorrRPhi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getEr (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getEz (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getEphi (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getDensity (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getPotential (const size_t iz, const size_t ir, const size_t iphi, const Side side) const |
DataT | getPhiVertex (const size_t indexPhi, const Side side) const |
DataT | getRVertex (const size_t indexR, const Side side) const |
DataT | getZVertex (const size_t indexZ, const Side side) const |
void | setOmegaTauT1T2 (const DataT omegaTau=0.32f, const DataT t1=1, const DataT t2=1) |
void | setC0C1 (const DataT c0, const DataT c1) |
void | setC2 (const DataT c2) |
void | initBField (const int field) |
void | setSimExBMisalignment (const bool simExBMisalignment) |
enable/disable calculation of distortions due to ExB misalignment | |
bool | getSimExBMisalignment () const |
void | setSimEDistortions (const bool simEDistortions) |
calculate distortions due to electric fields (space charge, boundary potential...) | |
bool | getSimEDistortions () const |
void | setUseInitialSCDensity (const bool useInitialSCDensity) |
void | dumpToFile (std::string_view file, const Side side, std::string_view option) const |
void | dumpToFile (std::string_view file) const |
void | dumpMetaData (std::string_view file, std::string_view option, const bool overwriteExisting=false) const |
void | readMetaData (std::string_view file) |
void | dumpToTree (const char *outFileName, const Side side, const int nZPoints=50, const int nRPoints=150, const int nPhiPoints=180, const bool randomize=false) const |
void | dumpToTree (const char *outFileName, const o2::tpc::Sector §or, const int nZPoints=50) const |
int | dumpElectricFields (std::string_view file, const Side side, std::string_view option) const |
void | setElectricFieldsFromFile (std::string_view file, const Side side) |
int | dumpPotential (std::string_view file, const Side side, std::string_view option) const |
void | setPotentialFromFile (std::string_view file, const Side side) |
void | fillPotential (const DataT potential, const size_t iz, const size_t ir, const size_t iphi, const Side side) |
int | dumpDensity (std::string_view file, const Side side, std::string_view option) const |
void | setDensityFromFile (std::string_view file, const Side side) |
void | fillDensity (const DataT density, const size_t iz, const size_t ir, const size_t iphi, const Side side) |
void | averageDensityPerSector (const Side side) |
average the sc density across all sectors per side | |
int | dumpGlobalDistortions (std::string_view file, const Side side, std::string_view option) const |
int | dumpAnalyticalCorrectionsDistortions (TFile &outf) const |
void | setAnalyticalCorrectionsDistortionsFromFile (std::string_view inpf) |
void | setGlobalDistortionsFromFile (std::string_view file, const Side side) |
template<typename DataTIn = DataT> | |
void | setGlobalDistortionsFromFile (TFile &inpf, const Side side) |
int | dumpGlobalCorrections (std::string_view file, const Side side, std::string_view option) const |
int | dumpGlobalCorrections (TFile &outf, const Side side) const |
void | setGlobalCorrectionsFromFile (std::string_view file, const Side side) |
template<typename DataTIn = DataT> | |
void | setGlobalCorrectionsFromFile (TFile &inpf, const Side side) |
int | dumpLocalCorrections (std::string_view file, const Side side, std::string_view option) const |
void | setLocalCorrectionsFromFile (std::string_view file, const Side side) |
int | dumpLocalDistortions (std::string_view file, const Side side, std::string_view option) const |
int | dumpLocalDistCorrVectors (std::string_view file, const Side side, std::string_view option) const |
void | setLocalDistortionsFromFile (std::string_view file, const Side side) |
void | setLocalDistCorrVectorsFromFile (std::string_view file, const Side side) |
DataT | regulateZ (const DataT posZ, const Side side) const |
DataT | regulateR (const DataT posR, const Side side) const |
set r coordinate between 'RMIN - 4 * GRIDSPACINGR' and 'RMAX + 2 * GRIDSPACINGR'. the r coordinate is not clamped to RMIN and RMAX to ensure correct interpolation at the borders of the grid. | |
DataT | regulatePhi (const DataT posPhi, const Side side) const |
set phi coordinate between min phi max phi | |
std::vector< std::pair< std::vector< o2::math_utils::Point3D< float > >, std::array< DataT, 3 > > > | calculateElectronDriftPath (const std::vector< GlobalPosition3D > &elePos, const int nSamplingPoints, const std::string_view outFile="electron_tracks.root") const |
bool | checkGridFromFile (std::string_view file, std::string_view tree) |
void | setDeltaVoltageCopperRodShiftIFC (const float deltaPot, const Side side, const int sector) |
void | setDeltaVoltageCopperRodShiftOFC (const float deltaPot, const Side side, const int sector) |
void | setDeltaVoltageStripsShiftIFC (const float deltaPot, const Side side, const int sector) |
void | setDeltaVoltageStripsShiftOFC (const float deltaPot, const Side side, const int sector) |
void | setDeltaVoltageRotatedClipIFC (const float deltaPot, const Side side, const int sector) |
void | setDeltaVoltageRotatedClipOFC (const float deltaPot, const Side side, const int sector) |
void | setIFCChargeUpRisingPot (const float deltaPot, const float zMaxDeltaPot, const int type, const float zStart, const float offs, const Side side) |
void | setIFCChargeUpFallingPot (const float deltaPot, const float zMaxDeltaPot, const int type, const float zEnd, const float offs, const Side side) |
void | setGlobalCorrections (const std::function< void(int sector, DataT gx, DataT gy, DataT gz, DataT &gCx, DataT &gCy, DataT &gCz)> &gCorr, const Side side) |
void | setROCMisalignmentShiftZ (const int sector, const int type, const float potential) |
void | setROCMisalignmentRotationAlongX (const int sector, const int type, const float potentialMin, const float potentialMax) |
void | setROCMisalignmentRotationAlongY (const int sector, const int type, const float potentialMin, const float potentialMax) |
void | subtractGlobalCorrections (const SpaceCharge< DataT > &otherSC, const Side side) |
void | subtractGlobalDistortions (const SpaceCharge< DataT > &otherSC, const Side side) |
void | scaleCorrections (const float scaleFac, const Side side) |
void | setMetaData (const SCMetaData &meta) |
setting meta data for this object | |
const auto & | getMetaData () const |
void | printMetaData () const |
float | getMeanLumi () const |
void | setMeanLumi (float lumi) |
void | initAfterReadingFromFile () |
float | getDCAr (float tgl, const int nPoints, const float phi, o2::utils::TreeStreamRedirector *pcstream=nullptr) const |
template<typename ElectricFields > | |
void | calcLocalDistortionsCorrections (const SpaceCharge< DataT >::Type type, const ElectricFields &formulaStruct) |
template<typename ElectricFields > | |
void | calcLocalDistortionsCorrectionsRK4 (const SpaceCharge< DataT >::Type type, const Side side) |
template<typename Formulas > | |
void | calcGlobalCorrections (const Formulas &formulaStruct, const int type) |
Static Public Member Functions | |
static void | convertIDCsToCharge (std::vector< CalDet< float > > &idcZero, const CalDet< float > &mapIBF, const float ionDriftTimeMS=200, const bool normToPadArea=true) |
static void | unsetSimOneSector () |
unsetting simulation of one sector | |
static DataT | getRadiusFromCartesian (const DataT x, const DataT y) |
convert x and y coordinates from cartesian to the radius in polar coordinates | |
static DataT | getPhiFromCartesian (const DataT x, const DataT y) |
convert x and y coordinates from cartesian to phi in polar coordinates | |
static DataT | getXFromPolar (const DataT r, const DataT phi) |
convert radius and phi coordinates from polar coordinates to x cartesian coordinates | |
static DataT | getYFromPolar (const DataT r, const DataT phi) |
convert radius and phi coordinates from polar coordinates to y cartesian coordinate | |
static constexpr DataT | getEzField (const Side side) |
Get constant electric field. | |
static Side | getSide (const DataT z) |
static int | getStepWidth () |
Get the step width which is used for the calculation of the correction/distortions in units of the z-bin. | |
static void | setNStep (const int nSteps) |
static int | getNStep () |
static int | getNThreads () |
get the number of threads used for some of the calculations | |
static void | setNThreads (const int nThreads) |
set the number of threads used for some of the calculations | |
static void | setNumericalIntegrationStrategy (const IntegrationStrategy strategy) |
static IntegrationStrategy | getNumericalIntegrationStrategy () |
static void | setGlobalDistType (const GlobalDistType globalDistType) |
static GlobalDistType | getGlobalDistType () |
static void | setGlobalDistCorrMethod (const GlobalDistCorrMethod globalDistCorrMethod) |
static GlobalDistCorrMethod | getGlobalDistCorrMethod () |
static void | setSimpsonNIteratives (const int nIter) |
static int | getSimpsonNIteratives () |
static void | setSCDistortionType (SCDistortionType distortionType) |
static SCDistortionType | getSCDistortionType () |
Get the space-charge distortions model. | |
static void | makeElectronDriftPathGif (const char *inpFile, TH2F &hBorder, const int type=0, const int gifSpeed=2, const int maxsamplingpoints=100, const char *outName="electron_drift_path") |
static void | normalizeHistoQVEps0 (TH3 &histoIonsPhiRZ) |
static int | getOMPMaxThreads () |
DataT | the data type which is used during the calculations |
Nz | number of vertices in z direction |
Nr | number of vertices in r direction |
Nphi | number of vertices in phi direction |
this class provides the algorithms for calculating the global distortions and corrections from the space charge density. The calculations can be done by a realistic space charge histogram as an input or by an analytical formula. An example of of the usage can be found in 'macro/calculateDistortionsCorrections.C'
Definition at line 72 of file SpaceCharge.h.
|
strong |
Enumerator | |
---|---|
IFC | inner field cage |
OFC | outer field cage |
Definition at line 128 of file SpaceCharge.h.
|
strong |
Enumerator | |
---|---|
LocalDistCorr | using local dis/corr interpolator for calculation of global distortions/corrections |
ElectricalField | using electric field for calculation of global distortions/corrections |
Definition at line 117 of file SpaceCharge.h.
|
strong |
Definition at line 111 of file SpaceCharge.h.
|
strong |
numerical integration strategys
Enumerator | |
---|---|
Trapezoidal | trapezoidal integration (https://en.wikipedia.org/wiki/Trapezoidal_rule). straight electron drift line assumed: z0->z1, r0->r0, phi0->phi0 |
Simpson | simpon integration. see: https://en.wikipedia.org/wiki/Simpson%27s_rule. straight electron drift line assumed: z0->z1, r0->r0, phi0->phi0 |
Root | Root integration. straight electron drift line assumed: z0->z1, r0->r0, phi0->phi0. |
SimpsonIterative | simpon integration, but using an iterative method to approximate the drift path. No straight electron drift line assumed: z0->z1, r0->r1, phi0->phi1 |
Definition at line 100 of file SpaceCharge.h.
|
strong |
Enumerator | |
---|---|
ShiftedClip | shifted copper rod clip |
RotatedClip | rotated mylar strips from FC |
RodShift | shifted copper rod |
Definition at line 122 of file SpaceCharge.h.
|
strong |
Enumerator | |
---|---|
Distortions | distortions |
Corrections | corrections |
Definition at line 106 of file SpaceCharge.h.
SpaceCharge::SpaceCharge | ( | const int | bfield, |
const unsigned short | nZVertices = 129 , |
||
const unsigned short | nRVertices = 129 , |
||
const unsigned short | nPhiVertices = 360 , |
||
const bool | initBuffers = false |
||
) |
constructor grid granularity has to set before constructing an object using the static function setGrid(nZVertices, nRVertices, nPhiVertices)!
bfield | magnetic field (-5, -2, 0, 2, 5) |
nZVertices | grid vertices in z direction |
nRVertices | grid vertices in r direction |
nPhiVertices | grid vertices in phi direction |
initBuffers | initialize all buffers |
Definition at line 61 of file SpaceCharge.cxx.
SpaceCharge::SpaceCharge | ( | ) |
constructor for empty object. Can be used when buffers are loaded from file
Definition at line 71 of file SpaceCharge.cxx.
|
default |
default move constructor
void SpaceCharge::addBoundaryPotential | ( | const SpaceCharge< DataT > & | other, |
const Side | side, | ||
const float | scaling = 1 |
||
) |
adding the boundary potential from other sc object which same number of vertices!
other | other SC object which boundary potential witll be added |
scaling | sclaing factor which used to scale the others boundary potential |
Definition at line 3584 of file SpaceCharge.cxx.
void SpaceCharge::addChargeDensity | ( | const SpaceCharge< DataT > & | otherSC | ) |
add space charge density from other object (this.mDensity = this.mDensity + other.mDensity)
otherSC | other space-charge object, which charge will be added to current object |
Definition at line 3415 of file SpaceCharge.cxx.
average the sc density across all sectors per side
Definition at line 3907 of file SpaceCharge.cxx.
step 2: calculate numerically the electric field from the potential
side | side of the TPC |
Definition at line 610 of file SpaceCharge.cxx.
void o2::tpc::SpaceCharge< DataT >::calcGlobalCorrections | ( | const Fields & | formulaStruct, |
const int | type = 3 |
||
) |
step 4: calculate global corrections by using the electric field or the local corrections
formulaStruct | struct containing a method to evaluate the electric field Er, Ez, Ephi or the local corrections |
type | how to treat the corrections at regions where the corrected value is out of the TPC volume: type=0: use last valid correction value, type=1 do linear extrapolation, type=2 do parabolic extrapolation, type=3 do NOT abort when reaching the CE or the IFC to get a smooth estimate of the corrections |
void o2::tpc::SpaceCharge< DataT >::calcGlobalCorrections | ( | const Formulas & | formulaStruct, |
const int | type | ||
) |
Definition at line 1587 of file SpaceCharge.cxx.
|
inline |
calculate the global corrections using the electric fields (interface for python)
Definition at line 340 of file SpaceCharge.h.
void SpaceCharge::calcGlobalCorrWithGlobalDistIterative | ( | const Side | side, |
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachZ = 1 , |
||
const DataT | approachR = 1 , |
||
const DataT | approachPhi = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
calculate global corrections from global distortions
scSCale | possible second sc object |
scale | scaling for second sc object |
Definition at line 722 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalCorrWithGlobalDistIterative | ( | const SpaceCharge< DataT > * | scSCale = nullptr , |
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachZ = 1 , |
||
const DataT | approachR = 1 , |
||
const DataT | approachPhi = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Definition at line 728 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalCorrWithGlobalDistIterativeLinearCartesian | ( | const Side | side, |
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachX = 1 , |
||
const DataT | approachY = 1 , |
||
const DataT | approachZ = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Calculate global corrections using the global distortions by scaling with second distortion object, which will be consistent with scaled corrections in cartesian coordinates
side | Side of the TPC |
scSCale | possible second sc object |
scale | scaling for second sc object |
Definition at line 878 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalCorrWithGlobalDistIterativeLinearCartesian | ( | const SpaceCharge< DataT > * | scSCale = nullptr , |
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachX = 1 , |
||
const DataT | approachY = 1 , |
||
const DataT | approachZ = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Definition at line 884 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistortions | ( | const Fields & | formulaStruct, |
const int | maxIterations = 3 * sSteps * 129 |
||
) |
step 5: calculate global distortions by using the electric field or the local distortions (SLOW)
formulaStruct | struct containing a method to evaluate the electric field Er, Ez, Ephi or the local distortions |
maxIterations | maximum steps which are are performed to reach the central electrode (in general this is not necessary, but in case of problems this value aborts the calculation) |
Definition at line 1503 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistWithGlobalCorrIterative | ( | const DistCorrInterpolator< DataT > & | globCorr, |
const int | maxIter = 100 , |
||
const DataT | approachZ = 1 , |
||
const DataT | approachR = 1 , |
||
const DataT | approachPhi = 1 , |
||
const DataT | diffCorr = 50e-6 , |
||
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 |
||
) |
step 5: calculate global distortions using the global corrections (FAST)
globCorr | interpolator for global corrections |
maxIter | maximum iterations per global distortion |
approachZ | when the difference between the desired z coordinate and the position of the global correction is deltaZ, approach the desired z coordinate by deltaZ * approachZ . |
approachR | when the difference between the desired r coordinate and the position of the global correction is deltaR, approach the desired r coordinate by deltaR * approachR . |
approachPhi | when the difference between the desired phi coordinate and the position of the global correction is deltaPhi, approach the desired phi coordinate by deltaPhi * approachPhi . |
diffCorr | if the absolute differences from the interpolated values for the global corrections from the last iteration compared to the current iteration is smaller than this value, set converged to true for current global distortion |
type | whether to calculate distortions or corrections |
Definition at line 700 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistWithGlobalCorrIterative | ( | const Side | side, |
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachZ = 1 , |
||
const DataT | approachR = 1 , |
||
const DataT | approachPhi = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
step 5: calculate global distortions using the global corrections (FAST)
scSCale | possible second sc object |
scale | scaling for second sc object |
Definition at line 706 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistWithGlobalCorrIterative | ( | const SpaceCharge< DataT > * | scSCale = nullptr , |
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachZ = 1 , |
||
const DataT | approachR = 1 , |
||
const DataT | approachPhi = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Definition at line 712 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistWithGlobalCorrIterativeLinearCartesian | ( | const Side | side, |
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachX = 1 , |
||
const DataT | approachY = 1 , |
||
const DataT | approachZ = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Calculate global distortions using the global corrections by scaling with second distortion object, which will be consistent with scaled corrections in cartesian coordinates (as done in the tracking)
side | Side of the TPC |
scSCale | possible second sc object |
scale | scaling for second sc object |
Definition at line 862 of file SpaceCharge.cxx.
void SpaceCharge::calcGlobalDistWithGlobalCorrIterativeLinearCartesian | ( | const SpaceCharge< DataT > * | scSCale = nullptr , |
float | scale = 0 , |
||
const int | maxIter = 100 , |
||
const DataT | approachX = 1 , |
||
const DataT | approachY = 1 , |
||
const DataT | approachZ = 1 , |
||
const DataT | diffCorr = 50e-6 |
||
) |
Definition at line 868 of file SpaceCharge.cxx.
void SpaceCharge::calcLocalDistortionCorrectionVector | ( | const ElectricFields & | formulaStruct | ) |
step 3b: calculate the local distortion and correction vectors with an electric field
formulaStruct | struct containing a method to evaluate the electric field Er, Ez, Ephi (analytical formula or by TriCubic interpolator) |
Definition at line 1336 of file SpaceCharge.cxx.
void o2::tpc::SpaceCharge< DataT >::calcLocalDistortionsCorrections | ( | const SpaceCharge< DataT >::Type | type, |
const ElectricFields & | formulaStruct | ||
) |
Definition at line 1252 of file SpaceCharge.cxx.
void o2::tpc::SpaceCharge< DataT >::calcLocalDistortionsCorrections | ( | const Type | type, |
const ElectricFields & | formulaStruct | ||
) |
step 3: calculate the local distortions and corrections with an electric field
type | calculate local corrections or local distortions: type = o2::tpc::SpaceCharge<>::Type::Distortions or o2::tpc::SpaceCharge<>::Type::Corrections |
formulaStruct | struct containing a method to evaluate the electric field Er, Ez, Ephi (analytical formula or by TriCubic interpolator) |
void o2::tpc::SpaceCharge< DataT >::calcLocalDistortionsCorrectionsRK4 | ( | const SpaceCharge< DataT >::Type | type, |
const Side | side | ||
) |
Definition at line 1372 of file SpaceCharge.cxx.
void o2::tpc::SpaceCharge< DataT >::calcLocalDistortionsCorrectionsRK4 | ( | const Type | type, |
const Side | side | ||
) |
step 3b: calculate the local distortions and corrections with the local distortion/correction vectors using Runge Kutta 4. calcLocalDistortionCorrectionVector() has to be called before this function
type | calculate local corrections or local distortions: type = o2::tpc::SpaceCharge<>::Type::Distortions or o2::tpc::SpaceCharge<>::Type::Corrections |
side | side of the TPC |
void SpaceCharge::calculateDistortionsCorrections | ( | const o2::tpc::Side | side, |
const bool | calcVectors = false |
||
) |
side | side of the TPC |
calcVectors | set to calculate also the local distortion and local correction vectors |
Definition at line 111 of file SpaceCharge.cxx.
std::vector< std::pair< std::vector< o2::math_utils::Point3D< float > >, std::array< DataT, 3 > > > SpaceCharge::calculateElectronDriftPath | ( | const std::vector< GlobalPosition3D > & | elePos, |
const int | nSamplingPoints, | ||
const std::string_view | outFile = "electron_tracks.root" |
||
) | const |
function to calculate the drift paths of the electron whose starting position is delivered. Electric fields must be set!
elePos | global position of the start position of the electron |
nSamplingPoints | number of output points of the electron drift path |
outFile | name of the output debug file (if empty no file is created) |
Definition at line 2195 of file SpaceCharge.cxx.
bool SpaceCharge::checkGridFromFile | ( | std::string_view | file, |
std::string_view | tree | ||
) |
compare currently set grid with stored grid (in case the grid definition differs this instance will be newly initalizes with correct grid)
file | input file |
tree | tree in input file |
Definition at line 3222 of file SpaceCharge.cxx.
|
static |
Convert the IDCs to the number of ions for the ion backflow (primary ionization is not considered)
idcZero | map containing the IDCs values which will be converted to the space-charge density |
mapIBF | map contains the pad-by-pad IBF in % |
ionDriftTimeMS | ion drift time in ms |
normToPadArea | normalize IDCs to pad area (should always be true as the normalization is performed in IDCFactorization::calcIDCZero |
Definition at line 3472 of file SpaceCharge.cxx.
void SpaceCharge::correctElectron | ( | GlobalPosition3D & | point | ) |
Correct electron position using correction lookup tables
point | 3D coordinates of the electron |
Definition at line 1741 of file SpaceCharge.cxx.
void SpaceCharge::distortElectron | ( | GlobalPosition3D & | point, |
const SpaceCharge< DataT > * | scSCale = nullptr , |
||
float | scale = 0 |
||
) | const |
Distort electron position using distortion lookup tables
point | 3D coordinates of the electron |
scSCale | other sc object which is used for scaling of the distortions |
scale | scaling value |
Definition at line 1756 of file SpaceCharge.cxx.
int SpaceCharge::dumpAnalyticalCorrectionsDistortions | ( | TFile & | outf | ) | const |
write analytical corrections and distortions to file
outf | output file where the analytical corrections and distortions will be written to |
Definition at line 2794 of file SpaceCharge.cxx.
int SpaceCharge::dumpDensity | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write density to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3210 of file SpaceCharge.cxx.
int SpaceCharge::dumpElectricFields | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write electric fields to file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 2967 of file SpaceCharge.cxx.
int SpaceCharge::dumpGlobalCorrections | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write global correction to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3044 of file SpaceCharge.cxx.
int SpaceCharge::dumpGlobalCorrections | ( | TFile & | outf, |
const Side | side | ||
) | const |
write global corrections to root file (deprecated)
outf | output file where the global corrections will be written to |
side | of the TPC |
Definition at line 3254 of file SpaceCharge.cxx.
int SpaceCharge::dumpGlobalDistortions | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write global distortions to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 2999 of file SpaceCharge.cxx.
int SpaceCharge::dumpLocalCorrections | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write local corrections to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3090 of file SpaceCharge.cxx.
int SpaceCharge::dumpLocalDistCorrVectors | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write local corrections to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3137 of file SpaceCharge.cxx.
int SpaceCharge::dumpLocalDistortions | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write local corrections to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3122 of file SpaceCharge.cxx.
void SpaceCharge::dumpMetaData | ( | std::string_view | file, |
std::string_view | option, | ||
const bool | overwriteExisting = false |
||
) | const |
dump meta data to file (mC0, mC1, mC2, RegularGrid)
file | output file |
option | "RECREATE" or "UPDATE" |
overwriteExisting | overwrite existing meta data in file |
Definition at line 3292 of file SpaceCharge.cxx.
int SpaceCharge::dumpPotential | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write potential to root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3186 of file SpaceCharge.cxx.
write all fields etc to a file for both sides
file | output file where the electrical fields will be written to |
Definition at line 3285 of file SpaceCharge.cxx.
void SpaceCharge::dumpToFile | ( | std::string_view | file, |
const Side | side, | ||
std::string_view | option | ||
) | const |
write all fields etc to a file
file | output file where the electrical fields will be written to |
side | of the TPC |
option | "RECREATE" or "UPDATE" |
Definition at line 3268 of file SpaceCharge.cxx.
void SpaceCharge::dumpToTree | ( | const char * | outFileName, |
const o2::tpc::Sector & | sector, | ||
const int | nZPoints = 50 |
||
) | const |
dump to tree evaluated on the pads for given sector
outFileName | name of the output file |
sector | TPC sector for which the TTree is created |
nZPoints | number of points in z to consider |
Definition at line 2653 of file SpaceCharge.cxx.
void SpaceCharge::dumpToTree | ( | const char * | outFileName, |
const Side | side, | ||
const int | nZPoints = 50 , |
||
const int | nRPoints = 150 , |
||
const int | nPhiPoints = 180 , |
||
const bool | randomize = false |
||
) | const |
dump sc density, potential, electric fields, global distortions/corrections to tree
outFileName | name of the output file |
side | of the TPC |
nZPoints | number of vertices of the output in z |
nRPoints | number of vertices of the output in r |
nPhiPoints | number of vertices of the output in phi |
randomize | randomize points |
Adding other TTree as friend: TFile f1("file_1.root"); TTree* tree1 = (TTree*)f1.Get("tree"); tree1->BuildIndex("globalIndex")
TFile f2("file_2.root"); TTree* tree2 = (TTree*)f2.Get("tree"); tree2->BuildIndex("globalIndex"); tree1->AddFriend(tree2, "t2");
Definition at line 2422 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeDensityFromCalDet | ( | const std::vector< CalDet< float > > & | calSCDensity3D | ) |
step 0: set the space charge density from std::vector<CalDet> containing the space charge density. Each entry in the object corresponds to one z slice
calSCDensity3D | pad-by-pad CalDet for the space charge density |
Definition at line 1091 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeDensityFromFile | ( | TFile & | fInp, |
const char * | name | ||
) |
step 0: set the charge density from TH3 histogram containing the space charge density
fInp | input file containing a histogram for the space charge density |
name | the name of the space charge density histogram in the file |
Definition at line 1083 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeDensityFromHisto | ( | const char * | file, |
const char * | nameA, | ||
const char * | nameC | ||
) |
step 0: set the space-charge density from two TH3 histograms containing the space charge density for A and C side separately which are stored in a ROOT file
file | path to root file containing the space-charge density |
nameA | name of the space-charge density histogram for the A-side |
nameC | name of the space-charge density histogram for the C-side |
Definition at line 3428 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeDensityFromHisto | ( | const TH3 & | hisSCDensity3D | ) |
step 0: set the charge density from TH3 histogram containing the space charge density
hisSCDensity3D | histogram for the space charge density |
Definition at line 1106 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeDensityFromHisto | ( | const TH3 & | hisSCDensity3D_A, |
const TH3 & | hisSCDensity3D_C | ||
) |
step 0: set the space-charge density from two TH3 histograms containing the space-charge density for A and C side seperately
hisSCDensity3D_A | histogram for the space charge density for A-side |
hisSCDensity3D_C | histogram for the space charge density for C-side |
Definition at line 3443 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeFromCalDet | ( | const std::vector< CalDet< float > > & | calCharge3D | ) |
step 0: set the charge (number of ions) from std::vector<CalDet> containing the charge. Each entry in the object corresponds to one z slice. Normalization to the space charge is also done automatically
calCharge3D | histogram for the charge |
Definition at line 1098 of file SpaceCharge.cxx.
void SpaceCharge::fillChargeFromIDCs | ( | std::vector< CalDet< float > > & | idcZero, |
const CalDet< float > & | mapIBF, | ||
const float | ionDriftTimeMS = 200 , |
||
const bool | normToPadArea = true |
||
) |
Convert the IDCs to the number of ions for the ion backflow (primary ionization is not considered)
idcZero | map containing the IDCs which will be converted to the space-charge density |
mapIBF | map contains the pad-by-pad IBF in % |
ionDriftTimeMS | ion drift time in ms |
normToPadArea | normalize IDCs to pad area (should always be true as the normalization is performed in IDCFactorization::calcIDCZero |
Definition at line 3522 of file SpaceCharge.cxx.
|
inline |
set the space charge density directly
density | space charege density which will be set |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | of the TPC |
Definition at line 1020 of file SpaceCharge.h.
|
inline |
set the potential directly
potential | potential which will be set |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | of the TPC |
Definition at line 1001 of file SpaceCharge.h.
|
inline |
Definition at line 405 of file SpaceCharge.h.
|
inline |
Definition at line 396 of file SpaceCharge.h.
|
inline |
Definition at line 399 of file SpaceCharge.h.
|
inline |
Definition at line 402 of file SpaceCharge.h.
void SpaceCharge::getCorrections | ( | const DataT | x, |
const DataT | y, | ||
const DataT | z, | ||
const Side | side, | ||
DataT & | corrX, | ||
DataT & | corrY, | ||
DataT & | corrZ | ||
) | const |
get the global corrections for given coordinate
x | global x coordinate |
y | global y coordinate |
z | global z coordinate |
corrX | returns corrections in x direction |
corrY | returns corrections in y direction |
corrZ | returns corrections in z direction |
Definition at line 1881 of file SpaceCharge.cxx.
|
inline |
get the analytical global corrections for given coordinate
x | global x coordinate |
y | global y coordinate |
z | global z coordinate |
corrX | returns corrections in x direction |
corrY | returns corrections in y direction |
corrZ | returns corrections in z direction |
Definition at line 488 of file SpaceCharge.h.
void SpaceCharge::getCorrectionsCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | corrZ, | ||
DataT & | corrR, | ||
DataT & | corrRPhi | ||
) | const |
get the global correction for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
corrZ | returns correction in z direction |
corrR | returns correction in r direction |
corrRPhi | returns correction in rphi direction |
Definition at line 1860 of file SpaceCharge.cxx.
void SpaceCharge::getCorrectionsCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | corrZ, | ||
std::vector< DataT > & | corrR, | ||
std::vector< DataT > & | corrRPhi | ||
) | const |
get the global correction for given coordinates
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
corrZ | returns corrections in z direction |
corrR | returns corrections in r direction |
corrRPhi | returns corrections in rphi direction |
Definition at line 1868 of file SpaceCharge.cxx.
float SpaceCharge::getDCAr | ( | float | tgl, |
const int | nPoints, | ||
const float | phi, | ||
o2::utils::TreeStreamRedirector * | pcstream = nullptr |
||
) | const |
get DCA in RPhi for high pt track
tgl | tgl of the track |
nPoints | number of points used to calculate the DCAr |
pcstream | if provided debug output is being created |
Definition at line 3985 of file SpaceCharge.cxx.
|
inline |
Get density for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 838 of file SpaceCharge.h.
DataT SpaceCharge::getDensityCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side | ||
) | const |
get the space charge density for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
Definition at line 1807 of file SpaceCharge.cxx.
void SpaceCharge::getDistortions | ( | const DataT | x, |
const DataT | y, | ||
const DataT | z, | ||
const Side | side, | ||
DataT & | distX, | ||
DataT & | distY, | ||
DataT & | distZ | ||
) | const |
get the global distortions for given coordinate
x | global x coordinate |
y | global y coordinate |
z | global z coordinate |
distX | returns distortion in x direction |
distY | returns distortion in y direction |
distZ | returns distortion in z direction |
Definition at line 1988 of file SpaceCharge.cxx.
|
inline |
get the global distortions for given coordinate for a possible analytical formula if it was specified
x | global x coordinate |
y | global y coordinate |
z | global z coordinate |
distX | returns distortion in x direction |
distY | returns distortion in y direction |
distZ | returns distortion in z direction |
Definition at line 578 of file SpaceCharge.h.
|
inline |
Definition at line 584 of file SpaceCharge.h.
void SpaceCharge::getDistortionsCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | distZ, | ||
DataT & | distR, | ||
DataT & | distRPhi | ||
) | const |
get the global distortions for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
distZ | returns distortion in z direction |
distR | returns distortion in r direction |
distRPhi | returns distortion in rphi direction |
Definition at line 1967 of file SpaceCharge.cxx.
void SpaceCharge::getDistortionsCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | distZ, | ||
std::vector< DataT > & | distR, | ||
std::vector< DataT > & | distRPhi | ||
) | const |
get the global distortions for given coordinate
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
distZ | returns distortions in z direction |
distR | returns distortions in r direction |
distRPhi | returns distortions in rphi direction |
Definition at line 1975 of file SpaceCharge.cxx.
void SpaceCharge::getElectricFieldsCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | eZ, | ||
DataT & | eR, | ||
DataT & | ePhi | ||
) | const |
get the electric field for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
eZ | returns correction in z direction |
eR | returns correction in r direction |
ePhi | returns correction in phi direction |
Definition at line 1831 of file SpaceCharge.cxx.
NumericalFields< DataT > SpaceCharge::getElectricFieldsInterpolator | ( | const Side | side | ) | const |
Get struct containing interpolators for the electrical fields
side | side of the TPC |
Definition at line 1033 of file SpaceCharge.cxx.
|
inline |
Get global electric Field Ephi for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 831 of file SpaceCharge.h.
|
inline |
Get global electric Field Er for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 817 of file SpaceCharge.h.
|
inline |
Get global electric Field Ez for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 824 of file SpaceCharge.h.
|
inlinestaticconstexpr |
Get constant electric field.
Definition at line 640 of file SpaceCharge.h.
DistCorrInterpolator< DataT > SpaceCharge::getGlobalCorrInterpolator | ( | const Side | side | ) | const |
Get struct containing interpolators for global corrections dR, dZ, dPhi
side | side of the TPC |
Definition at line 1073 of file SpaceCharge.cxx.
|
inline |
Get global correction dR for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 796 of file SpaceCharge.h.
|
inline |
Get global correction dRPhi for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 810 of file SpaceCharge.h.
|
inline |
Get global correction dZ for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 803 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 915 of file SpaceCharge.h.
DistCorrInterpolator< DataT > SpaceCharge::getGlobalDistInterpolator | ( | const Side | side | ) | const |
Get struct containing interpolators for global distortions dR, dZ, dPhi
side | side of the TPC |
Definition at line 1063 of file SpaceCharge.cxx.
|
inline |
Get global distortion dR for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 775 of file SpaceCharge.h.
|
inline |
Get global distortion dRPhi for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 789 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 912 of file SpaceCharge.h.
|
inline |
Get global distortion dZ for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 782 of file SpaceCharge.h.
|
inline |
Get the grid object.
Definition at line 664 of file SpaceCharge.h.
|
inline |
Get grid spacing in phi direction.
Definition at line 637 of file SpaceCharge.h.
|
inline |
Get grid spacing in r direction.
Definition at line 631 of file SpaceCharge.h.
|
inline |
Get grid spacing in z direction.
Definition at line 634 of file SpaceCharge.h.
void SpaceCharge::getLocalCorrectionsCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | lcorrZ, | ||
DataT & | lcorrR, | ||
DataT & | lcorrRPhi | ||
) | const |
get the local correction for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
lcorrZ | returns local correction in z direction |
lcorrR | returns local correction in r direction |
lcorrRPhi | returns local correction in rphi direction |
Definition at line 1839 of file SpaceCharge.cxx.
void SpaceCharge::getLocalCorrectionsCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | lcorrZ, | ||
std::vector< DataT > & | lcorrR, | ||
std::vector< DataT > & | lcorrRPhi | ||
) | const |
get the local correction for given coordinates
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
lcorrZ | returns local corrections in z direction |
lcorrR | returns local corrections in r direction |
lcorrRPhi | returns local corrections in rphi direction |
Definition at line 1847 of file SpaceCharge.cxx.
void SpaceCharge::getLocalCorrectionVectorCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | lveccorrZ, | ||
DataT & | lveccorrR, | ||
DataT & | lveccorrRPhi | ||
) | const |
get the local correction vector for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
ldistZ | returns local correction vector in z direction |
ldistR | returns local correction vector in r direction |
ldistRPhi | returns local correction vector in rphi direction |
Definition at line 1946 of file SpaceCharge.cxx.
void SpaceCharge::getLocalCorrectionVectorCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | lveccorrZ, | ||
std::vector< DataT > & | lveccorrR, | ||
std::vector< DataT > & | lveccorrRPhi | ||
) | const |
get the local correction vector for given coordinate
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
ldistZ | returns local correction vectors in z direction |
ldistR | returns local correction vectors in r direction |
ldistRPhi | returns local correction vectors in rphi direction |
Definition at line 1954 of file SpaceCharge.cxx.
DistCorrInterpolator< DataT > SpaceCharge::getLocalCorrInterpolator | ( | const Side | side | ) | const |
Get struct containing interpolators for local corrections dR, dZ, dPhi
side | side of the TPC |
Definition at line 1053 of file SpaceCharge.cxx.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 733 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 747 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 740 of file SpaceCharge.h.
DistCorrInterpolator< DataT > SpaceCharge::getLocalDistInterpolator | ( | const Side | side | ) | const |
Get struct containing interpolators for local distortions dR, dZ, dPhi
side | side of the TPC |
Definition at line 1043 of file SpaceCharge.cxx.
void SpaceCharge::getLocalDistortionsCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | ldistZ, | ||
DataT & | ldistR, | ||
DataT & | ldistRPhi | ||
) | const |
get the local distortions for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
ldistZ | returns local distortion in z direction |
ldistR | returns local distortion in r direction |
ldistRPhi | returns local distortion in rphi direction |
Definition at line 1904 of file SpaceCharge.cxx.
void SpaceCharge::getLocalDistortionsCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | ldistZ, | ||
std::vector< DataT > & | ldistR, | ||
std::vector< DataT > & | ldistRPhi | ||
) | const |
get the local distortions for given coordinates
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
ldistZ | returns local distortions in z direction |
ldistR | returns local distortions in r direction |
ldistRPhi | returns local distortions in rphi direction |
Definition at line 1912 of file SpaceCharge.cxx.
void SpaceCharge::getLocalDistortionVectorCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side, | ||
DataT & | lvecdistZ, | ||
DataT & | lvecdistR, | ||
DataT & | lvecdistRPhi | ||
) | const |
get the local distortion vector for given coordinates
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
lvecdistZ | returns local distortion vector in z direction |
lvecdistR | returns local distortion vector in r direction |
lvecdistRPhi | returns local distortion vector in rphi direction |
Definition at line 1925 of file SpaceCharge.cxx.
void SpaceCharge::getLocalDistortionVectorCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side, | ||
std::vector< DataT > & | lvecdistZ, | ||
std::vector< DataT > & | lvecdistR, | ||
std::vector< DataT > & | lvecdistRPhi | ||
) | const |
get the local distortion vector for given coordinate
z | global z coordinates |
r | global r coordinates |
phi | global phi coordinates |
lvecdistZ | returns local distortion vectors in z direction |
lvecdistR | returns local distortion vectors in r direction |
lvecdistRPhi | returns local distortion vectors in rphi direction |
Definition at line 1933 of file SpaceCharge.cxx.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 691 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 705 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 698 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 754 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 768 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 761 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 712 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 726 of file SpaceCharge.h.
|
inline |
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 719 of file SpaceCharge.h.
|
inline |
Definition at line 1224 of file SpaceCharge.h.
|
inline |
Definition at line 1222 of file SpaceCharge.h.
|
inline |
Definition at line 393 of file SpaceCharge.h.
|
inline |
Definition at line 390 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 898 of file SpaceCharge.h.
|
inlinestatic |
get the number of threads used for some of the calculations
Definition at line 901 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 909 of file SpaceCharge.h.
|
inline |
Definition at line 387 of file SpaceCharge.h.
Definition at line 105 of file SpaceCharge.cxx.
|
inlinestatic |
convert x and y coordinates from cartesian to phi in polar coordinates
Definition at line 596 of file SpaceCharge.h.
|
inline |
Get max phi.
Definition at line 658 of file SpaceCharge.h.
|
inline |
Get min phi.
Definition at line 649 of file SpaceCharge.h.
|
inline |
Get phi vertex position for index in phi direction
indexPhi | index in phi direction |
Definition at line 852 of file SpaceCharge.h.
|
inline |
Definition at line 407 of file SpaceCharge.h.
|
inline |
Get potential for vertex
iz | vertex in z dimension |
ir | vertex in r dimension |
iphi | vertex in phi dimension |
side | side of the TPC |
Definition at line 845 of file SpaceCharge.h.
DataT SpaceCharge::getPotentialCyl | ( | const DataT | z, |
const DataT | r, | ||
const DataT | phi, | ||
const Side | side | ||
) | const |
get the potential for given coordinate
z | global z coordinate |
r | global r coordinate |
phi | global phi coordinate |
Definition at line 1813 of file SpaceCharge.cxx.
std::vector< float > SpaceCharge::getPotentialCyl | ( | const std::vector< DataT > & | z, |
const std::vector< DataT > & | r, | ||
const std::vector< DataT > & | phi, | ||
const Side | side | ||
) | const |
get the potential for list of given coordinate
Definition at line 1819 of file SpaceCharge.cxx.
|
inlinestatic |
convert x and y coordinates from cartesian to the radius in polar coordinates
Definition at line 593 of file SpaceCharge.h.
|
inline |
Get max r.
Definition at line 652 of file SpaceCharge.h.
|
inline |
Get inner radius of tpc.
Definition at line 643 of file SpaceCharge.h.
|
inline |
Get r vertex position for index in r direction
indexR | index in r direction |
Definition at line 856 of file SpaceCharge.h.
|
inlinestatic |
Get the space-charge distortions model.
Definition at line 924 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 661 of file SpaceCharge.h.
|
inline |
Definition at line 892 of file SpaceCharge.h.
|
inline |
Definition at line 886 of file SpaceCharge.h.
|
inlinestatic |
Definition at line 918 of file SpaceCharge.h.
|
inlinestatic |
Get the step width which is used for the calculation of the correction/distortions in units of the z-bin.
Definition at line 848 of file SpaceCharge.h.
|
inline |
Definition at line 590 of file SpaceCharge.h.
|
inlinestatic |
convert radius and phi coordinates from polar coordinates to x cartesian coordinates
Definition at line 599 of file SpaceCharge.h.
|
inlinestatic |
convert radius and phi coordinates from polar coordinates to y cartesian coordinate
Definition at line 602 of file SpaceCharge.h.
|
inline |
Get max z.
Definition at line 655 of file SpaceCharge.h.
|
inline |
Get min z position which is used during the calaculations.
Definition at line 646 of file SpaceCharge.h.
|
inline |
Get z vertex position for index in z direction
indexZ | index in z direction |
Definition at line 860 of file SpaceCharge.h.
TODO is there a faster way to get the drift velocity
TODO fix hard coded values (ezField, t1, t2): export to Constants.h or get from somewhere?
TODO use this parameterization or fixed value(s) from Magboltz calculations?
Definition at line 2055 of file SpaceCharge.cxx.
Definition at line 3978 of file SpaceCharge.cxx.
set magnetic field which can be used for ExB misalignment distortions
field | magnetic field (-5, -2, 0, 2, 5) |
Definition at line 2838 of file SpaceCharge.cxx.
|
static |
inpFile | input file containing the electron tracks tree, which is the output file of calculateElectronDriftPath() |
hBorder | histogram which defines the borders for the drawing and also the axis titles |
type | setting dimensions: type=0: radius vs z, type=0: radius vs phi |
gifSpeed | speed of the output gif file (fastest is 2, slowest is 99) |
maxsamplingpoints | maximum number of frames which will be drawn (higher number increases processing time) |
outName | name of the output file |
Definition at line 2324 of file SpaceCharge.cxx.
void SpaceCharge::mirrorPotential | ( | const Side | sideRef, |
const Side | sideMirrored | ||
) |
mirror potential from one side to the other side
sideRef | side which contains the reference potential |
sideMirrored | side where the potential will be set from sideRef |
Definition at line 509 of file SpaceCharge.cxx.
normalize a histogram containing the charge to the space charge density
histoIonsPhiRZ | histogram which will be normalized to the sapce charge density |
Definition at line 2775 of file SpaceCharge.cxx.
|
default |
move assignment
void SpaceCharge::poissonSolver | ( | const DataT | stoppingConvergence = 1e-6 , |
const int | symmetry = 0 |
||
) |
step 1: use the O2TPCPoissonSolver class to numerically calculate the potential with set space charge density and boundary conditions from potential for A and C side in parallel
stoppingConvergence | stopping criterion used in the poisson solver |
symmetry | use symmetry or not in the poisson solver |
Definition at line 579 of file SpaceCharge.cxx.
void SpaceCharge::poissonSolver | ( | const Side | side, |
const DataT | stoppingConvergence = 1e-6 , |
||
const int | symmetry = 0 |
||
) |
step 1: use the O2TPCPoissonSolver class to numerically calculate the potential with set space charge density and boundary conditions from potential
side | side of the TPC |
stoppingConvergence | stopping criterion used in the poisson solver |
symmetry | use symmetry or not in the poisson solver |
Definition at line 569 of file SpaceCharge.cxx.
|
inline |
Definition at line 1223 of file SpaceCharge.h.
set meta data from file (mC0, mC1, mC2, RegularGrid)
file | input file to read from |
Definition at line 3321 of file SpaceCharge.cxx.
|
inline |
set phi coordinate between min phi max phi
Definition at line 1113 of file SpaceCharge.h.
set r coordinate between 'RMIN - 4 * GRIDSPACINGR' and 'RMAX + 2 * GRIDSPACINGR'. the r coordinate is not clamped to RMIN and RMAX to ensure correct interpolation at the borders of the grid.
Definition at line 198 of file SpaceCharge.cxx.
|
inline |
set z coordinate between min z max z
posZ | z position which will be regulated if needed |
Definition at line 1107 of file SpaceCharge.h.
void SpaceCharge::resetBoundaryPotentialToZeroInRangeZ | ( | float | zMin, |
float | zMax, | ||
const Side | side | ||
) |
setting the boundary potential to 0 for z<zMin and z>zMax
Definition at line 3624 of file SpaceCharge.cxx.
|
inline |
scale the space charge density by a scaling factor: sc_new = sc * scalingFactor
scalingFactor | factor to scale the space-charge density |
side | side for which the space-charge density will be scaled |
Definition at line 302 of file SpaceCharge.h.
void SpaceCharge::scaleChargeDensitySector | ( | const float | scalingFactor, |
const Sector | sector | ||
) |
scale the space-charge for one sector: space-charge density *= scalingFactor;
scalingFactor | scaling factor for the space-charge density |
Definition at line 3935 of file SpaceCharge.cxx.
void SpaceCharge::scaleChargeDensityStack | ( | const float | scalingFactor, |
const Sector | sector, | ||
const GEMstack | stack | ||
) |
scaling the space-charge density for given stack
Definition at line 3953 of file SpaceCharge.cxx.
void SpaceCharge::scaleCorrections | ( | const float | scaleFac, |
const Side | side | ||
) |
scale corrections by factor
scaleFac | global corrections are multiplied by this factor |
Definition at line 3899 of file SpaceCharge.cxx.
void SpaceCharge::setAnalyticalCorrectionsDistortionsFromFile | ( | std::string_view | inpf | ) |
set analytical corrections and distortions from file
inpf | input file where the analytical corrections and distortions are stored |
Definition at line 2805 of file SpaceCharge.cxx.
|
inline |
c0 | coefficient C0 (compare Jim Thomas's notes for definitions) |
c1 | coefficient C1 (compare Jim Thomas's notes for definitions) |
Definition at line 869 of file SpaceCharge.h.
|
inline |
c2 | coefficient C2 |
Definition at line 876 of file SpaceCharge.h.
void SpaceCharge::setChargeDensityFromFormula | ( | const AnalyticalFields< DataT > & | formulaStruct | ) |
step 0: this function fills the internal storage for the charge density using an analytical formula
formulaStruct | struct containing a method to evaluate the density |
Definition at line 212 of file SpaceCharge.cxx.
void SpaceCharge::setDefaultStaticDistortionsGEMFrameChargeUp | ( | const Side | side, |
const DataT | deltaPotential = 1000 |
||
) |
setting default potential (same potential for all GEM frames. The default value of 1000V are matched to distortions observed in laser data without X-Ray etc.
side | side of the TPC where the potential will be set |
deltaPotential | delta potential which will be set at the GEM frames |
Definition at line 229 of file SpaceCharge.cxx.
|
inline |
Define delta potential due to a possible shifted copper rod (delta potential spike at the copper rods) at the IFC
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1145 of file SpaceCharge.h.
|
inline |
Define delta potential due to a possible shifted copper rod (delta potential spike at the copper rods) at the OFC
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1149 of file SpaceCharge.h.
|
inline |
Define delta potential due to rotated clip at IFC. Only possible at sector 11 and 11+18. The delta potential increases linearly from neighbouring rod to the specified rod, the sign of the delta potential inverts and increases linearly to 0 to the other nerighbouring rod
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1162 of file SpaceCharge.h.
|
inline |
Define delta potential due to rotated clip at OFC. Only possible at sector 3 and 3+18. The delta potential increases linearly from neighbouring rod to the specified rod, the sign of the delta potential inverts and increases linearly to 0 to the other nerighbouring rod
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1167 of file SpaceCharge.h.
|
inline |
Define delta potential due to shifted copper rod and field cage strips at IFC (maximum of the delta potential at the copper rod and linear decreasing up to left and right neighbouring rods)
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1153 of file SpaceCharge.h.
|
inline |
Define delta potential due to shifted copper rod and field cage strips at OFC (maximum of the delta potential at the copper rod and linear decreasing up to left and right neighbouring rods)
deltaPot | delta potential which will be set at the copper rod |
Definition at line 1157 of file SpaceCharge.h.
void SpaceCharge::setDensityFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set density from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3242 of file SpaceCharge.cxx.
void SpaceCharge::setDistortionLookupTables | ( | const DataContainer & | distdZ, |
const DataContainer & | distdR, | ||
const DataContainer & | distdRPhi, | ||
const Side | side | ||
) |
set the distortions directly from a look up table
distdZ | distortions in z direction |
distdR | distortions in r direction |
distdRPhi | distortions in rphi direction |
side | side of the TPC |
Definition at line 2084 of file SpaceCharge.cxx.
|
inline |
set distortions and corrections by an analytical formula
Definition at line 581 of file SpaceCharge.h.
void SpaceCharge::setEFieldFromFormula | ( | const AnalyticalFields< DataT > & | formulaStruct | ) |
step 2a: set the electric field from an analytical formula
formulaStruct | struct containing a method to evaluate the electric fields |
Definition at line 589 of file SpaceCharge.cxx.
void SpaceCharge::setElectricFieldsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set electric field from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 2982 of file SpaceCharge.cxx.
set the density, potential, electric fields, local distortions/corrections, global distortions/corrections from a file for both sides. Missing objects in the file are ignored.
file | output file where the electrical fields will be written to |
Definition at line 3389 of file SpaceCharge.cxx.
set the density, potential, electric fields, local distortions/corrections, global distortions/corrections from a file. Missing objects in the file are ignored.
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3376 of file SpaceCharge.cxx.
void SpaceCharge::setGlobalCorrections | ( | const std::function< void(int sector, DataT gx, DataT gy, DataT gz, DataT &gCx, DataT &gCy, DataT &gCz)> & | gCorr, |
const Side | side | ||
) |
setting the global corrections directly from input function provided in global coordinates
gCorr | function returning global corrections for given global coordinate |
Definition at line 3740 of file SpaceCharge.cxx.
void SpaceCharge::setGlobalCorrectionsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set global corrections from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3059 of file SpaceCharge.cxx.
void SpaceCharge::setGlobalCorrectionsFromFile | ( | TFile & | inpf, |
const Side | side | ||
) |
set global corrections from root file (deprecated)
inpf | input file where the global corrections are stored |
side | of the TPC |
Definition at line 3078 of file SpaceCharge.cxx.
|
inlinestatic |
Definition at line 914 of file SpaceCharge.h.
void SpaceCharge::setGlobalDistortionsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set global distortions from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3014 of file SpaceCharge.cxx.
void SpaceCharge::setGlobalDistortionsFromFile | ( | TFile & | inpf, |
const Side | side | ||
) |
set global distortions from root file (deprecated)
inpf | input file where the global distortions are stored |
side | of the TPC |
Definition at line 3032 of file SpaceCharge.cxx.
|
inlinestatic |
Definition at line 911 of file SpaceCharge.h.
void SpaceCharge::setIFCChargeUpFallingPot | ( | const float | deltaPot, |
const float | zMaxDeltaPot, | ||
const int | type, | ||
const float | zEnd, | ||
const float | offs, | ||
const Side | side | ||
) |
IFC charge up: set a linear rising delta potential from the CE to given z position which falls linear down to 0 at the readout
deltaPot | maximum value of the delta potential in V |
zMaxDeltaPot | z position where the maximum delta potential of deltaPot will be set |
type | function which is used to set falling potential: 0=linear falling off, 1=1/x falling off, 2=1/x steeply falling, 3=linear with offset |
offs | if offs != 0 the potential doesnt fall to 0. E.g. deltaPot=400V and offs=-10V -> Potential falls from 400V at zMaxDeltaPot to -10V at z=250cm |
Definition at line 3696 of file SpaceCharge.cxx.
void SpaceCharge::setIFCChargeUpRisingPot | ( | const float | deltaPot, |
const float | zMaxDeltaPot, | ||
const int | type, | ||
const float | zStart, | ||
const float | offs, | ||
const Side | side | ||
) |
IFC charge up: set a linear rising delta potential from the CE to given z
deltaPot | maximum value of the delta potential in V |
zMaxDeltaPot | z position where the maximum delta potential of deltaPot will be set |
type | functional form of the delta potential: 0=linear, 1= 1/x, 2=flat, 3=linear falling, 4=flat no z dependence |
zStart | usually at 0 to start the rising of the potential at the IFC |
Definition at line 3658 of file SpaceCharge.cxx.
void SpaceCharge::setLocalCorrectionsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set local corrections from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3105 of file SpaceCharge.cxx.
void SpaceCharge::setLocalDistCorrVectorsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set local distortions from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3169 of file SpaceCharge.cxx.
void SpaceCharge::setLocalDistortionsFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set local distortions from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3152 of file SpaceCharge.cxx.
|
inline |
Definition at line 1225 of file SpaceCharge.h.
|
inline |
setting meta data for this object
Definition at line 1221 of file SpaceCharge.h.
|
inlinestatic |
set number of steps used for calculation of distortions/corrections per z bin
nSteps | number of steps per z bin |
Definition at line 896 of file SpaceCharge.h.
|
inlinestatic |
set the number of threads used for some of the calculations
Definition at line 904 of file SpaceCharge.h.
|
inlinestatic |
set which kind of numerical integration is used for calcution of the integrals int Er/Ez dz, int Ephi/Ez dz, int Ez dz
strategy | numerical integration strategy. see enum IntegrationStrategy for the different types |
Definition at line 908 of file SpaceCharge.h.
void SpaceCharge::setOmegaTauT1T2 | ( | const DataT | omegaTau = 0.32f , |
const DataT | t1 = 1 , |
||
const DataT | t2 = 1 |
||
) |
omegaTau | \omega \tau value |
t1 | value for t1 see: ??? |
t2 | value for t2 see: ??? |
Definition at line 3404 of file SpaceCharge.cxx.
void SpaceCharge::setPotential | ( | int | iz, |
int | ir, | ||
int | iphi, | ||
Side | side, | ||
float | val | ||
) |
setting the potential directly for given vertex
Definition at line 4044 of file SpaceCharge.cxx.
void SpaceCharge::setPotentialBoundaryFromFormula | ( | const AnalyticalFields< DataT > & | formulaStruct | ) |
step 0: this function fills the boundary of the potential using an analytical formula. The boundary is used in the PoissonSolver.
formulaStruct | struct containing a method to evaluate the potential |
Definition at line 523 of file SpaceCharge.cxx.
void SpaceCharge::setPotentialBoundaryGEMFrameAlongR | ( | const std::function< DataT(DataT)> & | potentialFunc, |
const Side | side | ||
) |
setting the boundary potential of the GEM stack along the radius
potentialFunc | potential funtion as a function of the radius |
Side | of the TPC |
Definition at line 263 of file SpaceCharge.cxx.
|
inline |
setting the boundary potential of the IROC on the bottom along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 225 of file SpaceCharge.h.
|
inline |
setting the potential from the IROC to the inner field cage
Definition at line 268 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the IROC on the top along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 230 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC1 on the bottom along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 235 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC1 on the top along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 240 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC2 on the bottom along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 245 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC2 on the top along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 250 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC3 on the bottom along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 255 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC3 on the top along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 265 of file SpaceCharge.h.
|
inline |
setting the boundary potential of the OROC3 on the top along phi
potentialFunc | potential funtion as a function of global phi |
Side | of the TPC |
Definition at line 260 of file SpaceCharge.h.
void SpaceCharge::setPotentialBoundaryInnerRadius | ( | const std::function< DataT(DataT)> & | potentialFunc, |
const Side | side | ||
) |
setting the boundary potential for the inner TPC radius along r
potentialFunc | potential funtion as a function of z |
Side | of the TPC |
Definition at line 464 of file SpaceCharge.cxx.
void SpaceCharge::setPotentialBoundaryOuterRadius | ( | const std::function< DataT(DataT)> & | potentialFunc, |
const Side | side | ||
) |
setting the boundary potential for the outer TPC radius along r
potentialFunc | potential funtion as a function of z |
Side | of the TPC |
Definition at line 478 of file SpaceCharge.cxx.
void SpaceCharge::setPotentialFromFile | ( | std::string_view | file, |
const Side | side | ||
) |
set potential from root file using RDataFrame
file | output file where the electrical fields will be written to |
side | of the TPC |
Definition at line 3198 of file SpaceCharge.cxx.
void SpaceCharge::setPotentialFromFormula | ( | const AnalyticalFields< DataT > & | formulaStruct | ) |
step 0: this function fills the potential using an analytical formula
formulaStruct | struct containing a method to evaluate the potential |
Definition at line 492 of file SpaceCharge.cxx.
void SpaceCharge::setROCMisalignmentRotationAlongX | ( | const int | sector, |
const int | type, | ||
const float | potentialMin, | ||
const float | potentialMax | ||
) |
set misalignment of ROC for rotation along local x
sector | sector for which the misalignment in z will be applied (if sector=-1 all sectors are shifted) |
type | 0=IROC, 1=OROC, 2=IROC+OROC |
potential | minimum delta potential |
potential | maximum delta potential |
Definition at line 3790 of file SpaceCharge.cxx.
void SpaceCharge::setROCMisalignmentRotationAlongY | ( | const int | sector, |
const int | type, | ||
const float | potentialMin, | ||
const float | potentialMax | ||
) |
set misalignment of ROC for rotation along local y
sector | sector for which the misalignment in z will be applied (if sector=-1 all sectors are shifted) |
type | 0=IROC, 1=OROC, 2=IROC+OROC |
potential | minimum delta potential |
potential | maximum delta potential |
Definition at line 3796 of file SpaceCharge.cxx.
void SpaceCharge::setROCMisalignmentShiftZ | ( | const int | sector, |
const int | type, | ||
const float | potential | ||
) |
set misalignment of ROC for shift in z
sector | sector for which the misalignment in z will be applied (if sector=-1 all sectors are shifted) |
type | 0=IROC, 1=OROC, 2=IROC+OROC |
potential | delta potential on which the ROCs are set |
Definition at line 3784 of file SpaceCharge.cxx.
|
inlinestatic |
Set the space-charge distortions model
distortionType | distortion type (constant or realistic) |
Definition at line 922 of file SpaceCharge.h.
|
inline |
calculate distortions due to electric fields (space charge, boundary potential...)
Definition at line 889 of file SpaceCharge.h.
|
inline |
enable/disable calculation of distortions due to ExB misalignment
Definition at line 883 of file SpaceCharge.h.
simulate only one sector instead of 18 per side. This makes currently only sense for the static distortions (ToDo: simplify usage) phi max will be restricted to 2Pi/18 for this instance and for global instance of poisson solver
Definition at line 3359 of file SpaceCharge.cxx.
|
inlinestatic |
Definition at line 917 of file SpaceCharge.h.
|
inline |
setting usage of the analytical formula for the distortions and corrections
Definition at line 587 of file SpaceCharge.h.
|
inline |
Definition at line 926 of file SpaceCharge.h.
void SpaceCharge::subtractGlobalCorrections | ( | const SpaceCharge< DataT > & | otherSC, |
const Side | side | ||
) |
substract global corrections from other sc object (global corrections -= other.global corrections) can be used to calculate the derivative: (this - other)/normalization for normalization see scaleCorrections()
Definition at line 3883 of file SpaceCharge.cxx.
void SpaceCharge::subtractGlobalDistortions | ( | const SpaceCharge< DataT > & | otherSC, |
const Side | side | ||
) |
substract global distortions from other sc object (global distortions -= other.global distortions) can be used to calculate the derivative: (this - other)/normalization
Definition at line 3891 of file SpaceCharge.cxx.
unsetting simulation of one sector
Definition at line 3370 of file SpaceCharge.cxx.