19#ifndef ALICEO2_TPC_SPACECHARGE_H_
20#define ALICEO2_TPC_SPACECHARGE_H_
45class TreeStreamRedirector;
71template <
typename DataT =
double>
78 using TH3DataT = std::conditional_t<std::is_same<DataT, double>::value, TH3D, TH3F>;
88 SpaceCharge(
const int bfield,
const unsigned short nZVertices = 129,
const unsigned short nRVertices = 129,
const unsigned short nPhiVertices = 360,
const bool initBuffers =
false);
335 template <
typename ElectricFields = AnalyticalFields<DataT>>
340 template <
typename ElectricFields = AnalyticalFields<DataT>>
347 template <
typename ElectricFields = AnalyticalFields<DataT>>
353 template <
typename Fields = AnalyticalFields<DataT>>
362 template <
typename Fields = AnalyticalFields<DataT>>
436 std::vector<float>
getDensityCyl(
const std::vector<DataT>&
z,
const std::vector<DataT>&
r,
const std::vector<DataT>& phi,
const Side side)
const;
445 std::vector<float>
getPotentialCyl(
const std::vector<DataT>&
z,
const std::vector<DataT>&
r,
const std::vector<DataT>& phi,
const Side side)
const;
472 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;
490 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;
508 void getCorrectionsAnalytical(
const DataT x,
const DataT y,
const DataT z,
const Side side,
DataT& corrX,
DataT& corrY,
DataT& corrZ)
const {
getDistortionsCorrectionsAnalytical(
x,
y,
z,
side, corrX, corrY, corrZ,
false); }
526 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;
544 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;
562 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;
580 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;
598 void getDistortionsAnalytical(
const DataT x,
const DataT y,
const DataT z,
const Side side,
DataT& distX,
DataT& distY,
DataT& distZ)
const {
getDistortionsCorrectionsAnalytical(
x,
y,
z,
side, distX, distY, distZ,
true); }
916 static void setNStep(
const int nSteps) { sSteps = nSteps; }
924 static void setNThreads(
const int nThreads) { sNThreads = nThreads; }
962 void dumpMetaData(std::string_view
file, std::string_view option,
const bool overwriteExisting =
false)
const;
985 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;
1067 template <
typename DataTIn = DataT>
1089 template <
typename DataTIn = DataT>
1140 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;
1148 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");
1271 inline static int sSimpsonNIteratives{3};
1272 inline static int sSteps{1};
1280 static constexpr int FNSIDES =
SIDES;
1281 bool mUseInitialSCDensity{
false};
1282 bool mInitLookUpTables{
false};
1283 bool mSimExBMisalignment{
false};
1284 bool mSimEDistortions{
true};
1285 bool mReadMetaData{
false};
1286 RegularGrid mGrid3D[FNSIDES]{{
GridProp::ZMIN,
GridProp::RMIN,
GridProp::PHIMIN, getSign(
Side::A) *
GridProp::getGridSpacingZ(mParamGrid.
NZVertices),
GridProp::getGridSpacingR(mParamGrid.
NRVertices),
GridProp::getGridSpacingPhi(mParamGrid.
NPhiVertices), mParamGrid}, {
GridProp::ZMIN,
GridProp::RMIN,
GridProp::PHIMIN, getSign(
Side::C) *
GridProp::getGridSpacingZ(mParamGrid.
NZVertices),
GridProp::getGridSpacingR(mParamGrid.
NRVertices),
GridProp::getGridSpacingPhi(mParamGrid.
NPhiVertices), mParamGrid}};
1288 DataContainer mLocalDistdR[FNSIDES]{};
1289 DataContainer mLocalDistdZ[FNSIDES]{};
1290 DataContainer mLocalDistdRPhi[FNSIDES]{};
1291 DataContainer mLocalVecDistdR[FNSIDES]{};
1292 DataContainer mLocalVecDistdZ[FNSIDES]{};
1293 DataContainer mLocalVecDistdRPhi[FNSIDES]{};
1294 DataContainer mLocalCorrdR[FNSIDES]{};
1295 DataContainer mLocalCorrdZ[FNSIDES]{};
1296 DataContainer mLocalCorrdRPhi[FNSIDES]{};
1297 DataContainer mGlobalDistdR[FNSIDES]{};
1298 DataContainer mGlobalDistdZ[FNSIDES]{};
1299 DataContainer mGlobalDistdRPhi[FNSIDES]{};
1300 DataContainer mGlobalCorrdR[FNSIDES]{};
1301 DataContainer mGlobalCorrdZ[FNSIDES]{};
1302 DataContainer mGlobalCorrdRPhi[FNSIDES]{};
1303 DataContainer mDensity[FNSIDES]{};
1304 DataContainer mPotential[FNSIDES]{};
1305 DataContainer mElectricFieldEr[FNSIDES]{};
1306 DataContainer mElectricFieldEz[FNSIDES]{};
1307 DataContainer mElectricFieldEphi[FNSIDES]{};
1317 AnalyticalDistCorr<DataT> mAnaDistCorr;
1318 bool mUseAnaDistCorr{
false};
1324 bool isCloseToZero(
const DataT valA,
const DataT valB)
const {
return std::abs(valA + valB) <
static_cast<DataT>(0.01); }
1346 template <
typename Fields = AnalyticalFields<DataT>>
1350 void langevinCylindricalE(
DataT& ddR,
DataT& ddPhi,
DataT& ddZ,
const DataT radius,
const DataT localIntErOverEz,
const DataT localIntEPhiOverEz,
const DataT localIntDeltaEz)
const;
1353 void langevinCylindricalB(
DataT& ddR,
DataT& ddPhi,
const DataT radius,
const DataT localIntBrOverBz,
const DataT localIntBPhiOverBz)
const;
1356 template <
typename Fields = AnalyticalFields<DataT>>
1357 void integrateEFieldsRoot(
const DataT p1r,
const DataT p1phi,
const DataT p1z,
const DataT p2z,
DataT& localIntErOverEz,
DataT& localIntEPhiOverEz,
DataT& localIntDeltaEz,
const Fields& formulaStruct,
const DataT ezField,
const Side side)
const;
1360 template <
typename Fields = AnalyticalFields<DataT>>
1361 void integrateEFieldsTrapezoidal(
const DataT p1r,
const DataT p1phi,
const DataT p1z,
const DataT p2z,
DataT& localIntErOverEz,
DataT& localIntEPhiOverEz,
DataT& localIntDeltaEz,
const Fields& formulaStruct,
const DataT ezField,
const Side side)
const;
1364 template <
typename Fields = AnalyticalFields<DataT>>
1365 void integrateEFieldsSimpson(
const DataT p1r,
const DataT p1phi,
const DataT p1z,
const DataT p2z,
DataT& localIntErOverEz,
DataT& localIntEPhiOverEz,
DataT& localIntDeltaEz,
const Fields& formulaStruct,
const DataT ezField,
const Side side)
const;
1368 template <
typename Fields = AnalyticalFields<DataT>>
1369 void integrateEFieldsSimpsonIterative(
const DataT p1r,
const DataT p2r,
const DataT p1phi,
const DataT p2phi,
const DataT p1z,
const DataT p2z,
DataT& localIntErOverEz,
DataT& localIntEPhiOverEz,
DataT& localIntDeltaEz,
const Fields& formulaStruct,
const DataT ezField,
const Side side)
const;
1372 void processGlobalDistCorr(
const DataT radius,
const DataT phi,
const DataT z0Tmp,
const DataT z1Tmp,
DataT& ddR,
DataT& ddPhi,
DataT& ddZ,
const AnalyticalFields<DataT>& formulaStruct)
const { calcDistCorr(radius, phi, z0Tmp, z1Tmp, ddR, ddPhi, ddZ, formulaStruct,
false, formulaStruct.getSide()); }
1375 void processGlobalDistCorr(
const DataT radius,
const DataT phi,
const DataT z0Tmp,
const DataT z1Tmp,
DataT& ddR,
DataT& ddPhi,
DataT& ddZ,
const NumericalFields<DataT>& formulaStruct)
const { calcDistCorr(radius, phi, z0Tmp, z1Tmp, ddR, ddPhi, ddZ, formulaStruct,
false, formulaStruct.getSide()); }
1378 void processGlobalDistCorr(
const DataT radius,
const DataT phi,
const DataT z0Tmp, [[maybe_unused]]
const DataT z1Tmp,
DataT& ddR,
DataT& ddPhi,
DataT& ddZ,
const DistCorrInterpolator<DataT>& localDistCorr)
const
1380 ddR = localDistCorr.evaldR(z0Tmp, radius, phi);
1381 ddZ = localDistCorr.evaldZ(z0Tmp, radius, phi);
1382 ddPhi = localDistCorr.evaldRPhi(z0Tmp, radius, phi) / radius;
1386 void dumpElectronTracksToTree(
const std::vector<std::pair<std::vector<
o2::math_utils::Point3D<float>>, std::array<DataT, 3>>>& electronTracks,
const int nSamplingPoints,
const char* outFile)
const;
1389 void setPotentialBoundaryGEMFrameAlongPhi(
const std::function<
DataT(
DataT)>& potentialFunc,
const GEMstack stack,
const bool bottom,
const Side side,
const bool outerFrame =
false);
1395 void initContainer(DataContainer&
data,
const bool initMem =
true);
1397 void initAllBuffers();
1399 void setBoundaryFromIndices(
const std::function<
DataT(
DataT)>& potentialFunc,
const std::vector<std::pair<size_t, float>>&
indices,
const Side side);
1402 std::vector<std::pair<size_t, float>> getPotentialBoundaryGEMFrameAlongRIndices(
const Side side)
const;
1405 std::vector<std::pair<size_t, float>> getPotentialBoundaryGEMFrameAlongPhiIndices(
const GEMstack stack,
const bool bottom,
const Side side,
const bool outerFrame,
const bool noGap =
false)
const;
1407 void setROCMisalignment(
int stackType,
int misalignmentType,
int sector,
const float potMin,
const float potMax);
1408 void fillROCMisalignment(
const std::vector<std::pair<size_t, float>>& indicesTop,
const std::vector<std::pair<size_t, float>>& indicesBottom,
int sector,
int misalignmentType,
const std::pair<float, float>& deltaPotPar);
1411 void initRodAlignmentVoltages(
const MisalignmentType misalignmentType,
const FCType fcType,
const int sector,
const Side side,
const float deltaPot);
1413 void calcGlobalDistCorrIterative(
const DistCorrInterpolator<DataT>& globCorr,
const int maxIter,
const DataT approachZ,
const DataT approachR,
const DataT approachPhi,
const DataT diffCorr,
const SpaceCharge<DataT>* scSCale,
float scale,
const Type type);
1414 void calcGlobalDistCorrIterativeLinearCartesian(
const DistCorrInterpolator<DataT>& globCorr,
const int maxIter,
const DataT approachX,
const DataT approachY,
const DataT approachZ,
const DataT diffCorr,
const SpaceCharge<DataT>* scSCale,
float scale,
const Type type);
This class provides a simple method to store values on a large 3-Dim grid with ROOT io functionality.
Definition of RegularGrid3D class.
This file provides all necesseray classes which are used during the calcution of the distortions and ...
Definition of TriCubic class.
struct for containing simple analytical distortions (as function of local coordinates) and the result...
DataT clampToGridCircular(DataT pos, const unsigned int dim) const
DataT getGridMinZ() const
DataT getRVertex(const size_t vertexY) const
DataT getGridMinPhi() const
DataT getGridMaxPhi() const
DataT getPhiVertex(const size_t vertexZ) const
DataT getGridMinR() const
DataT getInvSpacingR() const
DataT getInvSpacingZ() const
DataT getSpacingPhi() const
DataT getSpacingR() const
DataT getInvSpacingPhi() const
DataT getSpacingZ() const
DataT getGridMaxZ() const
DataT getZVertex(const size_t vertexX) const
DataT clampToGrid(const DataT pos, const unsigned int dim) const
DataT getGridMaxR() const
void setOmegaTauT1T2(const DataT omegaTau=0.32f, const DataT t1=1, const DataT t2=1)
void setC2(const DataT c2)
DataT getLocalDistRPhi(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
void setPotentialBoundaryGEMFrameOROC1TopAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
DataT getZMin(const Side side) const
Get min z position which is used during the calaculations.
void getDistortionsCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &distZ, DataT &distR, DataT &distRPhi) const
void addChargeDensity(const SpaceCharge< DataT > &otherSC)
void fillDensity(const DataT density, const size_t iz, const size_t ir, const size_t iphi, const Side side)
SpaceCharge & operator=(SpaceCharge &&)=default
move assignment
DataT getLocalCorrRPhi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
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
size_t getNearestPhiVertex(const DataT phi, const Side side) const
DataT getLocalCorrR(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void correctElectron(GlobalPosition3D &point)
void initBField(const int field)
static GlobalDistType getGlobalDistType()
void setPotentialBoundaryGEMFrameIROCBottomAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
DataT regulateR(const DataT posR, const Side side) const
set r coordinate between 'RMIN - 4 * GRIDSPACINGR' and 'RMAX + 2 * GRIDSPACINGR'. the r coordinate is...
size_t getPhiBinsGapFrame(const Side side) const
void setROCMisalignmentRotationAlongX(const int sector, const int type, const float potentialMin, const float potentialMax)
int dumpGlobalDistortions(std::string_view file, const Side side, std::string_view option) const
void addGlobalCorrections(const SpaceCharge< DataT > &otherSC, const Side side)
add global corrections from other space charge object
void setPotentialBoundaryInnerRadius(const std::function< DataT(DataT)> &potentialFunc, const Side side)
DataT getEz(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
DataT getRMin(const Side side) const
Get inner radius of tpc.
void setLocalCorrectionsFromFile(std::string_view file, const Side side)
static void setNumericalIntegrationStrategy(const IntegrationStrategy strategy)
static Side getSide(const DataT z)
DataT getEphi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setPotentialBoundaryGEMFrameAlongR(const std::function< DataT(DataT)> &potentialFunc, const Side side)
void setDistortionsCorrectionsAnalytical(const AnalyticalDistCorr< DataT > &formula)
set distortions and corrections by an analytical formula
static constexpr DataT getEzField(const Side side)
Get constant electric field.
DataT getGlobalDistR(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
static DataT getXFromPolar(const DataT r, const DataT phi)
convert radius and phi coordinates from polar coordinates to x cartesian coordinates
void setDeltaVoltageStripsShiftIFC(const float deltaPot, const Side side, const int sector)
DataT getPhiMin(const Side side) const
Get min phi.
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)
DataT getLocalCorrZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
int dumpLocalCorrections(std::string_view file, const Side side, std::string_view option) const
void setGlobalDistortionsFromFile(std::string_view file, const Side side)
bool checkGridFromFile(std::string_view file, std::string_view tree)
void setSimNSector(const int nSectors)
simulate N sectors
void setGlobalCorrectionsFromFile(std::string_view file, const Side side)
void setMetaData(const SCMetaData &meta)
setting meta data for this object
void getCorrectionsAnalytical(const DataT x, const DataT y, const DataT z, const Side side, DataT &corrX, DataT &corrY, DataT &corrZ) const
void setPotentialBoundaryGEMFrameIROCTopAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
bool getSimExBMisalignment() const
DataT getPotentialCyl(const DataT z, const DataT r, const DataT phi, const Side side) const
static int getNThreads()
get the number of threads used for some of the calculations
void setROCMisalignmentRotationAlongY(const int sector, const int type, const float potentialMin, const float potentialMax)
DataT regulatePhi(const DataT posPhi, const Side side) const
set phi coordinate between min phi max phi
int dumpDensity(std::string_view file, const Side side, std::string_view option) const
void setIFCChargeUpRisingPot(const float deltaPot, const float zMaxDeltaPot, const int type, const float zStart, const float offs, const Side side)
bool getUseAnalyticalDistCorr() const
bool getSimEDistortions() const
void getLocalDistortionVectorCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &lvecdistZ, DataT &lvecdistR, DataT &lvecdistRPhi) const
SpaceCharge(SpaceCharge &&)=default
default move constructor
DataT getLocalDistR(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
int dumpPotential(std::string_view file, const Side side, std::string_view option) const
void scalePotential(const DataT scalingFactor, const Side side)
void printMetaData() const
DataT regulateZ(const DataT posZ, const Side side) const
unsigned short getNPhiVertices() const
void setDefaultStaticDistortionsGEMFrameChargeUp(const Side side, const DataT deltaPotential=1000)
DataT getLocalVecCorrZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setAnalyticalCorrectionsDistortionsFromFile(std::string_view inpf)
void fillChargeFromCalDet(const std::vector< CalDet< float > > &calCharge3D)
const auto & getPotential(const Side side) const &
void setPotentialBoundaryGEMFrameOROC1BottomAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
DataT getGlobalCorrZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
const RegularGrid & getGrid3D(const Side side) const
Get the grid object.
void setPotentialBoundaryFromFormula(const AnalyticalFields< DataT > &formulaStruct)
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)
static void setGlobalDistType(const GlobalDistType globalDistType)
void setUseInitialSCDensity(const bool useInitialSCDensity)
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)
DataT getDensity(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setDeltaVoltageStripsShiftOFC(const float deltaPot, const Side side, const int sector)
DataT getGlobalDistRPhi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
@ ShiftedClip
shifted copper rod clip
@ RotatedClip
rotated mylar strips from FC
@ RodShift
shifted copper rod
void calcLocalDistortionsCorrections(const Type type, const ElectricFields &formulaStruct)
@ LocalDistCorr
using local dis/corr interpolator for calculation of global distortions/corrections
@ ElectricalField
using electric field for calculation of global distortions/corrections
int dumpGlobalCorrections(std::string_view file, const Side side, std::string_view option) const
IntegrationStrategy
numerical integration strategys
@ SimpsonIterative
simpon integration, but using an iterative method to approximate the drift path. No straight electron...
@ Trapezoidal
trapezoidal integration (https://en.wikipedia.org/wiki/Trapezoidal_rule). straight electron drift lin...
@ Root
Root integration. 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...
DataT getGridSpacingR(const Side side) const
Get grid spacing in r direction.
static int getStepWidth()
Get the step width which is used for the calculation of the correction/distortions in units of the z-...
DistCorrInterpolator< DataT > getGlobalCorrInterpolator(const Side side) const
DistCorrInterpolator< DataT > getLocalCorrInterpolator(const Side side) const
DataT getLocalDistZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setEFieldFromFormula(const AnalyticalFields< DataT > &formulaStruct)
static void setNThreads(const int nThreads)
set the number of threads used for some of the calculations
void subtractGlobalDistortions(const SpaceCharge< DataT > &otherSC, const Side side)
int dumpLocalDistCorrVectors(std::string_view file, const Side side, std::string_view option) const
DataT getGlobalDistZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void addBoundaryPotential(const SpaceCharge< DataT > &other, const Side side, const float scaling=1)
void getCorrections(const DataT x, const DataT y, const DataT z, const Side side, DataT &corrX, DataT &corrY, DataT &corrZ) const
void fillPotential(const DataT potential, const size_t iz, const size_t ir, const size_t iphi, const Side side)
DataT getGlobalCorrRPhi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void initAfterReadingFromFile()
static void normalizeHistoQVEps0(TH3 &histoIonsPhiRZ)
float getMeanLumi() const
void getDistortions(const DataT x, const DataT y, const DataT z, const Side side, DataT &distX, DataT &distY, DataT &distZ) const
void calcLocalDistortionsCorrectionsRK4(const Type type, const Side side)
void mirrorPotential(const Side sideRef, const Side sideMirrored)
static int getOMPMaxThreads()
int dumpAnalyticalCorrectionsDistortions(TFile &outf) const
static void unsetSimNSector()
unsetting simulation of one sector
DistCorrInterpolator< DataT > getLocalDistInterpolator(const Side side) const
void getLocalCorrectionsCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &lcorrZ, DataT &lcorrR, DataT &lcorrRPhi) const
void fillChargeFromIDCs(std::vector< CalDet< float > > &idcZero, const CalDet< float > &mapIBF, const float ionDriftTimeMS=200, const bool normToPadArea=true)
void scaleChargeDensityStack(const float scalingFactor, const Sector sector, const GEMstack stack)
scaling the space-charge density for given stack
void calcGlobalDistortions(const Fields &formulaStruct, const int maxIterations=3 *sSteps *129)
void setPotentialBoundaryGEMFrameOROC3ToOFCPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
void scaleChargeDensity(const DataT scalingFactor, const Side side)
DataT getGridSpacingPhi(const Side side) const
Get grid spacing in phi direction.
void dumpMetaData(std::string_view file, std::string_view option, const bool overwriteExisting=false) const
static void convertIDCsToCharge(std::vector< CalDet< float > > &idcZero, const CalDet< float > &mapIBF, const float ionDriftTimeMS=200, const bool normToPadArea=true)
@ Fast
interpolation of global corrections (use the global corrections to apply an iterative approach to obt...
@ Standard
classical method (start calculation of global distortion at each voxel in the tpc and follow electron...
size_t getNearestRVertex(const DataT r, const Side side) const
void getDistortionsAnalytical(const DataT x, const DataT y, const DataT z, const Side side, DataT &distX, DataT &distY, DataT &distZ) const
const auto & getDistortionsCorrectionsAnalytical() const &
NumericalFields< DataT > getElectricFieldsInterpolator(const Side side) const
void setDeltaVoltageCopperRodShiftOFC(const float deltaPot, const Side side, const int sector)
void subtractGlobalCorrections(const SpaceCharge< DataT > &otherSC, const Side side)
void getCorrectionsCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &corrZ, DataT &corrR, DataT &corrRPhi) const
void setDensityFromFile(std::string_view file, const Side side)
void setUseAnalyticalDistCorr(const bool useAna)
setting usage of the analytical formula for the distortions and corrections
void getLocalDistortionsCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &ldistZ, DataT &ldistR, DataT &ldistRPhi) const
void dumpToFile(std::string_view file, const Side side, std::string_view option) const
DataT getGridSpacingZ(const Side side) const
Get grid spacing in z direction.
void fillChargeDensityFromCalDet(const std::vector< CalDet< float > > &calSCDensity3D)
void setGlobalDistortions(const std::function< void(int sector, DataT gx, DataT gy, DataT gz, DataT &gCx, DataT &gCy, DataT &gCz)> &gDist, const Side side)
DataT getEr(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
static void setSCDistortionType(SCDistortionType distortionType)
void scaleCorrections(const float scaleFac, const Side side)
int dumpLocalDistortions(std::string_view file, const Side side, std::string_view option) const
unsigned short getNRVertices() const
void setMeanLumi(float lumi)
DataT getPhiMax(const Side side) const
Get max phi.
void setC0C1(const DataT c0, const DataT c1)
DataT getLocalVecDistZ(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setDistortionLookupTables(const DataContainer &distdZ, const DataContainer &distdR, const DataContainer &distdRPhi, const Side side)
static void setGlobalDistCorrMethod(const GlobalDistCorrMethod globalDistCorrMethod)
void setPotentialBoundaryOuterRadius(const std::function< DataT(DataT)> &potentialFunc, const Side side)
void setIFCChargeUpFallingPot(const float deltaPot, const float zMaxDeltaPot, const int type, const float zEnd, const float offs, const Side side)
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 setSimExBMisalignment(const bool simExBMisalignment)
enable/disable calculation of distortions due to ExB misalignment
void calcGlobalCorrectionsEField(const Side side, const int type=3)
calculate the global corrections using the electric fields (interface for python)
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")
void setElectricFieldsFromFile(std::string_view file, const Side side)
static DataT getYFromPolar(const DataT r, const DataT phi)
convert radius and phi coordinates from polar coordinates to y cartesian coordinate
static void setSimpsonNIteratives(const int nIter)
DataT getPhiVertex(const size_t indexPhi, const Side side) const
void calcEField(const Side side)
void setPotentialFromFile(std::string_view file, const Side side)
void fillChargeDensityFromHisto(const TH3 &hisSCDensity3D)
void setPotentialBoundaryGEMFrameOROC3TopAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
void getLocalCorrectionVectorCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &lveccorrZ, DataT &lveccorrR, DataT &lveccorrRPhi) const
void setPotentialBoundaryGEMFrameOROC2TopAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
void setPotential(int iz, int ir, int iphi, Side side, float val)
setting the potential directly for given vertex
void setFromFile(std::string_view file, const Side side)
void setSimEDistortions(const bool simEDistortions)
calculate distortions due to electric fields (space charge, boundary potential...)
DataT getLocalVecCorrR(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
static int getSimpsonNIteratives()
SpaceCharge()
constructor for empty object. Can be used when buffers are loaded from file
void setDeltaVoltageCopperRodShiftIFC(const float deltaPot, const Side side, const int sector)
void setChargeDensityFromFormula(const AnalyticalFields< DataT > &formulaStruct)
void setDeltaVoltageRotatedClipOFC(const float deltaPot, const Side side, const int sector)
void averageDensityPerSector(const Side side)
average the sc density across all sectors per side
void setDeltaVoltageRotatedClipIFC(const float deltaPot, const Side side, const int sector)
DataT getDensityCyl(const DataT z, const DataT r, const DataT phi, const Side side) const
DataT getLocalVecCorrRPhi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void getElectricFieldsCyl(const DataT z, const DataT r, const DataT phi, const Side side, DataT &eZ, DataT &eR, DataT &ePhi) const
float getDCAr(float tgl, const int nPoints, const float phi, float rStart=-1, o2::utils::TreeStreamRedirector *pcstream=nullptr) const
DistCorrInterpolator< DataT > getGlobalDistInterpolator(const Side side) const
DataT getZMax(const Side side) const
Get max z.
void setROCMisalignmentShiftZ(const int sector, const int type, const float potential)
void calculateDistortionsCorrections(const o2::tpc::Side side, const bool calcVectors=false)
void calcGlobalCorrections(const Fields &formulaStruct, const int type=3)
static IntegrationStrategy getNumericalIntegrationStrategy()
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 setPotentialBoundaryGEMFrameIROCToIFCPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
setting the potential from the IROC to the inner field cage
static SCDistortionType getSCDistortionType()
Get the space-charge distortions model.
const auto & getMetaData() const
DataT getLocalVecDistR(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
void setPotentialBoundaryGEMFrameOROC3BottomAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
static GlobalDistCorrMethod getGlobalDistCorrMethod()
void calcLocalDistortionCorrectionVector(const ElectricFields &formulaStruct)
DataT getRVertex(const size_t indexR, const Side side) const
static DataT getRadiusFromCartesian(const DataT x, const DataT y)
convert x and y coordinates from cartesian to the radius in polar coordinates
void setPotentialFromFormula(const AnalyticalFields< DataT > &formulaStruct)
void setLocalDistCorrVectorsFromFile(std::string_view file, const Side side)
void distortElectron(GlobalPosition3D &point, const SpaceCharge< DataT > *scSCale=nullptr, float scale=0) const
DataT getPotential(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
unsigned short getNZVertices() const
void setPotentialBoundaryGEMFrameOROC2BottomAlongPhi(const std::function< DataT(DataT)> &potentialFunc, const Side side)
static DataT getPhiFromCartesian(const DataT x, const DataT y)
convert x and y coordinates from cartesian to phi in polar coordinates
void readMetaData(std::string_view file)
void resetBoundaryPotentialToZeroInRangeZ(float zMin, float zMax, const Side side)
setting the boundary potential to 0 for z<zMin and z>zMax
void downSampleObject(const int nZNew, const int nRNew, const int nPhiNew)
DataT getLocalVecDistRPhi(const size_t iz, const size_t ir, const size_t iphi, const Side side) const
int dumpElectricFields(std::string_view file, const Side side, std::string_view option) const
DataT getRMax(const Side side) const
Get max r.
void setLocalDistortionsFromFile(std::string_view file, const Side side)
static void setNStep(const int nSteps)
void scaleChargeDensitySector(const float scalingFactor, const Sector sector)
DataT getZVertex(const size_t indexZ, const Side side) const
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
void fillChargeDensityFromFile(TFile &fInp, const char *name)
void poissonSolver(const Side side, const DataT stoppingConvergence=1e-6, const int symmetry=0)
GLuint const GLchar * name
GLint GLint GLsizei GLint GLenum GLenum type
GLsizei GLenum const void * indices
GLdouble GLdouble GLdouble z
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t1
GEMstack
TPC GEM stack types.
constexpr unsigned char SIDES
SCDistortionType
Enumerator for setting the space-charge distortion mode.
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Common utility functions.
Marks an empty item in the context.
static constexpr DataT PHIMIN
min phi coordinate
static constexpr DataT ZMIN
min z coordinate
static constexpr DataT getGridSpacingR(const unsigned int nR)
static constexpr DataT getGridSpacingZ(const unsigned int nZ)
static constexpr DataT getGridSpacingPhi(const unsigned int nPhi)
static constexpr DataT RMIN
min radius
unsigned short NZVertices
NRVertices number of vertices in z direction.
unsigned short NPhiVertices
NZVertices number of vertices in r direction.
unsigned short NRVertices
VectorOfTObjectPtrs other
o2::InteractionRecord ir(0, 0)
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))