7 #ifndef DOUBLE_SHIFT_QR_H 8 #define DOUBLE_SHIFT_QR_H 16 #include "../Util/TypeTraits.h" 20 template <
typename Scalar =
double>
55 unsigned char* nr = m_ref_nr.
data();
60 if (
abs(x3) < m_near_0)
63 if (
abs(x2) < m_near_0)
76 x2x3 = Eigen::numext::hypot(x2, x3);
81 Scalar x1_new = x1 - ((x1 <= 0) - (x1 > 0)) * Eigen::numext::hypot(x1, x2x3);
82 Scalar x_norm = Eigen::numext::hypot(x1_new, x2x3);
84 if (x_norm < m_near_0)
89 u[0] = x1_new / x_norm;
103 const Index bsize = iu - il + 1;
113 x01 = m_mat_H.
coeff(il, il + 1),
114 x10 = m_mat_H.
coeff(il + 1, il),
115 x11 = m_mat_H.
coeff(il + 1, il + 1);
127 apply_PX(m_mat_H.block(il, il, 2, m_n - il),
m_n, il);
128 apply_XP(m_mat_H.block(0, il, il + 2, 2),
m_n, il);
136 const Scalar m20 = m_mat_H.
coeff(il + 2, il + 1) * m_mat_H.
coeff(il + 1, il);
140 apply_PX(m_mat_H.block(il, il, 3, m_n - il),
m_n, il);
145 for (Index
i = 1;
i < bsize - 2;
i++)
149 apply_PX(m_mat_H.block(il +
i, il +
i - 1, 3, m_n - il -
i + 1),
m_n, il +
i);
157 apply_PX(m_mat_H.block(iu - 1, iu - 2, 2, m_n - iu + 2),
m_n, iu - 1);
158 apply_XP(m_mat_H.block(0, iu - 1, il + bsize, 2),
m_n, iu - 1);
165 void apply_PX(GenericMatrix
X, Index stride, Index u_ind)
const 167 const Index nr = m_ref_nr.
coeff(u_ind);
172 u1 = m_ref_u.
coeff(1, u_ind);
176 const Index nrow = X.rows();
177 const Index ncol = X.cols();
180 if (nr == 2 || nrow == 2)
182 for (Index
i = 0;
i < ncol;
i++, xptr += stride)
184 const Scalar tmp = u0_2 * xptr[0] + u1_2 * xptr[1];
193 for (Index
i = 0;
i < ncol;
i++, xptr += stride)
195 const Scalar tmp = u0_2 * xptr[0] + u1_2 * xptr[1] + u2_2 * xptr[2];
207 const Index nr = m_ref_nr.
coeff(u_ind);
212 u1 = m_ref_u.
coeff(1, u_ind),
216 const bool nr_is_2 = (nr == 2);
217 const Scalar dot2 =
Scalar(2) * (x[0] * u0 + x[1] * u1 + (nr_is_2 ? 0 : (x[2] *
u2)));
225 void apply_XP(GenericMatrix
X, Index stride, Index u_ind)
const 227 const Index nr = m_ref_nr.
coeff(u_ind);
232 u1 = m_ref_u.
coeff(1, u_ind);
236 const int nrow = X.rows();
237 const int ncol = X.cols();
238 Scalar *X0 = X.data(), *X1 = X0 + stride;
240 if (nr == 2 || ncol == 2)
245 for (Index
i = 0;
i < nrow;
i++)
247 const Scalar tmp = u0_2 * X0[
i] + u1_2 * X1[
i];
257 for (Index
i = 0;
i < nrow;
i++)
259 const Scalar tmp = u0_2 * X0[
i] + u1_2 * X1[
i] + u2_2 * X2[
i];
273 m_eps_abs(m_near_0 * (m_n / m_eps)),
287 m_eps_abs(m_near_0 * (m_n / m_eps)),
298 if (m_n != mat.cols())
299 throw std::invalid_argument(
"DoubleShiftQR: matrix must be square");
308 std::copy(mat.data(), mat.data() + mat.size(), m_mat_H.
data());
312 std::vector<int> zero_ind;
313 zero_ind.reserve(m_n - 1);
314 zero_ind.push_back(0);
316 for (Index
i = 0;
i < m_n - 2;
i++, Hii += (m_n + 1))
320 if (h <= 0 || h <= m_eps_rel * (
abs(Hii[0]) +
abs(Hii[m_n + 1])))
323 zero_ind.push_back(
i + 1);
327 std::fill(Hii + 2, Hii + m_n -
i,
Scalar(0));
329 zero_ind.push_back(m_n);
331 for (std::vector<int>::size_type
i = 0;
i < zero_ind.size() - 1;
i++)
333 const Index start = zero_ind[
i];
334 const Index
end = zero_ind[
i + 1] - 1;
345 throw std::logic_error(
"DoubleShiftQR: need to call compute() first");
355 throw std::logic_error(
"DoubleShiftQR: need to call compute() first");
358 const Index
n1 = m_n - 1;
359 for (Index
i = 0;
i <
n1;
i++, y_ptr++)
370 throw std::logic_error(
"DoubleShiftQR: need to call compute() first");
372 const Index nrow = Y.rows();
373 const Index
n2 = m_n - 2;
374 for (Index
i = 0;
i <
n2;
i++)
384 #endif // DOUBLE_SHIFT_QR_H Eigen::Array< unsigned char, Eigen::Dynamic, 1 > IntArray
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Namespace containing all symbols from the Eigen library.
Pose3 x2(Rot3::Ypr(0.0, 0.0, 0.0), l2)
void compute(ConstGenericMatrix &mat, const Scalar &s, const Scalar &t)
void apply_YQ(GenericMatrix Y) const
void compute_reflector(const Scalar *x, Index ind)
void matrix_QtHQ(Matrix &dest) const
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
void apply_XP(GenericMatrix X, Index stride, Index u_ind) const
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
void compute_reflector(const Scalar &x1, const Scalar &x2, const Scalar &x3, Index ind)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
void update_block(Index il, Index iu)
void apply_PX(GenericMatrix X, Index stride, Index u_ind) const
Pose3 x3(Rot3::Ypr(M_PI/4.0, 0.0, 0.0), l2)
A matrix or vector expression mapping an existing expression.
Eigen::Ref< Matrix > GenericMatrix
DoubleShiftQR(ConstGenericMatrix &mat, const Scalar &s, const Scalar &t)
static EIGEN_DEPRECATED const end_t end
void apply_PX(Scalar *x, Index u_ind) const
const Eigen::Ref< const Matrix > ConstGenericMatrix
Eigen::Matrix< Scalar, 3, Eigen::Dynamic > Matrix3X
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
DoubleShiftQR(Index size)
void apply_QtY(Vector &y) const