File indexing completed on 2024-04-06 12:28:37
0001 #ifndef RecoTracker_PixelTrackFitting_interface_alpaka_RiemannFit_h
0002 #define RecoTracker_PixelTrackFitting_interface_alpaka_RiemannFit_h
0003
0004 #include <alpaka/alpaka.hpp>
0005
0006 #include <Eigen/Core>
0007 #include <Eigen/Eigenvalues>
0008
0009 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0010 #include "RecoTracker/PixelTrackFitting/interface/alpaka/FitUtils.h"
0011
0012 namespace ALPAKA_ACCELERATOR_NAMESPACE::riemannFit {
0013 using namespace ::riemannFit;
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037 template <typename TAcc, typename VNd1, typename VNd2>
0038 ALPAKA_FN_ACC ALPAKA_FN_INLINE void computeRadLenUniformMaterial(const TAcc& acc,
0039 const VNd1& length_values,
0040 VNd2& rad_lengths) {
0041
0042
0043 constexpr double xx_0_inv = 0.06 / 16.;
0044 uint n = length_values.rows();
0045 rad_lengths(0) = length_values(0) * xx_0_inv;
0046 for (uint j = 1; j < n; ++j) {
0047 rad_lengths(j) = alpaka::math::abs(acc, length_values(j) - length_values(j - 1)) * xx_0_inv;
0048 }
0049 }
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070 template <typename TAcc, typename V4, typename VNd1, typename VNd2, int N>
0071 ALPAKA_FN_ACC ALPAKA_FN_INLINE auto scatterCovLine(const TAcc& acc,
0072 Matrix2d const* cov_sz,
0073 const V4& fast_fit,
0074 VNd1 const& s_arcs,
0075 VNd2 const& z_values,
0076 const double theta,
0077 const double bField,
0078 MatrixNd<N>& ret) {
0079 #ifdef RFIT_DEBUG
0080 riemannFit::printIt(&s_arcs, "Scatter_cov_line - s_arcs: ");
0081 #endif
0082 constexpr uint n = N;
0083 double p_t = alpaka::math::min(acc, 20., fast_fit(2) * bField);
0084 double p_2 = p_t * p_t * (1. + 1. / sqr(fast_fit(3)));
0085 VectorNd<N> rad_lengths_S;
0086
0087
0088
0089 VectorNd<N> s_values = s_arcs.array() * s_arcs.array() + z_values.array() * z_values.array();
0090 s_values = s_values.array().sqrt();
0091 computeRadLenUniformMaterial(acc, s_values, rad_lengths_S);
0092 VectorNd<N> sig2_S;
0093 sig2_S = .000225 / p_2 * (1. + 0.038 * rad_lengths_S.array().log()).abs2() * rad_lengths_S.array();
0094 #ifdef RFIT_DEBUG
0095 riemannFit::printIt(cov_sz, "Scatter_cov_line - cov_sz: ");
0096 #endif
0097 Matrix2Nd<N> tmp = Matrix2Nd<N>::Zero();
0098 for (uint k = 0; k < n; ++k) {
0099 tmp(k, k) = cov_sz[k](0, 0);
0100 tmp(k + n, k + n) = cov_sz[k](1, 1);
0101 tmp(k, k + n) = tmp(k + n, k) = cov_sz[k](0, 1);
0102 }
0103 for (uint k = 0; k < n; ++k) {
0104 for (uint l = k; l < n; ++l) {
0105 for (uint i = 0; i < uint(alpaka::math::min(acc, k, l)); ++i) {
0106 tmp(k + n, l + n) += alpaka::math::abs(acc, s_values(k) - s_values(i)) *
0107 alpaka::math::abs(acc, s_values(l) - s_values(i)) * sig2_S(i);
0108 }
0109 tmp(l + n, k + n) = tmp(k + n, l + n);
0110 }
0111 }
0112
0113
0114 #ifdef RFIT_DEBUG
0115 riemannFit::printIt(&tmp, "Scatter_cov_line - tmp: ");
0116 #endif
0117 ret = tmp.block(n, n, n, n);
0118 }
0119
0120
0121
0122
0123
0124
0125
0126
0127
0128
0129
0130
0131
0132
0133 template <typename TAcc, typename M2xN, typename V4, int N>
0134 ALPAKA_FN_ACC ALPAKA_FN_INLINE MatrixNd<N> scatter_cov_rad(
0135 const TAcc& acc, const M2xN& p2D, const V4& fast_fit, VectorNd<N> const& rad, double B) {
0136 constexpr uint n = N;
0137 double p_t = alpaka::math::min(acc, 20., fast_fit(2) * B);
0138 double p_2 = p_t * p_t * (1. + 1. / sqr(fast_fit(3)));
0139 double theta = atan(fast_fit(3));
0140 theta = theta < 0. ? theta + M_PI : theta;
0141 VectorNd<N> s_values;
0142 VectorNd<N> rad_lengths;
0143 const Vector2d oVec(fast_fit(0), fast_fit(1));
0144
0145
0146 for (uint i = 0; i < n; ++i) {
0147 Vector2d pVec = p2D.block(0, i, 2, 1) - oVec;
0148 const double cross = cross2D(acc, -oVec, pVec);
0149 const double dot = (-oVec).dot(pVec);
0150 const double tempAtan2 = atan2(cross, dot);
0151 s_values(i) = alpaka::math::abs(acc, tempAtan2 * fast_fit(2));
0152 }
0153 computeRadLenUniformMaterial(acc, s_values * sqrt(1. + 1. / sqr(fast_fit(3))), rad_lengths);
0154 MatrixNd<N> scatter_cov_rad = MatrixNd<N>::Zero();
0155 VectorNd<N> sig2 = (1. + 0.038 * rad_lengths.array().log()).abs2() * rad_lengths.array();
0156 sig2 *= 0.000225 / (p_2 * sqr(sin(theta)));
0157 for (uint k = 0; k < n; ++k) {
0158 for (uint l = k; l < n; ++l) {
0159 for (uint i = 0; i < uint(alpaka::math::min(acc, k, l)); ++i) {
0160 scatter_cov_rad(k, l) += (rad(k) - rad(i)) * (rad(l) - rad(i)) * sig2(i);
0161 }
0162 scatter_cov_rad(l, k) = scatter_cov_rad(k, l);
0163 }
0164 }
0165 #ifdef RFIT_DEBUG
0166 riemannFit::printIt(&scatter_cov_rad, "Scatter_cov_rad - scatter_cov_rad: ");
0167 #endif
0168 return scatter_cov_rad;
0169 }
0170
0171
0172
0173
0174
0175
0176
0177
0178
0179 template <typename TAcc, typename M2xN, int N>
0180 ALPAKA_FN_ACC ALPAKA_FN_INLINE Matrix2Nd<N> cov_radtocart(const TAcc& acc,
0181 const M2xN& p2D,
0182 const MatrixNd<N>& cov_rad,
0183 const VectorNd<N>& rad) {
0184 #ifdef RFIT_DEBUG
0185 printf("Address of p2D: %p\n", &p2D);
0186 #endif
0187 printIt(&p2D, "cov_radtocart - p2D:");
0188 constexpr uint n = N;
0189 Matrix2Nd<N> cov_cart = Matrix2Nd<N>::Zero();
0190 VectorNd<N> rad_inv = rad.cwiseInverse();
0191 printIt(&rad_inv, "cov_radtocart - rad_inv:");
0192 for (uint i = 0; i < n; ++i) {
0193 for (uint j = i; j < n; ++j) {
0194 cov_cart(i, j) = cov_rad(i, j) * p2D(1, i) * rad_inv(i) * p2D(1, j) * rad_inv(j);
0195 cov_cart(i + n, j + n) = cov_rad(i, j) * p2D(0, i) * rad_inv(i) * p2D(0, j) * rad_inv(j);
0196 cov_cart(i, j + n) = -cov_rad(i, j) * p2D(1, i) * rad_inv(i) * p2D(0, j) * rad_inv(j);
0197 cov_cart(i + n, j) = -cov_rad(i, j) * p2D(0, i) * rad_inv(i) * p2D(1, j) * rad_inv(j);
0198 cov_cart(j, i) = cov_cart(i, j);
0199 cov_cart(j + n, i + n) = cov_cart(i + n, j + n);
0200 cov_cart(j + n, i) = cov_cart(i, j + n);
0201 cov_cart(j, i + n) = cov_cart(i + n, j);
0202 }
0203 }
0204 return cov_cart;
0205 }
0206
0207
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217 template <typename TAcc, typename M2xN, int N>
0218 ALPAKA_FN_ACC ALPAKA_FN_INLINE VectorNd<N> cov_carttorad(const TAcc& acc,
0219 const M2xN& p2D,
0220 const Matrix2Nd<N>& cov_cart,
0221 const VectorNd<N>& rad) {
0222 constexpr uint n = N;
0223 VectorNd<N> cov_rad;
0224 const VectorNd<N> rad_inv2 = rad.cwiseInverse().array().square();
0225 for (uint i = 0; i < n; ++i) {
0226
0227 if (rad(i) < 1.e-4)
0228 cov_rad(i) = cov_cart(i, i);
0229 else {
0230 cov_rad(i) = rad_inv2(i) * (cov_cart(i, i) * sqr(p2D(1, i)) + cov_cart(i + n, i + n) * sqr(p2D(0, i)) -
0231 2. * cov_cart(i, i + n) * p2D(0, i) * p2D(1, i));
0232 }
0233 }
0234 return cov_rad;
0235 }
0236
0237
0238
0239
0240
0241
0242
0243
0244
0245
0246
0247
0248
0249 template <typename TAcc, typename M2xN, typename V4, int N>
0250 ALPAKA_FN_ACC ALPAKA_FN_INLINE VectorNd<N> cov_carttorad_prefit(
0251 const TAcc& acc, const M2xN& p2D, const Matrix2Nd<N>& cov_cart, V4& fast_fit, const VectorNd<N>& rad) {
0252 constexpr uint n = N;
0253 VectorNd<N> cov_rad;
0254 for (uint i = 0; i < n; ++i) {
0255
0256 if (rad(i) < 1.e-4)
0257 cov_rad(i) = cov_cart(i, i);
0258 else {
0259 Vector2d a = p2D.col(i);
0260 Vector2d b = p2D.col(i) - fast_fit.head(2);
0261 const double x2 = a.dot(b);
0262 const double y2 = cross2D(acc, a, b);
0263 const double tan_c = -y2 / x2;
0264 const double tan_c2 = sqr(tan_c);
0265 cov_rad(i) =
0266 1. / (1. + tan_c2) * (cov_cart(i, i) + cov_cart(i + n, i + n) * tan_c2 + 2 * cov_cart(i, i + n) * tan_c);
0267 }
0268 }
0269 return cov_rad;
0270 }
0271
0272
0273
0274
0275
0276
0277
0278
0279
0280
0281
0282
0283 template <typename TAcc, int N>
0284 ALPAKA_FN_ACC ALPAKA_FN_INLINE VectorNd<N> weightCircle(const TAcc& acc, const MatrixNd<N>& cov_rad_inv) {
0285 return cov_rad_inv.colwise().sum().transpose();
0286 }
0287
0288
0289
0290
0291
0292
0293
0294
0295
0296 template <typename TAcc, typename M2xN>
0297 ALPAKA_FN_ACC ALPAKA_FN_INLINE int32_t charge(const TAcc& acc, const M2xN& p2D, const Vector3d& par_uvr) {
0298 return ((p2D(0, 1) - p2D(0, 0)) * (par_uvr.y() - p2D(1, 0)) - (p2D(1, 1) - p2D(1, 0)) * (par_uvr.x() - p2D(0, 0)) >
0299 0)
0300 ? -1
0301 : 1;
0302 }
0303
0304
0305
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315
0316
0317
0318 template <typename TAcc>
0319 ALPAKA_FN_ACC ALPAKA_FN_INLINE Vector3d min_eigen3D(const TAcc& acc, const Matrix3d& A, double& chi2) {
0320 #ifdef RFIT_DEBUG
0321 printf("min_eigen3D - enter\n");
0322 #endif
0323 Eigen::SelfAdjointEigenSolver<Matrix3d> solver(3);
0324 solver.computeDirect(A);
0325 int min_index;
0326 chi2 = solver.eigenvalues().minCoeff(&min_index);
0327 #ifdef RFIT_DEBUG
0328 printf("min_eigen3D - exit\n");
0329 #endif
0330 return solver.eigenvectors().col(min_index);
0331 }
0332
0333
0334
0335
0336
0337
0338
0339
0340
0341
0342
0343
0344 template <typename TAcc>
0345 ALPAKA_FN_ACC ALPAKA_FN_INLINE Vector3d min_eigen3D_fast(const TAcc& acc, const Matrix3d& A) {
0346 Eigen::SelfAdjointEigenSolver<Matrix3f> solver(3);
0347 solver.computeDirect(A.cast<float>());
0348 int min_index;
0349 solver.eigenvalues().minCoeff(&min_index);
0350 return solver.eigenvectors().col(min_index).cast<double>();
0351 }
0352
0353
0354
0355
0356
0357
0358
0359
0360
0361
0362 template <typename TAcc>
0363 ALPAKA_FN_ACC ALPAKA_FN_INLINE Vector2d min_eigen2D(const TAcc& acc, const Matrix2d& aMat, double& chi2) {
0364 Eigen::SelfAdjointEigenSolver<Matrix2d> solver(2);
0365 solver.computeDirect(aMat);
0366 int min_index;
0367 chi2 = solver.eigenvalues().minCoeff(&min_index);
0368 return solver.eigenvectors().col(min_index);
0369 }
0370
0371
0372
0373
0374
0375
0376
0377
0378
0379
0380
0381
0382
0383
0384 template <typename TAcc, typename M3xN, typename V4>
0385 ALPAKA_FN_ACC ALPAKA_FN_INLINE void fastFit(const TAcc& acc, const M3xN& hits, V4& result) {
0386 constexpr uint32_t N = M3xN::ColsAtCompileTime;
0387 constexpr auto n = N;
0388 printIt(&hits, "Fast_fit - hits: ");
0389
0390
0391
0392 const Vector2d bVec = hits.block(0, n / 2, 2, 1) - hits.block(0, 0, 2, 1);
0393 const Vector2d cVec = hits.block(0, n - 1, 2, 1) - hits.block(0, 0, 2, 1);
0394 printIt(&bVec, "Fast_fit - b: ");
0395 printIt(&cVec, "Fast_fit - c: ");
0396
0397 auto b2 = bVec.squaredNorm();
0398 auto c2 = cVec.squaredNorm();
0399
0400
0401
0402
0403
0404
0405 bool flip = abs(bVec.x()) < abs(bVec.y());
0406 auto bx = flip ? bVec.y() : bVec.x();
0407 auto by = flip ? bVec.x() : bVec.y();
0408 auto cx = flip ? cVec.y() : cVec.x();
0409 auto cy = flip ? cVec.x() : cVec.y();
0410
0411 auto div = 2. * (cx * by - bx * cy);
0412
0413 auto y0 = (cx * b2 - bx * c2) / div;
0414 auto x0 = (0.5 * b2 - y0 * by) / bx;
0415 result(0) = hits(0, 0) + (flip ? y0 : x0);
0416 result(1) = hits(1, 0) + (flip ? x0 : y0);
0417 result(2) = sqrt(sqr(x0) + sqr(y0));
0418 printIt(&result, "Fast_fit - result: ");
0419
0420
0421 const Vector2d dVec = hits.block(0, 0, 2, 1) - result.head(2);
0422 const Vector2d eVec = hits.block(0, n - 1, 2, 1) - result.head(2);
0423 printIt(&eVec, "Fast_fit - e: ");
0424 printIt(&dVec, "Fast_fit - d: ");
0425
0426 auto dr = result(2) * atan2(cross2D(acc, dVec, eVec), dVec.dot(eVec));
0427
0428 auto dz = hits(2, n - 1) - hits(2, 0);
0429
0430 result(3) = (dr / dz);
0431
0432 #ifdef RFIT_DEBUG
0433 printf("Fast_fit: [%f, %f, %f, %f]\n", result(0), result(1), result(2), result(3));
0434 #endif
0435 }
0436
0437
0438
0439
0440
0441
0442
0443
0444
0445
0446
0447
0448
0449
0450
0451
0452
0453
0454
0455
0456
0457
0458
0459
0460
0461
0462
0463
0464 template <typename TAcc, typename M2xN, typename V4, int N>
0465 ALPAKA_FN_ACC ALPAKA_FN_INLINE CircleFit circleFit(const TAcc& acc,
0466 const M2xN& hits2D,
0467 const Matrix2Nd<N>& hits_cov2D,
0468 const V4& fast_fit,
0469 const VectorNd<N>& rad,
0470 const double bField,
0471 const bool error) {
0472 #ifdef RFIT_DEBUG
0473 printf("circle_fit - enter\n");
0474 #endif
0475
0476 Matrix2Nd<N> vMat = hits_cov2D;
0477 constexpr uint n = N;
0478 printIt(&hits2D, "circle_fit - hits2D:");
0479 printIt(&hits_cov2D, "circle_fit - hits_cov2D:");
0480
0481 #ifdef RFIT_DEBUG
0482 printf("circle_fit - WEIGHT COMPUTATION\n");
0483 #endif
0484
0485 VectorNd<N> weight;
0486 MatrixNd<N> gMat;
0487 double renorm;
0488 {
0489 MatrixNd<N> cov_rad = cov_carttorad_prefit(acc, hits2D, vMat, fast_fit, rad).asDiagonal();
0490 MatrixNd<N> scatterCovRadMat = scatter_cov_rad(acc, hits2D, fast_fit, rad, bField);
0491 printIt(&scatterCovRadMat, "circle_fit - scatter_cov_rad:");
0492 printIt(&hits2D, "circle_fit - hits2D bis:");
0493 #ifdef RFIT_DEBUG
0494 printf("Address of hits2D: a) %p\n", &hits2D);
0495 #endif
0496 vMat += cov_radtocart(acc, hits2D, scatterCovRadMat, rad);
0497 printIt(&vMat, "circle_fit - V:");
0498 cov_rad += scatterCovRadMat;
0499 printIt(&cov_rad, "circle_fit - cov_rad:");
0500 math::cholesky::invert(cov_rad, gMat);
0501
0502 renorm = gMat.sum();
0503 gMat *= 1. / renorm;
0504 weight = weightCircle(acc, gMat);
0505 }
0506 printIt(&weight, "circle_fit - weight:");
0507
0508
0509 #ifdef RFIT_DEBUG
0510 printf("circle_fit - SPACE TRANSFORMATION\n");
0511 #endif
0512
0513
0514 #ifdef RFIT_DEBUG
0515 printf("Address of hits2D: b) %p\n", &hits2D);
0516 #endif
0517 const Vector2d hCentroid = hits2D.rowwise().mean();
0518 printIt(&hCentroid, "circle_fit - h_:");
0519 Matrix3xNd<N> p3D;
0520 p3D.block(0, 0, 2, n) = hits2D.colwise() - hCentroid;
0521 printIt(&p3D, "circle_fit - p3D: a)");
0522 Vector2Nd<N> mc;
0523 mc << p3D.row(0).transpose(), p3D.row(1).transpose();
0524 printIt(&mc, "circle_fit - mc(centered hits):");
0525
0526
0527 const double tempQ = mc.squaredNorm();
0528 const double tempS = sqrt(n * 1. / tempQ);
0529 p3D.block(0, 0, 2, n) *= tempS;
0530
0531
0532 p3D.row(2) = p3D.block(0, 0, 2, n).colwise().squaredNorm();
0533 printIt(&p3D, "circle_fit - p3D: b)");
0534
0535 #ifdef RFIT_DEBUG
0536 printf("circle_fit - COST FUNCTION\n");
0537 #endif
0538
0539
0540
0541 Vector3d r0;
0542 r0.noalias() = p3D * weight;
0543 const Matrix3xNd<N> xMat = p3D.colwise() - r0;
0544 Matrix3d aMat = xMat * gMat * xMat.transpose();
0545 printIt(&aMat, "circle_fit - A:");
0546
0547 #ifdef RFIT_DEBUG
0548 printf("circle_fit - MINIMIZE\n");
0549 #endif
0550
0551 double chi2;
0552 Vector3d vVec = min_eigen3D(acc, aMat, chi2);
0553 #ifdef RFIT_DEBUG
0554 printf("circle_fit - AFTER MIN_EIGEN\n");
0555 #endif
0556 printIt(&vVec, "v BEFORE INVERSION");
0557 vVec *= (vVec(2) > 0) ? 1 : -1;
0558 printIt(&vVec, "v AFTER INVERSION");
0559
0560
0561 #ifdef RFIT_DEBUG
0562 printf("circle_fit - AFTER MIN_EIGEN 1\n");
0563 #endif
0564 Eigen::Matrix<double, 1, 1> cm;
0565 #ifdef RFIT_DEBUG
0566 printf("circle_fit - AFTER MIN_EIGEN 2\n");
0567 #endif
0568 cm = -vVec.transpose() * r0;
0569 #ifdef RFIT_DEBUG
0570 printf("circle_fit - AFTER MIN_EIGEN 3\n");
0571 #endif
0572 const double tempC = cm(0, 0);
0573
0574 #ifdef RFIT_DEBUG
0575 printf("circle_fit - COMPUTE CIRCLE PARAMETER\n");
0576 #endif
0577
0578
0579
0580 const double tempH = sqrt(1. - sqr(vVec(2)) - 4. * tempC * vVec(2));
0581 const double v2x2_inv = 1. / (2. * vVec(2));
0582 const double s_inv = 1. / tempS;
0583 Vector3d par_uvr;
0584 par_uvr << -vVec(0) * v2x2_inv, -vVec(1) * v2x2_inv, tempH * v2x2_inv;
0585
0586 CircleFit circle;
0587 circle.par << par_uvr(0) * s_inv + hCentroid(0), par_uvr(1) * s_inv + hCentroid(1), par_uvr(2) * s_inv;
0588 circle.qCharge = charge(acc, hits2D, circle.par);
0589 circle.chi2 = abs(chi2) * renorm / sqr(2 * vVec(2) * par_uvr(2) * tempS);
0590 printIt(&circle.par, "circle_fit - CIRCLE PARAMETERS:");
0591 printIt(&circle.cov, "circle_fit - CIRCLE COVARIANCE:");
0592 #ifdef RFIT_DEBUG
0593 printf("circle_fit - CIRCLE CHARGE: %d\n", circle.qCharge);
0594 #endif
0595
0596 #ifdef RFIT_DEBUG
0597 printf("circle_fit - ERROR PROPAGATION\n");
0598 #endif
0599
0600 if (error) {
0601 #ifdef RFIT_DEBUG
0602 printf("circle_fit - ERROR PRPAGATION ACTIVATED\n");
0603 #endif
0604 ArrayNd<N> vcsMat[2][2];
0605 MatrixNd<N> cMat[3][3];
0606 #ifdef RFIT_DEBUG
0607 printf("circle_fit - ERROR PRPAGATION ACTIVATED 2\n");
0608 #endif
0609 {
0610 Eigen::Matrix<double, 1, 1> cm;
0611 Eigen::Matrix<double, 1, 1> cm2;
0612 cm = mc.transpose() * vMat * mc;
0613 const double tempC2 = cm(0, 0);
0614 Matrix2Nd<N> tempVcsMat;
0615 tempVcsMat.template triangularView<Eigen::Upper>() =
0616 (sqr(tempS) * vMat + sqr(sqr(tempS)) * 1. / (4. * tempQ * n) *
0617 (2. * vMat.squaredNorm() + 4. * tempC2) *
0618 (mc * mc.transpose()));
0619
0620 printIt(&tempVcsMat, "circle_fit - Vcs:");
0621 cMat[0][0] = tempVcsMat.block(0, 0, n, n).template selfadjointView<Eigen::Upper>();
0622 vcsMat[0][1] = tempVcsMat.block(0, n, n, n);
0623 cMat[1][1] = tempVcsMat.block(n, n, n, n).template selfadjointView<Eigen::Upper>();
0624 vcsMat[1][0] = vcsMat[0][1].transpose();
0625 printIt(&tempVcsMat, "circle_fit - Vcs:");
0626 }
0627
0628 {
0629 const ArrayNd<N> t0 = (VectorXd::Constant(n, 1.) * p3D.row(0));
0630 const ArrayNd<N> t1 = (VectorXd::Constant(n, 1.) * p3D.row(1));
0631 const ArrayNd<N> t00 = p3D.row(0).transpose() * p3D.row(0);
0632 const ArrayNd<N> t01 = p3D.row(0).transpose() * p3D.row(1);
0633 const ArrayNd<N> t11 = p3D.row(1).transpose() * p3D.row(1);
0634 const ArrayNd<N> t10 = t01.transpose();
0635 vcsMat[0][0] = cMat[0][0];
0636 cMat[0][1] = vcsMat[0][1];
0637 cMat[0][2] = 2. * (vcsMat[0][0] * t0 + vcsMat[0][1] * t1);
0638 vcsMat[1][1] = cMat[1][1];
0639 cMat[1][2] = 2. * (vcsMat[1][0] * t0 + vcsMat[1][1] * t1);
0640 MatrixNd<N> tmp;
0641 tmp.template triangularView<Eigen::Upper>() =
0642 (2. * (vcsMat[0][0] * vcsMat[0][0] + vcsMat[0][0] * vcsMat[0][1] + vcsMat[1][1] * vcsMat[1][0] +
0643 vcsMat[1][1] * vcsMat[1][1]) +
0644 4. * (vcsMat[0][0] * t00 + vcsMat[0][1] * t01 + vcsMat[1][0] * t10 + vcsMat[1][1] * t11))
0645 .matrix();
0646 cMat[2][2] = tmp.template selfadjointView<Eigen::Upper>();
0647 }
0648 printIt(&cMat[0][0], "circle_fit - C[0][0]:");
0649
0650 Matrix3d c0Mat;
0651 for (uint i = 0; i < 3; ++i) {
0652 for (uint j = i; j < 3; ++j) {
0653 Eigen::Matrix<double, 1, 1> tmp;
0654 tmp = weight.transpose() * cMat[i][j] * weight;
0655
0656 const double tempC = tmp(0, 0);
0657 c0Mat(i, j) = tempC;
0658 c0Mat(j, i) = c0Mat(i, j);
0659 }
0660 }
0661 printIt(&c0Mat, "circle_fit - C0:");
0662
0663 const MatrixNd<N> wMat = weight * weight.transpose();
0664 const MatrixNd<N> hMat = MatrixNd<N>::Identity().rowwise() - weight.transpose();
0665 const MatrixNx3d<N> s_v = hMat * p3D.transpose();
0666 printIt(&wMat, "circle_fit - W:");
0667 printIt(&hMat, "circle_fit - H:");
0668 printIt(&s_v, "circle_fit - s_v:");
0669
0670 MatrixNd<N> dMat[3][3];
0671 dMat[0][0] = (hMat * cMat[0][0] * hMat.transpose()).cwiseProduct(wMat);
0672 dMat[0][1] = (hMat * cMat[0][1] * hMat.transpose()).cwiseProduct(wMat);
0673 dMat[0][2] = (hMat * cMat[0][2] * hMat.transpose()).cwiseProduct(wMat);
0674 dMat[1][1] = (hMat * cMat[1][1] * hMat.transpose()).cwiseProduct(wMat);
0675 dMat[1][2] = (hMat * cMat[1][2] * hMat.transpose()).cwiseProduct(wMat);
0676 dMat[2][2] = (hMat * cMat[2][2] * hMat.transpose()).cwiseProduct(wMat);
0677 dMat[1][0] = dMat[0][1].transpose();
0678 dMat[2][0] = dMat[0][2].transpose();
0679 dMat[2][1] = dMat[1][2].transpose();
0680 printIt(&dMat[0][0], "circle_fit - D_[0][0]:");
0681
0682 constexpr uint nu[6][2] = {{0, 0}, {0, 1}, {0, 2}, {1, 1}, {1, 2}, {2, 2}};
0683
0684 Matrix6d eMat;
0685 for (uint a = 0; a < 6; ++a) {
0686 const uint i = nu[a][0], j = nu[a][1];
0687 for (uint b = a; b < 6; ++b) {
0688 const uint k = nu[b][0], l = nu[b][1];
0689 VectorNd<N> t0(n);
0690 VectorNd<N> t1(n);
0691 if (l == k) {
0692 t0 = 2. * dMat[j][l] * s_v.col(l);
0693 if (i == j)
0694 t1 = t0;
0695 else
0696 t1 = 2. * dMat[i][l] * s_v.col(l);
0697 } else {
0698 t0 = dMat[j][l] * s_v.col(k) + dMat[j][k] * s_v.col(l);
0699 if (i == j)
0700 t1 = t0;
0701 else
0702 t1 = dMat[i][l] * s_v.col(k) + dMat[i][k] * s_v.col(l);
0703 }
0704
0705 if (i == j) {
0706 Eigen::Matrix<double, 1, 1> cm;
0707 cm = s_v.col(i).transpose() * (t0 + t1);
0708
0709 const double tempC = cm(0, 0);
0710 eMat(a, b) = 0. + tempC;
0711 } else {
0712 Eigen::Matrix<double, 1, 1> cm;
0713 cm = (s_v.col(i).transpose() * t0) + (s_v.col(j).transpose() * t1);
0714
0715 const double tempC = cm(0, 0);
0716 eMat(a, b) = 0. + tempC;
0717 }
0718 if (b != a)
0719 eMat(b, a) = eMat(a, b);
0720 }
0721 }
0722 printIt(&eMat, "circle_fit - E:");
0723
0724 Eigen::Matrix<double, 3, 6> j2Mat;
0725 for (uint a = 0; a < 6; ++a) {
0726 const uint i = nu[a][0], j = nu[a][1];
0727 Matrix3d delta = Matrix3d::Zero();
0728 delta(i, j) = delta(j, i) = abs(aMat(i, j) * epsilon);
0729 j2Mat.col(a) = min_eigen3D_fast(acc, aMat + delta);
0730 const int sign = (j2Mat.col(a)(2) > 0) ? 1 : -1;
0731 j2Mat.col(a) = (j2Mat.col(a) * sign - vVec) / delta(i, j);
0732 }
0733 printIt(&j2Mat, "circle_fit - J2:");
0734
0735 Matrix4d cvcMat;
0736 {
0737 Matrix3d t0 = j2Mat * eMat * j2Mat.transpose();
0738 Vector3d t1 = -t0 * r0;
0739 cvcMat.block(0, 0, 3, 3) = t0;
0740 cvcMat.block(0, 3, 3, 1) = t1;
0741 cvcMat.block(3, 0, 1, 3) = t1.transpose();
0742 Eigen::Matrix<double, 1, 1> cm1;
0743 Eigen::Matrix<double, 1, 1> cm3;
0744 cm1 = (vVec.transpose() * c0Mat * vVec);
0745
0746 cm3 = (r0.transpose() * t0 * r0);
0747
0748 const double tempC = cm1(0, 0) + (c0Mat.cwiseProduct(t0)).sum() + cm3(0, 0);
0749 cvcMat(3, 3) = tempC;
0750
0751 }
0752 printIt(&cvcMat, "circle_fit - Cvc:");
0753
0754 Eigen::Matrix<double, 3, 4> j3Mat;
0755 {
0756 const double t = 1. / tempH;
0757 j3Mat << -v2x2_inv, 0, vVec(0) * sqr(v2x2_inv) * 2., 0, 0, -v2x2_inv, vVec(1) * sqr(v2x2_inv) * 2., 0,
0758 vVec(0) * v2x2_inv * t, vVec(1) * v2x2_inv * t,
0759 -tempH * sqr(v2x2_inv) * 2. - (2. * tempC + vVec(2)) * v2x2_inv * t, -t;
0760 }
0761 printIt(&j3Mat, "circle_fit - J3:");
0762
0763 const RowVector2Nd<N> Jq = mc.transpose() * tempS * 1. / n;
0764 printIt(&Jq, "circle_fit - Jq:");
0765
0766 Matrix3d cov_uvr = j3Mat * cvcMat * j3Mat.transpose() * sqr(s_inv)
0767 + (par_uvr * par_uvr.transpose()) * (Jq * vMat * Jq.transpose());
0768
0769 circle.cov = cov_uvr;
0770 }
0771
0772 printIt(&circle.cov, "Circle cov:");
0773 #ifdef RFIT_DEBUG
0774 printf("circle_fit - exit\n");
0775 #endif
0776 return circle;
0777 }
0778
0779
0780
0781
0782
0783
0784
0785
0786
0787
0788
0789
0790
0791
0792
0793
0794
0795 template <typename TAcc, typename M3xN, typename M6xN, typename V4>
0796 ALPAKA_FN_ACC ALPAKA_FN_INLINE LineFit lineFit(const TAcc& acc,
0797 const M3xN& hits,
0798 const M6xN& hits_ge,
0799 const CircleFit& circle,
0800 const V4& fast_fit,
0801 const double bField,
0802 const bool error) {
0803 constexpr uint32_t N = M3xN::ColsAtCompileTime;
0804 constexpr auto n = N;
0805 double theta = -circle.qCharge * atan(fast_fit(3));
0806 theta = theta < 0. ? theta + M_PI : theta;
0807
0808
0809 Eigen::Matrix<double, 2, 2> rot;
0810 rot << sin(theta), cos(theta), -cos(theta), sin(theta);
0811
0812
0813
0814
0815
0816
0817
0818
0819
0820 Matrix2xNd<N> p2D = Matrix2xNd<N>::Zero();
0821 Eigen::Matrix<double, 2, 6> jxMat;
0822
0823 #ifdef RFIT_DEBUG
0824 printf("Line_fit - B: %g\n", bField);
0825 printIt(&hits, "Line_fit points: ");
0826 printIt(&hits_ge, "Line_fit covs: ");
0827 printIt(&rot, "Line_fit rot: ");
0828 #endif
0829
0830
0831
0832
0833
0834 const Vector2d oVec(circle.par(0), circle.par(1));
0835
0836
0837 Matrix6d covMat = Matrix6d::Zero();
0838 Matrix2d cov_sz[N];
0839 for (uint i = 0; i < n; ++i) {
0840 Vector2d pVec = hits.block(0, i, 2, 1) - oVec;
0841 const double cross = cross2D(acc, -oVec, pVec);
0842 const double dot = (-oVec).dot(pVec);
0843
0844
0845 const double tempQAtan2 = -circle.qCharge * atan2(cross, dot);
0846
0847 p2D(0, i) = tempQAtan2 * circle.par(2);
0848
0849
0850 const double temp0 = -circle.qCharge * circle.par(2) * 1. / (sqr(dot) + sqr(cross));
0851 double d_X0 = 0., d_Y0 = 0., d_R = 0.;
0852 if (error) {
0853 d_X0 = -temp0 * ((pVec(1) + oVec(1)) * dot - (pVec(0) - oVec(0)) * cross);
0854 d_Y0 = temp0 * ((pVec(0) + oVec(0)) * dot - (oVec(1) - pVec(1)) * cross);
0855 d_R = tempQAtan2;
0856 }
0857 const double d_x = temp0 * (oVec(1) * dot + oVec(0) * cross);
0858 const double d_y = temp0 * (-oVec(0) * dot + oVec(1) * cross);
0859 jxMat << d_X0, d_Y0, d_R, d_x, d_y, 0., 0., 0., 0., 0., 0., 1.;
0860
0861 covMat.block(0, 0, 3, 3) = circle.cov;
0862 covMat(3, 3) = hits_ge.col(i)[0];
0863 covMat(4, 4) = hits_ge.col(i)[2];
0864 covMat(5, 5) = hits_ge.col(i)[5];
0865 covMat(3, 4) = covMat(4, 3) = hits_ge.col(i)[1];
0866 covMat(3, 5) = covMat(5, 3) = hits_ge.col(i)[3];
0867 covMat(4, 5) = covMat(5, 4) = hits_ge.col(i)[4];
0868 Matrix2d tmp = jxMat * covMat * jxMat.transpose();
0869 cov_sz[i].noalias() = rot * tmp * rot.transpose();
0870 }
0871
0872 p2D.row(1) = hits.row(2);
0873
0874
0875
0876 MatrixNd<N> cov_with_ms;
0877 scatterCovLine(acc, cov_sz, fast_fit, p2D.row(0), p2D.row(1), theta, bField, cov_with_ms);
0878 #ifdef RFIT_DEBUG
0879 printIt(cov_sz, "line_fit - cov_sz:");
0880 printIt(&cov_with_ms, "line_fit - cov_with_ms: ");
0881 #endif
0882
0883
0884 Matrix2xNd<N> p2D_rot = rot * p2D;
0885
0886 #ifdef RFIT_DEBUG
0887 printf("Fast fit Tan(theta): %g\n", fast_fit(3));
0888 printf("Rotation angle: %g\n", theta);
0889 printIt(&rot, "Rotation Matrix:");
0890 printIt(&p2D, "Original Hits(s,z):");
0891 printIt(&p2D_rot, "Rotated hits(S3D, Z'):");
0892 printIt(&rot, "Rotation Matrix:");
0893 #endif
0894
0895
0896 Matrix2xNd<N> aMat;
0897 aMat << MatrixXd::Ones(1, n), p2D_rot.row(0);
0898
0899 #ifdef RFIT_DEBUG
0900 printIt(&aMat, "A Matrix:");
0901 #endif
0902
0903
0904 MatrixNd<N> vyInvMat;
0905 math::cholesky::invert(cov_with_ms, vyInvMat);
0906
0907 Eigen::Matrix<double, 2, 2> covParamsMat = aMat * vyInvMat * aMat.transpose();
0908
0909 math::cholesky::invert(covParamsMat, covParamsMat);
0910
0911
0912
0913
0914 Eigen::Matrix<double, 2, 1> sol = covParamsMat * aMat * vyInvMat * p2D_rot.row(1).transpose();
0915
0916 #ifdef RFIT_DEBUG
0917 printIt(&sol, "Rotated solutions:");
0918 #endif
0919
0920
0921 const auto sinTheta = sin(theta);
0922 const auto cosTheta = cos(theta);
0923 auto common_factor = 1. / (sinTheta - sol(1, 0) * cosTheta);
0924 Eigen::Matrix<double, 2, 2> jMat;
0925 jMat << 0., common_factor * common_factor, common_factor, sol(0, 0) * cosTheta * common_factor * common_factor;
0926
0927 double tempM = common_factor * (sol(1, 0) * sinTheta + cosTheta);
0928 double tempQ = common_factor * sol(0, 0);
0929 auto cov_mq = jMat * covParamsMat * jMat.transpose();
0930
0931 VectorNd<N> res = p2D_rot.row(1).transpose() - aMat.transpose() * sol;
0932 double chi2 = res.transpose() * vyInvMat * res;
0933
0934 LineFit line;
0935 line.par << tempM, tempQ;
0936 line.cov << cov_mq;
0937 line.chi2 = chi2;
0938
0939 #ifdef RFIT_DEBUG
0940 printf("Common_factor: %g\n", common_factor);
0941 printIt(&jMat, "Jacobian:");
0942 printIt(&sol, "Rotated solutions:");
0943 printIt(&covParamsMat, "Cov_params:");
0944 printIt(&cov_mq, "Rotated Covariance Matrix:");
0945 printIt(&(line.par), "Real Parameters:");
0946 printIt(&(line.cov), "Real Covariance Matrix:");
0947 printf("Chi2: %g\n", chi2);
0948 #endif
0949
0950 return line;
0951 }
0952
0953 }
0954
0955 namespace riemannFit {
0956
0957
0958
0959
0960
0961
0962
0963
0964
0965
0966
0967
0968
0969
0970
0971
0972
0973
0974
0975
0976
0977
0978
0979
0980
0981
0982
0983
0984
0985
0986
0987
0988
0989 template <int N>
0990 class helixFit {
0991 public:
0992 template <typename TAcc>
0993 ALPAKA_FN_ACC ALPAKA_FN_INLINE void operator()(const TAcc& acc,
0994 const Matrix3xNd<N>* hits,
0995 const Eigen::Matrix<float, 6, N>* hits_ge,
0996 const double bField,
0997 const bool error,
0998 HelixFit* helix) const {
0999 constexpr uint n = N;
1000 VectorNd<4> rad = (hits->block(0, 0, 2, n).colwise().norm());
1001
1002
1003 Vector4d fast_fit;
1004 ALPAKA_ACCELERATOR_NAMESPACE::riemannFit::fastFit(acc, *hits, fast_fit);
1005 riemannFit::Matrix2Nd<N> hits_cov = MatrixXd::Zero(2 * n, 2 * n);
1006 ALPAKA_ACCELERATOR_NAMESPACE::riemannFit::loadCovariance2D(acc, *hits_ge, hits_cov);
1007 CircleFit circle = ALPAKA_ACCELERATOR_NAMESPACE::riemannFit::circleFit(
1008 acc, hits->block(0, 0, 2, n), hits_cov, fast_fit, rad, bField, error);
1009 LineFit line =
1010 ALPAKA_ACCELERATOR_NAMESPACE::riemannFit::lineFit(acc, *hits, *hits_ge, circle, fast_fit, bField, error);
1011
1012 ALPAKA_ACCELERATOR_NAMESPACE::riemannFit::par_uvrtopak(acc, circle, bField, error);
1013
1014 helix->par << circle.par, line.par;
1015 if (error) {
1016 helix->cov = MatrixXd::Zero(5, 5);
1017 helix->cov.block(0, 0, 3, 3) = circle.cov;
1018 helix->cov.block(3, 3, 2, 2) = line.cov;
1019 }
1020 helix->qCharge = circle.qCharge;
1021 helix->chi2_circle = circle.chi2;
1022 helix->chi2_line = line.chi2;
1023 }
1024 };
1025
1026 }
1027
1028 #endif