1 #ifndef UPPER_HESSENBERG_QR_H
2 #define UPPER_HESSENBERG_QR_H
20 template <
typename Scalar =
double>
24 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Matrix;
25 typedef Eigen::Matrix<Scalar, 2, 2> Matrix22;
26 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> Vector;
27 typedef Eigen::Matrix<Scalar, 1, Eigen::Dynamic> RowVector;
28 typedef Eigen::Array<Scalar, Eigen::Dynamic, 1> Array;
30 typedef Eigen::Ref<Matrix> GenericMatrix;
31 typedef const Eigen::Ref<const Matrix> ConstGenericMatrix;
80 virtual void compute(ConstGenericMatrix &mat)
84 rot_cos.resize(n - 1);
85 rot_sin.resize(n - 1);
87 mat_T = mat.template triangularView<Eigen::Upper>();
88 mat_T.diagonal(-1) = mat.diagonal(-1);
90 Scalar xi, xj, r, c, s;
92 for(
int i = 0; i < n - 2; i++)
96 r = std::sqrt(xi * xi + xj * xj);
97 rot_cos[i] = c = xi / r;
98 rot_sin[i] = s = -xj / r;
106 mat_T.block(i, i, 2, n - i) = Gt * mat_T.block(i, i, 2, n - i);
115 xi = mat_T(n - 2, n - 2);
116 xj = mat_T(n - 1, n - 2);
117 r = std::sqrt(xi * xi + xj * xj);
118 rot_cos[n - 2] = c = xi / r;
119 rot_sin[n - 2] = s = -xj / r;
121 mat_T.template block<2, 2>(n - 2, n - 2) = Gt * mat_T.template block<2, 2>(n - 2, n - 2);
136 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
151 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
154 Matrix RQ = mat_T.template triangularView<Eigen::Upper>();
156 Scalar *c = rot_cos.data(),
158 for(
int i = 0; i < n - 1; i++)
163 Vector Yi = RQ.block(0, i, i + 2, 1);
164 RQ.block(0, i, i + 2, 1) = (*c) * Yi - (*s) * RQ.block(0, i + 1, i + 2, 1);
165 RQ.block(0, i + 1, i + 2, 1) = (*s) * Yi + (*c) * RQ.block(0, i + 1, i + 2, 1);
185 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
188 for(
int i = n - 2; i >= 0; i--)
194 Y[i] = rot_cos[i] * tmp + rot_sin[i] * Y[i + 1];
195 Y[i + 1] = -rot_sin[i] * tmp + rot_cos[i] * Y[i + 1];
211 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
214 for(
int i = 0; i < n - 1; i++)
220 Y[i] = rot_cos[i] * tmp - rot_sin[i] * Y[i + 1];
221 Y[i + 1] = rot_sin[i] * tmp + rot_cos[i] * Y[i + 1];
238 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
240 Scalar *c = rot_cos.data() + n - 2,
241 *s = rot_sin.data() + n - 2;
242 RowVector Yi(Y.cols()), Yi1(Y.cols());
243 for(
int i = n - 2; i >= 0; i--)
250 Y.row(i) = (*c) * Yi + (*s) * Yi1;
251 Y.row(i + 1) = -(*s) * Yi + (*c) * Yi1;
270 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
272 Scalar *c = rot_cos.data(),
274 RowVector Yi(Y.cols()), Yi1(Y.cols());
275 for(
int i = 0; i < n - 1; i++)
282 Y.row(i) = (*c) * Yi - (*s) * Yi1;
283 Y.row(i + 1) = (*s) * Yi + (*c) * Yi1;
302 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
304 Scalar *c = rot_cos.data(),
307 for(
int i = 0; i < n - 1; i++)
313 Y.col(i) = (*c) * Yi - (*s) * Y.col(i + 1);
314 Y.col(i + 1) = (*s) * Yi + (*c) * Y.col(i + 1);
333 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
335 Scalar *c = rot_cos.data() + n - 2,
336 *s = rot_sin.data() + n - 2;
338 for(
int i = n - 2; i >= 0; i--)
344 Y.col(i) = (*c) * Yi + (*s) * Y.col(i + 1);
345 Y.col(i + 1) = -(*s) * Yi + (*c) * Y.col(i + 1);
361 template <
typename Scalar =
double>
365 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Matrix;
366 typedef const Eigen::Ref<const Matrix> ConstGenericMatrix;
404 this->n = mat.rows();
405 this->mat_T.resize(this->n, this->n);
406 this->rot_cos.resize(this->n - 1);
407 this->rot_sin.resize(this->n - 1);
409 this->mat_T.setZero();
410 this->mat_T.diagonal() = mat.diagonal();
411 this->mat_T.diagonal(1) = mat.diagonal(-1);
412 this->mat_T.diagonal(-1) = mat.diagonal(-1);
415 Scalar *Tii = this->mat_T.data(),
417 *c = this->rot_cos.data(),
418 *s = this->rot_sin.data(),
420 for(
int i = 0; i < this->n - 2; i++)
424 r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
444 ptr[0] = (*c) * tmp - (*s) * ptr[1];
445 ptr[1] = (*s) * tmp + (*c) * ptr[1];
450 ptr[0] = -(*s) * ptr[1];
467 r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
474 ptr[0] = (*c) * tmp - (*s) * ptr[1];
475 ptr[1] = (*s) * tmp + (*c) * ptr[1];
477 this->computed =
true;
490 throw std::logic_error(
"TridiagQR: need to call compute() first");
493 Matrix RQ(this->n, this->n);
495 RQ.diagonal() = this->mat_T.diagonal();
496 RQ.diagonal(1) = this->mat_T.diagonal(1);
500 Scalar *m11 = RQ.data(), *m12, *m21, *m22,
501 *c = this->rot_cos.data(),
502 *s = this->rot_sin.data(),
504 for(
int i = 0; i < this->n - 1; i++)
512 *m11 = (*c) * (*m11) - (*s) * (*m12);
513 *m21 = (*c) * tmp - (*s) * (*m22);
514 *m22 = (*s) * tmp + (*c) * (*m22);
523 RQ.diagonal(1) = RQ.diagonal(-1);
531 #endif // UPPER_HESSENBERG_QR_H
void compute(ConstGenericMatrix &mat)
void apply_YQt(GenericMatrix Y)
virtual Matrix matrix_RQ()
virtual void compute(ConstGenericMatrix &mat)
void apply_YQ(GenericMatrix Y)
UpperHessenbergQR(ConstGenericMatrix &mat)
TridiagQR(ConstGenericMatrix &mat)
void apply_QtY(GenericMatrix Y)
void apply_QY(GenericMatrix Y)
void apply_QtY(Vector &Y)