7 #ifndef UPPER_HESSENBERG_QR_H
8 #define UPPER_HESSENBERG_QR_H
24 template <
typename Scalar =
double>
28 typedef arma::Mat<Scalar> Matrix;
29 typedef arma::Col<Scalar> Vector;
30 typedef arma::Row<Scalar> RowVector;
81 rot_cos.set_size(n - 1);
82 rot_sin.set_size(n - 1);
87 Scalar xi, xj, r, c, s, eps = std::numeric_limits<Scalar>::epsilon();
89 for(
int i = 0; i < n - 1; i++)
91 Tii = mat_T.colptr(i) + i;
95 std::fill(Tii + 2, Tii + n - i, Scalar(0));
99 r = std::sqrt(xi * xi + xj * xj);
105 rot_cos[i] = c = xi / r;
106 rot_sin[i] = s = -xj / r;
119 for(
int j = i + 1; j < n; j++, ptr += n)
122 ptr[0] = c * tmp - s * ptr[1];
123 ptr[1] = s * tmp + c * ptr[1];
146 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
161 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
166 Scalar *c = rot_cos.memptr(),
167 *s = rot_sin.memptr();
168 for(
int i = 0; i < n - 1; i++)
177 for(
int j = 0; j < i + 2; j++)
180 Yi[j] = (*c) * tmp - (*s) * Yi1[j];
181 Yi1[j] = (*s) * tmp + (*c) * Yi1[j];
206 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
209 for(
int i = n - 2; i >= 0; i--)
215 Y[i] = rot_cos[i] * tmp + rot_sin[i] * Y[i + 1];
216 Y[i + 1] = -rot_sin[i] * tmp + rot_cos[i] * Y[i + 1];
232 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
235 for(
int i = 0; i < n - 1; i++)
241 Y[i] = rot_cos[i] * tmp - rot_sin[i] * Y[i + 1];
242 Y[i + 1] = rot_sin[i] * tmp + rot_cos[i] * Y[i + 1];
258 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
260 Scalar *c = rot_cos.memptr() + n - 2,
261 *s = rot_sin.memptr() + n - 2;
262 RowVector Yi(Y.n_cols), Yi1(Y.n_cols);
263 for(
int i = n - 2; i >= 0; i--)
270 Y.row(i) = (*c) * Yi + (*s) * Yi1;
271 Y.row(i + 1) = -(*s) * Yi + (*c) * Yi1;
289 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
291 Scalar *c = rot_cos.memptr(),
292 *s = rot_sin.memptr();
293 RowVector Yi(Y.n_cols), Yi1(Y.n_cols);
294 for(
int i = 0; i < n - 1; i++)
301 Y.row(i) = (*c) * Yi - (*s) * Yi1;
302 Y.row(i + 1) = (*s) * Yi + (*c) * Yi1;
320 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
322 Scalar *c = rot_cos.memptr(),
323 *s = rot_sin.memptr();
336 Scalar *Y_col_i, *Y_col_i1;
338 for(
int i = 0; i < n - 1; i++)
340 Y_col_i = Y.colptr(i);
341 Y_col_i1 = Y.colptr(i + 1);
342 for(
int j = 0; j < nrow; j++)
344 Scalar tmp = Y_col_i[j];
345 Y_col_i[j] = (*c) * tmp - (*s) * Y_col_i1[j];
346 Y_col_i1[j] = (*s) * tmp + (*c) * Y_col_i1[j];
365 throw std::logic_error(
"UpperHessenbergQR: need to call compute() first");
367 Scalar *c = rot_cos.memptr() + n - 2,
368 *s = rot_sin.memptr() + n - 2;
370 for(
int i = n - 2; i >= 0; i--)
376 Y.col(i) = (*c) * Yi + (*s) * Y.col(i + 1);
377 Y.col(i + 1) = -(*s) * Yi + (*c) * Y.col(i + 1);
395 template <
typename Scalar =
double>
399 typedef arma::Mat<Scalar> Matrix;
400 typedef arma::Col<Scalar> Vector;
436 this->n = mat.n_rows;
437 this->mat_T.set_size(this->n, this->n);
438 this->rot_cos.set_size(this->n - 1);
439 this->rot_sin.set_size(this->n - 1);
442 this->mat_T.diag() = mat.diag();
443 this->mat_T.diag(1) = mat.diag(-1);
444 this->mat_T.diag(-1) = mat.diag(-1);
447 Scalar *Tii = this->mat_T.memptr(),
449 *c = this->rot_cos.memptr(),
450 *s = this->rot_sin.memptr(),
452 eps = std::numeric_limits<Scalar>::epsilon();
453 for(
int i = 0; i < this->n - 2; i++)
457 r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
483 ptr[0] = (*c) * tmp - (*s) * ptr[1];
484 ptr[1] = (*s) * tmp + (*c) * ptr[1];
489 ptr[0] = -(*s) * ptr[1];
506 r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
519 ptr[0] = (*c) * tmp - (*s) * ptr[1];
520 ptr[1] = (*s) * tmp + (*c) * ptr[1];
522 this->computed =
true;
535 throw std::logic_error(
"TridiagQR: need to call compute() first");
538 Matrix RQ(this->n, this->n, arma::fill::zeros);
539 RQ.diag() = this->mat_T.diag();
540 RQ.diag(1) = this->mat_T.diag(1);
544 Scalar *m11 = RQ.memptr(), *m12, *m21, *m22,
545 *c = this->rot_cos.memptr(),
546 *s = this->rot_sin.memptr(),
548 for(
int i = 0; i < this->n - 1; i++)
556 *m11 = (*c) * (*m11) - (*s) * (*m12);
557 *m21 = (*c) * tmp - (*s) * (*m22);
558 *m22 = (*s) * tmp + (*c) * (*m22);
567 RQ.diag(1) = RQ.diag(-1);
575 #endif // UPPER_HESSENBERG_QR_H
virtual void compute(const Matrix &mat)
virtual Matrix matrix_RQ()
void apply_YQt(Matrix &Y)
void compute(const Matrix &mat)
void apply_QtY(Matrix &Y)
UpperHessenbergQR(const Matrix &mat)
TridiagQR(const Matrix &mat)
void apply_QtY(Vector &Y)