CMS 3D CMS Logo

FitUtils.h
Go to the documentation of this file.
1 #ifndef RecoPixelVertexing_PixelTrackFitting_alpaka_FitUtils_h
2 #define RecoPixelVertexing_PixelTrackFitting_alpaka_FitUtils_h
3 #include <alpaka/alpaka.hpp>
6 namespace riemannFit {
7 
8  constexpr double epsilon = 1.e-4;
9 
12  template <int N>
13  using MatrixNd = Eigen::Matrix<double, N, N>;
14  template <int N>
15  using MatrixNplusONEd = Eigen::Matrix<double, N + 1, N + 1>;
16  template <int N>
17  using ArrayNd = Eigen::Array<double, N, N>;
18  template <int N>
19  using Matrix2Nd = Eigen::Matrix<double, 2 * N, 2 * N>;
20  template <int N>
21  using Matrix3Nd = Eigen::Matrix<double, 3 * N, 3 * N>;
22  template <int N>
23  using Matrix2xNd = Eigen::Matrix<double, 2, N>;
24  template <int N>
25  using Array2xNd = Eigen::Array<double, 2, N>;
26  template <int N>
27  using MatrixNx3d = Eigen::Matrix<double, N, 3>;
28  template <int N>
29  using MatrixNx5d = Eigen::Matrix<double, N, 5>;
30  template <int N>
31  using VectorNd = Eigen::Matrix<double, N, 1>;
32  template <int N>
33  using VectorNplusONEd = Eigen::Matrix<double, N + 1, 1>;
34  template <int N>
35  using Vector2Nd = Eigen::Matrix<double, 2 * N, 1>;
36  template <int N>
37  using Vector3Nd = Eigen::Matrix<double, 3 * N, 1>;
38  template <int N>
39  using RowVectorNd = Eigen::Matrix<double, 1, 1, N>;
40  template <int N>
41  using RowVector2Nd = Eigen::Matrix<double, 1, 2 * N>;
42 
43  using Matrix2x3d = Eigen::Matrix<double, 2, 3>;
44 
48  using Vector6f = Eigen::Matrix<double, 6, 1>;
49  // transformation between the "perigee" to cmssw localcoord frame
50  // the plane of the latter is the perigee plane...
51  // from //!<(phi,Tip,q/pt,cotan(theta)),Zip)
52  // to q/p,dx/dz,dy/dz,x,z
53  template <typename VI5, typename MI5, typename VO5, typename MO5>
54  inline void transformToPerigeePlane(VI5 const& ip, MI5 const& icov, VO5& op, MO5& ocov) {
55  auto sinTheta2 = 1. / (1. + ip(3) * ip(3));
56  auto sinTheta = std::sqrt(sinTheta2);
57  auto cosTheta = ip(3) * sinTheta;
58 
59  op(0) = sinTheta * ip(2);
60  op(1) = 0.;
61  op(2) = -ip(3);
62  op(3) = ip(1);
63  op(4) = -ip(4);
64 
65  Matrix5d jMat = Matrix5d::Zero();
66 
67  jMat(0, 2) = sinTheta;
68  jMat(0, 3) = -sinTheta2 * cosTheta * ip(2);
69  jMat(1, 0) = 1.;
70  jMat(2, 3) = -1.;
71  jMat(3, 1) = 1.;
72  jMat(4, 4) = -1;
73 
74  ocov = jMat * icov * jMat.transpose();
75  }
76 
77 } // namespace riemannFit
78 
80  namespace riemannFit {
81  using namespace ::riemannFit;
82 
83  template <typename TAcc, class C>
84  ALPAKA_FN_ACC void printIt(const TAcc& acc, C* m, const char* prefix = "") {
85 #ifdef RFIT_DEBUG
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));
89  }
90  }
91 #endif
92  }
93 
97  template <typename T>
98  constexpr T sqr(const T a) {
99  return a * a;
100  }
101 
110  template <typename TAcc>
111  ALPAKA_FN_ACC ALPAKA_FN_INLINE double cross2D(const TAcc& acc, const Vector2d& a, const Vector2d& b) {
112  return a.x() * b.y() - a.y() * b.x();
113  }
114 
119  template <typename TAcc, typename M6xNf, typename M2Nd>
120  ALPAKA_FN_ACC void loadCovariance2D(const TAcc& acc, M6xNf const& ge, M2Nd& hits_cov) {
121  // Index numerology:
122  // i: index of the hits/point (0,..,3)
123  // j: index of space component (x,y,z)
124  // l: index of space components (x,y,z)
125  // ge is always in sync with the index i and is formatted as:
126  // ge[] ==> [xx, xy, yy, xz, yz, zz]
127  // in (j,l) notation, we have:
128  // ge[] ==> [(0,0), (0,1), (1,1), (0,2), (1,2), (2,2)]
129  // so the index ge_idx corresponds to the matrix elements:
130  // | 0 1 3 |
131  // | 1 2 4 |
132  // | 3 4 5 |
133  constexpr uint32_t hits_in_fit = M6xNf::ColsAtCompileTime;
134  for (uint32_t i = 0; i < hits_in_fit; ++i) {
135  {
136  constexpr uint32_t ge_idx = 0, j = 0, l = 0;
137  hits_cov(i + j * hits_in_fit, i + l * hits_in_fit) = ge.col(i)[ge_idx];
138  }
139  {
140  constexpr uint32_t ge_idx = 2, j = 1, l = 1;
141  hits_cov(i + j * hits_in_fit, i + l * hits_in_fit) = ge.col(i)[ge_idx];
142  }
143  {
144  constexpr uint32_t ge_idx = 1, j = 1, l = 0;
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) =
146  ge.col(i)[ge_idx];
147  }
148  }
149  }
150 
151  template <typename TAcc, typename M6xNf, typename M3xNd>
152  ALPAKA_FN_ACC void loadCovariance(const TAcc& acc, M6xNf const& ge, M3xNd& hits_cov) {
153  // Index numerology:
154  // i: index of the hits/point (0,..,3)
155  // j: index of space component (x,y,z)
156  // l: index of space components (x,y,z)
157  // ge is always in sync with the index i and is formatted as:
158  // ge[] ==> [xx, xy, yy, xz, yz, zz]
159  // in (j,l) notation, we have:
160  // ge[] ==> [(0,0), (0,1), (1,1), (0,2), (1,2), (2,2)]
161  // so the index ge_idx corresponds to the matrix elements:
162  // | 0 1 3 |
163  // | 1 2 4 |
164  // | 3 4 5 |
165  constexpr uint32_t hits_in_fit = M6xNf::ColsAtCompileTime;
166  for (uint32_t i = 0; i < hits_in_fit; ++i) {
167  {
168  constexpr uint32_t ge_idx = 0, j = 0, l = 0;
169  hits_cov(i + j * hits_in_fit, i + l * hits_in_fit) = ge.col(i)[ge_idx];
170  }
171  {
172  constexpr uint32_t ge_idx = 2, j = 1, l = 1;
173  hits_cov(i + j * hits_in_fit, i + l * hits_in_fit) = ge.col(i)[ge_idx];
174  }
175  {
176  constexpr uint32_t ge_idx = 5, j = 2, l = 2;
177  hits_cov(i + j * hits_in_fit, i + l * hits_in_fit) = ge.col(i)[ge_idx];
178  }
179  {
180  constexpr uint32_t ge_idx = 1, j = 1, l = 0;
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) =
182  ge.col(i)[ge_idx];
183  }
184  {
185  constexpr uint32_t ge_idx = 3, j = 2, l = 0;
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) =
187  ge.col(i)[ge_idx];
188  }
189  {
190  constexpr uint32_t ge_idx = 4, j = 2, l = 1;
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) =
192  ge.col(i)[ge_idx];
193  }
194  }
195  }
196 
205  template <typename TAcc>
206  ALPAKA_FN_ACC ALPAKA_FN_INLINE void par_uvrtopak(const TAcc& acc,
207  CircleFit& circle,
208  const double B,
209  const bool error) {
210  Vector3d par_pak;
211  const double temp0 = circle.par.head(2).squaredNorm();
212  const double temp1 = alpaka::math::sqrt(acc, temp0);
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;
215  if (error) {
216  const double temp2 = sqr(circle.par(0)) * 1. / temp0;
217  const double temp3 = 1. / temp1 * circle.qCharge;
218  Matrix3d j4Mat;
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();
222  }
223  circle.par = par_pak;
224  }
225 
232  template <typename TAcc>
233  ALPAKA_FN_ACC ALPAKA_FN_INLINE void fromCircleToPerigee(const TAcc& acc, CircleFit& circle) {
234  Vector3d par_pak;
235  const double temp0 = circle.par.head(2).squaredNorm();
236  const double temp1 = alpaka::math::sqrt(acc, temp0);
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);
239 
240  const double temp2 = sqr(circle.par(0)) * 1. / temp0;
241  const double temp3 = 1. / temp1 * circle.qCharge;
242  Matrix3d j4Mat;
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();
246 
247  circle.par = par_pak;
248  }
249 
250  } // namespace riemannFit
251 
252 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
253 #endif // RecoPixelVertexing_PixelTrackFitting_interface_FitUtils_h
Eigen::Matrix< double, 2, N > Matrix2xNd
Definition: FitUtils.h:23
Eigen::Matrix< double, 6, 1 > Vector6f
Definition: FitUtils.h:48
Eigen::Vector3d Vector3d
Definition: FitResult.h:13
Eigen::Matrix< double, 3 *N, 1 > Vector3Nd
Definition: FitUtils.h:37
Eigen::Matrix3d Matrix3d
Definition: FitResult.h:17
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: FitResult.h:19
Eigen::Matrix< double, 1, 2 *N > RowVector2Nd
Definition: FitUtils.h:41
Definition: APVGainStruct.h:7
ALPAKA_FN_ACC void loadCovariance2D(const TAcc &acc, M6xNf const &ge, M2Nd &hits_cov)
Definition: FitUtils.h:120
Eigen::Matrix< double, 2 *N, 1 > Vector2Nd
Definition: FitUtils.h:35
Eigen::Matrix< double, 2, 3 > Matrix2x3d
Definition: FitUtils.h:43
Eigen::Matrix< double, 3 *N, 3 *N > Matrix3Nd
Definition: FitUtils.h:21
Eigen::VectorXd VectorXd
Definition: FitUtils.h:10
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...
Definition: FitUtils.h:111
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...
Definition: FitUtils.h:206
Eigen::Array< double, N, N > ArrayNd
Definition: FitUtils.h:17
Eigen::Matrix< double, N+1, 1 > VectorNplusONEd
Definition: FitUtils.h:33
constexpr T sqr(const T a)
raise to square.
Definition: FitUtils.h:98
void transformToPerigeePlane(VI5 const &ip, MI5 const &icov, VO5 &op, MO5 &ocov)
Definition: FitUtils.h:54
Eigen::Matrix< double, N, 1 > VectorNd
Definition: FitUtils.h:31
T sqrt(T t)
Definition: SSEVec.h:19
Eigen::Matrix< double, N, N > MatrixNd
Definition: FitUtils.h:13
constexpr double epsilon
used in numerical derivative (J2 in Circle_fit())
Definition: FitUtils.h:8
Eigen::Matrix3f Matrix3f
Definition: FitUtils.h:45
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...
Definition: FitUtils.h:233
Eigen::MatrixXd MatrixXd
Definition: FitUtils.h:11
Eigen::Matrix< double, N, 5 > MatrixNx5d
Definition: FitUtils.h:29
Eigen::Matrix< double, 2 *N, 2 *N > Matrix2Nd
Definition: FitUtils.h:19
ALPAKA_FN_ACC void loadCovariance(const TAcc &acc, M6xNf const &ge, M3xNd &hits_cov)
Definition: FitUtils.h:152
Eigen::Array< double, 2, N > Array2xNd
Definition: FitUtils.h:25
Eigen::Vector2d Vector2d
Definition: FitResult.h:12
ALPAKA_FN_ACC void printIt(const TAcc &acc, C *m, const char *prefix="")
Definition: FitUtils.h:84
double b
Definition: hdecay.h:120
Eigen::Matrix< double, N, 3 > MatrixNx3d
Definition: FitUtils.h:27
Eigen::Matrix< double, N+1, N+1 > MatrixNplusONEd
Definition: FitUtils.h:15
double a
Definition: hdecay.h:121
Eigen::Matrix< double, 1, 1, N > RowVectorNd
Definition: FitUtils.h:39
Eigen::Vector3f Vector3f
Definition: FitUtils.h:46
long double T
Eigen::Vector4f Vector4f
Definition: FitUtils.h:47