Project
Loading...
Searching...
No Matches
TrackResiduals.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
17
18#ifndef ALICEO2_TPC_TRACKRESIDUALS_H_
19#define ALICEO2_TPC_TRACKRESIDUALS_H_
20
21#include <memory>
22#include <vector>
23#include <array>
24#include <bitset>
25#include <string>
26
27#include "DataFormatsTPC/Defs.h"
31
32#include "TTree.h"
33#include "TFile.h"
34
35namespace o2
36{
37namespace tpc
38{
39
46{
47 public:
49 TrackResiduals() = default;
50
54
56 enum { VoxZ,
60 VoxDim = 3,
61 VoxHDim = 4 };
62
64 enum { ResX,
69
71 enum { DistDone = 1 << 0,
72 DispDone = 1 << 1,
73 SmoothDone = 1 << 2,
74 Masked = 1 << 7 };
75
77 Gaussian };
78
80 struct VoxRes {
81 VoxRes() = default;
82 std::array<float, ResDim> D{};
83 std::array<float, ResDim> E{};
84 std::array<float, ResDim> DS{};
85 std::array<float, ResDim> DC{};
86 float EXYCorr{0.f};
87 float dYSigMAD{0.f};
88 float dZSigLTM{0.f};
89 std::array<float, VoxHDim> stat{};
90 std::array<unsigned char, VoxDim> bvox{};
91 unsigned char bsec{0};
92 unsigned char flags{0};
94 };
95
98 struct LocalResid {
99 LocalResid() = default;
100 LocalResid(short dyIn, short dzIn, short tgSlpIn, std::array<unsigned char, VoxDim> bvoxIn) : dy(dyIn), dz(dzIn), tgSlp(tgSlpIn), bvox(bvoxIn) {}
101 short dy{0};
102 short dz{0};
103 short tgSlp{0};
104 std::array<unsigned char, VoxDim> bvox{};
106 };
107
109 struct VoxStats {
110 VoxStats() = default;
111 std::array<float, VoxDim> meanPos{};
112 float nEntries{0.f};
114 };
115
116 // -------------------------------------- initialization --------------------------------------------------
119 void init(bool doBinning = true);
121 void initBinning();
125 void initResultsContainer(int iSec);
129 void reset();
130
131 // -------------------------------------- settings --------------------------------------------------
133 void setPrintMemoryUsage() { mPrintMem = true; }
142 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);
143
146 void setSmoothPol2(int dim, bool flag) { mSmoothPol2[dim] = flag; }
147
148 void setVdriftCorr(float corr) { mEffVdriftCorr = corr; }
149
150 void setT0Corr(float corr) { mEffT0Corr = corr; }
151
152 // -------------------------------------- I/O --------------------------------------------------
153
154 std::vector<LocalResid>& getLocalResVec() { return mLocalResidualsIn; }
155 std::vector<VoxStats>** getVoxStatPtr() { return &mVoxStatsInPtr; }
156
157 const std::array<std::vector<VoxRes>, SECTORSPERSIDE * SIDES>& getVoxelResults() const { return mVoxelResults; }
158
159 // -------------------------------------- steering functions --------------------------------------------------
160
163 void processSectorResiduals(Int_t iSec);
164
170 void processVoxelResiduals(std::vector<float>& dy, std::vector<float>& dz, std::vector<float>& tg, VoxRes& resVox);
171
176 void processVoxelDispersions(std::vector<float>& tg, std::vector<float>& dy, VoxRes& resVox);
177
182 int validateVoxels(int iSec);
183
186 void smooth(int iSec);
187
188 // -------------------------------------- statistics --------------------------------------------------
189
198 float fitPoly1Robust(std::vector<float>& x, std::vector<float>& y, std::array<float, 2>& res, std::array<float, 3>& err, float cutLTM) const;
199
204 float getMAD2Sigma(std::vector<float> data) const;
205
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;
216
226 float roFunc(int nPoints, int offset, const std::vector<float>& x, const std::vector<float>& y, float b, float& aa) const;
227
235 float selectKthMin(const int k, std::vector<float>& data) const;
236
245 bool getSmoothEstimate(int iSec, float x, float p, float z, std::array<float, ResDim>& res, int whichDim = 0);
246
251 double getKernelWeight(std::array<double, 3> u2vec) const;
252
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);
262
267 static bool fitPoly1(int nCl, std::array<float, param::NPadRows>& x, std::array<float, param::NPadRows>& y, std::array<float, 2>& res);
268
269 // -------------------------------------- binning / geometry --------------------------------------------------
270
273 void setNXBins(int nBins) { mNXBins = nBins; }
274 int getNXBins() const { return mNXBins; }
275
278 void setNY2XBins(int nBins) { mNY2XBins = nBins; }
279 int getNY2XBins() const { return mNY2XBins; }
280
283 void setNZ2XBins(int nBins) { mNZ2XBins = nBins; }
284 int getNZ2XBins() const { return mNZ2XBins; }
285
287 int getNVoxelsPerSector() const { return mNVoxPerSector; }
288
291 void setY2XBinning(const std::vector<float>& binning);
292
295 void setZ2XBinning(const std::vector<float>& binning);
296
302 size_t getGlbVoxBin(int ix, int ip, int iz) const;
303
307 size_t getGlbVoxBin(const std::array<unsigned char, VoxDim>& bvox) const;
308
318 void getVoxelCoordinates(int isec, int ix, int ip, int iz, float& x, float& p, float& z) const;
319
323 float getX(int i) const;
324
329 float getY2X(int ix, int ip) const;
330
334 float getZ2X(int iz) const;
335
340 bool getXBinIgnored(int iSec, int bin) const { return mXBinsIgnore[iSec].test(bin); }
341
349 void findVoxel(float x, float y2x, float z2x, int& ix, int& ip, int& iz) const;
350
352 bool findVoxelBin(int secID, float x, float y, float z, std::array<unsigned char, VoxDim>& bvox) const;
353
357 int getXBin(float x) const;
358
362 int getXBinExact(float x) const;
363
367 int getRowID(float x) const;
368
374 int getY2XBinExact(float y2x, int ix) const;
375
381 int getY2XBin(float y2x, int ix) const;
382
387 int getZ2XBinExact(float z2x) const;
388
393 int getZ2XBin(float z2x) const;
394
398 float getDXI(int ix) const;
399
404 float getDY2XI(int ix, int iy = 0) const;
405
409 float getDZ2X(int iz = 0) const;
410
414 float getDZ2XI(int iz = 0) const;
415
416 // -------------------------------------- debugging --------------------------------------------------
417
419 void printMem() const;
420
423 void dumpResults(int iSec);
424
426 void createOutputFile(const char* filename = "debugVoxRes.root");
427
429 void closeOutputFile();
430
432 TFile* getOutputFilePtr() { return mFileOut.get(); }
433
435 void setStats(const std::vector<TrackResiduals::VoxStats>& statsIn, int iSec);
436
438 void fillStats(int iSec);
439
441 void clear();
442
444 TTree* getOutputTree() { return mTreeOut.get(); }
445
447 void setAdhocScalingFactorX(const std::array<float, 2>& scaling) { mAdhocScalingX = scaling; }
448
450 void doAdhocCorrectionZ2X(bool corr) { mDoAdhocCorrectionZ2X = corr; }
451
452 private:
453 std::bitset<SECTORSPERSIDE * SIDES> mInitResultsContainer{};
454
455 // some constants
456 static constexpr float sFloatEps{1.e-7f};
457 static constexpr float sDeadZone{1.5f};
458 static constexpr int sSmtLinDim{4};
459 static constexpr int sMaxSmtDim{7};
460
461 // settings
462 const SpacePointsCalibConfParam* mParams = nullptr;
463
464 // input data
465 std::vector<LocalResid> mLocalResidualsIn;
466 std::vector<VoxStats> mVoxStatsIn, *mVoxStatsInPtr{&mVoxStatsIn};
467 // output data
468 std::unique_ptr<TFile> mFileOut;
469 std::unique_ptr<TTree> mTreeOut;
470 // status flags
471 bool mIsInitialized{};
472 bool mPrintMem{};
473 // binning
474 int mNXBins{param::NPadRows};
475 int mNY2XBins{param::NY2XBins};
476 int mNZ2XBins{param::NZ2XBins};
477 int mNVoxPerSector{};
478 float mDX{};
479 float mDXI{};
480 std::vector<float> mMaxY2X{};
481 std::vector<float> mDY2X{};
482 std::vector<float> mDY2XI{};
483 std::vector<float> mY2XBinsDH{};
484 std::vector<float> mY2XBinsDI{};
485 std::vector<float> mY2XBinsCenter{};
486 float mDZ2X{};
487 float mDZ2XI{};
488 std::vector<float> mZ2XBinsDH{};
489 std::vector<float> mZ2XBinsDI{};
490 std::vector<float> mZ2XBinsCenter{};
491 float mMaxZ2X{1.f};
492 std::array<bool, VoxDim> mUniformBins{true, true, true};
493 // smoothing
495 bool mUseErrInSmoothing{true};
496 std::array<bool, VoxDim> mSmoothPol2{};
497 std::array<int, SECTORSPERSIDE * SIDES> mNSmoothingFailedBins{};
498 std::array<int, VoxDim> mStepKern{};
499 std::array<float, VoxDim> mKernelScaleEdge{};
500 std::array<float, VoxDim> mKernelWInv{};
501 std::array<double, ResDim * sMaxSmtDim> mLastSmoothingRes{};
502 // calibrated parameters
503 float mEffVdriftCorr{0.f};
504 float mEffT0Corr{0.f};
505 // (intermediate) results
506 std::array<std::bitset<param::NPadRows>, SECTORSPERSIDE * SIDES> mXBinsIgnore{};
507 std::array<std::array<float, param::NPadRows>, SECTORSPERSIDE * SIDES> mValidFracXBins{};
508 std::array<std::vector<VoxRes>, SECTORSPERSIDE * SIDES> mVoxelResults{};
509 VoxRes mVoxelResultsOut{};
510 VoxRes* mVoxelResultsOutPtr{&mVoxelResultsOut};
511 std::array<float, 2> mAdhocScalingX{0, 0};
512 bool mDoAdhocCorrectionZ2X{false};
513
514 ClassDefNV(TrackResiduals, 3);
515};
516
517//_____________________________________________________
518inline size_t TrackResiduals::getGlbVoxBin(const std::array<unsigned char, VoxDim>& bvox) const
519{
520 return bvox[VoxX] + (bvox[VoxF] + bvox[VoxZ] * mNY2XBins) * mNXBins;
521}
522
523//_____________________________________________________
524inline size_t TrackResiduals::getGlbVoxBin(int ix, int ip, int iz) const
525{
526 return ix + (ip + iz * mNY2XBins) * mNXBins;
527}
528
529//_____________________________________________________
530inline void TrackResiduals::getVoxelCoordinates(int isec, int ix, int ip, int iz, float& x, float& p, float& z) const
531{
532 x = getX(ix);
533 p = getY2X(ix, ip);
534 z = getZ2X(iz);
535 if (isec >= SECTORSPERSIDE) {
536 z = -z;
537 }
538}
539
540//_____________________________________________________
541inline float TrackResiduals::getDXI(int ix) const
542{
543 if (mUniformBins[VoxX]) {
544 return mDXI;
545 } else {
546 if (ix < param::NRowsPerROC[0]) {
547 // we are in the IROC
548 return 1.f / param::RowDX[0];
549 } else if (ix > param::NRowsAccumulated[param::NROCTypes - 2]) {
550 // we are in the last OROC
551 return 1.f / param::RowDX[param::NROCTypes - 1];
552 } else if (ix < param::NRowsAccumulated[1]) {
553 // OROC1
554 return 1.f / param::RowDX[1];
555 } else {
556 // OROC2
557 return 1.f / param::RowDX[2];
558 }
559 }
560}
561
562//_____________________________________________________
563inline float TrackResiduals::getX(int i) const
564{
565 return mUniformBins[VoxX] ? param::MinX + (i + 0.5) * mDX : param::RowX[i];
566}
567
568//_____________________________________________________
569inline float TrackResiduals::getY2X(int ix, int ip) const
570{
571 if (mUniformBins[VoxF]) {
572 return (0.5f + ip) * mDY2X[ix] - mMaxY2X[ix];
573 }
574 return mMaxY2X[ix] * mY2XBinsCenter[ip];
575}
576
577//_____________________________________________________
578inline float TrackResiduals::getZ2X(int iz) const
579{
580 if (mUniformBins[VoxZ]) {
581 return (0.5f + iz) * getDZ2X();
582 }
583 return mZ2XBinsCenter[iz];
584}
585
586//_____________________________________________________
587inline float TrackResiduals::getDY2XI(int ix, int iy) const
588{
589 if (mUniformBins[VoxF]) {
590 return mDY2XI[ix];
591 }
592 return mY2XBinsDI[iy] / mMaxY2X[ix];
593}
594
595//_____________________________________________________
596inline float TrackResiduals::getDZ2X(int iz) const
597{
598 if (mUniformBins[VoxZ]) {
599 return mDZ2X;
600 }
601 return 2.f * mZ2XBinsDH[iz];
602}
603
604//_____________________________________________________
605inline float TrackResiduals::getDZ2XI(int iz) const
606{
607 if (mUniformBins[VoxZ]) {
608 return mDZ2XI;
609 }
610 return mZ2XBinsDI[iz];
611}
612
613//_____________________________________________________
614inline void TrackResiduals::findVoxel(float x, float y2x, float z2x, int& ix, int& ip, int& iz) const
615{
616 ix = getXBin(x);
617 ip = getY2XBin(y2x, ix);
618 iz = getZ2XBin(z2x);
619}
620
621//_____________________________________________________
622inline int TrackResiduals::getXBinExact(float x) const
623{
624 if (mUniformBins[VoxX]) {
625 int ix = (x - param::MinX) * mDXI;
626 return (ix < 0 || ix >= mNXBins) ? -1 : ix;
627 }
628 return getRowID(x);
629}
630
631//_____________________________________________________
632inline int TrackResiduals::getXBin(float x) const
633{
634 int bx = getXBinExact(x);
635 return bx > -1 ? (bx < mNXBins ? bx : mNXBins - 1) : 0;
636}
637
638//_____________________________________________________
639inline int TrackResiduals::getY2XBinExact(float y2x, int ix) const
640{
641 if (y2x < -mMaxY2X[ix]) {
642 return -1;
643 }
644 if (y2x > mMaxY2X[ix]) {
645 return mNY2XBins;
646 }
647 if (mUniformBins[VoxF]) {
648 return static_cast<int>((y2x + mMaxY2X[ix]) * getDY2XI(ix));
649 }
650 y2x /= mMaxY2X[ix];
651 for (int iBin = 0; iBin < mNY2XBins; ++iBin) {
652 if (y2x < mY2XBinsCenter[iBin] + mY2XBinsDH[iBin]) {
653 return iBin;
654 }
655 }
656 return mNY2XBins;
657}
658
659//_____________________________________________________
660inline int TrackResiduals::getY2XBin(float y2x, int ix) const
661{
662 int bp = getY2XBinExact(y2x, ix);
663 return bp > -1 ? (bp < mNY2XBins ? bp : mNY2XBins - 1) : 0;
664}
665
666//_____________________________________________________
667inline int TrackResiduals::getZ2XBinExact(float z2x) const
668{
669 if (mUniformBins[VoxZ]) {
670 float bz = z2x * getDZ2XI();
671 if (bz >= mNZ2XBins) {
672 return -1;
673 }
674 if (bz < 0) {
675 // accounting for clusters which were moved to the wrong side
676 // TODO: how can this happen?
677 bz = 0;
678 }
679 return static_cast<int>(bz);
680 }
681 for (int iBin = 0; iBin < mNZ2XBins; ++iBin) {
682 if (z2x < mZ2XBinsCenter[iBin] + mZ2XBinsDH[iBin]) {
683 return iBin;
684 }
685 }
686 return -1;
687}
688
689//_____________________________________________________
690inline int TrackResiduals::getZ2XBin(float z2x) const
691{
692 int iz = getZ2XBinExact(z2x);
693 return iz < 0 ? mNZ2XBins - 1 : iz;
694}
695
696} // namespace tpc
697
698} // namespace o2
699
700#endif
int32_t i
uint32_t res
Definition RawData.h:0
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)
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
@ VoxDim
dimensionality of the voxels
@ VoxV
voxel dispersions
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
float getX(int i) 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.
void setAdhocScalingFactorX(const std::array< float, 2 > &scaling)
Ad-hoc radial scaling factor A/C-Side.
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)
@ ResD
index for dispersions
@ DispDone
voxel dispersions have been processed
@ DistDone
voxel residuals have been processed
@ SmoothDone
voxel has been smoothed
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
float getDZ2X(int iz=0) const
void processVoxelResiduals(std::vector< float > &dy, std::vector< float > &dz, std::vector< float > &tg, VoxRes &resVox)
void doAdhocCorrectionZ2X(bool corr)
Ad-hoc correction of Z/X.
const std::array< std::vector< VoxRes >, SECTORSPERSIDE *SIDES > & getVoxelResults() const
std::vector< VoxStats > ** getVoxStatPtr()
GLint GLenum GLint x
Definition glcorearb.h:403
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLboolean * data
Definition glcorearb.h:298
GLintptr offset
Definition glcorearb.h:660
GLbitfield flags
Definition glcorearb.h:1570
GLboolean r
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
constexpr unsigned char SIDES
Definition Defs.h:41
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string filename()
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
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