ARPACK-Armadillo
UpperHessenbergQR.h
1 // Copyright (C) 2015 Yixuan Qiu
2 //
3 // This Source Code Form is subject to the terms of the Mozilla Public
4 // License, v. 2.0. If a copy of the MPL was not distributed with this
5 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
6 
7 #ifndef UPPER_HESSENBERG_QR_H
8 #define UPPER_HESSENBERG_QR_H
9 
10 #include <armadillo>
11 #include <cmath> // std::sqrt
12 #include <algorithm> // std::fill
13 #include <limits> // std::numeric_limits
14 #include <stdexcept> // std::logic_error
15 
24 template <typename Scalar = double>
26 {
27 private:
28  typedef arma::Mat<Scalar> Matrix;
29  typedef arma::Col<Scalar> Vector;
30  typedef arma::Row<Scalar> RowVector;
31 
32 protected:
33  int n;
34  Matrix mat_T;
35  // Gi = [ cos[i] sin[i]]
36  // [-sin[i] cos[i]]
37  // Q = G1 * G2 * ... * G_{n-1}
38  Vector rot_cos;
39  Vector rot_sin;
40  bool computed;
41 public:
47  n(0), computed(false)
48  {}
49 
59  UpperHessenbergQR(const Matrix &mat) :
60  n(mat.n_rows),
61  mat_T(n, n),
62  rot_cos(n - 1),
63  rot_sin(n - 1),
64  computed(false)
65  {
66  compute(mat);
67  }
68 
77  virtual void compute(const Matrix &mat)
78  {
79  n = mat.n_rows;
80  mat_T.set_size(n, n);
81  rot_cos.set_size(n - 1);
82  rot_sin.set_size(n - 1);
83 
84  // Make a copy of mat
85  mat_T = mat;
86 
87  Scalar xi, xj, r, c, s, eps = std::numeric_limits<Scalar>::epsilon();
88  Scalar *Tii, *ptr;
89  for(int i = 0; i < n - 1; i++)
90  {
91  Tii = mat_T.colptr(i) + i;
92 
93  // Make sure mat_T is upper Hessenberg
94  // Zero the elements below mat_T(i + 1, i)
95  std::fill(Tii + 2, Tii + n - i, Scalar(0));
96 
97  xi = Tii[0]; // mat_T(i, i)
98  xj = Tii[1]; // mat_T(i + 1, i)
99  r = std::sqrt(xi * xi + xj * xj);
100  if(r <= eps)
101  {
102  rot_cos[i] = c = 1;
103  rot_sin[i] = s = 0;
104  } else {
105  rot_cos[i] = c = xi / r;
106  rot_sin[i] = s = -xj / r;
107  }
108  // For a complete QR decomposition,
109  // we first obtain the rotation matrix
110  // G = [ cos sin]
111  // [-sin cos]
112  // and then do T[i:(i + 1), i:(n - 1)] = G' * T[i:(i + 1), i:(n - 1)]
113 
114  // Gt << c << -s << arma::endr << s << c << arma::endr;
115  // mat_T.submat(i, i, i + 1, n - 1) = Gt * mat_T.submat(i, i, i + 1, n - 1);
116  Tii[0] = r; // mat_T(i, i) => r
117  Tii[1] = 0; // mat_T(i + 1, i) => 0
118  ptr = Tii + n; // mat_T(i, k), k = i+1, i+2, ..., n-1
119  for(int j = i + 1; j < n; j++, ptr += n)
120  {
121  Scalar tmp = ptr[0];
122  ptr[0] = c * tmp - s * ptr[1];
123  ptr[1] = s * tmp + c * ptr[1];
124  }
125 
126  // If we do not need to calculate the R matrix, then
127  // only the cos and sin sequences are required.
128  // In such case we only update T[i + 1, (i + 1):(n - 1)]
129  // mat_T(i + 1, arma::span(i + 1, n - 1)) *= c;
130  // mat_T(i + 1, arma::span(i + 1, n - 1)) += s * mat_T(i, arma::span(i + 1, n - 1));
131  }
132 
133  computed = true;
134  }
135 
143  Matrix matrix_R()
144  {
145  if(!computed)
146  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
147 
148  return mat_T;
149  }
150 
158  virtual Matrix matrix_RQ()
159  {
160  if(!computed)
161  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
162 
163  // Make a copy of the R matrix
164  Matrix RQ = mat_T;
165 
166  Scalar *c = rot_cos.memptr(),
167  *s = rot_sin.memptr();
168  for(int i = 0; i < n - 1; i++)
169  {
170  // RQ[, i:(i + 1)] = RQ[, i:(i + 1)] * Gi
171  // Gi = [ cos[i] sin[i]]
172  // [-sin[i] cos[i]]
173 
174  Scalar *Yi, *Yi1;
175  Yi = RQ.colptr(i);
176  Yi1 = Yi + n; // RQ(0, i + 1)
177  for(int j = 0; j < i + 2; j++)
178  {
179  Scalar tmp = Yi[j];
180  Yi[j] = (*c) * tmp - (*s) * Yi1[j];
181  Yi1[j] = (*s) * tmp + (*c) * Yi1[j];
182  }
183 
184  /* Vector Yi = RQ(arma::span(0, i + 1), i);
185  RQ(arma::span(0, i + 1), i) = (*c) * Yi - (*s) * RQ(arma::span(0, i + 1), i + 1);
186  RQ(arma::span(0, i + 1), i + 1) = (*s) * Yi + (*c) * RQ(arma::span(0, i + 1), i + 1); */
187  c++;
188  s++;
189  }
190 
191  return RQ;
192  }
193 
202  // Y -> QY = G1 * G2 * ... * Y
203  void apply_QY(Vector &Y)
204  {
205  if(!computed)
206  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
207 
208  Scalar tmp;
209  for(int i = n - 2; i >= 0; i--)
210  {
211  // Y[i:(i + 1)] = Gi * Y[i:(i + 1)]
212  // Gi = [ cos[i] sin[i]]
213  // [-sin[i] cos[i]]
214  tmp = Y[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];
217  }
218  }
219 
228  // Y -> Q'Y = G_{n-1}' * ... * G2' * G1' * Y
229  void apply_QtY(Vector &Y)
230  {
231  if(!computed)
232  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
233 
234  Scalar tmp;
235  for(int i = 0; i < n - 1; i++)
236  {
237  // Y[i:(i + 1)] = Gi' * Y[i:(i + 1)]
238  // Gi = [ cos[i] sin[i]]
239  // [-sin[i] cos[i]]
240  tmp = Y[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];
243  }
244  }
245 
254  // Y -> QY = G1 * G2 * ... * Y
255  void apply_QY(Matrix &Y)
256  {
257  if(!computed)
258  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
259 
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--)
264  {
265  // Y[i:(i + 1), ] = Gi * Y[i:(i + 1), ]
266  // Gi = [ cos[i] sin[i]]
267  // [-sin[i] cos[i]]
268  Yi = Y.row(i);
269  Yi1 = Y.row(i + 1);
270  Y.row(i) = (*c) * Yi + (*s) * Yi1;
271  Y.row(i + 1) = -(*s) * Yi + (*c) * Yi1;
272  c--;
273  s--;
274  }
275  }
276 
285  // Y -> Q'Y = G_{n-1}' * ... * G2' * G1' * Y
286  void apply_QtY(Matrix &Y)
287  {
288  if(!computed)
289  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
290 
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++)
295  {
296  // Y[i:(i + 1), ] = Gi' * Y[i:(i + 1), ]
297  // Gi = [ cos[i] sin[i]]
298  // [-sin[i] cos[i]]
299  Yi = Y.row(i);
300  Yi1 = Y.row(i + 1);
301  Y.row(i) = (*c) * Yi - (*s) * Yi1;
302  Y.row(i + 1) = (*s) * Yi + (*c) * Yi1;
303  c++;
304  s++;
305  }
306  }
307 
316  // Y -> YQ = Y * G1 * G2 * ...
317  void apply_YQ(Matrix &Y)
318  {
319  if(!computed)
320  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
321 
322  Scalar *c = rot_cos.memptr(),
323  *s = rot_sin.memptr();
324  /*Vector Yi(Y.n_rows);
325  for(int i = 0; i < n - 1; i++)
326  {
327  // Y[, i:(i + 1)] = Y[, i:(i + 1)] * Gi
328  // Gi = [ cos[i] sin[i]]
329  // [-sin[i] cos[i]]
330  Yi = Y.col(i);
331  Y.col(i) = (*c) * Yi - (*s) * Y.col(i + 1);
332  Y.col(i + 1) = (*s) * Yi + (*c) * Y.col(i + 1);
333  c++;
334  s++;
335  }*/
336  Scalar *Y_col_i, *Y_col_i1;
337  int nrow = Y.n_rows;
338  for(int i = 0; i < n - 1; i++)
339  {
340  Y_col_i = Y.colptr(i);
341  Y_col_i1 = Y.colptr(i + 1);
342  for(int j = 0; j < nrow; j++)
343  {
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];
347  }
348  c++;
349  s++;
350  }
351  }
352 
361  // Y -> YQ' = Y * G_{n-1}' * ... * G2' * G1'
362  void apply_YQt(Matrix &Y)
363  {
364  if(!computed)
365  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
366 
367  Scalar *c = rot_cos.memptr() + n - 2,
368  *s = rot_sin.memptr() + n - 2;
369  Vector Yi(Y.n_rows);
370  for(int i = n - 2; i >= 0; i--)
371  {
372  // Y[, i:(i + 1)] = Y[, i:(i + 1)] * Gi'
373  // Gi = [ cos[i] sin[i]]
374  // [-sin[i] cos[i]]
375  Yi = Y.col(i);
376  Y.col(i) = (*c) * Yi + (*s) * Y.col(i + 1);
377  Y.col(i + 1) = -(*s) * Yi + (*c) * Y.col(i + 1);
378  c--;
379  s--;
380  }
381  }
382 };
383 
384 
385 
395 template <typename Scalar = double>
396 class TridiagQR: public UpperHessenbergQR<Scalar>
397 {
398 private:
399  typedef arma::Mat<Scalar> Matrix;
400  typedef arma::Col<Scalar> Vector;
401 
402 public:
408  UpperHessenbergQR<Scalar>()
409  {}
410 
420  TridiagQR(const Matrix &mat) :
421  UpperHessenbergQR<Scalar>()
422  {
423  this->compute(mat);
424  }
425 
434  void compute(const Matrix &mat)
435  {
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);
440 
441  this->mat_T.zeros();
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);
445 
446  // A number of pointers to avoid repeated address calculation
447  Scalar *Tii = this->mat_T.memptr(), // pointer to T[i, i]
448  *ptr, // some location relative to Tii
449  *c = this->rot_cos.memptr(), // pointer to the cosine vector
450  *s = this->rot_sin.memptr(), // pointer to the sine vector
451  r, tmp,
452  eps = std::numeric_limits<Scalar>::epsilon();
453  for(int i = 0; i < this->n - 2; i++)
454  {
455  // Tii[0] == T[i, i]
456  // Tii[1] == T[i + 1, i]
457  r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
458  if(r <= eps)
459  {
460  *c = 1;
461  *s = 0;
462  } else {
463  *c = Tii[0] / r;
464  *s = -Tii[1] / r;
465  }
466 
467  // For a complete QR decomposition,
468  // we first obtain the rotation matrix
469  // G = [ cos sin]
470  // [-sin cos]
471  // and then do T[i:(i + 1), i:(i + 2)] = G' * T[i:(i + 1), i:(i + 2)]
472 
473  // Update T[i, i] and T[i + 1, i]
474  // The updated value of T[i, i] is known to be r
475  // The updated value of T[i + 1, i] is known to be 0
476  Tii[0] = r;
477  Tii[1] = 0;
478  // Update T[i, i + 1] and T[i + 1, i + 1]
479  // ptr[0] == T[i, i + 1]
480  // ptr[1] == T[i + 1, i + 1]
481  ptr = Tii + this->n;
482  tmp = *ptr;
483  ptr[0] = (*c) * tmp - (*s) * ptr[1];
484  ptr[1] = (*s) * tmp + (*c) * ptr[1];
485  // Update T[i, i + 2] and T[i + 1, i + 2]
486  // ptr[0] == T[i, i + 2] == 0
487  // ptr[1] == T[i + 1, i + 2]
488  ptr += this->n;
489  ptr[0] = -(*s) * ptr[1];
490  ptr[1] *= (*c);
491 
492  // Move from T[i, i] to T[i + 1, i + 1]
493  Tii += this->n + 1;
494  // Increase c and s by 1
495  c++;
496  s++;
497 
498 
499  // If we do not need to calculate the R matrix, then
500  // only the cos and sin sequences are required.
501  // In such case we only update T[i + 1, (i + 1):(i + 2)]
502  // this->mat_T(i + 1, i + 1) = (*c) * this->mat_T(i + 1, i + 1) + (*s) * this->mat_T(i, i + 1);
503  // this->mat_T(i + 1, i + 2) *= (*c);
504  }
505  // For i = n - 2
506  r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
507  if(r <= eps)
508  {
509  *c = 1;
510  *s = 0;
511  } else {
512  *c = Tii[0] / r;
513  *s = -Tii[1] / r;
514  }
515  Tii[0] = r;
516  Tii[1] = 0;
517  ptr = Tii + this->n; // points to T[i - 2, i - 1]
518  tmp = *ptr;
519  ptr[0] = (*c) * tmp - (*s) * ptr[1];
520  ptr[1] = (*s) * tmp + (*c) * ptr[1];
521 
522  this->computed = true;
523  }
524 
532  Matrix matrix_RQ()
533  {
534  if(!this->computed)
535  throw std::logic_error("TridiagQR: need to call compute() first");
536 
537  // Make a copy of the R matrix
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);
541 
542  // [m11 m12] will point to RQ[i:(i+1), i:(i+1)]
543  // [m21 m22]
544  Scalar *m11 = RQ.memptr(), *m12, *m21, *m22,
545  *c = this->rot_cos.memptr(),
546  *s = this->rot_sin.memptr(),
547  tmp;
548  for(int i = 0; i < this->n - 1; i++)
549  {
550  m21 = m11 + 1;
551  m12 = m11 + this->n;
552  m22 = m12 + 1;
553  tmp = *m21;
554 
555  // Update diagonal and the below-subdiagonal
556  *m11 = (*c) * (*m11) - (*s) * (*m12);
557  *m21 = (*c) * tmp - (*s) * (*m22);
558  *m22 = (*s) * tmp + (*c) * (*m22);
559 
560  // Move m11 to RQ[i+1, i+1]
561  m11 = m22;
562  c++;
563  s++;
564  }
565 
566  // Copy the below-subdiagonal to above-subdiagonal
567  RQ.diag(1) = RQ.diag(-1);
568 
569  return RQ;
570  }
571 };
572 
573 
574 
575 #endif // UPPER_HESSENBERG_QR_H
virtual void compute(const Matrix &mat)
virtual Matrix matrix_RQ()
Matrix matrix_RQ()
void apply_YQt(Matrix &Y)
void apply_YQ(Matrix &Y)
void apply_QY(Matrix &Y)
void compute(const Matrix &mat)
void apply_QY(Vector &Y)
void apply_QtY(Matrix &Y)
UpperHessenbergQR(const Matrix &mat)
TridiagQR(const Matrix &mat)
void apply_QtY(Vector &Y)