Project
Loading...
Searching...
No Matches
TrackResiduals.cxx
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
18
22#include "MathUtils/fit.h"
23
24#include "TMatrixDSym.h"
25#include "TDecompChol.h"
26#include "TVectorD.h"
27
28#include <cmath>
29#include <cstring>
30#include <algorithm>
31
32// for debugging
33#include "TStopwatch.h"
34#include "TSystem.h"
35#include <iostream>
36#include <fstream>
37#include <limits>
38#include <iomanip>
39
40#include <fairlogger/Logger.h>
41
42using namespace o2::tpc;
43
49
50//______________________________________________________________________________
51void TrackResiduals::init(bool doBinning)
52{
53 mSmoothPol2[VoxX] = true;
54 mSmoothPol2[VoxF] = true;
57 mMaxZ2X = mParams->maxZ2X;
58 mIsInitialized = true;
59
60 if (doBinning) {
61 // initialize binning
63 }
64
65 LOG(info) << "Initialization complete";
66}
67
68//______________________________________________________________________________
69void TrackResiduals::setY2XBinning(const std::vector<float>& binning)
70{
71 if (mIsInitialized) {
72 LOG(error) << "Binning already initialized, not changing y/x binning";
73 return;
74 }
75
76 if (binning.size() == 0) {
77 LOGP(info, "Empty binning provided, will use default uniform y/x binning with {} bins", mNY2XBins);
78 return;
79 } else if (binning.size() == 1) {
80 const int bins = static_cast<int>(binning.at(0));
82 LOGP(info, "Setting uniform binning for y/x with {} bins", bins - 1);
83 return;
84 }
85
86 const int nBins = binning.size() - 1;
87 if (fabsf(binning[0] + 1.f) > param::sEps || fabsf(binning[nBins] - 1.f) > param::sEps) {
88 LOG(error) << "Provided binning for y/x not in range -1 to 1: " << binning[0] << " - " << binning[nBins] << ". Not changing y/x binning";
89 return;
90 }
91
92 LOGP(info, "Setting custom binning for y/x with {} bins", nBins);
93 setNY2XBins(nBins);
94 mUniformBins[VoxF] = false;
95 mY2XBinsDH.clear();
96 mY2XBinsDI.clear();
97 mY2XBinsCenter.clear();
98 for (int iBin = 0; iBin < nBins; ++iBin) {
99 mY2XBinsDH.push_back(.5f * (binning[iBin + 1] - binning[iBin]));
100 mY2XBinsDI.push_back(.5f / mY2XBinsDH[iBin]);
101 mY2XBinsCenter.push_back(binning[iBin] + mY2XBinsDH[iBin]);
102 LOGF(info, "Bin %i: center (%.3f), half bin width (%.3f)", iBin, mY2XBinsCenter.back(), mY2XBinsDH.back());
103 }
104}
105
106//______________________________________________________________________________
107void TrackResiduals::setZ2XBinning(const std::vector<float>& binning)
108{
109 if (mIsInitialized) {
110 LOG(error) << "Binning already initialized, not changing z/x binning";
111 return;
112 }
113
114 if (binning.size() == 0) {
115 LOGP(info, "Empty binning provided, will use default uniform z/x binning with {} bins", mNZ2XBins);
116 return;
117 } else if (binning.size() == 1) {
118 const int bins = static_cast<int>(binning.at(0));
120 LOGP(info, "Setting uniform binning for z/x with {} bins", bins);
121 return;
122 }
123
124 int nBins = binning.size() - 1;
125 if (fabsf(binning[0]) > param::sEps || fabsf(binning[nBins] - 1.f) > param::sEps) {
126 LOG(error) << "Provided binning for z/x not in range 0 to 1: " << binning[0] << " - " << binning[nBins] << ". Not changing z/x binning";
127 return;
128 }
129
130 LOGP(info, "Setting custom binning for z/x with {} bins", nBins);
131 setNZ2XBins(nBins);
132 mUniformBins[VoxZ] = false;
133 mZ2XBinsDH.clear();
134 mZ2XBinsDI.clear();
135 mZ2XBinsCenter.clear();
136 for (int iBin = 0; iBin < nBins; ++iBin) {
137 mZ2XBinsDH.push_back(.5f * (binning[iBin + 1] - binning[iBin]) * mMaxZ2X);
138 mZ2XBinsDI.push_back(.5f / mZ2XBinsDH[iBin]);
139 mZ2XBinsCenter.push_back(binning[iBin] * mMaxZ2X + mZ2XBinsDH[iBin]);
140 LOGF(info, "Bin %i: center (%.3f), half bin width (%.3f)", iBin, mZ2XBinsCenter.back(), mZ2XBinsDH.back());
141 }
142}
143
144//______________________________________________________________________________
146{
147 // initialize binning structures
148 //
149 // X binning
150 if (mNXBins > 0 && mNXBins < param::NPadRows) {
151 // uniform binning in X
152 LOGF(info, "X-binning is uniform with %i bins from %.2f to %.2f", mNXBins, param::MinX, param::MaxX);
153 mDXI = mNXBins / (param::MaxX - param::MinX);
154 mDX = 1.0f / mDXI;
155 mUniformBins[VoxX] = true;
156 } else {
157 // binning per pad row
158 LOGF(info, "X-binning is per pad-row");
159 mNXBins = param::NPadRows;
160 mUniformBins[VoxX] = false;
161 mDX = param::RowDX[0];
162 mDXI = 1.f / mDX; // should not be used
163 }
164 //
165 // Y/X binning
166 mMaxY2X.resize(mNXBins);
167 mDY2XI.resize(mNXBins);
168 mDY2X.resize(mNXBins);
169 //
170 for (int ix = 0; ix < mNXBins; ++ix) {
171 float x = getX(ix);
172 mMaxY2X[ix] = tan(.5f * SECPHIWIDTH) - sDeadZone / x;
173 mDY2XI[ix] = mNY2XBins / (2.f * mMaxY2X[ix]);
174 mDY2X[ix] = 1.f / mDY2XI[ix];
175 }
176 if (mUniformBins[VoxF]) {
177 LOGF(info, "Y/X-binning is uniform with %i bins from -MaxY2X to +MaxY2X (values depend on X-bin)", mNY2XBins);
178 for (int ip = 0; ip < mNY2XBins; ++ip) {
179 mY2XBinsDH.push_back(1.f / mNY2XBins);
180 mY2XBinsDI.push_back(.5f / mY2XBinsDH[ip]);
181 mY2XBinsCenter.push_back(-1.f + (ip + 0.5f) * 2.f * mY2XBinsDH[ip]);
182 LOGF(info, "Bin %i: center (%.3f), half bin width (%.3f)", ip, mY2XBinsCenter.back(), mY2XBinsDH.back());
183 }
184 }
185 //
186 // Z/X binning
187 mDZ2XI = mNZ2XBins / mMaxZ2X;
188 mDZ2X = 1.0f / mDZ2XI; // for uniform case only
189 if (mUniformBins[VoxZ]) {
190 LOGF(info, "Z/X-binning is uniform with %i bins from 0 to %f", mNZ2XBins, mMaxZ2X);
191 for (int iz = 0; iz < mNZ2XBins; ++iz) {
192 mZ2XBinsDH.push_back(.5f * mDZ2X);
193 mZ2XBinsDI.push_back(mDZ2XI);
194 mZ2XBinsCenter.push_back((iz + 0.5f) * mDZ2X);
195 LOGF(info, "Bin %i: center (%.3f), half bin width (%.3f)", iz, mZ2XBinsCenter.back(), mZ2XBinsDH.back());
196 }
197 }
198 //
199 mNVoxPerSector = mNY2XBins * mNZ2XBins * mNXBins;
200 LOGF(info, "Each TPC sector is divided into %i voxels", mNVoxPerSector);
201}
202
203//______________________________________________________________________________
205{
206 if (mInitResultsContainer.test(iSec)) {
207 return;
208 }
209 mInitResultsContainer.set(iSec);
210 mVoxelResults[iSec].resize(mNVoxPerSector);
211 for (int ix = 0; ix < mNXBins; ++ix) {
212 for (int ip = 0; ip < mNY2XBins; ++ip) {
213 for (int iz = 0; iz < mNZ2XBins; ++iz) {
214 const size_t binGlb = getGlbVoxBin(ix, ip, iz);
215 VoxRes& resVox = mVoxelResults[iSec][binGlb];
216 resVox.bvox[VoxX] = ix;
217 resVox.bvox[VoxF] = ip;
218 resVox.bvox[VoxZ] = iz;
219 resVox.bsec = iSec;
220 }
221 }
222 }
223 LOG(debug) << "initialized the container for the main results";
224}
225
226//______________________________________________________________________________
228{
229 for (int iSec = 0; iSec < SECTORSPERSIDE * SIDES; ++iSec) {
230 mXBinsIgnore[iSec].reset();
231 std::fill(mVoxelResults[iSec].begin(), mVoxelResults[iSec].end(), VoxRes());
232 std::fill(mValidFracXBins[iSec].begin(), mValidFracXBins[iSec].end(), 0);
233 }
234}
235
236//______________________________________________________________________________
238{
239 int ix;
240 if (x < param::RowX[param::NRowsAccumulated[0] - 1] + param::RowDX[0]) {
241 // we are in the IROC
242 ix = (x - (param::RowX[0] - .5f * param::RowDX[0])) / param::RowDX[0];
243 if (ix < 0) {
244 // x is smaller than the inner radius of the first pad row
245 ix = -1;
246 }
247 } else if (x >= param::RowX[param::NRowsAccumulated[param::NROCTypes - 2]] - .5f * param::RowDX[param::NROCTypes - 1]) {
248 // we are in the OROC3
249 ix = (x - (param::RowX[param::NRowsAccumulated[param::NROCTypes - 2]] - .5f * param::RowDX[param::NROCTypes - 1])) / param::RowDX[param::NROCTypes - 1] + param::NRowsAccumulated[param::NROCTypes - 2];
250 if (ix >= param::NPadRows) {
251 // x is larger than the outer radius of the last OROC pad row
252 ix = -1;
253 }
254 } else if (x > param::RowX[param::NRowsAccumulated[0]] - .5f * param::RowDX[1] && x < param::RowX[param::NRowsAccumulated[1] - 1] + .5f * param::RowDX[1]) {
255 // we are in the OROC1
256 ix = (x - (param::RowX[param::NRowsAccumulated[0]] - .5f * param::RowDX[1])) / param::RowDX[1] + param::NRowsAccumulated[0];
257 } else if (x > param::RowX[param::NRowsAccumulated[1]] - .5f * param::RowDX[2] && x < param::RowX[param::NRowsAccumulated[2] - 1] + .5f * param::RowDX[2]) {
258 // we are in the OROC2
259 ix = (x - (param::RowX[param::NRowsAccumulated[1]] - .5f * param::RowDX[2])) / param::RowDX[2] + param::NRowsAccumulated[1];
260 } else {
261 // x is in one of the gaps between the ROCs
262 ix = -1;
263 }
264 return ix;
265}
266
267bool TrackResiduals::findVoxelBin(int secID, float x, float y, float z, std::array<unsigned char, VoxDim>& bvox) const
268{
269 // Z/X bin
270 if (fabs(z / x) > mMaxZ2X) {
271 return false;
272 }
273 int bz = getZ2XBinExact(secID < SECTORSPERSIDE ? z / x : -z / x);
274 if (bz < 0) {
275 return false;
276 }
277 // X bin
278 int bx = getXBinExact(x);
279 if (bx < 0 || bx >= mNXBins) {
280 return false;
281 }
282 // Y/X bin
283 int bp = getY2XBinExact(y / x, bx);
284 if (bp < 0 || bp >= mNY2XBins) {
285 return false;
286 }
287
288 bvox[VoxZ] = bz;
289 bvox[VoxX] = bx;
290 bvox[VoxF] = bp;
291 return true;
292}
293
294void TrackResiduals::setKernelType(KernelType kernel, float bwX, float bwP, float bwZ, float scX, float scP, float scZ)
295{
296 // set kernel type and widths in terms of binning in x, y/x, z/x and define aux variables
297 mKernelType = kernel;
298
299 mKernelScaleEdge[VoxX] = scX;
300 mKernelScaleEdge[VoxF] = scP;
301 mKernelScaleEdge[VoxZ] = scZ;
302
303 mKernelWInv[VoxX] = (bwX > 0) ? 1. / bwX : 1.;
304 mKernelWInv[VoxF] = (bwP > 0) ? 1. / bwP : 1.;
305 mKernelWInv[VoxZ] = (bwZ > 0) ? 1. / bwZ : 1.;
306
307 if (mKernelType == KernelType::Epanechnikov) {
308 // bandwidth 1
309 mStepKern[VoxX] = static_cast<int>(nearbyint(bwX + 0.5));
310 mStepKern[VoxF] = static_cast<int>(nearbyint(bwP + 0.5));
311 mStepKern[VoxZ] = static_cast<int>(nearbyint(bwZ + 0.5));
312 } else if (mKernelType == KernelType::Gaussian) {
313 // look in ~5 sigma
314 mStepKern[VoxX] = static_cast<int>(nearbyint(bwX * 5. + 0.5));
315 mStepKern[VoxF] = static_cast<int>(nearbyint(bwP * 5. + 0.5));
316 mStepKern[VoxZ] = static_cast<int>(nearbyint(bwZ * 5. + 0.5));
317 } else {
318 LOG(error) << "given kernel type is not defined";
319 }
320 for (int i = VoxDim; i--;) {
321 if (mStepKern[i] < 1) {
322 mStepKern[i] = 1;
323 }
324 }
325}
326
332
333void TrackResiduals::setStats(const std::vector<TrackResiduals::VoxStats>& statsIn, int iSec)
334{
336 std::vector<VoxRes>& secDataTmp = mVoxelResults[iSec];
337 for (int iVox = 0; iVox < mNVoxPerSector; ++iVox) {
338 secDataTmp[iVox].stat[VoxX] = statsIn[iVox].meanPos[VoxX];
339 secDataTmp[iVox].stat[VoxF] = statsIn[iVox].meanPos[VoxF];
340 secDataTmp[iVox].stat[VoxZ] = statsIn[iVox].meanPos[VoxZ];
341 secDataTmp[iVox].stat[VoxDim] = statsIn[iVox].nEntries;
342 }
343}
344
346{
348 std::vector<VoxRes>& secDataTmp = mVoxelResults[iSec];
349 for (int iVox = 0; iVox < mNVoxPerSector; ++iVox) {
350 const auto& voxStat = mVoxStatsIn[iVox];
351 VoxRes& resVox = secDataTmp[iVox];
352 for (int iDim = VoxDim; iDim--;) {
353 const auto sumStat = (resVox.stat[VoxDim] + voxStat.nEntries);
354 if (sumStat == 0) {
355 continue;
356 }
357 double norm = 1. / sumStat;
358 resVox.stat[iDim] = (resVox.stat[iDim] * resVox.stat[VoxDim] + voxStat.meanPos[iDim] * voxStat.nEntries) * norm;
359 }
360 resVox.stat[VoxDim] += voxStat.nEntries;
361 }
362}
363
364//______________________________________________________________________________
366{
367 LOGP(info, "Processing {} voxel residuals for sector {}", mLocalResidualsIn.size(), iSec);
369 // effective t0 correction changes sign between A-/C-side
370 float effT0corr = (iSec < SECTORSPERSIDE) ? mEffT0Corr : -1. * mEffT0Corr;
371 std::vector<size_t> binData;
372 for (const auto& res : mLocalResidualsIn) {
373 binData.push_back(getGlbVoxBin(res.bvox));
374 }
375 // sort in voxel increasing order
376 std::vector<size_t> binIndices(binData.size());
377 o2::math_utils::SortData(binData, binIndices);
378 // fill the voxel statistics into the results container
379 std::vector<VoxRes>& secData = mVoxelResults[iSec];
380
381 // vectors holding the data for one voxel at a time
382 std::vector<float> dyVec;
383 std::vector<float> dzVec;
384 std::vector<float> tgVec;
385 // assuming we will always have around 1000 entries per voxel
386 dyVec.reserve(1e3);
387 dzVec.reserve(1e3);
388 tgVec.reserve(1e3);
389 size_t currVoxBin = -1;
390 unsigned int nPointsInVox = 0;
391 unsigned int nProcessed = 0;
392 while (nProcessed < binData.size()) {
393 // read all points, voxel by voxel
394 int idx = binIndices[nProcessed];
395 if (currVoxBin != binData[idx]) {
396 if (nPointsInVox) {
397 VoxRes& resVox = secData[currVoxBin];
398 processVoxelResiduals(dyVec, dzVec, tgVec, resVox);
399 }
400 currVoxBin = binData[idx];
401 nPointsInVox = 0;
402 dyVec.clear();
403 dzVec.clear();
404 tgVec.clear();
405 }
406 dyVec.push_back(mLocalResidualsIn[idx].dy * param::MaxResid / 0x7fff);
407 dzVec.push_back(mLocalResidualsIn[idx].dz * param::MaxResid / 0x7fff -
408 mEffVdriftCorr * secData[currVoxBin].stat[VoxZ] * secData[currVoxBin].stat[VoxX] -
409 effT0corr);
410 tgVec.push_back(mLocalResidualsIn[idx].tgSlp * param::MaxTgSlp / 0x7fff);
411
412 ++nPointsInVox;
413 ++nProcessed;
414 }
415 if (nPointsInVox) {
416 // process last voxel
417 VoxRes& resVox = secData[currVoxBin];
418 processVoxelResiduals(dyVec, dzVec, tgVec, resVox);
419 }
420 LOG(info) << "extracted residuals for sector " << iSec;
421
422 int nRowsOK = validateVoxels(iSec);
423 LOG(info) << "number of validated X rows: " << nRowsOK;
424 if (!nRowsOK) {
425 LOG(warning) << "sector " << iSec << ": all X-bins disabled, abandon smoothing";
426 return;
427 } else {
428 smooth(iSec);
429 }
430
431 // process dispersions
432 dyVec.clear();
433 tgVec.clear();
434 currVoxBin = -1;
435 nProcessed = 0;
436 nPointsInVox = 0;
437 while (nProcessed < binData.size()) {
438 int idx = binIndices[nProcessed];
439 if (currVoxBin != binData[idx]) {
440 if (nPointsInVox) {
441 VoxRes& resVox = secData[currVoxBin];
442 if (!getXBinIgnored(iSec, resVox.bvox[VoxX])) {
443 processVoxelDispersions(tgVec, dyVec, resVox);
444 }
445 }
446 currVoxBin = binData[idx];
447 nPointsInVox = 0;
448 dyVec.clear();
449 tgVec.clear();
450 }
451 dyVec.push_back(mLocalResidualsIn[idx].dy * param::MaxResid / 0x7fff);
452 tgVec.push_back(mLocalResidualsIn[idx].tgSlp * param::MaxTgSlp / 0x7fff);
453 ++nPointsInVox;
454 ++nProcessed;
455 }
456 if (nPointsInVox) {
457 // process last voxel
458 VoxRes& resVox = secData[currVoxBin];
459 if (!getXBinIgnored(iSec, resVox.bvox[VoxX])) {
460 processVoxelDispersions(tgVec, dyVec, resVox);
461 }
462 }
463 // smooth dispersions
464 for (int ix = 0; ix < mNXBins; ++ix) {
465 if (getXBinIgnored(iSec, ix)) {
466 continue;
467 }
468 for (int iz = 0; iz < mNZ2XBins; ++iz) {
469 for (int ip = 0; ip < mNY2XBins; ++ip) {
470 int voxBin = getGlbVoxBin(ix, ip, iz);
471 VoxRes& resVox = secData[voxBin];
472 getSmoothEstimate(iSec, resVox.stat[VoxX], resVox.stat[VoxF], resVox.stat[VoxZ], resVox.DS, 0x1 << VoxV);
473 }
474 }
475 }
476 LOG(info) << "Done processing residuals for sector " << iSec;
477 dumpResults(iSec);
478}
479
480//______________________________________________________________________________
481void TrackResiduals::processVoxelResiduals(std::vector<float>& dy, std::vector<float>& dz, std::vector<float>& tg, VoxRes& resVox)
482{
483 int nPoints = dy.size();
484 if (nPoints < mParams->minEntriesPerVoxel) {
485 LOG(debug) << "voxel " << getGlbVoxBin(resVox.bvox) << " is skipped due to too few entries (" << nPoints << " < " << mParams->minEntriesPerVoxel << ")";
486 return;
487 } else {
488 LOGF(debug, "Processing voxel %i with %i entries", getGlbVoxBin(resVox.bvox), nPoints);
489 }
490 std::array<float, 7> zResults;
491 resVox.flags = 0;
492 std::vector<size_t> indices(dz.size());
493 if (!o2::math_utils::LTMUnbinned(dz, indices, zResults, mParams->LTMCut)) {
494 LOG(debug) << "failed trimming input array for voxel " << getGlbVoxBin(resVox.bvox);
495 return;
496 }
497 if (!mParams->isBfieldZero) {
498 std::array<float, 2> res{0.f};
499 std::array<float, 3> err{0.f};
500 float sigMAD = fitPoly1Robust(tg, dy, res, err, mParams->LTMCut);
501 if (sigMAD < 0) {
502 LOG(debug) << "failed robust linear fit, sigMAD = " << sigMAD;
503 return;
504 }
505 float corrErr = err[0] * err[2];
506 corrErr = corrErr > 0 ? err[1] / std::sqrt(corrErr) : -999;
507 //
508 resVox.D[ResX] = -res[1];
509 resVox.D[ResY] = res[0];
510 resVox.D[ResZ] = zResults[1];
511 resVox.E[ResX] = std::sqrt(err[2]);
512 resVox.E[ResY] = std::sqrt(err[0]);
513 resVox.E[ResZ] = zResults[4];
514 resVox.EXYCorr = corrErr;
515 resVox.D[ResD] = resVox.dYSigMAD = sigMAD; // later will be overwritten by real dispersion
516 resVox.dZSigLTM = zResults[2];
517 } else {
518 // for B=0 we cannot disentangle radial distortions from distortions in y,
519 // so simply use average for dy as well and set distortion in X to zero
520 std::array<float, 7> yResults;
521 std::vector<size_t> indicesY(dy.size());
522 if (!o2::math_utils::LTMUnbinned(dy, indicesY, yResults, mParams->LTMCut)) {
523 LOG(debug) << "failed trimming input array for voxel " << getGlbVoxBin(resVox.bvox);
524 return;
525 }
526 resVox.D[ResX] = 0; // force to zero
527 resVox.D[ResY] = yResults[1];
528 resVox.D[ResZ] = zResults[1];
529 resVox.E[ResX] = 0;
530 resVox.E[ResY] = yResults[4];
531 resVox.E[ResZ] = zResults[4];
532 resVox.EXYCorr = 0;
533 resVox.D[ResD] = resVox.dYSigMAD = yResults[2];
534 resVox.dZSigLTM = zResults[2];
535 }
536
537 LOGF(debug, "D[0]=%.2f, D[1]=%.2f, D[2]=%.2f, E[0]=%.2f, E[1]=%.2f, E[2]=%.2f, EXYCorr=%.4f, dYSigMAD=%.3f, dZSigLTM=%.3f",
538 resVox.D[0], resVox.D[1], resVox.D[2], resVox.E[0], resVox.E[1], resVox.E[2], resVox.EXYCorr, resVox.dYSigMAD, resVox.dZSigLTM);
539
540 resVox.flags |= DistDone;
541
542 return;
543}
544
545void TrackResiduals::processVoxelDispersions(std::vector<float>& tg, std::vector<float>& dy, VoxRes& resVox)
546{
547 size_t nPoints = tg.size();
548 LOG(debug) << "processing voxel dispersions for vox " << getGlbVoxBin(resVox.bvox) << " with " << nPoints << " points";
549 if (nPoints < 2) {
550 return;
551 }
552 for (size_t i = nPoints; i--;) {
553 dy[i] -= resVox.DS[ResY] - resVox.DS[ResX] * tg[i];
554 }
555 resVox.D[ResD] = getMAD2Sigma(dy);
556 resVox.E[ResD] = resVox.D[ResD] / sqrt(2.f * nPoints); // a la gaussioan RMS error (very crude)
557 resVox.flags |= DispDone;
558}
559
560//______________________________________________________________________________
562{
563 // apply voxel validation cuts
564 // return number of good voxels for given sector
565 int cntMasked = 0; // number of voxels masked for any reason (either low statistics or bad fit)
566 int cntLowStat = 0; // number of voxels which were not processed due to too low statistics
567 mXBinsIgnore[iSec].reset();
568 std::vector<VoxRes>& secData = mVoxelResults[iSec];
569
570 int cntMaskedFit = 0;
571 int cntMaskedSigma = 0;
572
573 // find bad voxels in sector
574 for (int ix = 0; ix < mNXBins; ++ix) {
575 int cntValid = 0;
576 for (int ip = 0; ip < mNY2XBins; ++ip) {
577 for (int iz = 0; iz < mNZ2XBins; ++iz) {
578 int binGlb = getGlbVoxBin(ix, ip, iz);
579 VoxRes& resVox = secData[binGlb];
580 if ((resVox.flags & DistDone) == 0) {
581 ++cntLowStat;
582 }
583 bool voxelOK = (resVox.flags & DistDone) && !(resVox.flags & Masked);
584 if (voxelOK) {
585 // check fit errors
586 if (resVox.E[ResY] * resVox.E[ResY] > mParams->maxFitErrY2 ||
587 resVox.E[ResX] * resVox.E[ResX] > mParams->maxFitErrX2 ||
588 fabs(resVox.EXYCorr) > mParams->maxFitCorrXY) {
589 voxelOK = false;
590 ++cntMaskedFit;
591 }
592 // check raw distribution sigmas
593 if (resVox.dYSigMAD > mParams->maxSigY ||
594 resVox.dZSigLTM > mParams->maxSigZ) {
595 voxelOK = false;
596 ++cntMaskedSigma;
597 }
598 }
599 if (voxelOK) {
600 ++cntValid;
601 } else {
602 ++cntMasked;
603 resVox.flags |= Masked;
604 }
605 } // loop over Z
606 } // loop over Y/X
607 mValidFracXBins[iSec][ix] = static_cast<float>(cntValid) / (mNY2XBins * mNZ2XBins);
608 LOGP(debug, "Sector {}: xBin {} has {} % of voxels valid. Total masked due to fit: {} ,and sigma: {}",
609 iSec, ix, mValidFracXBins[iSec][ix] * 100., cntMaskedFit, cntMaskedSigma);
610 } // loop over X
611
612 // mask X-bins which cannot be smoothed
613
614 short nBadReg = 0; // count bad regions (one or more consecutive bad X-bins)
615 std::array<short, param::NPadRows> badStart; // to store indices to the beginnings of the bad regions
616 std::array<short, param::NPadRows> badEnd; // to store indices to the end of the bad regions
617 bool prevBad = false;
618 float fracBadRows = 0.f;
619 for (int ix = 0; ix < mNXBins; ++ix) {
620 if (mValidFracXBins[iSec][ix] < mParams->minValidVoxFracDrift) {
621 LOG(debug) << "row " << ix << " is bad";
622 ++fracBadRows;
623 if (prevBad) {
624 badEnd[nBadReg] = ix;
625 } else {
626 badStart[nBadReg] = ix;
627 badEnd[nBadReg] = ix;
628 prevBad = true;
629 }
630 } else {
631 if (prevBad) {
632 ++nBadReg;
633 prevBad = false;
634 }
635 }
636 }
637 if (prevBad) {
638 ++nBadReg;
639 }
640 fracBadRows /= mNXBins;
641 if (fracBadRows > mParams->maxFracBadRowsPerSector) {
642 LOG(warning) << "sector " << iSec << ": Fraction of bad X-bins: " << fracBadRows << " -> masking whole sector";
643 mXBinsIgnore[iSec].set();
644 } else {
645 for (int iBad = 0; iBad < nBadReg; ++iBad) {
646 LOG(debug) << "masking bad region " << iBad;
647 short badInReg = badEnd[iBad] - badStart[iBad] + 1;
648 short badInNextReg = iBad < (nBadReg - 1) ? badEnd[iBad] - badStart[iBad] + 1 : 0;
649 if (badInReg > mParams->maxBadXBinsToCover) {
650 // disable too large bad patches
651 for (int i = 0; i < badInReg; ++i) {
652 LOG(debug) << "disabling too large patch in bad region " << iBad << ", badStart(" << badStart[iBad] << "), i(" << i << ")";
653 mXBinsIgnore[iSec].set(badStart[iBad] + i);
654 }
655 }
656 if (badInNextReg > mParams->maxBadXBinsToCover && (badStart[iBad + 1] - badEnd[iBad] - 1) < mParams->minGoodXBinsToCover) {
657 // disable too small isolated good patches`
658 for (int i = badEnd[iBad] + 1; i < badStart[iBad + 1]; ++i) {
659 LOG(debug) << "disabling too small good patch before bad region " << iBad + 1 << ", badStart(" << badEnd[iBad] << "), badEnd(" << badStart[iBad + 1] << ")";
660 mXBinsIgnore[iSec].set(i);
661 }
662 }
663 }
664 if (nBadReg) {
665 if (mXBinsIgnore[iSec].test(badStart[0]) && badStart[0] < mParams->minGoodXBinsToCover) {
666 // 1st good patch is too small
667 for (int i = 0; i < badStart[0]; ++i) {
668 LOG(debug) << "disabling too small first good patch badStart(0), badEnd(" << badStart[0] << ")";
669 mXBinsIgnore[iSec].set(i);
670 }
671 }
672 if (mXBinsIgnore[iSec].test(badStart[nBadReg - 1]) && (mNXBins - badEnd[nBadReg - 1] - 1) < mParams->minGoodXBinsToCover) {
673 // last good patch is too small
674 for (int i = badEnd[nBadReg - 1] + 1; i < mNXBins; ++i) {
675 LOG(debug) << "disabling too small last good patch badStart(" << badEnd[nBadReg - 1] << "), badEnd(" << mNXBins << ")";
676 mXBinsIgnore[iSec].set(i);
677 }
678 }
679 }
680 }
681 //
682 int nMaskedRows = mXBinsIgnore[iSec].count();
683 LOGP(info, "Sector {}: out of {} voxels {} are masked. {} (low stat), {} (invalid fit) and {} (raw distrib sigma)",
684 iSec, mNVoxPerSector, cntMasked, cntLowStat, cntMaskedFit, cntMaskedSigma);
685 //
686 return mNXBins - nMaskedRows;
687}
688
690{
691 std::vector<VoxRes>& secData = mVoxelResults[iSec];
692 for (int ix = 0; ix < mNXBins; ++ix) {
693 if (getXBinIgnored(iSec, ix)) {
694 continue;
695 }
696 for (int ip = 0; ip < mNY2XBins; ++ip) {
697 for (int iz = 0; iz < mNZ2XBins; ++iz) {
698 int voxBin = getGlbVoxBin(ix, ip, iz);
699 VoxRes& resVox = secData[voxBin];
700 resVox.flags &= ~SmoothDone;
701 bool res = getSmoothEstimate(resVox.bsec, resVox.stat[VoxX], resVox.stat[VoxF], resVox.stat[VoxZ], resVox.DS, (0x1 << VoxX | 0x1 << VoxF | 0x1 << VoxZ));
702 if (!res) {
703 mNSmoothingFailedBins[iSec]++;
704 } else {
705 resVox.flags |= SmoothDone;
706 }
707 }
708 }
709 }
710 // substract dX contribution to dZ
711 for (int ix = 0; ix < mNXBins; ++ix) {
712 if (getXBinIgnored(iSec, ix)) {
713 continue;
714 }
715 for (int ip = 0; ip < mNY2XBins; ++ip) {
716 for (int iz = 0; iz < mNZ2XBins; ++iz) {
717 int voxBin = getGlbVoxBin(ix, ip, iz);
718 VoxRes& resVox = secData[voxBin];
719 if (!(resVox.flags & SmoothDone)) {
720 continue;
721 }
722 // TODO: Usage of Z/X is bug???
723 float z2x = resVox.stat[VoxZ];
724 if (mDoAdhocCorrectionZ2X) {
725 //
726 const float z = z2x * resVox.stat[VoxX] - resVox.DS[ResZ];
727 const float x = resVox.stat[VoxX] - resVox.DS[ResX]; // is subration of DS[ResX] correct?
728 z2x = z / x;
729 }
730 resVox.DS[ResZ] += z2x * resVox.DS[ResX]; // remove slope*dX contribution from dZ
731 resVox.D[ResZ] += z2x * resVox.DS[ResX]; // remove slope*dX contribution from dZ
732 //
733 if (mAdhocScalingX[iSec >= 18] != 0) {
734 const float aDX = resVox.DS[ResX] * mAdhocScalingX[iSec >= 18];
735 resVox.D[ResX] += aDX;
736 resVox.DS[ResX] += aDX;
737 resVox.D[ResY] += aDX * resVox.stat[VoxF];
738 resVox.DS[ResY] += aDX * resVox.stat[VoxF];
739 resVox.D[ResZ] += aDX * z2x;
740 resVox.DS[ResZ] += aDX * z2x;
741 }
742 }
743 }
744 }
745}
746
747bool TrackResiduals::getSmoothEstimate(int iSec, float x, float p, float z, std::array<float, ResDim>& res, int whichDim)
748{
749 // get smooth estimate for distortions for point in sector coordinates
751
752 std::array<int, VoxDim> minPointsDir{0}; // min number of points per direction
753 const float kTrialStep = 0.5;
754 std::array<bool, ResDim> doDim{false};
755 for (int i = 0; i < ResDim; ++i) {
756 doDim[i] = (whichDim & (0x1 << i)) > 0;
757 if (doDim[i]) {
758 res[i] = 0.f;
759 }
760 }
761
762 int matSize = sSmtLinDim;
763 for (int i = 0; i < VoxDim; ++i) {
764 minPointsDir[i] = 3; // for pol1 smoothing require at least 3 points
765 if (mSmoothPol2[i]) {
766 ++minPointsDir[i];
767 ++matSize;
768 }
769 }
770
771 int ix0, ip0, iz0;
772 findVoxel(x, p, iSec < SECTORSPERSIDE ? z : -z, ix0, ip0, iz0); // find nearest voxel
773 std::vector<VoxRes>& secData = mVoxelResults[iSec];
774 int binCenter = getGlbVoxBin(ix0, ip0, iz0); // global bin of nearest voxel
775 VoxRes& voxCenter = secData[binCenter]; // nearest voxel
776 LOG(debug) << "getting smooth estimate around voxel " << binCenter;
777
778 // cache
779 // \todo maybe a 1-D cache would be more efficient?
780 std::array<std::array<double, sMaxSmtDim*(sMaxSmtDim + 1) / 2>, ResDim> cmat;
781 int maxNeighb = 10 * 10 * 10;
782 std::vector<VoxRes*> currVox;
783 currVox.reserve(maxNeighb);
784 std::vector<float> currCache;
785 currCache.reserve(maxNeighb * VoxHDim);
786
787 std::array<int, VoxDim> maxTrials;
788 maxTrials[VoxZ] = mNZ2XBins / 2;
789 maxTrials[VoxF] = mNY2XBins / 2;
790 maxTrials[VoxX] = mParams->maxBadXBinsToCover * 2;
791
792 std::array<int, VoxDim> trial{0};
793
794 while (true) {
795 std::fill(mLastSmoothingRes.begin(), mLastSmoothingRes.end(), 0);
796 memset(&cmat[0][0], 0, sizeof(cmat));
797
798 int nbOK = 0; // accounted neighbours
799
800 float stepX = mStepKern[VoxX] * (1. + kTrialStep * trial[VoxX]);
801 float stepF = mStepKern[VoxF] * (1. + kTrialStep * trial[VoxF]);
802 float stepZ = mStepKern[VoxZ] * (1. + kTrialStep * trial[VoxZ]);
803
804 if (!(voxCenter.flags & DistDone) || (voxCenter.flags & Masked) || getXBinIgnored(iSec, ix0)) {
805 // closest voxel has no data -> increase smoothing step
806 stepX += kTrialStep * mStepKern[VoxX];
807 stepF += kTrialStep * mStepKern[VoxF];
808 stepZ += kTrialStep * mStepKern[VoxZ];
809 }
810
811 // effective kernel widths accounting for the increased bandwidth at the edges and missing data
812 float kWXI = getDXI(ix0) * mKernelWInv[VoxX] * mStepKern[VoxX] / stepX;
813 float kWFI = getDY2XI(ix0, ip0) * mKernelWInv[VoxF] * mStepKern[VoxF] / stepF;
814 float kWZI = getDZ2XI(iz0) * mKernelWInv[VoxZ] * mStepKern[VoxZ] / stepZ;
815 int iStepX = static_cast<int>(nearbyint(stepX + 0.5));
816 int iStepF = static_cast<int>(nearbyint(stepF + 0.5));
817 int iStepZ = static_cast<int>(nearbyint(stepZ + 0.5));
818 int ixMin = ix0 - iStepX;
819 int ixMax = ix0 + iStepX;
820 if (ixMin < 0) {
821 ixMin = 0;
822 ixMax = std::min(static_cast<int>(nearbyint(ix0 + stepX * mKernelScaleEdge[VoxX])), mNXBins - 1);
823 kWXI /= mKernelScaleEdge[VoxX];
824 }
825 if (ixMax >= mNXBins) {
826 ixMax = mNXBins - 1;
827 ixMin = std::max(static_cast<int>(nearbyint(ix0 - stepX * mKernelScaleEdge[VoxX])), 0);
828 kWXI /= mKernelScaleEdge[VoxX];
829 }
830
831 int ipMin = ip0 - iStepF;
832 int ipMax = ip0 + iStepF;
833 if (ipMin < 0) {
834 ipMin = 0;
835 ipMax = std::min(static_cast<int>(nearbyint(ip0 + stepF * mKernelScaleEdge[VoxF])), mNY2XBins - 1);
836 kWFI /= mKernelScaleEdge[VoxF];
837 }
838 if (ipMax >= mNY2XBins) {
839 ipMax = mNY2XBins - 1;
840 ipMin = std::max(static_cast<int>(nearbyint(ip0 - stepF * mKernelScaleEdge[VoxF])), 0);
841 kWFI /= mKernelScaleEdge[VoxF];
842 }
843
844 int izMin = iz0 - iStepZ;
845 int izMax = iz0 + iStepZ;
846 if (izMin < 0) {
847 izMin = 0;
848 izMax = std::min(static_cast<int>(nearbyint(iz0 + stepZ * mKernelScaleEdge[VoxZ])), mNZ2XBins - 1);
849 kWZI /= mKernelScaleEdge[VoxZ];
850 }
851 if (izMax >= mNZ2XBins) {
852 izMax = mNZ2XBins - 1;
853 izMin = std::max(static_cast<int>(nearbyint(iz0 - stepZ * mKernelScaleEdge[VoxZ])), 0);
854 kWZI /= mKernelScaleEdge[VoxZ];
855 }
856
857 std::vector<unsigned short> nOccX(ixMax - ixMin + 1, 0);
858 std::vector<unsigned short> nOccF(ipMax - ipMin + 1, 0);
859 std::vector<unsigned short> nOccZ(izMax - izMin + 1, 0);
860
861 int nbCheck = (ixMax - ixMin + 1) * (ipMax - ipMin + 1) * (izMax - izMin + 1);
862 if (nbCheck >= maxNeighb) {
863 maxNeighb = nbCheck + 100;
864 currCache.reserve(maxNeighb * VoxHDim);
865 currVox.reserve(maxNeighb);
866 }
867 std::array<double, 3> u2Vec;
868
869 // first loop, check presence of enough points
870 for (int ix = ixMin; ix <= ixMax; ++ix) {
871 for (int ip = ipMin; ip <= ipMax; ++ip) {
872 for (int iz = izMin; iz <= izMax; ++iz) {
873 int binNb = getGlbVoxBin(ix, ip, iz);
874 VoxRes& voxNb = secData[binNb];
875 if (!(voxNb.flags & DistDone) ||
876 (voxNb.flags & Masked) ||
877 getXBinIgnored(iSec, ix)) {
878 // skip voxels w/o data
879 continue;
880 }
881 // estimate weighted distance
882 float dx = voxNb.stat[VoxX] - x;
883 float df = voxNb.stat[VoxF] - p;
884 float dz = voxNb.stat[VoxZ] - z;
885 float dxw = dx * kWXI;
886 float dfw = df * kWFI;
887 float dzw = dz * kWZI;
888 u2Vec[0] = dxw * dxw;
889 u2Vec[1] = dfw * dfw;
890 u2Vec[2] = dzw * dzw;
891 double kernelWeight = getKernelWeight(u2Vec);
892 if (kernelWeight < 1e-6) {
893 continue;
894 }
895 // new point is validated
896 ++nOccX[ix - ixMin];
897 ++nOccF[ip - ipMin];
898 ++nOccZ[iz - izMin];
899 currVox[nbOK] = &voxNb;
900 currCache[nbOK * VoxHDim + VoxX] = dx;
901 currCache[nbOK * VoxHDim + VoxF] = df;
902 currCache[nbOK * VoxHDim + VoxZ] = dz;
903 currCache[nbOK * VoxHDim + VoxV] = kernelWeight;
904 ++nbOK;
905 }
906 }
907 }
908
909 // check if we have enough points in every dimension
910 std::array<int, VoxDim> nPoints{0};
911 for (int i = ixMax - ixMin + 1; i--;) {
912 if (nOccX[i]) {
913 ++nPoints[VoxX];
914 }
915 }
916 for (int i = ipMax - ipMin + 1; i--;) {
917 if (nOccF[i]) {
918 ++nPoints[VoxF];
919 }
920 }
921 for (int i = izMax - izMin + 1; i--;) {
922 if (nOccZ[i]) {
923 ++nPoints[VoxZ];
924 }
925 }
926 bool enoughPoints = true;
927 std::array<bool, VoxDim> incrDone{false};
928 for (int i = 0; i < VoxDim; ++i) {
929 if (nPoints[i] < minPointsDir[i]) {
930 // need to extend smoothing neighbourhood
931 enoughPoints = false;
932 if (trial[i] < maxTrials[i] && !incrDone[i]) {
933 // try to increment only missing direction
934 ++trial[i];
935 incrDone[i] = true;
936 } else if (trial[i] == maxTrials[i]) {
937 // cannot increment missing direction, try others
938 for (int j = VoxDim; j--;) {
939 if (i != j && trial[j] < maxTrials[j] && !incrDone[j]) {
940 ++trial[j];
941 incrDone[j] = true;
942 }
943 }
944 }
945 }
946 }
947
948 if (!enoughPoints) {
949 if (!(incrDone[VoxX] || incrDone[VoxF] || incrDone[VoxZ])) {
950 LOG(error) << fmt::format("trial limit reached, skipping this voxel: incrDone[VoxX] {}, incrDone[VoxF] {}, incrDone[VoxZ] {}", incrDone[VoxX], incrDone[VoxF], incrDone[VoxZ]);
951 return false;
952 }
953 LOG(debug) << "sector " << iSec << ": increasing filter bandwidth around voxel " << binCenter;
954 // printf("Sector:%2d x=%.2f y/x=%.2f z/x=%.2f (iX: %d iY2X:%d iZ2X:%d)\n", iSec, x, p, z, ix0, ip0, iz0);
955 // printf("not enough neighbours (need min %d) %d %d %d (tot: %d) | Steps: %.1f %.1f %.1f\n", 2, nPoints[VoxX], nPoints[VoxF], nPoints[VoxZ], nbOK, stepX, stepF, stepZ);
956 // printf("trying to increase filter bandwidth (trialXFZ: %d %d %d)\n", trial[VoxX], trial[VoxF], trial[VoxZ]);
957 continue;
958 }
959
960 // now fill matrices and solve
961 for (int iNb = 0; iNb < nbOK; ++iNb) {
962 double wiCache = currCache[iNb * VoxHDim + VoxV];
963 double dxi = currCache[iNb * VoxHDim + VoxX];
964 double dfi = currCache[iNb * VoxHDim + VoxF];
965 double dzi = currCache[iNb * VoxHDim + VoxZ];
966 double dxi2 = dxi * dxi;
967 double dfi2 = dfi * dfi;
968 double dzi2 = dzi * dzi;
969 const VoxRes* voxNb = currVox[iNb];
970 for (int iDim = 0; iDim < ResDim; ++iDim) {
971 if (!doDim[iDim]) {
972 continue;
973 }
974 double vi = voxNb->D[iDim];
975 double wi = wiCache;
976 if (mUseErrInSmoothing && fabs(voxNb->E[iDim]) > 1e-6) {
977 // account for point error apart from kernel value
978 wi /= (voxNb->E[iDim] * voxNb->E[iDim]);
979 }
980 std::array<double, sMaxSmtDim*(sMaxSmtDim + 1) / 2>& cmatD = cmat[iDim];
981 double* rhsD = &mLastSmoothingRes[iDim * sMaxSmtDim];
982 unsigned short iMat = 0;
983 unsigned short iRhs = 0;
984 // linear part
985 cmatD[iMat++] += wi;
986 rhsD[iRhs++] += wi * vi;
987 //
988 cmatD[iMat++] += wi * dxi;
989 cmatD[iMat++] += wi * dxi2;
990 rhsD[iRhs++] += wi * dxi * vi;
991 //
992 cmatD[iMat++] += wi * dfi;
993 cmatD[iMat++] += wi * dxi * dfi;
994 cmatD[iMat++] += wi * dfi2;
995 rhsD[iRhs++] += wi * dfi * vi;
996 //
997 cmatD[iMat++] += wi * dzi;
998 cmatD[iMat++] += wi * dxi * dzi;
999 cmatD[iMat++] += wi * dfi * dzi;
1000 cmatD[iMat++] += wi * dzi2;
1001 rhsD[iRhs++] += wi * dzi * vi;
1002 //
1003 // check if quadratic part is needed
1004 if (mSmoothPol2[VoxX]) {
1005 cmatD[iMat++] += wi * dxi2;
1006 cmatD[iMat++] += wi * dxi * dxi2;
1007 cmatD[iMat++] += wi * dfi * dxi2;
1008 cmatD[iMat++] += wi * dzi * dxi2;
1009 cmatD[iMat++] += wi * dxi2 * dxi2;
1010 rhsD[iRhs++] += wi * dxi2 * vi;
1011 }
1012 if (mSmoothPol2[VoxF]) {
1013 cmatD[iMat++] += wi * dfi2;
1014 cmatD[iMat++] += wi * dxi * dfi2;
1015 cmatD[iMat++] += wi * dfi * dfi2;
1016 cmatD[iMat++] += wi * dzi * dfi2;
1017 cmatD[iMat++] += wi * dxi2 * dfi2;
1018 cmatD[iMat++] += wi * dfi2 * dfi2;
1019 rhsD[iRhs++] += wi * dfi2 * vi;
1020 }
1021 if (mSmoothPol2[VoxZ]) {
1022 cmatD[iMat++] += wi * dzi2;
1023 cmatD[iMat++] += wi * dxi * dzi2;
1024 cmatD[iMat++] += wi * dfi * dzi2;
1025 cmatD[iMat++] += wi * dzi * dzi2;
1026 cmatD[iMat++] += wi * dxi2 * dzi2;
1027 cmatD[iMat++] += wi * dfi2 * dzi2;
1028 cmatD[iMat++] += wi * dzi2 * dzi2;
1029 rhsD[iRhs++] += wi * dzi2 * vi;
1030 }
1031 }
1032 }
1033
1034 bool fitRes = true;
1035
1036 // solve system of linear equations
1037
1038 TMatrixDSym matrix(matSize);
1039 TDecompChol chol(matSize);
1040 TVectorD rhsVec(matSize);
1041 for (int iDim = 0; iDim < ResDim; ++iDim) {
1042 if (!doDim[iDim]) {
1043 continue;
1044 }
1045 matrix.Zero(); // reset matrix
1046 std::array<double, sMaxSmtDim*(sMaxSmtDim + 1) / 2>& cmatD = cmat[iDim];
1047 double* rhsD = &mLastSmoothingRes[iDim * sMaxSmtDim];
1048 short iMat = -1;
1049 short row = -1;
1050
1051 // with the studid implementation of TMatrixDSym we need to set all elements of the matrix explicitly (or maybe only upper triangle?)
1052 matrix(++row, 0) = cmatD[++iMat];
1053 matrix(++row, 0) = cmatD[++iMat];
1054 matrix(row, 1) = cmatD[++iMat];
1055 matrix(0, row) = matrix(row, 0);
1056 matrix(++row, 0) = cmatD[++iMat];
1057 matrix(row, 1) = cmatD[++iMat];
1058 matrix(row, 2) = cmatD[++iMat];
1059 matrix(0, row) = matrix(row, 0);
1060 matrix(1, row) = matrix(row, 1);
1061 matrix(++row, 0) = cmatD[++iMat];
1062 matrix(row, 1) = cmatD[++iMat];
1063 matrix(row, 2) = cmatD[++iMat];
1064 matrix(row, 3) = cmatD[++iMat];
1065 matrix(0, row) = matrix(row, 0);
1066 matrix(1, row) = matrix(row, 1);
1067 matrix(2, row) = matrix(row, 2);
1068 // add pol2 elements if needed
1069 if (mSmoothPol2[VoxX]) {
1070 const int colLim = (++row) + 1;
1071 for (int iCol = 0; iCol < colLim; ++iCol) {
1072 matrix(row, iCol) = cmatD[++iMat];
1073 matrix(iCol, row) = matrix(row, iCol);
1074 }
1075 }
1076 if (mSmoothPol2[VoxF]) {
1077 const int colLim = (++row) + 1;
1078 for (int iCol = 0; iCol < colLim; ++iCol) {
1079 matrix(row, iCol) = cmatD[++iMat];
1080 matrix(iCol, row) = matrix(row, iCol);
1081 }
1082 }
1083 if (mSmoothPol2[VoxZ]) {
1084 const int colLim = (++row) + 1;
1085 for (int iCol = 0; iCol < colLim; ++iCol) {
1086 matrix(row, iCol) = cmatD[++iMat];
1087 matrix(iCol, row) = matrix(row, iCol);
1088 }
1089 }
1090 rhsVec.SetElements(rhsD);
1091 chol.SetMatrix(matrix);
1092 chol.Decompose();
1093 fitRes = chol.Solve(rhsVec);
1094 if (!fitRes) {
1095 for (int i = VoxDim; i--;) {
1096 trial[i]++;
1097 }
1098 LOG(error) << "solution for smoothing failed, trying to increase filter bandwidth";
1099 continue;
1100 }
1101 res[iDim] = rhsVec[0];
1102 }
1103
1104 break;
1105 }
1106
1107 return true;
1108}
1109
1110double TrackResiduals::getKernelWeight(std::array<double, 3> u2vec) const
1111{
1112 double w = 1.;
1113 if (mKernelType == KernelType::Epanechnikov) {
1114 for (size_t i = u2vec.size(); i--;) {
1115 if (u2vec[i] > 1) {
1116 return 0.;
1117 }
1118 w *= 3. / 4. * (1. - u2vec[i]);
1119 }
1120 } else if (mKernelType == KernelType::Gaussian) {
1121 double u2 = 0.;
1122 for (size_t i = u2vec.size(); i--;) {
1123 u2 += u2vec[i];
1124 }
1125 w = u2 < mParams->maxGaussStdDev * mParams->maxGaussStdDev * u2vec.size() ? std::exp(-u2) / std::sqrt(2. * M_PI) : 0;
1126 }
1127 return w;
1128}
1129
1135
1136float TrackResiduals::fitPoly1Robust(std::vector<float>& x, std::vector<float>& y, std::array<float, 2>& res, std::array<float, 3>& err, float cutLTM) const
1137{
1138 // robust pol1 fit, modifies input arrays order
1139 if (x.size() != y.size()) {
1140 LOG(error) << "x and y must not have different sizes for fitPoly1Robust (" << x.size() << " != " << y.size() << ")";
1141 }
1142 size_t nPoints = x.size();
1143 res[0] = res[1] = 0.f;
1144 if (nPoints < 2) {
1145 return -1;
1146 }
1147 std::array<float, 7> yResults;
1148 std::vector<size_t> indY(nPoints);
1149 if (!o2::math_utils::LTMUnbinned(y, indY, yResults, cutLTM)) {
1150 return -1;
1151 }
1152 // rearrange used events in increasing order
1155 //
1156 // 1st fit to get crude slope
1157 int nPointsUsed = std::lrint(yResults[0]);
1158 int vecOffset = std::lrint(yResults[5]);
1159 // use only entries selected by LTM for the fit
1160 float a, b;
1161 medFit(nPointsUsed, vecOffset, x, y, a, b, err);
1162 //
1163 std::vector<float> ycm(nPoints);
1164 for (size_t i = nPoints; i-- > 0;) {
1165 ycm[i] = y[i] - (a + b * x[i]);
1166 }
1167 std::vector<size_t> indices(nPoints);
1172 //
1173 // robust estimate of sigma after crude slope correction
1174 float sigMAD = getMAD2Sigma({ycm.begin() + vecOffset, ycm.begin() + vecOffset + nPointsUsed});
1175 // find LTM estimate matching to sigMAD, keaping at least given fraction
1176 if (!o2::math_utils::LTMUnbinnedSig(ycm, indY, yResults, mParams->minFracLTM, sigMAD, true)) {
1177 return -1;
1178 }
1179 // final fit
1180 nPointsUsed = std::lrint(yResults[0]);
1181 vecOffset = std::lrint(yResults[5]);
1182 medFit(nPointsUsed, vecOffset, x, y, a, b, err);
1183 res[0] = a;
1184 res[1] = b;
1185 return sigMAD;
1186}
1187
1188//___________________________________________________________________
1189void TrackResiduals::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
1190{
1191 // fitting a straight line y(x|a, b) = a + b * x
1192 // to given x and y data minimizing the absolute deviation
1193 float aa, bb, chi2 = 0.f;
1194 if (nPoints < 2) {
1195 a = b = 0.f;
1196 err[0] = err[1] = err[2] = 999.f;
1197 return;
1198 }
1199 // do least squares minimization as first guess
1200 float sx = 0.f, sxx = 0.f, sy = 0.f, sxy = 0.f;
1201 for (int j = nPoints + offset; j-- > offset;) { // same order as in AliRoot version such that resulting sums are identical
1202 sx += x[j];
1203 sxx += x[j] * x[j];
1204 sy += y[j];
1205 sxy += x[j] * y[j];
1206 }
1207 float del = nPoints * sxx - sx * sx;
1208 float delI = 1. / del;
1209 aa = (sxx * sy - sx * sxy) * delI;
1210 bb = (nPoints * sxy - sx * sy) * delI;
1211 err[0] = sxx * delI;
1212 err[1] = sx * delI;
1213 err[2] = nPoints * delI;
1214
1215 for (int j = nPoints + offset; j-- > offset;) {
1216 float tmp = y[j] - (aa + bb * x[j]);
1217 chi2 += tmp * tmp;
1218 }
1219 float sigb = std::sqrt(chi2 * delI); // expected sigma for b
1220 float b1 = bb;
1221 float f1 = roFunc(nPoints, offset, x, y, b1, aa);
1222 if (sigb > 0) {
1223 float b2 = bb + std::copysign(3.f * sigb, f1);
1224 float f2 = roFunc(nPoints, offset, x, y, b2, aa);
1225 if (fabs(f1 - f2) < sFloatEps) {
1226 a = aa;
1227 b = bb;
1228 return;
1229 }
1230 while (f1 * f2 > 0.f) {
1231 bb = b2 + 1.6f * (b2 - b1);
1232 b1 = b2;
1233 f1 = f2;
1234 b2 = bb;
1235 f2 = roFunc(nPoints, offset, x, y, b2, aa);
1236 }
1237 sigb = .01f * sigb;
1238 while (fabs(b2 - b1) > sigb) {
1239 bb = b1 + .5f * (b2 - b1);
1240 if (bb == b1 || bb == b2) {
1241 break;
1242 }
1243 float f = roFunc(nPoints, offset, x, y, bb, aa);
1244 if (f * f1 >= .0f) {
1245 f1 = f;
1246 b1 = bb;
1247 } else {
1248 f2 = f;
1249 b2 = bb;
1250 }
1251 }
1252 }
1253 a = aa;
1254 b = bb;
1255}
1256
1257float TrackResiduals::roFunc(int nPoints, int offset, const std::vector<float>& x, const std::vector<float>& y, float b, float& aa) const
1258{
1259 // calculate sum(x_i * sgn(y_i - a - b * x_i)) for given b
1260 // see numberical recipies paragraph 15.7.3
1261 std::vector<float> vecTmp(nPoints);
1262 float sum = 0.f;
1263 for (int j = nPoints; j-- > 0;) {
1264 vecTmp[j] = y[j + offset] - b * x[j + offset];
1265 }
1266 int nPointsHalf = nPoints / 2;
1267 if (nPoints < 20) { // it is faster to do insertion sort
1268 for (int i = 1; i < nPoints; i++) {
1269 float v = vecTmp[i];
1270 int j;
1271 for (j = i; j--;) {
1272 if (vecTmp[j] > v) {
1273 vecTmp[j + 1] = vecTmp[j];
1274 } else {
1275 break;
1276 }
1277 }
1278 vecTmp[j + 1] = v;
1279 }
1280 aa = (nPoints & 0x1) ? vecTmp[nPointsHalf] : .5f * (vecTmp[nPointsHalf - 1] + vecTmp[nPointsHalf]);
1281 } else {
1282 std::vector<float>::iterator nth = vecTmp.begin() + vecTmp.size() / 2;
1283 if (nPoints & 0x1) {
1284 std::nth_element(vecTmp.begin(), nth, vecTmp.end());
1285 aa = *nth;
1286 } else {
1287 std::nth_element(vecTmp.begin(), nth - 1, vecTmp.end());
1288 std::nth_element(nth, nth, vecTmp.end());
1289 aa = 0.5 * (*(nth - 1) + *(nth));
1290 }
1291 // aa = (nPoints & 0x1) ? selectKthMin(nPointsHalf, vecTmp) : .5f * (selectKthMin(nPointsHalf - 1, vecTmp) + selectKthMin(nPointsHalf, vecTmp));
1292 }
1293 for (int j = nPoints; j-- > 0;) {
1294 float d = y[j + offset] - (b * x[j + offset] + aa);
1295 if (y[j + offset] != 0.f) {
1296 d /= fabs(y[j + offset]);
1297 }
1298 if (fabs(d) > sFloatEps) {
1299 sum += (d >= 0.f ? x[j + offset] : -x[j + offset]);
1300 }
1301 }
1302 return sum;
1303}
1304
1305//______________________________________________________________________________
1306float TrackResiduals::selectKthMin(const int k, std::vector<float>& data) const
1307{
1308 // Returns the k th smallest value in the vector. The input vector will be rearranged
1309 // to have this value in location data[k] , with all smaller elements moved before it
1310 // (in arbitrary order) and all larger elements after (also in arbitrary order).
1311 // From Numerical Recipes in C++ (paragraph 8.5)
1312 // Probably it is not needed anymore, since std::nth_element() can also be used
1313
1314 int i, ir, j, l, mid, n = data.size();
1315 float a; // partitioning element
1316 l = 0; // left hand side of active partition
1317 ir = n - 1; // right hand side of active partition
1318
1319 while (true) {
1320 if (ir <= l + 1) { // active partition with 1 or 2 elements
1321 if (ir == l + 1 && data[ir] < data[l]) { // case of 2 elements
1322 std::swap(data[l], data[ir]);
1323 }
1324 return data[k];
1325 } else {
1326 mid = (l + ir) >> 1;
1327 std::swap(data[mid], data[l + 1]);
1328 if (data[l] > data[ir]) {
1329 std::swap(data[l], data[ir]);
1330 }
1331 if (data[l + 1] > data[ir]) {
1332 std::swap(data[l + 1], data[ir]);
1333 }
1334 if (data[l] > data[l + 1]) {
1335 std::swap(data[l], data[l + 1]);
1336 }
1337 i = l + 1; // initialize pointers for partitioning
1338 j = ir;
1339 a = data[l + 1];
1340 while (true) {
1341 // innermost loop used for partitioning
1342 do {
1343 i++;
1344 } while (data[i] < a);
1345 do {
1346 j--;
1347 } while (data[j] > a);
1348 if (j < i) {
1349 // pointers crossed -> partitioning complete
1350 break;
1351 }
1352 std::swap(data[i], data[j]);
1353 }
1354 data[l + 1] = data[j];
1355 data[j] = a;
1356 // keep the partition which contains kth element active
1357 if (j >= k) {
1358 ir = j - 1;
1359 }
1360 if (j <= k) {
1361 l = i;
1362 }
1363 }
1364 }
1365}
1366
1367//___________________________________________________________________
1368float TrackResiduals::getMAD2Sigma(std::vector<float> data) const
1369{
1370 // Sigma calculated from median absolute deviations
1371 // see: https://en.wikipedia.org/wiki/Median_absolute_deviation
1372 // the data is passed by value (copied!), such that the original vector
1373 // is not rearranged
1374
1375 int nPoints = data.size();
1376 if (nPoints < 2) {
1377 return 0;
1378 }
1379
1380 // calculate median of the input data
1381 float medianOfData;
1382 std::vector<float>::iterator nth = data.begin() + data.size() / 2;
1383 if (nPoints & 0x1) {
1384 std::nth_element(data.begin(), nth, data.end());
1385 medianOfData = *nth;
1386 } else {
1387 std::nth_element(data.begin(), nth - 1, data.end());
1388 std::nth_element(nth, nth, data.end());
1389 medianOfData = .5f * (*(nth - 1) + (*nth));
1390 }
1391
1392 // fill vector with absolute deviations to median
1393 for (auto& entry : data) {
1394 entry = fabs(entry - medianOfData);
1395 }
1396
1397 // calculate median of abs deviations
1398 float medianOfAbsDeviations;
1399 if (nPoints & 0x1) {
1400 std::nth_element(data.begin(), nth, data.end());
1401 medianOfAbsDeviations = *nth;
1402 } else {
1403 std::nth_element(data.begin(), nth - 1, data.end());
1404 std::nth_element(nth, nth, data.end());
1405 medianOfAbsDeviations = .5f * (*(nth - 1) + (*nth));
1406 }
1407
1408 float k = 1.4826f; // scale factor for normally distributed data
1409 return k * medianOfAbsDeviations;
1410}
1411
1412void TrackResiduals::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)
1413{
1414 // this fast algebraic circle fit is described here:
1415 // https://dtcenter.org/met/users/docs/write_ups/circle_fit.pdf
1416 float xMean = 0.f;
1417 float yMean = 0.f;
1418
1419 for (int i = nCl; i--;) {
1420 xMean += x[i];
1421 yMean += y[i];
1422 }
1423 xMean /= nCl;
1424 yMean /= nCl;
1425 // define sums needed for circular fit
1426 float su2 = 0.f, sv2 = 0.f, suv = 0.f, su3 = 0.f, sv3 = 0.f, su2v = 0.f, suv2 = 0.f;
1427 for (int i = nCl; i--;) {
1428 float ui = x[i] - xMean;
1429 float vi = y[i] - yMean;
1430 float ui2 = ui * ui;
1431 float vi2 = vi * vi;
1432 suv += ui * vi;
1433 su2 += ui2;
1434 sv2 += vi2;
1435 su3 += ui2 * ui;
1436 sv3 += vi2 * vi;
1437 su2v += ui2 * vi;
1438 suv2 += ui * vi2;
1439 }
1440 float rhsU = .5f * (su3 + suv2);
1441 float rhsV = .5f * (sv3 + su2v);
1442 float det = su2 * sv2 - suv * suv;
1443 float uc = (rhsU * sv2 - rhsV * suv) / det;
1444 float vc = (su2 * rhsV - suv * rhsU) / det;
1445 float r2 = uc * uc + vc * vc + (su2 + sv2) / nCl;
1446 xc = uc + xMean;
1447 yc = vc + yMean;
1448 r = sqrt(r2);
1449 // write residuals to residHelixY
1450 for (int i = nCl; i--;) {
1451 float dx = x[i] - xc;
1452 float dxr = r2 - dx * dx;
1453 float ys = dxr > 0 ? sqrt(dxr) : 0.f; // distance of point in y from the circle center (using fit results for r and xc)
1454 float dy = y[i] - yc; // distance of point in y from the circle center (using fit result for yc)
1455 float dysp = dy - ys;
1456 float dysm = dy + ys;
1457 residHelixY[i] = fabsf(dysp) < fabsf(dysm) ? dysp : dysm;
1458 }
1459 // printf("r = %.4f m, xc = %.4f, yc = %.4f\n", r/100.f, xc, yc);
1460}
1461
1462bool TrackResiduals::fitPoly1(int nCl, std::array<float, param::NPadRows>& x, std::array<float, param::NPadRows>& y, std::array<float, 2>& res)
1463{
1464 // fit a straight line y = ax + b to a given set of points (x,y)
1465 // no measurement errors assumed, no fit errors calculated
1466 // res[0] = a (slope)
1467 // res[1] = b (offset)
1468 if (nCl < 2) {
1469 // not enough points
1470 return false;
1471 }
1472 float sumX = 0.f, sumY = 0.f, sumXY = 0.f, sumX2 = 0.f, nInv = 1.f / nCl;
1473 for (int i = nCl; i--;) {
1474 sumX += x[i];
1475 sumY += y[i];
1476 sumXY += x[i] * y[i];
1477 sumX2 += x[i] * x[i];
1478 }
1479 float det = sumX2 - nInv * sumX * sumX;
1480 if (fabsf(det) < 1e-12f) {
1481 return false;
1482 }
1483 res[0] = (sumXY - nInv * sumX * sumY) / det;
1484 res[1] = nInv * sumY - nInv * res[0] * sumX;
1485 return true;
1486}
1487
1493
1495{
1496 if (getNVoxelsPerSector() == 0) {
1497 LOG(warn) << "For the tree aliases to work you must initialize the binning before calling createOutputFile()";
1498 }
1499 mFileOut = std::make_unique<TFile>(filename, "recreate");
1500 mTreeOut = std::make_unique<TTree>("voxResTree", "Voxel results and statistics");
1501 mTreeOut->SetAlias("z2xBin", "bvox[0]");
1502 mTreeOut->SetAlias("y2xBin", "bvox[1]");
1503 mTreeOut->SetAlias("xBin", "bvox[2]");
1504 mTreeOut->SetAlias("z2xAV", "stat[0]");
1505 mTreeOut->SetAlias("y2xAV", "stat[1]");
1506 mTreeOut->SetAlias("xAV", "stat[2]");
1507 mTreeOut->SetAlias("fsector", "bsec+0.5+9.*(y2xAV)/pi");
1508 mTreeOut->SetAlias("phi", "(bsec%18+0.5+9.*(stat[1])/pi)/9*pi");
1509 mTreeOut->SetAlias("r", "stat[2]");
1510 mTreeOut->SetAlias("z", "z2xAV*xAV");
1511 mTreeOut->SetAlias("dX", "D[0]");
1512 mTreeOut->SetAlias("dY", "D[1]");
1513 mTreeOut->SetAlias("dZ", "D[2]");
1514 mTreeOut->SetAlias("dXS", "DS[0]");
1515 mTreeOut->SetAlias("dYS", "DS[1]");
1516 mTreeOut->SetAlias("dZS", "DS[2]");
1517 mTreeOut->SetAlias("dXE", "E[0]");
1518 mTreeOut->SetAlias("dYE", "E[1]");
1519 mTreeOut->SetAlias("dZE", "E[2]");
1520 mTreeOut->SetAlias("voxelIndex", Form("xBin + %i * (y2xBin + %i * z2xBin) + %i * bsec", getNXBins(), getNY2XBins(), getNVoxelsPerSector()));
1521 mTreeOut->SetAlias("entries", "stat[3]");
1522 mTreeOut->SetAlias("fitOK", Form("(flags & %u) == %u", DistDone, DistDone));
1523 mTreeOut->SetAlias("dispOK", Form("(flags & %u) == %u", DispDone, DispDone));
1524 mTreeOut->SetAlias("smtOK", Form("(flags & %u) == %u", SmoothDone, SmoothDone));
1525 mTreeOut->SetAlias("masked", Form("(flags & %u) == %u", Masked, Masked));
1526 mTreeOut->Branch("voxRes", &mVoxelResultsOutPtr);
1527}
1528
1530{
1531 mFileOut->cd();
1532 mTreeOut->Write();
1533 mTreeOut.reset();
1534 mFileOut->Close();
1535 mFileOut.reset();
1536}
1537
1539{
1540 if (mTreeOut) {
1541 printf("Dumping results for sector %i. Don't forget the call to closeOutputFile() in the end...\n", iSec);
1542 for (int i = 0; i < mNVoxPerSector; ++i) {
1543 mVoxelResultsOut = mVoxelResults[iSec][i];
1544 mTreeOut->Fill();
1545 }
1546 }
1547}
1548
1550{
1551 static float mres = 0, mvir = 0, mres0 = 0, mvir0 = 0;
1552 static ProcInfo_t procInfo;
1553 static TStopwatch sw;
1554 const Long_t kMB = 1024;
1555 gSystem->GetProcInfo(&procInfo);
1556 mres = float(procInfo.fMemResident) / kMB;
1557 mvir = float(procInfo.fMemVirtual) / kMB;
1558 sw.Stop();
1559 printf("RSS: %.3f(%.3f) VMEM: %.3f(%.3f) MB | CpuTime:%.3f RealTime:%.3f s\n",
1560 mres, mres - mres0, mvir, mvir - mvir0, sw.CpuTime(), sw.RealTime());
1561 mres0 = mres;
1562 mvir0 = mvir;
1563 sw.Start();
1564}
1565
1567{
1568 getLocalResVec().clear();
1569 mInitResultsContainer.reset();
1570}
Base track model for the Barrel, params only, w/o covariance.
const auto bins
Definition PID.cxx:49
std::ostringstream debug
int32_t i
const int16_t bb
useful math constants
uint32_t j
Definition RawData.h:0
uint32_t res
Definition RawData.h:0
Definition of the TrackResiduals class.
std::vector< LocalResid > & getLocalResVec()
float getDZ2XI(int iz=0) const
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)
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)
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)
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 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)
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.
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)
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
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)
@ 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 setNZ2XBins(int nBins)
void processVoxelResiduals(std::vector< float > &dy, std::vector< float > &dz, std::vector< float > &tg, VoxRes &resVox)
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLuint entry
Definition glcorearb.h:5735
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
GLdouble f
Definition glcorearb.h:310
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLboolean * data
Definition glcorearb.h:298
GLintptr offset
Definition glcorearb.h:660
GLsizei GLenum const void * indices
Definition glcorearb.h:400
GLboolean r
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
void Reorder(std::vector< T > &data, const std::vector< size_t > &index)
Definition fit.h:626
bool LTMUnbinnedSig(const std::vector< T > &data, std::vector< size_t > &index, std::array< float, 7 > &params, float fracKeepMin, float sigTgt, bool sorted=false)
Definition fit.h:649
void SortData(std::vector< T > const &values, std::vector< size_t > &index)
Definition fit.h:540
bool LTMUnbinned(const std::vector< T > &data, std::vector< size_t > &index, std::array< float, 7 > &params, float fracKeep)
Definition fit.h:570
Global TPC definitions and constants.
Definition SimTraits.h:168
constexpr double SECPHIWIDTH
Definition Defs.h:45
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
Enum< T >::Iterator begin(Enum< T >)
Definition Defs.h:156
constexpr unsigned char SIDES
Definition Defs.h:41
DataT sum(const Vector< DataT, N > &a)
Definition Vector.h:107
FIXME: do not use data model tables.
std::string filename()
float minFracLTM
minimum fraction of points to keep when trimming data to fit expected sigma
float maxFitErrX2
maximum fit error for X2
float minValidVoxFracDrift
if more than this fraction of bins are bad for one pad row the whole pad row is declared bad
float maxGaussStdDev
maximum number of sigmas to be considered for gaussian kernel smoothing
int minEntriesPerVoxel
minimum number of points in voxel for processing
float LTMCut
fraction op points to keep when trimming input data
int maxBadXBinsToCover
a lower number of consecutive bad X bins will not be declared bad
float maxFitErrY2
maximum fit error for Y2
float maxSigY
maximum sigma for y of the voxel
float maxFitCorrXY
maximum fit correlation for x and y
float maxSigZ
maximum sigma for z of the voxel
int minGoodXBinsToCover
minimum number of consecutive good bins, otherwise bins are declared bad
float maxFracBadRowsPerSector
maximum fraction of bad rows before whole sector is masked
bool isBfieldZero
for B=0 we set the radial distortions to zero and don't fit dy vs tgSlp
Structure which gets filled with the results for each voxel.
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.
unsigned char flags
status flag
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
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
o2::InteractionRecord ir(0, 0)
TStopwatch sw
std::vector< int > row