18#ifndef ALICEO2_TPC_TRACKRESIDUALS_H_
19#define ALICEO2_TPC_TRACKRESIDUALS_H_
82 std::array<float, ResDim>
D{};
83 std::array<float, ResDim>
E{};
84 std::array<float, ResDim>
DS{};
85 std::array<float, ResDim>
DC{};
89 std::array<float, VoxHDim>
stat{};
90 std::array<unsigned char, VoxDim>
bvox{};
100 LocalResid(
short dyIn,
short dzIn,
short tgSlpIn, std::array<unsigned char, VoxDim> bvoxIn) :
dy(dyIn),
dz(dzIn),
tgSlp(tgSlpIn),
bvox(bvoxIn) {}
104 std::array<unsigned char, VoxDim>
bvox{};
119 void init(
bool doBinning =
true);
170 void processVoxelResiduals(std::vector<float>& dy, std::vector<float>& dz, std::vector<float>& tg, VoxRes& resVox);
198 float fitPoly1Robust(std::vector<float>&
x, std::vector<float>&
y, std::array<float, 2>&
res, std::array<float, 3>& err,
float cutLTM)
const;
215 void medFit(
int nPoints,
int offset,
const std::vector<float>&
x,
const std::vector<float>&
y,
float&
a,
float&
b, std::array<float, 3>& err)
const;
226 float roFunc(
int nPoints,
int offset,
const std::vector<float>&
x,
const std::vector<float>&
y,
float b,
float&
aa)
const;
245 bool getSmoothEstimate(
int iSec,
float x,
float p,
float z, std::array<float, ResDim>&
res,
int whichDim = 0);
261 static void fitCircle(
int nCl, std::array<float, param::NPadRows>&
x, std::array<float, param::NPadRows>&
y,
float& xc,
float& yc,
float&
r, std::array<float, param::NPadRows>& residHelixY);
267 static bool fitPoly1(
int nCl, std::array<float, param::NPadRows>&
x, std::array<float, param::NPadRows>&
y, std::array<float, 2>&
res);
307 size_t getGlbVoxBin(
const std::array<unsigned char, VoxDim>& bvox)
const;
323 float getX(
int i)
const;
329 float getY2X(
int ix,
int ip)
const;
334 float getZ2X(
int iz)
const;
340 bool getXBinIgnored(
int iSec,
int bin)
const {
return mXBinsIgnore[iSec].test(bin); }
349 void findVoxel(
float x,
float y2x,
float z2x,
int& ix,
int& ip,
int& iz)
const;
352 bool findVoxelBin(
int secID,
float x,
float y,
float z, std::array<unsigned char, VoxDim>& bvox)
const;
398 float getDXI(
int ix)
const;
404 float getDY2XI(
int ix,
int iy = 0)
const;
409 float getDZ2X(
int iz = 0)
const;
435 void setStats(
const std::vector<TrackResiduals::VoxStats>& statsIn,
int iSec);
447 std::bitset<SECTORSPERSIDE * SIDES> mInitResultsContainer{};
450 static constexpr float sFloatEps{1.e-7f};
451 static constexpr float sDeadZone{1.5f};
452 static constexpr int sSmtLinDim{4};
453 static constexpr int sMaxSmtDim{7};
456 const SpacePointsCalibConfParam* mParams =
nullptr;
459 std::vector<LocalResid> mLocalResidualsIn;
460 std::vector<VoxStats> mVoxStatsIn, *mVoxStatsInPtr{&mVoxStatsIn};
462 std::unique_ptr<TFile> mFileOut;
463 std::unique_ptr<TTree> mTreeOut;
465 bool mIsInitialized{};
468 int mNXBins{param::NPadRows};
469 int mNY2XBins{param::NY2XBins};
470 int mNZ2XBins{param::NZ2XBins};
471 int mNVoxPerSector{};
474 std::vector<float> mMaxY2X{};
475 std::vector<float> mDY2X{};
476 std::vector<float> mDY2XI{};
477 std::vector<float> mY2XBinsDH{};
478 std::vector<float> mY2XBinsDI{};
479 std::vector<float> mY2XBinsCenter{};
482 std::vector<float> mZ2XBinsDH{};
483 std::vector<float> mZ2XBinsDI{};
484 std::vector<float> mZ2XBinsCenter{};
486 std::array<bool, VoxDim> mUniformBins{
true,
true,
true};
489 bool mUseErrInSmoothing{
true};
490 std::array<bool, VoxDim> mSmoothPol2{};
491 std::array<int, SECTORSPERSIDE * SIDES> mNSmoothingFailedBins{};
492 std::array<int, VoxDim> mStepKern{};
493 std::array<float, VoxDim> mKernelScaleEdge{};
494 std::array<float, VoxDim> mKernelWInv{};
495 std::array<double, ResDim * sMaxSmtDim> mLastSmoothingRes{};
497 float mEffVdriftCorr{0.f};
498 float mEffT0Corr{0.f};
503 VoxRes mVoxelResultsOut{};
504 VoxRes* mVoxelResultsOutPtr{&mVoxelResultsOut};
512 return bvox[
VoxX] + (bvox[
VoxF] + bvox[
VoxZ] * mNY2XBins) * mNXBins;
518 return ix + (ip + iz * mNY2XBins) * mNXBins;
535 if (mUniformBins[
VoxX]) {
538 if (ix < param::NRowsPerROC[0]) {
540 return 1.f / param::RowDX[0];
541 }
else if (ix > param::NRowsAccumulated[param::NROCTypes - 2]) {
543 return 1.f / param::RowDX[param::NROCTypes - 1];
544 }
else if (ix < param::NRowsAccumulated[1]) {
546 return 1.f / param::RowDX[1];
549 return 1.f / param::RowDX[2];
557 return mUniformBins[
VoxX] ? param::MinX + (
i + 0.5) * mDX : param::RowX[
i];
563 if (mUniformBins[
VoxF]) {
564 return (0.5f + ip) * mDY2X[ix] - mMaxY2X[ix];
566 return mMaxY2X[ix] * mY2XBinsCenter[ip];
572 if (mUniformBins[
VoxZ]) {
573 return (0.5f + iz) *
getDZ2X();
575 return mZ2XBinsCenter[iz];
581 if (mUniformBins[
VoxF]) {
584 return mY2XBinsDI[iy] / mMaxY2X[ix];
590 if (mUniformBins[
VoxZ]) {
593 return 2.f * mZ2XBinsDH[iz];
599 if (mUniformBins[
VoxZ]) {
602 return mZ2XBinsDI[iz];
616 if (mUniformBins[
VoxX]) {
617 int ix = (
x - param::MinX) * mDXI;
618 return (ix < 0 || ix >= mNXBins) ? -1 : ix;
627 return bx > -1 ? (bx < mNXBins ? bx : mNXBins - 1) : 0;
633 if (y2x < -mMaxY2X[ix]) {
636 if (y2x > mMaxY2X[ix]) {
639 if (mUniformBins[
VoxF]) {
640 return static_cast<int>((y2x + mMaxY2X[ix]) *
getDY2XI(ix));
643 for (
int iBin = 0; iBin < mNY2XBins; ++iBin) {
644 if (y2x < mY2XBinsCenter[iBin] + mY2XBinsDH[iBin]) {
655 return bp > -1 ? (bp < mNY2XBins ? bp : mNY2XBins - 1) : 0;
661 if (mUniformBins[
VoxZ]) {
663 if (bz >= mNZ2XBins) {
671 return static_cast<int>(bz);
673 for (
int iBin = 0; iBin < mNZ2XBins; ++iBin) {
674 if (z2x < mZ2XBinsCenter[iBin] + mZ2XBinsDH[iBin]) {
685 return iz < 0 ? mNZ2XBins - 1 : iz;
Parameters used for TPC space point calibration.
Definition of the TrackInterpolation class.
std::vector< LocalResid > & getLocalResVec()
TrackResiduals()=default
Default constructor.
float getDZ2XI(int iz=0) const
void setNXBins(int nBins)
void setVdriftCorr(float corr)
void setKernelType(KernelType kernel=KernelType::Epanechnikov, float bwX=2.1f, float bwP=2.1f, float bwZ=1.7f, float scX=1.f, float scP=1.f, float scZ=1.f)
int getZ2XBin(float z2x) const
bool getSmoothEstimate(int iSec, float x, float p, float z, std::array< float, ResDim > &res, int whichDim=0)
int getXBinExact(float x) const
void setY2XBinning(const std::vector< float > &binning)
@ DispDone
voxel dispersions have been processed
@ DistDone
voxel residuals have been processed
@ SmoothDone
voxel has been smoothed
void dumpResults(int iSec)
int getXBin(float x) const
void processSectorResiduals(Int_t iSec)
static bool fitPoly1(int nCl, std::array< float, param::NPadRows > &x, std::array< float, param::NPadRows > &y, std::array< float, 2 > &res)
float getZ2X(int iz) const
size_t getGlbVoxBin(int ix, int ip, int iz) const
void setStats(const std::vector< TrackResiduals::VoxStats > &statsIn, int iSec)
Set the voxel statistics directly from outside.
void closeOutputFile()
Closes the file with the debug output.
void clear()
clear member to be able to process new sector or new input files
void fillStats(int iSec)
Fill statistics from TTree.
void setSmoothPol2(int dim, bool flag)
TrackResiduals(const TrackResiduals &)=delete
Copying and assigning is forbidden.
TFile * getOutputFilePtr()
Allow to access the output file from outside.
void printMem() const
Prints the current memory usage.
void setZ2XBinning(const std::vector< float > &binning)
void initResultsContainer(int iSec)
int getNVoxelsPerSector() const
Get the total number of voxels per TPC sector (mNXBins * mNY2XBins * mNZ2XBins)
float getY2X(int ix, int ip) const
int getY2XBin(float y2x, int ix) const
int getRowID(float x) const
void createOutputFile(const char *filename="debugVoxRes.root")
Creates a file for the debug output.
void setPrintMemoryUsage()
Sets a flag to print the memory usage at certain points in the program for performance studies.
float getDY2XI(int ix, int iy=0) const
void findVoxel(float x, float y2x, float z2x, int &ix, int &ip, int &iz) const
float roFunc(int nPoints, int offset, const std::vector< float > &x, const std::vector< float > &y, float b, float &aa) const
void medFit(int nPoints, int offset, const std::vector< float > &x, const std::vector< float > &y, float &a, float &b, std::array< float, 3 > &err) const
static void fitCircle(int nCl, std::array< float, param::NPadRows > &x, std::array< float, param::NPadRows > &y, float &xc, float &yc, float &r, std::array< float, param::NPadRows > &residHelixY)
TrackResiduals & operator=(const TrackResiduals &)=delete
float getDXI(int ix) const
float fitPoly1Robust(std::vector< float > &x, std::vector< float > &y, std::array< float, 2 > &res, std::array< float, 3 > &err, float cutLTM) const
void getVoxelCoordinates(int isec, int ix, int ip, int iz, float &x, float &p, float &z) const
int getY2XBinExact(float y2x, int ix) const
bool findVoxelBin(int secID, float x, float y, float z, std::array< unsigned char, VoxDim > &bvox) const
Calculates the bin indices for given x, y, z in sector coordinates.
bool getXBinIgnored(int iSec, int bin) const
void initBinning()
Initializes the binning in X, Y/X and Z.
float getMAD2Sigma(std::vector< float > data) const
void reset()
Resets all (also intermediate) results.
int getZ2XBinExact(float z2x) const
double getKernelWeight(std::array< double, 3 > u2vec) const
void processVoxelDispersions(std::vector< float > &tg, std::vector< float > &dy, VoxRes &resVox)
void init(bool doBinning=true)
@ VoxDim
dimensionality of the voxels
@ ResD
index for dispersions
float selectKthMin(const int k, std::vector< float > &data) const
void setNY2XBins(int nBins)
void initVoxelStats()
Initialize the statistics for the local residuals container.
void setNZ2XBins(int nBins)
void setT0Corr(float corr)
TTree * getOutputTree()
output tree
int validateVoxels(int iSec)
float getDZ2X(int iz=0) const
void processVoxelResiduals(std::vector< float > &dy, std::vector< float > &dz, std::vector< float > &tg, VoxRes &resVox)
const std::array< std::vector< VoxRes >, SECTORSPERSIDE *SIDES > & getVoxelResults() const
std::vector< VoxStats > ** getVoxStatPtr()
GLboolean GLboolean GLboolean b
GLboolean GLboolean GLboolean GLboolean a
GLdouble GLdouble GLdouble z
constexpr unsigned char SECTORSPERSIDE
constexpr unsigned char SIDES
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LocalResid(short dyIn, short dzIn, short tgSlpIn, std::array< unsigned char, VoxDim > bvoxIn)
std::array< unsigned char, VoxDim > bvox
voxel identifier: VoxZ, VoxF, VoxX
short dy
residual in y, ranges from -param::sMaxResid to +param::sMaxResid
ClassDefNV(LocalResid, 1)
short tgSlp
tangens of the phi angle between padrow and track, ranges from -param::MaxTgSlp to +param::MaxTgSlp
short dz
residual in z, ranges from -param::sMaxResid to +param::sMaxResid
Structure which gets filled with the results for each voxel.
std::array< float, ResDim > DC
Cheb parameterized residual.
std::array< unsigned char, VoxDim > bvox
voxel identifier: VoxZ, VoxF, VoxX
std::array< float, ResDim > E
their errors
float dYSigMAD
MAD estimator of dY sigma (dispersion after slope removal)
std::array< float, ResDim > D
values of extracted distortions
float dZSigLTM
Z sigma from unbinned LTM estimator.
std::array< float, ResDim > DS
smoothed residual
unsigned char bsec
sector ID (0-35)
float EXYCorr
correlation between extracted X and Y
std::array< float, VoxHDim > stat
statistics: averages of each voxel dimension + entries
Structure which holds the statistics for each voxel.
std::array< float, VoxDim > meanPos