ARPACK-Eigen
UpperHessenbergQR.h
1 #ifndef UPPER_HESSENBERG_QR_H
2 #define UPPER_HESSENBERG_QR_H
3 
4 #include <Eigen/Core>
5 #include <stdexcept>
6 
11 
20 template <typename Scalar = double>
22 {
23 private:
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;
29 
30  typedef Eigen::Ref<Matrix> GenericMatrix;
31  typedef const Eigen::Ref<const Matrix> ConstGenericMatrix;
32 
33 protected:
34  int n;
35  Matrix mat_T;
36  // Gi = [ cos[i] sin[i]]
37  // [-sin[i] cos[i]]
38  // Q = G1 * G2 * ... * G_{n-1}
39  Array rot_cos;
40  Array rot_sin;
41  bool computed;
42 public:
48  n(0), computed(false)
49  {}
50 
61  UpperHessenbergQR(ConstGenericMatrix &mat) :
62  n(mat.rows()),
63  mat_T(n, n),
64  rot_cos(n - 1),
65  rot_sin(n - 1),
66  computed(false)
67  {
68  compute(mat);
69  }
70 
80  virtual void compute(ConstGenericMatrix &mat)
81  {
82  n = mat.rows();
83  mat_T.resize(n, n);
84  rot_cos.resize(n - 1);
85  rot_sin.resize(n - 1);
86 
87  mat_T = mat.template triangularView<Eigen::Upper>();
88  mat_T.diagonal(-1) = mat.diagonal(-1);
89 
90  Scalar xi, xj, r, c, s;
91  Matrix22 Gt;
92  for(int i = 0; i < n - 2; i++)
93  {
94  xi = mat_T(i, i);
95  xj = mat_T(i + 1, i);
96  r = std::sqrt(xi * xi + xj * xj);
97  rot_cos[i] = c = xi / r;
98  rot_sin[i] = s = -xj / r;
99  // For a complete QR decomposition,
100  // we first obtain the rotation matrix
101  // G = [ cos sin]
102  // [-sin cos]
103  // and then do T[i:(i + 1), i:(n - 1)] = G' * T[i:(i + 1), i:(n - 1)]
104 
105  Gt << c, -s, s, c;
106  mat_T.block(i, i, 2, n - i) = Gt * mat_T.block(i, i, 2, n - i);
107 
108  // If we do not need to calculate the R matrix, then
109  // only the cos and sin sequences are required.
110  // In such case we only update T[i + 1, (i + 1):(n - 1)]
111  // mat_T.block(i + 1, i + 1, 1, n - i - 1) *= c;
112  // mat_T.block(i + 1, i + 1, 1, n - i - 1) += s * mat_T.block(i, i + 1, 1, n - i - 1);
113  }
114  // For i = n - 2
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;
120  Gt << c, -s, s, c;
121  mat_T.template block<2, 2>(n - 2, n - 2) = Gt * mat_T.template block<2, 2>(n - 2, n - 2);
122 
123  computed = true;
124  }
125 
133  Matrix matrix_R()
134  {
135  if(!computed)
136  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
137 
138  return mat_T;
139  }
140 
148  virtual Matrix matrix_RQ()
149  {
150  if(!computed)
151  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
152 
153  // Make a copy of the R matrix
154  Matrix RQ = mat_T.template triangularView<Eigen::Upper>();
155 
156  Scalar *c = rot_cos.data(),
157  *s = rot_sin.data();
158  for(int i = 0; i < n - 1; i++)
159  {
160  // RQ[, i:(i + 1)] = RQ[, i:(i + 1)] * Gi
161  // Gi = [ cos[i] sin[i]]
162  // [-sin[i] cos[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);
166  c++;
167  s++;
168  }
169 
170  return RQ;
171  }
172 
181  // Y -> QY = G1 * G2 * ... * Y
182  void apply_QY(Vector &Y)
183  {
184  if(!computed)
185  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
186 
187  Scalar tmp;
188  for(int i = n - 2; i >= 0; i--)
189  {
190  // Y[i:(i + 1)] = Gi * Y[i:(i + 1)]
191  // Gi = [ cos[i] sin[i]]
192  // [-sin[i] cos[i]]
193  tmp = Y[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];
196  }
197  }
198 
207  // Y -> Q'Y = G_{n-1}' * ... * G2' * G1' * Y
208  void apply_QtY(Vector &Y)
209  {
210  if(!computed)
211  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
212 
213  Scalar tmp;
214  for(int i = 0; i < n - 1; i++)
215  {
216  // Y[i:(i + 1)] = Gi' * Y[i:(i + 1)]
217  // Gi = [ cos[i] sin[i]]
218  // [-sin[i] cos[i]]
219  tmp = Y[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];
222  }
223  }
224 
234  // Y -> QY = G1 * G2 * ... * Y
235  void apply_QY(GenericMatrix Y)
236  {
237  if(!computed)
238  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
239 
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--)
244  {
245  // Y[i:(i + 1), ] = Gi * Y[i:(i + 1), ]
246  // Gi = [ cos[i] sin[i]]
247  // [-sin[i] cos[i]]
248  Yi = Y.row(i);
249  Yi1 = Y.row(i + 1);
250  Y.row(i) = (*c) * Yi + (*s) * Yi1;
251  Y.row(i + 1) = -(*s) * Yi + (*c) * Yi1;
252  c--;
253  s--;
254  }
255  }
256 
266  // Y -> Q'Y = G_{n-1}' * ... * G2' * G1' * Y
267  void apply_QtY(GenericMatrix Y)
268  {
269  if(!computed)
270  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
271 
272  Scalar *c = rot_cos.data(),
273  *s = rot_sin.data();
274  RowVector Yi(Y.cols()), Yi1(Y.cols());
275  for(int i = 0; i < n - 1; i++)
276  {
277  // Y[i:(i + 1), ] = Gi' * Y[i:(i + 1), ]
278  // Gi = [ cos[i] sin[i]]
279  // [-sin[i] cos[i]]
280  Yi = Y.row(i);
281  Yi1 = Y.row(i + 1);
282  Y.row(i) = (*c) * Yi - (*s) * Yi1;
283  Y.row(i + 1) = (*s) * Yi + (*c) * Yi1;
284  c++;
285  s++;
286  }
287  }
288 
298  // Y -> YQ = Y * G1 * G2 * ...
299  void apply_YQ(GenericMatrix Y)
300  {
301  if(!computed)
302  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
303 
304  Scalar *c = rot_cos.data(),
305  *s = rot_sin.data();
306  Vector Yi(Y.rows());
307  for(int i = 0; i < n - 1; i++)
308  {
309  // Y[, i:(i + 1)] = Y[, i:(i + 1)] * Gi
310  // Gi = [ cos[i] sin[i]]
311  // [-sin[i] cos[i]]
312  Yi = Y.col(i);
313  Y.col(i) = (*c) * Yi - (*s) * Y.col(i + 1);
314  Y.col(i + 1) = (*s) * Yi + (*c) * Y.col(i + 1);
315  c++;
316  s++;
317  }
318  }
319 
329  // Y -> YQ' = Y * G_{n-1}' * ... * G2' * G1'
330  void apply_YQt(GenericMatrix Y)
331  {
332  if(!computed)
333  throw std::logic_error("UpperHessenbergQR: need to call compute() first");
334 
335  Scalar *c = rot_cos.data() + n - 2,
336  *s = rot_sin.data() + n - 2;
337  Vector Yi(Y.rows());
338  for(int i = n - 2; i >= 0; i--)
339  {
340  // Y[, i:(i + 1)] = Y[, i:(i + 1)] * Gi'
341  // Gi = [ cos[i] sin[i]]
342  // [-sin[i] cos[i]]
343  Yi = Y.col(i);
344  Y.col(i) = (*c) * Yi + (*s) * Y.col(i + 1);
345  Y.col(i + 1) = -(*s) * Yi + (*c) * Y.col(i + 1);
346  c--;
347  s--;
348  }
349  }
350 };
351 
352 
353 
361 template <typename Scalar = double>
362 class TridiagQR: public UpperHessenbergQR<Scalar>
363 {
364 private:
365  typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Matrix;
366  typedef const Eigen::Ref<const Matrix> ConstGenericMatrix;
367 
368 public:
374  UpperHessenbergQR<Scalar>()
375  {}
376 
387  TridiagQR(ConstGenericMatrix &mat) :
388  UpperHessenbergQR<Scalar>()
389  {
390  this->compute(mat);
391  }
392 
402  void compute(ConstGenericMatrix &mat)
403  {
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);
408 
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);
413 
414  // A number of pointers to avoid repeated address calculation
415  Scalar *Tii = this->mat_T.data(), // pointer to T[i, i]
416  *ptr, // some location relative to Tii
417  *c = this->rot_cos.data(), // pointer to the cosine vector
418  *s = this->rot_sin.data(), // pointer to the sine vector
419  r, tmp;
420  for(int i = 0; i < this->n - 2; i++)
421  {
422  // Tii[0] == T[i, i]
423  // Tii[1] == T[i + 1, i]
424  r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
425  *c = Tii[0] / r;
426  *s = -Tii[1] / r;
427 
428  // For a complete QR decomposition,
429  // we first obtain the rotation matrix
430  // G = [ cos sin]
431  // [-sin cos]
432  // and then do T[i:(i + 1), i:(i + 2)] = G' * T[i:(i + 1), i:(i + 2)]
433 
434  // Update T[i, i] and T[i + 1, i]
435  // The updated value of T[i, i] is known to be r
436  // The updated value of T[i + 1, i] is known to be 0
437  Tii[0] = r;
438  Tii[1] = 0;
439  // Update T[i, i + 1] and T[i + 1, i + 1]
440  // ptr[0] == T[i, i + 1]
441  // ptr[1] == T[i + 1, i + 1]
442  ptr = Tii + this->n;
443  tmp = *ptr;
444  ptr[0] = (*c) * tmp - (*s) * ptr[1];
445  ptr[1] = (*s) * tmp + (*c) * ptr[1];
446  // Update T[i, i + 2] and T[i + 1, i + 2]
447  // ptr[0] == T[i, i + 2] == 0
448  // ptr[1] == T[i + 1, i + 2]
449  ptr += this->n;
450  ptr[0] = -(*s) * ptr[1];
451  ptr[1] *= (*c);
452 
453  // Move from T[i, i] to T[i + 1, i + 1]
454  Tii += this->n + 1;
455  // Increase c and s by 1
456  c++;
457  s++;
458 
459 
460  // If we do not need to calculate the R matrix, then
461  // only the cos and sin sequences are required.
462  // In such case we only update T[i + 1, (i + 1):(i + 2)]
463  // this->mat_T(i + 1, i + 1) = (*c) * this->mat_T(i + 1, i + 1) + (*s) * this->mat_T(i, i + 1);
464  // this->mat_T(i + 1, i + 2) *= (*c);
465  }
466  // For i = n - 2
467  r = std::sqrt(Tii[0] * Tii[0] + Tii[1] * Tii[1]);
468  *c = Tii[0] / r;
469  *s = -Tii[1] / r;
470  Tii[0] = r;
471  Tii[1] = 0;
472  ptr = Tii + this->n; // points to T[i - 2, i - 1]
473  tmp = *ptr;
474  ptr[0] = (*c) * tmp - (*s) * ptr[1];
475  ptr[1] = (*s) * tmp + (*c) * ptr[1];
476 
477  this->computed = true;
478  }
479 
487  Matrix matrix_RQ()
488  {
489  if(!this->computed)
490  throw std::logic_error("TridiagQR: need to call compute() first");
491 
492  // Make a copy of the R matrix
493  Matrix RQ(this->n, this->n);
494  RQ.setZero();
495  RQ.diagonal() = this->mat_T.diagonal();
496  RQ.diagonal(1) = this->mat_T.diagonal(1);
497 
498  // [m11 m12] will point to RQ[i:(i+1), i:(i+1)]
499  // [m21 m22]
500  Scalar *m11 = RQ.data(), *m12, *m21, *m22,
501  *c = this->rot_cos.data(),
502  *s = this->rot_sin.data(),
503  tmp;
504  for(int i = 0; i < this->n - 1; i++)
505  {
506  m21 = m11 + 1;
507  m12 = m11 + this->n;
508  m22 = m12 + 1;
509  tmp = *m21;
510 
511  // Update diagonal and the below-subdiagonal
512  *m11 = (*c) * (*m11) - (*s) * (*m12);
513  *m21 = (*c) * tmp - (*s) * (*m22);
514  *m22 = (*s) * tmp + (*c) * (*m22);
515 
516  // Move m11 to RQ[i+1, i+1]
517  m11 = m22;
518  c++;
519  s++;
520  }
521 
522  // Copy the below-subdiagonal to above-subdiagonal
523  RQ.diagonal(1) = RQ.diagonal(-1);
524 
525  return RQ;
526  }
527 };
528 
529 
530 
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)
Matrix matrix_RQ()
UpperHessenbergQR(ConstGenericMatrix &mat)
TridiagQR(ConstGenericMatrix &mat)
void apply_QtY(GenericMatrix Y)
void apply_QY(Vector &Y)
void apply_QY(GenericMatrix Y)
void apply_QtY(Vector &Y)