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