192 GPUd()
bool recalculatePCAWithErrors(
int cand = 0);
198 auto m = calcPCACovMatrix(
cand);
199 return {
static_cast<float>(
m(0, 0)),
static_cast<float>(
m(1, 0)),
static_cast<float>(
m(1, 1)),
static_cast<float>(
m(2, 0)),
static_cast<float>(
m(2, 1)),
static_cast<float>(
m(2, 2))};
206 GPUdi()
void setPropagateToPCA(
bool v = true) { mPropagateToPCA =
v; }
207 GPUdi()
void setMaxIter(
int n = 20) { mMaxIter =
n > 2 ?
n : 2; }
208 GPUdi()
void setMaxR(
float r = 200.) { mMaxR2 =
r *
r; }
209 GPUdi()
void setMaxDZIni(
float d = 4.) { mMaxDZIni = d; }
210 GPUdi()
void setMaxDXYIni(
float d = 4.) { mMaxDXYIni = d > 0 ? d : 1e9; }
211 GPUdi()
void setMaxChi2(
float chi2 = 999.) { mMaxChi2 = chi2; }
213 GPUdi()
void setMinParamChange(
float x = 1e-3) { mMinParamChange =
x > 1e-4 ?
x : 1.e-4; }
214 GPUdi()
void setMinRelChi2Change(
float r = 0.9) { mMinRelChi2Change =
r > 0.1 ?
r : 999.; }
215 GPUdi()
void setUseAbsDCA(
bool v) { mUseAbsDCA =
v; }
216 GPUdi()
void setWeightedFinalPCA(
bool v) { mWeightedFinalPCA =
v; }
217 GPUdi()
void setMaxDistance2ToMerge(
float v) { mMaxDist2ToMergeSeeds =
v; }
219 GPUdi()
void setUsePropagator(
bool v) { mUsePropagator =
v; }
220 GPUdi()
void setRefitWithMatCorr(
bool v) { mRefitWithMatCorr =
v; }
221 GPUdi()
void setMaxSnp(
float s) { mMaxSnp =
s; }
222 GPUdi()
void setMaxStep(
float s) { mMaxStep =
s; }
223 GPUdi()
void setMinXSeed(
float x) { mMinXSeed =
x; }
224 GPUdi()
void setCollinear(
bool isCollinear) { mIsCollinear = isCollinear; }
226 GPUdi()
int getNCandidates()
const {
return mCurHyp; }
227 GPUdi()
int getMaxIter()
const {
return mMaxIter; }
228 GPUdi() float getMaxR()
const {
return o2::gpu::GPUCommonMath::Sqrt(mMaxR2); }
229 GPUdi() float getMaxDZIni()
const {
return mMaxDZIni; }
230 GPUdi() float getMaxDXYIni()
const {
return mMaxDXYIni; }
231 GPUdi() float getMaxChi2()
const {
return mMaxChi2; }
232 GPUdi() float getMinParamChange()
const {
return mMinParamChange; }
233 GPUdi() float getBz()
const {
return mBz; }
234 GPUdi() float getMaxDistance2ToMerge()
const {
return mMaxDist2ToMergeSeeds; }
235 GPUdi() bool getUseAbsDCA()
const {
return mUseAbsDCA; }
236 GPUdi() bool getWeightedFinalPCA()
const {
return mWeightedFinalPCA; }
237 GPUdi() bool getPropagateToPCA()
const {
return mPropagateToPCA; }
239 GPUdi() bool getUsePropagator()
const {
return mUsePropagator; }
240 GPUdi() bool getRefitWithMatCorr()
const {
return mRefitWithMatCorr; }
241 GPUdi() float getMaxSnp()
const {
return mMaxSnp; }
242 GPUdi() float getMasStep()
const {
return mMaxStep; }
243 GPUdi() float getMinXSeed()
const {
return mMinXSeed; }
245 template <
class... Tr>
249 GPUdi()
int getFitterID()
const {
return mFitterID; }
251 GPUdi() size_t getCallID()
const {
return mCallID; }
255 GPUd()
bool calcInverseWeight();
256 GPUd()
void calcResidDerivatives();
257 GPUd()
void calcResidDerivativesNoErr();
259 GPUd()
void calcChi2Derivatives();
260 GPUd()
void calcChi2DerivativesNoErr();
263 GPUd()
void calcTrackResiduals();
264 GPUd()
void calcTrackDerivatives();
265 GPUd()
double calcChi2() const;
266 GPUd()
double calcChi2NoErr() const;
269 GPUd()
bool minimizeChi2NoErr();
270 GPUd()
bool roughDZCut() const;
271 GPUd()
bool closerToAlternative() const;
273 GPUd()
bool propagateParamToX(
o2::track::TrackPar&
t,
float x);
280 GPUd() float getTrackX(
int i,
int cand = 0)
const {
return getTrackPos(
i,
cand)[0]; }
282 GPUd() MatStd3D getTrackRotMatrix(
int i) const
286 mat(0, 0) = mat(1, 1) = mTrAux[
i].c;
287 mat(0, 1) = -mTrAux[
i].s;
288 mat(1, 0) = mTrAux[
i].s;
292 GPUd() MatSym3D getTrackCovMatrix(
int i,
int cand = 0) const
294 const auto& trc = mCandTr[mOrder[
cand]][
i];
296 mat(0, 0) = trc.getSigmaY2() * XerrFactor;
297 mat(1, 1) = trc.getSigmaY2();
298 mat(2, 2) = trc.getSigmaZ2();
299 mat(2, 1) = trc.getSigmaZY();
304 template <
class T,
class... Tr>
307#ifndef GPUCA_GPUCODE_DEVICE
308 static_assert(std::is_convertible<T, Track>(),
"Wrong track type");
317 mAllowAltPreference =
true;
329 mLoggerBadCov.clear();
330 mLoggerBadInv.clear();
331 mLoggerBadProp.clear();
349 CrossInfo mCrossings;
362 LogLogThrottler mLoggerBadCov{};
363 LogLogThrottler mLoggerBadInv{};
364 LogLogThrottler mLoggerBadProp{};
369 int mCrossIDAlt = -1;
371 bool mAllowAltPreference =
true;
372 bool mUseAbsDCA =
false;
373 bool mWeightedFinalPCA =
false;
374 bool mPropagateToPCA =
true;
375 bool mUsePropagator =
false;
376 bool mRefitWithMatCorr =
false;
377 bool mIsCollinear =
false;
381 float mMaxR2 = 200. * 200.;
382 float mMinXSeed = -50.;
383 float mMaxDZIni = 4.;
384 float mMaxDXYIni = 4.;
385 float mMinParamChange = 1e-3;
386 float mMinRelChi2Change = 0.9;
387 float mMaxChi2 = 100;
388 float mMaxDist2ToMergeSeeds = 1.;
389 float mMaxSnp = 0.95;
390 float mMaxStep = 2.0;
397template <
int N,
typename... Args>
398template <
class... Tr>
403 static_assert(
sizeof...(args) == N,
"incorrect number of input tracks");
406 for (
int i = 0;
i < N;
i++) {
407 mTrAux[
i].set(*mOrigTrPtr[
i], mBz);
409 if (!mCrossings.set(mTrAux[0], *mOrigTrPtr[0], mTrAux[1], *mOrigTrPtr[1], mMaxDXYIni, mIsCollinear)) {
412 for (
int ih = 0; ih < MAXHYP; ih++) {
413 mPropFailed[ih] =
false;
418 if (mCrossings.nDCA == MAXHYP) {
419 auto dst2 = (mCrossings.xDCA[0] - mCrossings.xDCA[1]) * (mCrossings.xDCA[0] - mCrossings.xDCA[1]) +
420 (mCrossings.yDCA[0] - mCrossings.yDCA[1]) * (mCrossings.yDCA[0] - mCrossings.yDCA[1]);
421 if (dst2 < mMaxDist2ToMergeSeeds) {
423 mCrossings.xDCA[0] = 0.5 * (mCrossings.xDCA[0] + mCrossings.xDCA[1]);
424 mCrossings.yDCA[0] = 0.5 * (mCrossings.yDCA[0] + mCrossings.yDCA[1]);
428 for (
int ic = 0; ic < mCrossings.nDCA; ic++) {
430 if (mCrossings.xDCA[ic] * mCrossings.xDCA[ic] + mCrossings.yDCA[ic] * mCrossings.yDCA[ic] > mMaxR2) {
434 mCrossIDAlt = (mCrossings.nDCA == 2 && mAllowAltPreference) ? 1 - ic : -1;
435 mNIters[mCurHyp] = 0;
436 mTrPropDone[mCurHyp] =
false;
437 mChi2[mCurHyp] = -1.;
438 mPCA[mCurHyp][0] = mCrossings.xDCA[ic];
439 mPCA[mCurHyp][1] = mCrossings.yDCA[ic];
441 if (mUseAbsDCA ? minimizeChi2NoErr() : minimizeChi2()) {
442 mOrder[mCurHyp] = mCurHyp;
443 if (mPropagateToPCA && !propagateTracksToVertex(mCurHyp)) {
450 for (
int i = mCurHyp;
i--;) {
451 for (
int j =
i;
j--;) {
452 if (mChi2[mOrder[
i]] < mChi2[mOrder[
j]]) {
453 o2::gpu::GPUCommonMath::Swap(mOrder[
i], mOrder[
j]);
457 if (mUseAbsDCA && mWeightedFinalPCA) {
458 for (
int i = mCurHyp;
i--;) {
459 recalculatePCAWithErrors(
i);
466template <
int N,
typename... Args>
467GPUd() bool DCAFitterN<N, Args...>::calcPCACoefs()
470 if (!calcInverseWeight()) {
473 for (
int i = N;
i--;) {
474 const auto& taux = mTrAux[
i];
475 const auto& tcov = mTrcEInv[mCurHyp][
i];
477 miei[0][0] = taux.c * tcov.sxx;
478 miei[0][1] = -taux.s * tcov.syy;
479 miei[0][2] = -taux.s * tcov.syz;
480 miei[1][0] = taux.s * tcov.sxx;
481 miei[1][1] = taux.c * tcov.syy;
482 miei[1][2] = taux.c * tcov.syz;
484 miei[2][1] = tcov.syz;
485 miei[2][2] = tcov.szz;
486 mTrCFVT[mCurHyp][
i] = mWeightInv * miei;
492template <
int N,
typename... Args>
493GPUd() bool DCAFitterN<N, Args...>::calcInverseWeight()
496 auto* arrmat = mWeightInv.Array();
497 memset(arrmat, 0,
sizeof(mWeightInv));
504 for (
int i = N;
i--;) {
505 const auto& taux = mTrAux[
i];
506 const auto& tcov = mTrcEInv[mCurHyp][
i];
507 arrmat[XX] += taux.cc * tcov.sxx + taux.ss * tcov.syy;
508 arrmat[XY] += taux.cs * (tcov.sxx - tcov.syy);
509 arrmat[XZ] += -taux.s * tcov.syz;
510 arrmat[YY] += taux.cc * tcov.syy + taux.ss * tcov.sxx;
511 arrmat[YZ] += taux.c * tcov.syz;
512 arrmat[ZZ] += tcov.szz;
515 return mWeightInv.Invert();
519template <
int N,
typename... Args>
520GPUd()
void DCAFitterN<N, Args...>::calcResidDerivatives()
524 for (
int i = N;
i--;) {
525 const auto& taux = mTrAux[
i];
526 for (
int j = N;
j--;) {
527 const auto& matT = mTrCFVT[mCurHyp][
j];
528 const auto& trDx = mTrDer[mCurHyp][
j];
529 auto& dr1 = mDResidDx[
i][
j];
530 auto& dr2 = mD2ResidDx2[
i][
j];
532 matMT[0][0] = taux.c * matT[0][0] + taux.s * matT[1][0];
533 matMT[0][1] = taux.c * matT[0][1] + taux.s * matT[1][1];
534 matMT[0][2] = taux.c * matT[0][2] + taux.s * matT[1][2];
535 matMT[1][0] = -taux.s * matT[0][0] + taux.c * matT[1][0];
536 matMT[1][1] = -taux.s * matT[0][1] + taux.c * matT[1][1];
537 matMT[1][2] = -taux.s * matT[0][2] + taux.c * matT[1][2];
538 matMT[2][0] = matT[2][0];
539 matMT[2][1] = matT[2][1];
540 matMT[2][2] = matT[2][2];
543 dr1[0] = -(matMT[0][0] + matMT[0][1] * trDx.dydx + matMT[0][2] * trDx.dzdx);
544 dr1[1] = -(matMT[1][0] + matMT[1][1] * trDx.dydx + matMT[1][2] * trDx.dzdx);
545 dr1[2] = -(matMT[2][0] + matMT[2][1] * trDx.dydx + matMT[2][2] * trDx.dzdx);
548 dr2[0] = -(matMT[0][1] * trDx.d2ydx2 + matMT[0][2] * trDx.d2zdx2);
549 dr2[1] = -(matMT[1][1] * trDx.d2ydx2 + matMT[1][2] * trDx.d2zdx2);
550 dr2[2] = -(matMT[2][1] * trDx.d2ydx2 + matMT[2][2] * trDx.d2zdx2);
557 dr2[1] += trDx.d2ydx2;
558 dr2[2] += trDx.d2zdx2;
565template <
int N,
typename... Args>
566GPUd()
void DCAFitterN<N, Args...>::calcResidDerivativesNoErr()
569 constexpr double NInv1 = 1. - NInv;
570 for (
int i = N;
i--;) {
571 const auto& trDxi = mTrDer[mCurHyp][
i];
572 auto& dr1ii = mDResidDx[
i][
i];
573 auto& dr2ii = mD2ResidDx2[
i][
i];
575 dr1ii[1] = NInv1 * trDxi.dydx;
576 dr1ii[2] = NInv1 * trDxi.dzdx;
579 dr2ii[1] = NInv1 * trDxi.d2ydx2;
580 dr2ii[2] = NInv1 * trDxi.d2zdx2;
582 for (
int j =
i;
j--;) {
583 auto& dr1ij = mDResidDx[
i][
j];
584 auto& dr1ji = mDResidDx[
j][
i];
585 const auto& trDxj = mTrDer[mCurHyp][
j];
586 auto cij = mCosDif[
i][
j], sij = mSinDif[
i][
j];
589 dr1ij[0] = -(cij + sij * trDxj.dydx);
590 dr1ij[1] = -(-sij + cij * trDxj.dydx);
591 dr1ij[2] = -trDxj.dzdx * NInv;
594 dr1ji[0] = -(cij - sij * trDxi.dydx);
595 dr1ji[1] = -(sij + cij * trDxi.dydx);
596 dr1ji[2] = -trDxi.dzdx * NInv;
598 auto& dr2ij = mD2ResidDx2[
i][
j];
599 auto& dr2ji = mD2ResidDx2[
j][
i];
601 dr2ij[0] = -sij * trDxj.d2ydx2;
602 dr2ij[1] = -cij * trDxj.d2ydx2;
603 dr2ij[2] = -trDxj.d2zdx2 * NInv;
606 dr2ji[0] = sij * trDxi.d2ydx2;
607 dr2ji[1] = -cij * trDxi.d2ydx2;
608 dr2ji[2] = -trDxi.d2zdx2 * NInv;
615template <
int N,
typename... Args>
616GPUd()
void DCAFitterN<N, Args...>::calcRMatrices()
619 for (
int i = N;
i--;) {
620 const auto& mi = mTrAux[
i];
621 for (
int j =
i;
j--;) {
622 const auto& mj = mTrAux[
j];
623 mCosDif[
i][
j] = (mi.c * mj.c + mi.s * mj.s) * NInv;
624 mSinDif[
i][
j] = (mi.s * mj.c - mi.c * mj.s) * NInv;
630template <
int N,
typename... Args>
631GPUd()
void DCAFitterN<N, Args...>::calcChi2Derivatives()
637 for (
int i = N;
i--;) {
638 auto& dchi1 = mDChi2Dx[
i];
640 for (
int j = N;
j--;) {
641 const auto&
res = mTrRes[mCurHyp][
j];
642 const auto& covI = mTrcEInv[mCurHyp][
j];
643 const auto& dr1 = mDResidDx[
j][
i];
644 auto& cidr = covIDrDx[
i][
j];
645 cidr[0] = covI.sxx * dr1[0];
646 cidr[1] = covI.syy * dr1[1] + covI.syz * dr1[2];
647 cidr[2] = covI.syz * dr1[1] + covI.szz * dr1[2];
653 for (
int i = N;
i--;) {
654 for (
int j =
i + 1;
j--;) {
655 auto& dchi2 = mD2Chi2Dx2[
i][
j];
657 for (
int k = N; k--;) {
658 const auto& dr1j = mDResidDx[k][
j];
659 const auto& cidrkj = covIDrDx[
i][k];
662 const auto&
res = mTrRes[mCurHyp][k];
663 const auto& covI = mTrcEInv[mCurHyp][k];
664 const auto& dr2ij = mD2ResidDx2[k][
j];
665 dchi2 +=
res[0] * covI.sxx * dr2ij[0] +
res[1] * (covI.syy * dr2ij[1] + covI.syz * dr2ij[2]) +
res[2] * (covI.syz * dr2ij[1] + covI.szz * dr2ij[2]);
673template <
int N,
typename... Args>
674GPUd()
void DCAFitterN<N, Args...>::calcChi2DerivativesNoErr()
677 for (
int i = N;
i--;) {
678 auto& dchi1 = mDChi2Dx[
i];
680 for (
int j = N;
j--;) {
681 const auto&
res = mTrRes[mCurHyp][
j];
682 const auto& dr1 = mDResidDx[
j][
i];
686 auto& dchi2 = mD2Chi2Dx2[
i][
j];
688 for (
int k = N; k--;) {
697template <
int N,
typename... Args>
698GPUd()
void DCAFitterN<N, Args...>::calcPCA()
701 mPCA[mCurHyp] = mTrCFVT[mCurHyp][N - 1] * mTrPos[mCurHyp][N - 1];
702 for (
int i = N - 1;
i--;) {
703 mPCA[mCurHyp] += mTrCFVT[mCurHyp][
i] * mTrPos[mCurHyp][
i];
708template <
int N,
typename... Args>
709GPUd() bool DCAFitterN<N, Args...>::recalculatePCAWithErrors(
int cand)
712 if (isPropagateTracksToVertexDone(cand) && !propagateTracksToVertex(cand)) {
715 int saveCurHyp = mCurHyp;
716 mCurHyp = mOrder[cand];
718 for (
int i = N;
i--;) {
719 if (!mTrcEInv[mCurHyp][
i].
set(mCandTr[mCurHyp][
i], XerrFactor)) {
720 if (mLoggerBadCov.needToLog()) {
722 printf(
"fitter %d: error (%ld muted): overrode invalid track covariance from %s\n",
723 mFitterID, mLoggerBadCov.getNMuted(), mCandTr[mCurHyp][
i].asString().c_str());
725 printf(
"fitter %d: error (%ld muted): overrode invalid track covariance cyy:%e czz:%e cyz:%e\n",
726 mFitterID, mLoggerBadCov.getNMuted(), mCandTr[mCurHyp][
i].getSigmaY2(), mCandTr[mCurHyp][
i].getSigmaZ2(), mCandTr[mCurHyp][
i].getSigmaZY());
728 mLoggerBadCov.evCountPrev = mLoggerBadCov.evCount;
730 if (mBadCovPolicy == Discard) {
732 }
else if (mBadCovPolicy == OverrideAndFlag) {
733 mPropFailed[mCurHyp] =
true;
737 if (!calcPCACoefs()) {
738 mCurHyp = saveCurHyp;
742 auto oldPCA = mPCA[mOrder[cand]];
744 mCurHyp = saveCurHyp;
749template <
int N,
typename... Args>
750GPUd()
void DCAFitterN<N, Args...>::calcPCANoErr()
753 auto& pca = mPCA[mCurHyp];
754 o2::math_utils::rotateZd(mTrPos[mCurHyp][N - 1][0], mTrPos[mCurHyp][N - 1][1], pca[0], pca[1], mTrAux[N - 1].s, mTrAux[N - 1].
c);
756 pca[2] = mTrPos[mCurHyp][N - 1][2];
757 for (
int i = N - 1;
i--;) {
763 pca[2] += mTrPos[mCurHyp][
i][2];
771template <
int N,
typename... Args>
772GPUd()
o2::math_utils::SMatrix<
double, 3, 3,
o2::math_utils::MatRepSym<
double, 3>>
DCAFitterN<N, Args...>::calcPCACovMatrix(
int cand)
const
777 for (
int i = N;
i--;) {
781 if (covTr.Invert()) {
786 if (nAdded && covm.Invert()) {
791 for (
int i = N;
i--;) {
798template <
int N,
typename... Args>
799GPUd()
void DCAFitterN<N, Args...>::calcTrackResiduals()
803 for (
int i = N;
i--;) {
804 mTrRes[mCurHyp][
i] = mTrPos[mCurHyp][
i];
805 vtxLoc = mPCA[mCurHyp];
807 mTrRes[mCurHyp][
i] -= vtxLoc;
812template <
int N,
typename... Args>
813GPUdi()
void DCAFitterN<N, Args...>::calcTrackDerivatives()
816 for (
int i = N;
i--;) {
817 mTrDer[mCurHyp][
i].set(mCandTr[mCurHyp][
i], mBz);
822template <
int N,
typename... Args>
823GPUdi() double DCAFitterN<N, Args...>::calcChi2()
const
827 for (
int i = N;
i--;) {
828 const auto&
res = mTrRes[mCurHyp][
i];
829 const auto& covI = mTrcEInv[mCurHyp][
i];
830 chi2 +=
res[0] *
res[0] * covI.sxx +
res[1] *
res[1] * covI.syy +
res[2] *
res[2] * covI.szz + 2. *
res[1] *
res[2] * covI.syz;
836template <
int N,
typename... Args>
837GPUdi() double DCAFitterN<N, Args...>::calcChi2NoErr()
const
841 for (
int i = N;
i--;) {
842 const auto&
res = mTrRes[mCurHyp][
i];
849template <
int N,
typename... Args>
850GPUd() bool DCAFitterN<N, Args...>::correctTracks(const VecND& corrX)
853 for (
int i = N;
i--;) {
854 const auto& trDer = mTrDer[mCurHyp][
i];
855 auto dx2h = 0.5 * corrX[
i] * corrX[
i];
856 mTrPos[mCurHyp][
i][0] -= corrX[
i];
857 mTrPos[mCurHyp][
i][1] -= trDer.dydx * corrX[
i] - dx2h * trDer.d2ydx2;
858 mTrPos[mCurHyp][
i][2] -= trDer.dzdx * corrX[
i] - dx2h * trDer.d2zdx2;
864template <
int N,
typename... Args>
865GPUd() bool DCAFitterN<N, Args...>::propagateTracksToVertex(
int icand)
868 int ord = mOrder[icand];
869 if (mTrPropDone[ord]) {
874 if (mRefitWithMatCorr) {
875 int curHypSav = mCurHyp, curCrosIDAlt = mCrossIDAlt;
878 auto restore = [
this, curHypSav, curCrosIDAlt]() { this->mCurHyp = curHypSav; this->mCrossIDAlt = curCrosIDAlt; };
879 if (!(mUseAbsDCA ? minimizeChi2NoErr() : minimizeChi2())) {
886 for (
int i = N;
i--;) {
887 if (mUseAbsDCA || mUsePropagator || mMatCorr != o2::base::Propagator::MatCorrType::USEMatCorrNONE) {
888 mCandTr[ord][
i] = *mOrigTrPtr[
i];
890 auto x = mTrAux[
i].c * mPCA[ord][0] + mTrAux[
i].s * mPCA[ord][1];
891 if (!propagateToX(mCandTr[ord][
i],
x)) {
896 mTrPropDone[ord] =
true;
901template <
int N,
typename... Args>
905 int ord = mOrder[icand];
907 if (!mTrPropDone[ord]) {
908 auto x = mTrAux[
i].c * mPCA[ord][0] + mTrAux[
i].s * mPCA[ord][1];
909 if (!propagateParamToX(trc,
x)) {
917template <
int N,
typename... Args>
918GPUdi() double DCAFitterN<N, Args...>::getAbsMax(const VecND&
v)
921 for (
int i = N;
i--;) {
922 auto vai = o2::gpu::GPUCommonMath::Abs(
v[
i]);
931template <
int N,
typename... Args>
932GPUd() bool DCAFitterN<N, Args...>::minimizeChi2()
935 for (
int i = N;
i--;) {
936 mCandTr[mCurHyp][
i] = *mOrigTrPtr[
i];
937 auto x = mTrAux[
i].c * mPCA[mCurHyp][0] + mTrAux[
i].s * mPCA[mCurHyp][1];
938 if (
x < mMinXSeed || !propagateToX(mCandTr[mCurHyp][
i],
x)) {
941 setTrackPos(mTrPos[mCurHyp][
i], mCandTr[mCurHyp][
i]);
942 if (!mTrcEInv[mCurHyp][
i].
set(mCandTr[mCurHyp][
i], XerrFactor)) {
943 if (mLoggerBadCov.needToLog()) {
945 printf(
"fitter %d: error (%ld muted): overrode invalid track covariance from %s\n",
946 mFitterID, mLoggerBadCov.getNMuted(), mCandTr[mCurHyp][
i].asString().c_str());
948 printf(
"fitter %d: error (%ld muted): overrode invalid track covariance cyy:%e czz:%e cyz:%e\n",
949 mFitterID, mLoggerBadCov.getNMuted(), mCandTr[mCurHyp][
i].getSigmaY2(), mCandTr[mCurHyp][
i].getSigmaZ2(), mCandTr[mCurHyp][
i].getSigmaZY());
951 mLoggerBadCov.evCountPrev = mLoggerBadCov.evCount;
953 if (mBadCovPolicy == Discard) {
955 }
else if (mBadCovPolicy == OverrideAndFlag) {
956 mPropFailed[mCurHyp] =
true;
961 if (mMaxDZIni > 0 && !roughDZCut()) {
965 if (!calcPCACoefs()) {
969 calcTrackResiduals();
970 float chi2Upd, chi2 = calcChi2();
972 calcTrackDerivatives();
973 calcResidDerivatives();
974 calcChi2Derivatives();
977 if (!mD2Chi2Dx2.Invert()) {
978 if (mLoggerBadInv.needToLog()) {
979 printf(
"fitter %d: error (%ld muted): Inversion failed\n", mFitterID, mLoggerBadCov.getNMuted());
980 mLoggerBadInv.evCountPrev = mLoggerBadInv.evCount;
984 VecND dx = mD2Chi2Dx2 * mDChi2Dx;
985 if (!correctTracks(dx)) {
989 if (mCrossIDAlt >= 0 && closerToAlternative()) {
990 mAllowAltPreference =
false;
993 calcTrackResiduals();
994 chi2Upd = calcChi2();
995 if (getAbsMax(dx) < mMinParamChange || chi2Upd > chi2 * mMinRelChi2Change) {
1000 }
while (++mNIters[mCurHyp] < mMaxIter);
1002 mChi2[mCurHyp] = chi2 * NInv;
1003 return mChi2[mCurHyp] < mMaxChi2;
1007template <
int N,
typename... Args>
1008GPUd() bool DCAFitterN<N, Args...>::minimizeChi2NoErr()
1012 for (
int i = N;
i--;) {
1013 mCandTr[mCurHyp][
i] = *mOrigTrPtr[
i];
1014 auto x = mTrAux[
i].c * mPCA[mCurHyp][0] + mTrAux[
i].s * mPCA[mCurHyp][1];
1015 if (
x < mMinXSeed || !propagateParamToX(mCandTr[mCurHyp][
i],
x)) {
1018 setTrackPos(mTrPos[mCurHyp][
i], mCandTr[mCurHyp][
i]);
1020 if (mMaxDZIni > 0 && !roughDZCut()) {
1025 calcTrackResiduals();
1026 float chi2Upd, chi2 = calcChi2NoErr();
1028 calcTrackDerivatives();
1029 calcResidDerivativesNoErr();
1030 calcChi2DerivativesNoErr();
1033 if (!mD2Chi2Dx2.Invert()) {
1034 if (mLoggerBadInv.needToLog()) {
1035 printf(
"itter %d: error (%ld muted): Inversion failed\n", mFitterID, mLoggerBadCov.getNMuted());
1036 mLoggerBadInv.evCountPrev = mLoggerBadInv.evCount;
1040 VecND dx = mD2Chi2Dx2 * mDChi2Dx;
1041 if (!correctTracks(dx)) {
1045 if (mCrossIDAlt >= 0 && closerToAlternative()) {
1046 mAllowAltPreference =
false;
1049 calcTrackResiduals();
1050 chi2Upd = calcChi2NoErr();
1051 if (getAbsMax(dx) < mMinParamChange || chi2Upd > chi2 * mMinRelChi2Change) {
1056 }
while (++mNIters[mCurHyp] < mMaxIter);
1058 mChi2[mCurHyp] = chi2 * NInv;
1059 return mChi2[mCurHyp] < mMaxChi2;
1063template <
int N,
typename... Args>
1064GPUd() bool DCAFitterN<N, Args...>::roughDZCut()
const
1068 for (
int i = N; accept &&
i--;) {
1069 for (
int j =
i;
j--;) {
1070 if (o2::gpu::GPUCommonMath::Abs(mCandTr[mCurHyp][
i].getZ() - mCandTr[mCurHyp][
j].getZ()) > mMaxDZIni) {
1080template <
int N,
typename... Args>
1081GPUd() bool DCAFitterN<N, Args...>::closerToAlternative()
const
1084 auto dxCur = mPCA[mCurHyp][0] - mCrossings.xDCA[mCrossIDCur], dyCur = mPCA[mCurHyp][1] - mCrossings.yDCA[mCrossIDCur];
1085 auto dxAlt = mPCA[mCurHyp][0] - mCrossings.xDCA[mCrossIDAlt], dyAlt = mPCA[mCurHyp][1] - mCrossings.yDCA[mCrossIDAlt];
1086 return dxCur * dxCur + dyCur * dyCur > dxAlt * dxAlt + dyAlt * dyAlt;
1090template <
int N,
typename... Args>
1093#ifndef GPUCA_GPUCODE_DEVICE
1094 LOG(info) << N <<
"-prong vertex fitter in " << (mUseAbsDCA ?
"abs." :
"weighted") <<
" distance minimization mode, collinear tracks mode: " << (mIsCollinear ?
"ON" :
"OFF");
1095 LOG(info) <<
"Bz: " << mBz <<
" MaxIter: " << mMaxIter <<
" MaxChi2: " << mMaxChi2 <<
" MatCorrType: " <<
int(mMatCorr);
1096 LOG(info) <<
"Stopping condition: Max.param change < " << mMinParamChange <<
" Rel.Chi2 change > " << mMinRelChi2Change;
1097 LOG(info) <<
"Discard candidates for : Rvtx > " << getMaxR() <<
" DZ between tracks > " << mMaxDZIni;
1098 LOG(info) <<
"PropagateToPCA:" << mPropagateToPCA <<
" WeightedFinalPCA:" << mWeightedFinalPCA <<
" UsePropagator:" << mUsePropagator <<
" RefitWithMatCorr:" << mRefitWithMatCorr;
1100 for (
int i = 0;
i < mCrossings.nDCA;
i++) {
1101 rep += fmt::format(
"seed{}:{}/{} ",
i, mTrPropDone[
i], mPropFailed[
i]);
1103 LOG(info) <<
"Last call: NCand:" << mCurHyp <<
" from " << mCrossings.nDCA <<
" seeds, prop.done/failed: " << rep;
1106 printf(
"%d-prong vertex fitter in abs. distance minimization mode\n", N);
1108 printf(
"%d-prong vertex fitter in weighted distance minimization mode\n", N);
1110 printf(
"Bz: %1.f MaxIter: %3.d MaxChi2: %2.3f\n", mBz, mMaxIter, mMaxChi2);
1111 printf(
"Stopping condition: Max.param change < %2.3f Rel.Chi2 change > %2.3f\n", mMinParamChange, mMinRelChi2Change);
1112 printf(
"Discard candidates for : Rvtx > %2.3f DZ between tracks > %2.3f\n", getMaxR(), mMaxDZIni);
1117template <
int N,
typename... Args>
1118GPUd()
o2::track::
TrackParCov DCAFitterN<N, Args...>::createParentTrackParCov(
int cand,
bool sectorAlpha)
const
1120 const auto& trP = getTrack(0, cand);
1121 const auto& trN = getTrack(1, cand);
1125 for (
int it = 0; it < N; it++) {
1126 const auto& trc = getTrack(it, cand);
1129 trc.getPxPyPzGlo(pvecT);
1130 trc.getCovXYZPxPyPzGlo(covT);
1131 constexpr int MomInd[6] = {9, 13, 14, 18, 19, 20};
1132 for (
int i = 0;
i < 6;
i++) {
1133 covV[MomInd[
i]] += covT[MomInd[
i]];
1135 for (
int i = 0;
i < 3;
i++) {
1136 pvecV[
i] += pvecT[
i];
1138 q += trc.getCharge();
1140 auto covVtxV = calcPCACovMatrix(cand);
1141 covV[0] = covVtxV(0, 0);
1142 covV[1] = covVtxV(1, 0);
1143 covV[2] = covVtxV(1, 1);
1144 covV[3] = covVtxV(2, 0);
1145 covV[4] = covVtxV(2, 1);
1146 covV[5] = covVtxV(2, 2);
1151template <
int N,
typename... Args>
1152GPUd()
o2::track::TrackPar DCAFitterN<N, Args...>::createParentTrackPar(
int cand,
bool sectorAlpha)
const
1154 const auto& trP = getTrack(0, cand);
1155 const auto& trN = getTrack(1, cand);
1156 const auto& wvtx = getPCACandidate(cand);
1159 for (
int it = 0; it < N; it++) {
1160 const auto& trc = getTrack(it, cand);
1162 trc.getPxPyPzGlo(pvecT);
1163 for (
int i = 0;
i < 3;
i++) {
1164 pvecV[
i] += pvecT[
i];
1166 q += trc.getCharge();
1173template <
int N,
typename... Args>
1174GPUdi() bool DCAFitterN<N, Args...>::propagateParamToX(
o2::track::TrackPar& t,
float x)
1177 if (mUsePropagator || mMatCorr != o2::base::Propagator::MatCorrType::USEMatCorrNONE) {
1178#ifndef GPUCA_GPUCODE
1182 res = t.propagateParamTo(
x, mBz);
1185 mPropFailed[mCurHyp] =
true;
1186 if (mLoggerBadProp.needToLog()) {
1187#ifndef GPUCA_GPUCODE
1188 printf(
"fitter %d: error (%ld muted): propagation failed for %s\n", mFitterID, mLoggerBadProp.getNMuted(), t.asString().c_str());
1190 printf(
"fitter %d: error (%ld muted): propagation failed\n", mFitterID, mLoggerBadProp.getNMuted());
1192 mLoggerBadProp.evCountPrev = mLoggerBadProp.evCount;
1199template <
int N,
typename... Args>
1203 if (mUsePropagator || mMatCorr != o2::base::Propagator::MatCorrType::USEMatCorrNONE) {
1204#ifndef GPUCA_GPUCODE
1208 res = t.propagateTo(
x, mBz);
1211 mPropFailed[mCurHyp] =
true;
1212 if (mLoggerBadProp.needToLog()) {
1213#ifndef GPUCA_GPUCODE
1214 printf(
"fitter %d: error (%ld muted): propagation failed for %s\n", mFitterID, mLoggerBadProp.getNMuted(), t.asString().c_str());
1216 printf(
"fitter %d: error (%ld muted): propagation failed\n", mFitterID, mLoggerBadProp.getNMuted());
1218 mLoggerBadProp.evCountPrev = mLoggerBadProp.evCount;
1229template <
typename Fitter>
1230void print(
const int nBlocks,
const int nThreads, Fitter& ft);
1232template <
typename Fitter,
class... Tr>
1233int process(
const int nBlocks,
const int nThreads, Fitter&, Tr&... args);
1235template <
class Fitter,
class... Tr>
1236void processBulk(
const int nBlocks,
const int nThreads,
const int nBatches, std::vector<Fitter>& fitters, std::vector<int>& results, std::vector<Tr>&... args);