1 #ifndef RecoPixelVertexing_PixelTrackFitting_alpaka_FitUtils_h 2 #define RecoPixelVertexing_PixelTrackFitting_alpaka_FitUtils_h 3 #include <alpaka/alpaka.hpp> 17 using ArrayNd = Eigen::Array<double, N, N>;
19 using Matrix2Nd = Eigen::Matrix<double, 2 * N, 2 * N>;
21 using Matrix3Nd = Eigen::Matrix<double, 3 * N, 3 * N>;
53 template <
typename VI5,
typename MI5,
typename VO5,
typename MO5>
55 auto sinTheta2 = 1. / (1. + ip(3) * ip(3));
57 auto cosTheta = ip(3) * sinTheta;
59 op(0) = sinTheta * ip(2);
67 jMat(0, 2) = sinTheta;
68 jMat(0, 3) = -sinTheta2 * cosTheta * ip(2);
74 ocov = jMat * icov * jMat.transpose();
83 template <
typename TAcc,
class C>
86 for (
uint r = 0; r <
m->rows(); ++r) {
87 for (
uint c = 0;
c <
m->cols(); ++
c) {
88 printf(
"%s Matrix(%d,%d) = %g\n",
prefix, r,
c, (*
m)(r,
c));
110 template <
typename TAcc>
112 return a.x() *
b.y() -
a.y() *
b.x();
119 template <
typename TAcc,
typename M6xNf,
typename M2Nd>
133 constexpr uint32_t hits_in_fit = M6xNf::ColsAtCompileTime;
134 for (uint32_t
i = 0;
i < hits_in_fit; ++
i) {
137 hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) = ge.col(
i)[ge_idx];
141 hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) = ge.col(
i)[ge_idx];
145 hits_cov(
i +
l * hits_in_fit,
i +
j * hits_in_fit) = hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) =
151 template <
typename TAcc,
typename M6xNf,
typename M3xNd>
152 ALPAKA_FN_ACC
void loadCovariance(
const TAcc& acc, M6xNf
const& ge, M3xNd& hits_cov) {
165 constexpr uint32_t hits_in_fit = M6xNf::ColsAtCompileTime;
166 for (uint32_t
i = 0;
i < hits_in_fit; ++
i) {
169 hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) = ge.col(
i)[ge_idx];
173 hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) = ge.col(
i)[ge_idx];
177 hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) = ge.col(
i)[ge_idx];
181 hits_cov(
i +
l * hits_in_fit,
i +
j * hits_in_fit) = hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) =
186 hits_cov(
i +
l * hits_in_fit,
i +
j * hits_in_fit) = hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) =
191 hits_cov(
i +
l * hits_in_fit,
i +
j * hits_in_fit) = hits_cov(
i +
j * hits_in_fit,
i +
l * hits_in_fit) =
205 template <
typename TAcc>
211 const double temp0 = circle.par.head(2).squaredNorm();
213 par_pak << alpaka::math::atan2(acc, circle.qCharge * circle.par(0), -circle.qCharge * circle.par(1)),
214 circle.qCharge * (temp1 - circle.par(2)), circle.par(2) *
B;
216 const double temp2 =
sqr(circle.par(0)) * 1. / temp0;
217 const double temp3 = 1. / temp1 * circle.qCharge;
219 j4Mat << -circle.par(1) * temp2 * 1. /
sqr(circle.par(0)), temp2 * 1. / circle.par(0), 0.,
220 circle.par(0) * temp3, circle.par(1) * temp3, -circle.qCharge, 0., 0.,
B;
221 circle.cov = j4Mat * circle.cov * j4Mat.transpose();
223 circle.par = par_pak;
232 template <
typename TAcc>
235 const double temp0 = circle.par.head(2).squaredNorm();
237 par_pak << alpaka::math::atan2(acc, circle.qCharge * circle.par(0), -circle.qCharge * circle.par(1)),
238 circle.qCharge * (temp1 - circle.par(2)), circle.qCharge / circle.par(2);
240 const double temp2 =
sqr(circle.par(0)) * 1. / temp0;
241 const double temp3 = 1. / temp1 * circle.qCharge;
243 j4Mat << -circle.par(1) * temp2 * 1. /
sqr(circle.par(0)), temp2 * 1. / circle.par(0), 0., circle.par(0) * temp3,
244 circle.par(1) * temp3, -circle.qCharge, 0., 0., -circle.qCharge / (circle.par(2) * circle.par(2));
245 circle.cov = j4Mat * circle.cov * j4Mat.transpose();
247 circle.par = par_pak;
253 #endif // RecoPixelVertexing_PixelTrackFitting_interface_FitUtils_h Eigen::Matrix< double, 2, N > Matrix2xNd
Eigen::Matrix< double, 6, 1 > Vector6f
Eigen::Matrix< double, 3 *N, 1 > Vector3Nd
Eigen::Matrix< double, 5, 5 > Matrix5d
Eigen::Matrix< double, 1, 2 *N > RowVector2Nd
ALPAKA_FN_ACC void loadCovariance2D(const TAcc &acc, M6xNf const &ge, M2Nd &hits_cov)
Eigen::Matrix< double, 2 *N, 1 > Vector2Nd
Eigen::Matrix< double, 2, 3 > Matrix2x3d
Eigen::Matrix< double, 3 *N, 3 *N > Matrix3Nd
ALPAKA_FN_ACC ALPAKA_FN_INLINE double cross2D(const TAcc &acc, const Vector2d &a, const Vector2d &b)
Compute cross product of two 2D vector (assuming z component 0), returning z component of the result...
ALPAKA_FN_ACC ALPAKA_FN_INLINE void par_uvrtopak(const TAcc &acc, CircleFit &circle, const double B, const bool error)
Transform circle parameter from (X0,Y0,R) to (phi,Tip,p_t) and consequently covariance matrix...
Eigen::Array< double, N, N > ArrayNd
Eigen::Matrix< double, N+1, 1 > VectorNplusONEd
constexpr T sqr(const T a)
raise to square.
void transformToPerigeePlane(VI5 const &ip, MI5 const &icov, VO5 &op, MO5 &ocov)
Eigen::Matrix< double, N, 1 > VectorNd
Eigen::Matrix< double, N, N > MatrixNd
constexpr double epsilon
used in numerical derivative (J2 in Circle_fit())
ALPAKA_FN_ACC ALPAKA_FN_INLINE void fromCircleToPerigee(const TAcc &acc, CircleFit &circle)
Transform circle parameter from (X0,Y0,R) to (phi,Tip,q/R) and consequently covariance matrix...
Eigen::Matrix< double, N, 5 > MatrixNx5d
Eigen::Matrix< double, 2 *N, 2 *N > Matrix2Nd
ALPAKA_FN_ACC void loadCovariance(const TAcc &acc, M6xNf const &ge, M3xNd &hits_cov)
Eigen::Array< double, 2, N > Array2xNd
ALPAKA_FN_ACC void printIt(const TAcc &acc, C *m, const char *prefix="")
Eigen::Matrix< double, N, 3 > MatrixNx3d
Eigen::Matrix< double, N+1, N+1 > MatrixNplusONEd
Eigen::Matrix< double, 1, 1, N > RowVectorNd