7 #ifndef UPPER_HESSENBERG_QR_H 8 #define UPPER_HESSENBERG_QR_H 42 template <
typename Scalar =
double>
76 const Scalar xabs = x * xsign;
77 const Scalar yabs = y * ysign;
156 if (m_n != mat.cols())
157 throw std::invalid_argument(
"UpperHessenbergQR: matrix must be square");
161 m_rot_cos.
resize(m_n - 1);
162 m_rot_sin.
resize(m_n - 1);
165 std::copy(mat.data(), mat.data() + mat.size(), m_mat_T.
data());
166 m_mat_T.diagonal().array() -=
m_shift;
170 const Index
n1 = m_n - 1;
171 for (Index
i = 0;
i <
n1;
i++)
177 std::fill(Tii + 2, Tii + m_n -
i,
Scalar(0));
196 for (Index
j =
i + 1;
j <
m_n;
j++, ptr +=
m_n)
199 ptr[0] = c * tmp - s * ptr[1];
200 ptr[1] = s * tmp + c * ptr[1];
223 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
238 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
242 std::copy(m_mat_T.
data(), m_mat_T.
data() + m_mat_T.size(), dest.
data());
245 const Index
n1 = m_n - 1;
246 for (Index
i = 0;
i <
n1;
i++)
256 const Index i2 =
i + 2;
257 for (Index
j = 0;
j < i2;
j++)
260 Yi[
j] = c * tmp - s * Yi1[
j];
261 Yi1[
j] = s * tmp + c * Yi1[
j];
270 dest.diagonal().array() +=
m_shift;
285 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
287 for (Index
i = m_n - 2;
i >= 0;
i--)
295 Y[
i] = c * tmp + s * Y[
i + 1];
296 Y[
i + 1] = -s * tmp + c * Y[
i + 1];
312 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
314 const Index
n1 = m_n - 1;
315 for (Index
i = 0;
i <
n1;
i++)
323 Y[
i] = c * tmp - s * Y[
i + 1];
324 Y[
i + 1] = s * tmp + c * Y[
i + 1];
341 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
343 RowVector Yi(Y.cols()), Yi1(Y.cols());
344 for (Index
i = m_n - 2;
i >= 0;
i--)
351 Yi.noalias() = Y.row(
i);
352 Yi1.noalias() = Y.row(
i + 1);
353 Y.row(
i) = c * Yi + s * Yi1;
354 Y.row(
i + 1) = -s * Yi + c * Yi1;
371 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
373 RowVector Yi(Y.cols()), Yi1(Y.cols());
374 const Index
n1 = m_n - 1;
375 for (Index
i = 0;
i <
n1;
i++)
382 Yi.noalias() = Y.row(
i);
383 Yi1.noalias() = Y.row(
i + 1);
384 Y.row(
i) = c * Yi - s * Yi1;
385 Y.row(
i + 1) = s * Yi + c * Yi1;
402 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
416 Scalar *Y_col_i, *Y_col_i1;
417 const Index
n1 = m_n - 1;
418 const Index nrow = Y.rows();
419 for (Index
i = 0;
i <
n1;
i++)
424 Y_col_i = &Y.coeffRef(0,
i);
425 Y_col_i1 = &Y.coeffRef(0,
i + 1);
426 for (Index
j = 0;
j < nrow;
j++)
429 Y_col_i[
j] = c * tmp - s * Y_col_i1[
j];
430 Y_col_i1[
j] = s * tmp + c * Y_col_i1[
j];
448 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
451 for (Index
i = m_n - 2;
i >= 0;
i--)
458 Yi.noalias() = Y.col(
i);
459 Y.col(
i) = c * Yi + s * Y.col(
i + 1);
460 Y.col(
i + 1) = -s * Yi + c * Y.col(
i + 1);
474 template <
typename Scalar =
double>
528 this->m_n = mat.rows();
529 if (this->m_n != mat.cols())
530 throw std::invalid_argument(
"TridiagQR: matrix must be square");
532 this->m_shift = shift;
533 m_T_diag.
resize(this->m_n);
534 m_T_lsub.
resize(this->m_n - 1);
535 m_T_usub.
resize(this->m_n - 1);
536 m_T_usub2.
resize(this->m_n - 2);
537 this->m_rot_cos.
resize(this->m_n - 1);
538 this->m_rot_sin.
resize(this->m_n - 1);
540 m_T_diag.array() = mat.diagonal().array() - this->
m_shift;
541 m_T_lsub.noalias() = mat.diagonal(-1);
542 m_T_usub.noalias() = m_T_lsub;
546 *
s = this->m_rot_sin.
data(),
548 const Index
n1 = this->m_n - 1;
549 for (Index
i = 0;
i <
n1;
i++)
593 this->m_computed =
true;
605 if (!this->m_computed)
606 throw std::logic_error(
"TridiagQR: need to call compute() first");
608 Matrix
R = Matrix::Zero(this->m_n, this->m_n);
609 R.diagonal().noalias() = m_T_diag;
610 R.diagonal(1).noalias() = m_T_usub;
611 R.diagonal(2).noalias() = m_T_usub2;
625 if (!this->m_computed)
626 throw std::logic_error(
"TridiagQR: need to call compute() first");
629 dest.
resize(this->m_n, this->m_n);
631 dest.diagonal().noalias() = m_T_diag;
641 const Index
n1 = this->m_n - 1;
642 for (Index
i = 0;
i <
n1;
i++)
648 m22 = m_T_diag.
coeff(
i + 1);
657 dest.diagonal(1).noalias() = dest.diagonal(-1);
660 dest.diagonal().array() += this->
m_shift;
670 #endif // UPPER_HESSENBERG_QR_H virtual void compute(ConstGenericMatrix &mat, const Scalar &shift=Scalar(0))
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
void compute(ConstGenericMatrix &mat, const Scalar &shift=Scalar(0))
void apply_YQt(GenericMatrix Y) const
UpperHessenbergQR(Index size)
Rot2 R(Rot2::fromAngle(0.1))
void apply_QtY(Vector &Y) const
const Eigen::Ref< const Matrix > ConstGenericMatrix
UpperHessenbergQR(ConstGenericMatrix &mat, const Scalar &shift=Scalar(0))
Eigen::Array< Scalar, Eigen::Dynamic, 1 > Array
virtual Matrix matrix_R() const
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 y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size ratio
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Eigen::Ref< Matrix > GenericMatrix
virtual void matrix_QtHQ(Matrix &dest) const
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
void apply_QY(GenericMatrix Y) const
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 apply_QY(Vector &Y) const
virtual ~UpperHessenbergQR()
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
A matrix or vector expression mapping an existing expression.
const Eigen::Ref< const Matrix > ConstGenericMatrix
void matrix_QtHQ(Matrix &dest) const
Eigen::Matrix< Scalar, 1, Eigen::Dynamic > RowVector
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Jet< T, N > sqrt(const Jet< T, N > &f)
void apply_YQ(GenericMatrix Y) const
TridiagQR(ConstGenericMatrix &mat, const Scalar &shift=Scalar(0))
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
void apply_QtY(GenericMatrix Y) const
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
static void compute_rotation(const Scalar &x, const Scalar &y, Scalar &r, Scalar &c, Scalar &s)