Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
GenEigsBase.h
1 // Copyright (C) 2018-2022 Yixuan Qiu <yixuan.qiu@cos.name>
2 //
3 // This Source Code Form is subject to the terms of the Mozilla
4 // Public License v. 2.0. If a copy of the MPL was not distributed
5 // with this file, You can obtain one at https://mozilla.org/MPL/2.0/.
6 
7 #ifndef SPECTRA_GEN_EIGS_BASE_H
8 #define SPECTRA_GEN_EIGS_BASE_H
9 
10 #include <Eigen/Core>
11 #include <vector> // std::vector
12 #include <cmath> // std::abs, std::pow, std::sqrt
13 #include <algorithm> // std::min, std::copy
14 #include <complex> // std::complex, std::conj, std::norm, std::abs
15 #include <stdexcept> // std::invalid_argument
16 
17 #include "Util/Version.h"
18 #include "Util/TypeTraits.h"
19 #include "Util/SelectionRule.h"
20 #include "Util/CompInfo.h"
21 #include "Util/SimpleRandom.h"
22 #include "MatOp/internal/ArnoldiOp.h"
23 #include "LinAlg/UpperHessenbergQR.h"
24 #include "LinAlg/DoubleShiftQR.h"
25 #include "LinAlg/UpperHessenbergEigen.h"
26 #include "LinAlg/Arnoldi.h"
27 
28 namespace Spectra {
29 
37 template <typename OpType, typename BOpType>
39 {
40 private:
41  using Scalar = typename OpType::Scalar;
42  using Index = Eigen::Index;
43  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
44  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
45  using Array = Eigen::Array<Scalar, Eigen::Dynamic, 1>;
46  using BoolArray = Eigen::Array<bool, Eigen::Dynamic, 1>;
47  using MapMat = Eigen::Map<Matrix>;
48  using MapVec = Eigen::Map<Vector>;
49  using MapConstVec = Eigen::Map<const Vector>;
50 
51  using Complex = std::complex<Scalar>;
52  using ComplexMatrix = Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic>;
53  using ComplexVector = Eigen::Matrix<Complex, Eigen::Dynamic, 1>;
54 
55  using ArnoldiOpType = ArnoldiOp<Scalar, OpType, BOpType>;
56  using ArnoldiFac = Arnoldi<Scalar, ArnoldiOpType>;
57 
58 protected:
59  // clang-format off
60  OpType& m_op; // object to conduct matrix operation,
61  // e.g. matrix-vector product
62  const Index m_n; // dimension of matrix A
63  const Index m_nev; // number of eigenvalues requested
64  const Index m_ncv; // dimension of Krylov subspace in the Arnoldi method
65  Index m_nmatop; // number of matrix operations called
66  Index m_niter; // number of restarting iterations
67 
68  ArnoldiFac m_fac; // Arnoldi factorization
69 
70  ComplexVector m_ritz_val; // Ritz values
71  ComplexMatrix m_ritz_vec; // Ritz vectors
72  ComplexVector m_ritz_est; // last row of m_ritz_vec, also called the Ritz estimates
73 
74 private:
75  BoolArray m_ritz_conv; // indicator of the convergence of Ritz values
76  CompInfo m_info; // status of the computation
77  // clang-format on
78 
79  // Real Ritz values calculated from UpperHessenbergEigen have exact zero imaginary part
80  // Complex Ritz values have exact conjugate pairs
81  // So we use exact tests here
82  static bool is_complex(const Complex& v) { return v.imag() != Scalar(0); }
83  static bool is_conj(const Complex& v1, const Complex& v2) { return v1 == Eigen::numext::conj(v2); }
84 
85  // Implicitly restarted Arnoldi factorization
86  void restart(Index k, SortRule selection)
87  {
88  using std::norm;
89 
90  if (k >= m_ncv)
91  return;
92 
93  DoubleShiftQR<Scalar> decomp_ds(m_ncv);
94  UpperHessenbergQR<Scalar> decomp_hb(m_ncv);
95  Matrix Q = Matrix::Identity(m_ncv, m_ncv);
96 
97  for (Index i = k; i < m_ncv; i++)
98  {
99  if (is_complex(m_ritz_val[i]) && is_conj(m_ritz_val[i], m_ritz_val[i + 1]))
100  {
101  // H - mu * I = Q1 * R1
102  // H <- R1 * Q1 + mu * I = Q1' * H * Q1
103  // H - conj(mu) * I = Q2 * R2
104  // H <- R2 * Q2 + conj(mu) * I = Q2' * H * Q2
105  //
106  // (H - mu * I) * (H - conj(mu) * I) = Q1 * Q2 * R2 * R1 = Q * R
107  const Scalar s = Scalar(2) * m_ritz_val[i].real();
108  const Scalar t = norm(m_ritz_val[i]);
109 
110  decomp_ds.compute(m_fac.matrix_H(), s, t);
111 
112  // Q -> Q * Qi
113  decomp_ds.apply_YQ(Q);
114  // H -> Q'HQ
115  // Matrix Q = Matrix::Identity(m_ncv, m_ncv);
116  // decomp_ds.apply_YQ(Q);
117  // m_fac_H = Q.transpose() * m_fac_H * Q;
118  m_fac.compress_H(decomp_ds);
119 
120  i++;
121  }
122  else
123  {
124  // QR decomposition of H - mu * I, mu is real
125  decomp_hb.compute(m_fac.matrix_H(), m_ritz_val[i].real());
126 
127  // Q -> Q * Qi
128  decomp_hb.apply_YQ(Q);
129  // H -> Q'HQ = RQ + mu * I
130  m_fac.compress_H(decomp_hb);
131  }
132  }
133 
134  m_fac.compress_V(Q);
135  m_fac.factorize_from(k, m_ncv, m_nmatop);
136 
137  retrieve_ritzpair(selection);
138  }
139 
140  // Calculates the number of converged Ritz values
141  Index num_converged(const Scalar& tol)
142  {
143  using std::pow;
144 
145  // The machine precision, ~= 1e-16 for the "double" type
146  constexpr Scalar eps = TypeTraits<Scalar>::epsilon();
147  // std::pow() is not constexpr, so we do not declare eps23 to be constexpr
148  // But most compilers should be able to compute eps23 at compile time
149  const Scalar eps23 = pow(eps, Scalar(2) / 3);
150 
151  // thresh = tol * max(eps23, abs(theta)), theta for Ritz value
152  Array thresh = tol * m_ritz_val.head(m_nev).array().abs().max(eps23);
153  Array resid = m_ritz_est.head(m_nev).array().abs() * m_fac.f_norm();
154  // Converged "wanted" Ritz values
155  m_ritz_conv = (resid < thresh);
156 
157  return m_ritz_conv.count();
158  }
159 
160  // Returns the adjusted nev for restarting
161  Index nev_adjusted(Index nconv)
162  {
163  using std::abs;
164 
165  // A very small value, but 1.0 / near_0 does not overflow
166  // ~= 1e-307 for the "double" type
167  constexpr Scalar near_0 = TypeTraits<Scalar>::min() * Scalar(10);
168 
169  Index nev_new = m_nev;
170  for (Index i = m_nev; i < m_ncv; i++)
171  if (abs(m_ritz_est[i]) < near_0)
172  nev_new++;
173 
174  // Adjust nev_new, according to dnaup2.f line 660~674 in ARPACK
175  nev_new += (std::min)(nconv, (m_ncv - nev_new) / 2);
176  if (nev_new == 1 && m_ncv >= 6)
177  nev_new = m_ncv / 2;
178  else if (nev_new == 1 && m_ncv > 3)
179  nev_new = 2;
180 
181  if (nev_new > m_ncv - 2)
182  nev_new = m_ncv - 2;
183 
184  // Increase nev by one if ritz_val[nev - 1] and
185  // ritz_val[nev] are conjugate pairs
186  if (is_complex(m_ritz_val[nev_new - 1]) &&
187  is_conj(m_ritz_val[nev_new - 1], m_ritz_val[nev_new]))
188  {
189  nev_new++;
190  }
191 
192  return nev_new;
193  }
194 
195  // Retrieves and sorts Ritz values and Ritz vectors
196  void retrieve_ritzpair(SortRule selection)
197  {
198  UpperHessenbergEigen<Scalar> decomp(m_fac.matrix_H());
199  const ComplexVector& evals = decomp.eigenvalues();
200  ComplexMatrix evecs = decomp.eigenvectors();
201 
202  // Sort Ritz values and put the wanted ones at the beginning
203  std::vector<Index> ind;
204  switch (selection)
205  {
207  {
208  SortEigenvalue<Complex, SortRule::LargestMagn> sorting(evals.data(), m_ncv);
209  sorting.swap(ind);
210  break;
211  }
213  {
214  SortEigenvalue<Complex, SortRule::LargestReal> sorting(evals.data(), m_ncv);
215  sorting.swap(ind);
216  break;
217  }
219  {
220  SortEigenvalue<Complex, SortRule::LargestImag> sorting(evals.data(), m_ncv);
221  sorting.swap(ind);
222  break;
223  }
225  {
226  SortEigenvalue<Complex, SortRule::SmallestMagn> sorting(evals.data(), m_ncv);
227  sorting.swap(ind);
228  break;
229  }
231  {
232  SortEigenvalue<Complex, SortRule::SmallestReal> sorting(evals.data(), m_ncv);
233  sorting.swap(ind);
234  break;
235  }
237  {
238  SortEigenvalue<Complex, SortRule::SmallestImag> sorting(evals.data(), m_ncv);
239  sorting.swap(ind);
240  break;
241  }
242  default:
243  throw std::invalid_argument("unsupported selection rule");
244  }
245 
246  // Copy the Ritz values and vectors to m_ritz_val and m_ritz_vec, respectively
247  for (Index i = 0; i < m_ncv; i++)
248  {
249  m_ritz_val[i] = evals[ind[i]];
250  m_ritz_est[i] = evecs(m_ncv - 1, ind[i]);
251  }
252  for (Index i = 0; i < m_nev; i++)
253  {
254  m_ritz_vec.col(i).noalias() = evecs.col(ind[i]);
255  }
256  }
257 
258 protected:
259  // Sorts the first nev Ritz pairs in the specified order
260  // This is used to return the final results
261  virtual void sort_ritzpair(SortRule sort_rule)
262  {
263  std::vector<Index> ind;
264  switch (sort_rule)
265  {
267  {
268  SortEigenvalue<Complex, SortRule::LargestMagn> sorting(m_ritz_val.data(), m_nev);
269  sorting.swap(ind);
270  break;
271  }
273  {
274  SortEigenvalue<Complex, SortRule::LargestReal> sorting(m_ritz_val.data(), m_nev);
275  sorting.swap(ind);
276  break;
277  }
279  {
280  SortEigenvalue<Complex, SortRule::LargestImag> sorting(m_ritz_val.data(), m_nev);
281  sorting.swap(ind);
282  break;
283  }
285  {
286  SortEigenvalue<Complex, SortRule::SmallestMagn> sorting(m_ritz_val.data(), m_nev);
287  sorting.swap(ind);
288  break;
289  }
291  {
292  SortEigenvalue<Complex, SortRule::SmallestReal> sorting(m_ritz_val.data(), m_nev);
293  sorting.swap(ind);
294  break;
295  }
297  {
298  SortEigenvalue<Complex, SortRule::SmallestImag> sorting(m_ritz_val.data(), m_nev);
299  sorting.swap(ind);
300  break;
301  }
302  default:
303  throw std::invalid_argument("unsupported sorting rule");
304  }
305 
306  ComplexVector new_ritz_val(m_ncv);
307  ComplexMatrix new_ritz_vec(m_ncv, m_nev);
308  BoolArray new_ritz_conv(m_nev);
309 
310  for (Index i = 0; i < m_nev; i++)
311  {
312  new_ritz_val[i] = m_ritz_val[ind[i]];
313  new_ritz_vec.col(i).noalias() = m_ritz_vec.col(ind[i]);
314  new_ritz_conv[i] = m_ritz_conv[ind[i]];
315  }
316 
317  m_ritz_val.swap(new_ritz_val);
318  m_ritz_vec.swap(new_ritz_vec);
319  m_ritz_conv.swap(new_ritz_conv);
320  }
321 
322 public:
324 
325  GenEigsBase(OpType& op, const BOpType& Bop, Index nev, Index ncv) :
326  m_op(op),
327  m_n(m_op.rows()),
328  m_nev(nev),
329  m_ncv(ncv > m_n ? m_n : ncv),
330  m_nmatop(0),
331  m_niter(0),
332  m_fac(ArnoldiOpType(op, Bop), m_ncv),
333  m_info(CompInfo::NotComputed)
334  {
335  if (nev < 1 || nev > m_n - 2)
336  throw std::invalid_argument("nev must satisfy 1 <= nev <= n - 2, n is the size of matrix");
337 
338  if (ncv < nev + 2 || ncv > m_n)
339  throw std::invalid_argument("ncv must satisfy nev + 2 <= ncv <= n, n is the size of matrix");
340  }
341 
345  virtual ~GenEigsBase() {}
346 
348 
358  void init(const Scalar* init_resid)
359  {
360  // Reset all matrices/vectors to zero
361  m_ritz_val.resize(m_ncv);
362  m_ritz_vec.resize(m_ncv, m_nev);
363  m_ritz_est.resize(m_ncv);
364  m_ritz_conv.resize(m_nev);
365 
366  m_ritz_val.setZero();
367  m_ritz_vec.setZero();
368  m_ritz_est.setZero();
369  m_ritz_conv.setZero();
370 
371  m_nmatop = 0;
372  m_niter = 0;
373 
374  // Initialize the Arnoldi factorization
375  MapConstVec v0(init_resid, m_n);
376  m_fac.init(v0, m_nmatop);
377  }
378 
386  void init()
387  {
388  SimpleRandom<Scalar> rng(0);
389  Vector init_resid = rng.random_vec(m_n);
390  init(init_resid.data());
391  }
392 
417  Index compute(SortRule selection = SortRule::LargestMagn, Index maxit = 1000,
418  Scalar tol = 1e-10, SortRule sorting = SortRule::LargestMagn)
419  {
420  // The m-step Arnoldi factorization
421  m_fac.factorize_from(1, m_ncv, m_nmatop);
422  retrieve_ritzpair(selection);
423  // Restarting
424  Index i, nconv = 0, nev_adj;
425  for (i = 0; i < maxit; i++)
426  {
427  nconv = num_converged(tol);
428  if (nconv >= m_nev)
429  break;
430 
431  nev_adj = nev_adjusted(nconv);
432  restart(nev_adj, selection);
433  }
434  // Sorting results
435  sort_ritzpair(sorting);
436 
437  m_niter += i + 1;
438  m_info = (nconv >= m_nev) ? CompInfo::Successful : CompInfo::NotConverging;
439 
440  return (std::min)(m_nev, nconv);
441  }
442 
447  CompInfo info() const { return m_info; }
448 
452  Index num_iterations() const { return m_niter; }
453 
457  Index num_operations() const { return m_nmatop; }
458 
466  ComplexVector eigenvalues() const
467  {
468  const Index nconv = m_ritz_conv.cast<Index>().sum();
469  ComplexVector res(nconv);
470 
471  if (!nconv)
472  return res;
473 
474  Index j = 0;
475  for (Index i = 0; i < m_nev; i++)
476  {
477  if (m_ritz_conv[i])
478  {
479  res[j] = m_ritz_val[i];
480  j++;
481  }
482  }
483 
484  return res;
485  }
486 
496  ComplexMatrix eigenvectors(Index nvec) const
497  {
498  const Index nconv = m_ritz_conv.cast<Index>().sum();
499  nvec = (std::min)(nvec, nconv);
500  ComplexMatrix res(m_n, nvec);
501 
502  if (!nvec)
503  return res;
504 
505  ComplexMatrix ritz_vec_conv(m_ncv, nvec);
506  Index j = 0;
507  for (Index i = 0; i < m_nev && j < nvec; i++)
508  {
509  if (m_ritz_conv[i])
510  {
511  ritz_vec_conv.col(j).noalias() = m_ritz_vec.col(i);
512  j++;
513  }
514  }
515 
516  res.noalias() = m_fac.matrix_V() * ritz_vec_conv;
517 
518  return res;
519  }
520 
524  ComplexMatrix eigenvectors() const
525  {
526  return eigenvectors(m_nev);
527  }
528 };
529 
530 } // namespace Spectra
531 
532 #endif // SPECTRA_GEN_EIGS_BASE_H
ComplexVector eigenvalues() const
Definition: GenEigsBase.h:466
Index num_operations() const
Definition: GenEigsBase.h:457
ComplexMatrix eigenvectors() const
Definition: GenEigsBase.h:524
ComplexMatrix eigenvectors(Index nvec) const
Definition: GenEigsBase.h:496
Index compute(SortRule selection=SortRule::LargestMagn, Index maxit=1000, Scalar tol=1e-10, SortRule sorting=SortRule::LargestMagn)
Definition: GenEigsBase.h:417
Index num_iterations() const
Definition: GenEigsBase.h:452
CompInfo info() const
Definition: GenEigsBase.h:447
void init(const Scalar *init_resid)
Definition: GenEigsBase.h:358
CompInfo
Definition: CompInfo.h:18
@ SmallestImag
Select eigenvalues with smallest imaginary part (in magnitude). Only for general eigen solvers.
@ SmallestReal
Select eigenvalues with smallest real part. Only for general eigen solvers.
@ LargestImag
Select eigenvalues with largest imaginary part (in magnitude). Only for general eigen solvers.
@ LargestReal
Select eigenvalues with largest real part. Only for general eigen solvers.
@ Successful
Computation was successful.