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);
267 static bool fitPoly1(int nCl, std::array<float, param::NPadRows>& x, std::array<float, param::NPadRows>& y, std::array<float, 2>& res);
269
270 // -------------------------------------- binning / geometry --------------------------------------------------
271
274 void setNXBins(int nBins) { mNXBins = nBins; }
275 int getNXBins() const { return mNXBins; }
276
279 void setNY2XBins(int nBins) { mNY2XBins = nBins; }
280 int getNY2XBins() const { return mNY2XBins; }
281
284 void setNZ2XBins(int nBins) { mNZ2XBins = nBins; }
285 int getNZ2XBins() const { return mNZ2XBins; }
286
288 int getNVoxelsPerSector() const { return mNVoxPerSector; }
289
292 void setY2XBinning(const std::vector<float>& binning);
293
296 void setZ2XBinning(const std::vector<float>& binning);
297
303 size_t getGlbVoxBin(int ix, int ip, int iz) const;
304
308 size_t getGlbVoxBin(const std::array<unsigned char, VoxDim>& bvox) const;
309
319 void getVoxelCoordinates(int isec, int ix, int ip, int iz, float& x, float& p, float& z) const;
320
324 float getX(int ix) const;
325
329 float getMaxY2X(int ix) const;
330
335 float getY2X(int ix, int ip) const;
336
340 float getZ2X(int iz) const;
341
346 bool getXBinIgnored(int iSec, int bin) const { return mXBinsIgnore[iSec].test(bin); }
347
355 void findVoxel(float x, float y2x, float z2x, int& ix, int& ip, int& iz) const;
356
358 bool findVoxelBin(int secID, float x, float y, float z, std::array<unsigned char, VoxDim>& bvox) const;
359
363 int getXBin(float x) const;
364
368 int getXBinExact(float x) const;
369
373 int getRowID(float x) const;
374
380 int getY2XBinExact(float y2x, int ix) const;
381
387 int getY2XBin(float y2x, int ix) const;
388
393 int getZ2XBinExact(float z2x) const;
394
399 int getZ2XBin(float z2x) const;
400
404 float getDXI(int ix) const;
405
410 float getDY2XI(int ix, int iy = 0) const;
411
415 float getDZ2X(int iz = 0) const;
416
420 float getDZ2XI(int iz = 0) const;
421
422 // -------------------------------------- debugging --------------------------------------------------
423
425 void printMem() const;
426
429 void dumpResults(int iSec);
430
432 void createOutputFile(const char* filename = "debugVoxRes.root");
433
435 void closeOutputFile();
436
438 TFile* getOutputFilePtr() { return mFileOut.get(); }
439
441 void setStats(const std::vector<TrackResiduals::VoxStats>& statsIn, int iSec);
442
444 void fillStats(int iSec);
445
447 void clear();
448
450 TTree* getOutputTree() { return mTreeOut.get(); }
451
453 void setAdhocScalingFactorX(const std::array<float, 2>& scaling) { mAdhocScalingX = scaling; }
454
456 void doAdhocCorrectionZ2X(bool corr) { mDoAdhocCorrectionZ2X = corr; }
457
458 private:
459 std::bitset<SECTORSPERSIDE * SIDES> mInitResultsContainer{};
460
461 // some constants
462 static constexpr float sFloatEps{1.e-7f};
463 static constexpr float sDeadZone{1.5f};
464 static constexpr int sSmtLinDim{4};
465 static constexpr int sMaxSmtDim{7};
466
467 // settings
468 const SpacePointsCalibConfParam* mParams = nullptr;
469
470 // input data
471 std::vector<LocalResid> mLocalResidualsIn;
472 std::vector<VoxStats> mVoxStatsIn, *mVoxStatsInPtr{&mVoxStatsIn};
473 // output data
474 std::unique_ptr<TFile> mFileOut;
475 std::unique_ptr<TTree> mTreeOut;
476 // status flags
477 bool mIsInitialized{};
478 bool mPrintMem{};
479 // binning
480 int mNXBins{param::NPadRows};
481 int mNY2XBins{param::NY2XBins};
482 int mNZ2XBins{param::NZ2XBins};
483 int mNVoxPerSector{};
484 float mDX{};
485 float mDXI{};
486 std::vector<float> mMaxY2X{};
487 std::vector<float> mDY2X{};
488 std::vector<float> mDY2XI{};
489 std::vector<float> mY2XBinsDH{};
490 std::vector<float> mY2XBinsDI{};
491 std::vector<float> mY2XBinsCenter{};
492 float mDZ2X{};
493 float mDZ2XI{};
494 std::vector<float> mZ2XBinsDH{};
495 std::vector<float> mZ2XBinsDI{};
496 std::vector<float> mZ2XBinsCenter{};
497 float mMaxZ2X{1.f};
498 std::array<bool, VoxDim> mUniformBins{true, true, true};
499 // smoothing
501 bool mUseErrInSmoothing{true};
502 std::array<bool, VoxDim> mSmoothPol2{};
503 std::array<int, SECTORSPERSIDE * SIDES> mNSmoothingFailedBins{};
504 std::array<int, VoxDim> mStepKern{};
505 std::array<float, VoxDim> mKernelScaleEdge{};
506 std::array<float, VoxDim> mKernelWInv{};
507 std::array<double, ResDim * sMaxSmtDim> mLastSmoothingRes{};
508 // calibrated parameters
509 float mEffVdriftCorr{0.f};
510 float mEffT0Corr{0.f};
511 // (intermediate) results
512 std::array<std::bitset<param::NPadRows>, SECTORSPERSIDE * SIDES> mXBinsIgnore{};
513 std::array<std::array<float, param::NPadRows>, SECTORSPERSIDE * SIDES> mValidFracXBins{};
514 std::array<std::vector<VoxRes>, SECTORSPERSIDE * SIDES> mVoxelResults{};
515 VoxRes mVoxelResultsOut{};
516 VoxRes* mVoxelResultsOutPtr{&mVoxelResultsOut};
517 std::array<float, 2> mAdhocScalingX{0, 0};
518 bool mDoAdhocCorrectionZ2X{false};
519
520 ClassDefNV(TrackResiduals, 3);
521};
522
523//_____________________________________________________
524inline size_t TrackResiduals::getGlbVoxBin(const std::array<unsigned char, VoxDim>& bvox) const
525{
526 return bvox[VoxX] + (bvox[VoxF] + bvox[VoxZ] * mNY2XBins) * mNXBins;
527}
528
529//_____________________________________________________
530inline size_t TrackResiduals::getGlbVoxBin(int ix, int ip, int iz) const
531{
532 return ix + (ip + iz * mNY2XBins) * mNXBins;
533}
534
535//_____________________________________________________
536inline void TrackResiduals::getVoxelCoordinates(int isec, int ix, int ip, int iz, float& x, float& p, float& z) const
537{
538 x = getX(ix);
539 p = getY2X(ix, ip);
540 z = getZ2X(iz);
541 if (isec >= SECTORSPERSIDE) {
542 z = -z;
543 }
544}
545
546//_____________________________________________________
547inline float TrackResiduals::getDXI(int ix) const
548{
549 if (mUniformBins[VoxX]) {
550 return mDXI;
551 } else {
552 if (ix < param::NRowsPerROC[0]) {
553 // we are in the IROC
554 return 1.f / param::RowDX[0];
555 } else if (ix > param::NRowsAccumulated[param::NROCTypes - 2]) {
556 // we are in the last OROC
557 return 1.f / param::RowDX[param::NROCTypes - 1];
558 } else if (ix < param::NRowsAccumulated[1]) {
559 // OROC1
560 return 1.f / param::RowDX[1];
561 } else {
562 // OROC2
563 return 1.f / param::RowDX[2];
564 }
565 }
566}
567
568//_____________________________________________________
569inline float TrackResiduals::getX(int ix) const
570{
571 return mUniformBins[VoxX] ? param::MinX + (ix + 0.5) * mDX : param::RowX[ix];
572}
573
574//_____________________________________________________
575inline float TrackResiduals::getMaxY2X(int ix) const
576{
577 return mMaxY2X[ix];
578}
579
580//_____________________________________________________
581inline float TrackResiduals::getY2X(int ix, int ip) const
582{
583 if (mUniformBins[VoxF]) {
584 return (0.5f + ip) * mDY2X[ix] - mMaxY2X[ix];
585 }
586 return mMaxY2X[ix] * mY2XBinsCenter[ip];
587}
588
589//_____________________________________________________
590inline float TrackResiduals::getZ2X(int iz) const
591{
592 if (mUniformBins[VoxZ]) {
593 return (0.5f + iz) * getDZ2X();
594 }
595 return mZ2XBinsCenter[iz];
596}
597
598//_____________________________________________________
599inline float TrackResiduals::getDY2XI(int ix, int iy) const
600{
601 if (mUniformBins[VoxF]) {
602 return mDY2XI[ix];
603 }
604 return mY2XBinsDI[iy] / mMaxY2X[ix];
605}
606
607//_____________________________________________________
608inline float TrackResiduals::getDZ2X(int iz) const
609{
610 if (mUniformBins[VoxZ]) {
611 return mDZ2X;
612 }
613 return 2.f * mZ2XBinsDH[iz];
614}
615
616//_____________________________________________________
617inline float TrackResiduals::getDZ2XI(int iz) const
618{
619 if (mUniformBins[VoxZ]) {
620 return mDZ2XI;
621 }
622 return mZ2XBinsDI[iz];
623}
624
625//_____________________________________________________
626inline void TrackResiduals::findVoxel(float x, float y2x, float z2x, int& ix, int& ip, int& iz) const
627{
628 ix = getXBin(x);
629 ip = getY2XBin(y2x, ix);
630 iz = getZ2XBin(z2x);
631}
632
633//_____________________________________________________
634inline int TrackResiduals::getXBinExact(float x) const
635{
636 if (mUniformBins[VoxX]) {
637 int ix = (x - param::MinX) * mDXI;
638 return (ix < 0 || ix >= mNXBins) ? -1 : ix;
639 }
640 return getRowID(x);
641}
642
643//_____________________________________________________
644inline int TrackResiduals::getXBin(float x) const
645{
646 int bx = getXBinExact(x);
647 return bx > -1 ? (bx < mNXBins ? bx : mNXBins - 1) : 0;
648}
649
650//_____________________________________________________
651inline int TrackResiduals::getY2XBinExact(float y2x, int ix) const
652{
653 if (y2x < -mMaxY2X[ix]) {
654 return -1;
655 }
656 if (y2x > mMaxY2X[ix]) {
657 return mNY2XBins;
658 }
659 if (mUniformBins[VoxF]) {
660 return static_cast<int>((y2x + mMaxY2X[ix]) * getDY2XI(ix));
661 }
662 y2x /= mMaxY2X[ix];
663 for (int iBin = 0; iBin < mNY2XBins; ++iBin) {
664 if (y2x < mY2XBinsCenter[iBin] + mY2XBinsDH[iBin]) {
665 return iBin;
666 }
667 }
668 return mNY2XBins;
669}
670
671//_____________________________________________________
672inline int TrackResiduals::getY2XBin(float y2x, int ix) const
673{
674 int bp = getY2XBinExact(y2x, ix);
675 return bp > -1 ? (bp < mNY2XBins ? bp : mNY2XBins - 1) : 0;
676}
677
678//_____________________________________________________
679inline int TrackResiduals::getZ2XBinExact(float z2x) const
680{
681 if (mUniformBins[VoxZ]) {
682 float bz = z2x * getDZ2XI();
683 if (bz >= mNZ2XBins) {
684 return -1;
685 }
686 if (bz < 0) {
687 // accounting for clusters which were moved to the wrong side
688 // TODO: how can this happen?
689 bz = 0;
690 }
691 return static_cast<int>(bz);
692 }
693 for (int iBin = 0; iBin < mNZ2XBins; ++iBin) {
694 if (z2x < mZ2XBinsCenter[iBin] + mZ2XBinsDH[iBin]) {
695 return iBin;
696 }
697 }
698 return -1;
699}
700
701//_____________________________________________________
702inline int TrackResiduals::getZ2XBin(float z2x) const
703{
704 int iz = getZ2XBinExact(z2x);
705 return iz < 0 ? mNZ2XBins - 1 : iz;
706}
707
708} // namespace tpc
709
710} // namespace o2
711
712#endif
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
float getX(int ix) const
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
@ DispDone
voxel dispersions have been processed
@ DistDone
voxel residuals have been processed
@ SmoothDone
voxel has been smoothed
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
@ ResD
index for dispersions
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)
float selectKthMin(const int k, std::vector< float > &data) const
float getMaxY2X(int ix) 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
GLenum const GLfloat * params
Definition glcorearb.h:272
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