Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
SymEigsBase.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_SYM_EIGS_BASE_H
8 #define SPECTRA_SYM_EIGS_BASE_H
9 
10 #include <Eigen/Core>
11 #include <vector> // std::vector
12 #include <cmath> // std::abs, std::pow
13 #include <algorithm> // std::min
14 #include <stdexcept> // std::invalid_argument
15 #include <utility> // std::move
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/TridiagEigen.h"
25 #include "LinAlg/Lanczos.h"
26 
27 namespace Spectra {
28 
34 
42 template <typename OpType, typename BOpType>
44 {
45 private:
46  using Scalar = typename OpType::Scalar;
47  using Index = Eigen::Index;
48  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
49  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
50  using Array = Eigen::Array<Scalar, Eigen::Dynamic, 1>;
51  using BoolArray = Eigen::Array<bool, Eigen::Dynamic, 1>;
52  using MapMat = Eigen::Map<Matrix>;
53  using MapVec = Eigen::Map<Vector>;
54  using MapConstVec = Eigen::Map<const Vector>;
55 
56  using ArnoldiOpType = ArnoldiOp<Scalar, OpType, BOpType>;
57  using LanczosFac = Lanczos<Scalar, ArnoldiOpType>;
58 
59 protected:
60  // clang-format off
61 
62  // In SymEigsSolver and SymEigsShiftSolver, the A operator is an lvalue provided by
63  // the user. In SymGEigsSolver, the A operator is an rvalue. To avoid copying objects,
64  // we use the following scheme:
65  // 1. If the op parameter in the constructor is an lvalue, make m_op a const reference to op
66  // 2. If op is an rvalue, move op to m_op_container, and then make m_op a const
67  // reference to m_op_container[0]
68  std::vector<OpType> m_op_container;
69  const OpType& m_op; // matrix operator for A
70  const Index m_n; // dimension of matrix A
71  const Index m_nev; // number of eigenvalues requested
72  const Index m_ncv; // dimension of Krylov subspace in the Lanczos method
73  Index m_nmatop; // number of matrix operations called
74  Index m_niter; // number of restarting iterations
75 
76  LanczosFac m_fac; // Lanczos factorization
77  Vector m_ritz_val; // Ritz values
78 
79 private:
80  Matrix m_ritz_vec; // Ritz vectors
81  Vector m_ritz_est; // last row of m_ritz_vec, also called the Ritz estimates
82  BoolArray m_ritz_conv; // indicator of the convergence of Ritz values
83  CompInfo m_info; // status of the computation
84  // clang-format on
85 
86  // Move rvalue object to the container
87  static std::vector<OpType> create_op_container(OpType&& rval)
88  {
89  std::vector<OpType> container;
90  container.emplace_back(std::move(rval));
91  return container;
92  }
93 
94  // Implicitly restarted Lanczos factorization
95  void restart(Index k, SortRule selection)
96  {
97  using std::abs;
98 
99  if (k >= m_ncv)
100  return;
101 
102  TridiagQR<Scalar> decomp(m_ncv);
103  Matrix Q = Matrix::Identity(m_ncv, m_ncv);
104 
105  // Apply large shifts first
106  const int nshift = m_ncv - k;
107  Vector shifts = m_ritz_val.tail(nshift);
108  std::sort(shifts.data(), shifts.data() + nshift, [](const Scalar& v1, const Scalar& v2) { return abs(v1) > abs(v2); });
109 
110  for (Index i = 0; i < nshift; i++)
111  {
112  // QR decomposition of H-mu*I, mu is the shift
113  decomp.compute(m_fac.matrix_H(), shifts[i]);
114 
115  // Q -> Q * Qi
116  decomp.apply_YQ(Q);
117  // H -> Q'HQ
118  // Since QR = H - mu * I, we have H = QR + mu * I
119  // and therefore Q'HQ = RQ + mu * I
120  m_fac.compress_H(decomp);
121  }
122 
123  m_fac.compress_V(Q);
124  m_fac.factorize_from(k, m_ncv, m_nmatop);
125 
126  retrieve_ritzpair(selection);
127  }
128 
129  // Calculates the number of converged Ritz values
130  Index num_converged(const Scalar& tol)
131  {
132  using std::pow;
133 
134  // The machine precision, ~= 1e-16 for the "double" type
135  constexpr Scalar eps = TypeTraits<Scalar>::epsilon();
136  // std::pow() is not constexpr, so we do not declare eps23 to be constexpr
137  // But most compilers should be able to compute eps23 at compile time
138  const Scalar eps23 = pow(eps, Scalar(2) / 3);
139 
140  // thresh = tol * max(eps23, abs(theta)), theta for Ritz value
141  Array thresh = tol * m_ritz_val.head(m_nev).array().abs().max(eps23);
142  Array resid = m_ritz_est.head(m_nev).array().abs() * m_fac.f_norm();
143  // Converged "wanted" Ritz values
144  m_ritz_conv = (resid < thresh);
145 
146  // std::cout << "================ num_converged() ================" << std::endl;
147  // std::cout << "resid = " << resid.matrix().transpose() << std::endl;
148  // std::cout << "thresh = " << thresh.matrix().transpose() << std::endl;
149  // std::cout << "================ num_converged() ================" << std::endl
150  // << std::endl
151  // << std::endl
152  // << std::endl;
153 
154  return m_ritz_conv.count();
155  }
156 
157  // Returns the adjusted nev for restarting
158  Index nev_adjusted(Index nconv)
159  {
160  using std::abs;
161 
162  // A very small value, but 1.0 / near_0 does not overflow
163  // ~= 1e-307 for the "double" type
164  constexpr Scalar near_0 = TypeTraits<Scalar>::min() * Scalar(10);
165 
166  Index nev_new = m_nev;
167  for (Index i = m_nev; i < m_ncv; i++)
168  if (abs(m_ritz_est[i]) < near_0)
169  nev_new++;
170 
171  // Adjust nev_new, according to dsaup2.f line 677~684 in ARPACK
172  nev_new += (std::min)(nconv, (m_ncv - nev_new) / 2);
173  if (nev_new == 1 && m_ncv >= 6)
174  nev_new = m_ncv / 2;
175  else if (nev_new == 1 && m_ncv > 2)
176  nev_new = 2;
177 
178  if (nev_new > m_ncv - 1)
179  nev_new = m_ncv - 1;
180 
181  return nev_new;
182  }
183 
184  // Retrieves and sorts Ritz values and Ritz vectors
185  void retrieve_ritzpair(SortRule selection)
186  {
187  TridiagEigen<Scalar> decomp(m_fac.matrix_H());
188  const Vector& evals = decomp.eigenvalues();
189  const Matrix& evecs = decomp.eigenvectors();
190 
191  // std::cout << "================ retrieve_ritzpair() ================" << std::endl;
192  // std::cout << "H = \n"
193  // << m_fac.matrix_H().topLeftCorner(5, 5) << std::endl;
194  // std::cout << "evals = " << evals.transpose() << std::endl;
195  // std::cout << "================ retrieve_ritzpair() ================" << std::endl
196  // << std::endl;
197 
198  // Sort Ritz values and put the wanted ones at the beginning
199  std::vector<Index> ind = argsort(selection, evals, m_ncv);
200 
201  // Copy the Ritz values and vectors to m_ritz_val and m_ritz_vec, respectively
202  for (Index i = 0; i < m_ncv; i++)
203  {
204  m_ritz_val[i] = evals[ind[i]];
205  m_ritz_est[i] = evecs(m_ncv - 1, ind[i]);
206  }
207  for (Index i = 0; i < m_nev; i++)
208  {
209  m_ritz_vec.col(i).noalias() = evecs.col(ind[i]);
210  }
211  }
212 
213 protected:
214  // Sorts the first nev Ritz pairs in the specified order
215  // This is used to return the final results
216  virtual void sort_ritzpair(SortRule sort_rule)
217  {
218  if ((sort_rule != SortRule::LargestAlge) && (sort_rule != SortRule::LargestMagn) &&
219  (sort_rule != SortRule::SmallestAlge) && (sort_rule != SortRule::SmallestMagn))
220  throw std::invalid_argument("unsupported sorting rule");
221 
222  std::vector<Index> ind = argsort(sort_rule, m_ritz_val, m_nev);
223 
224  Vector new_ritz_val(m_ncv);
225  Matrix new_ritz_vec(m_ncv, m_nev);
226  BoolArray new_ritz_conv(m_nev);
227 
228  for (Index i = 0; i < m_nev; i++)
229  {
230  new_ritz_val[i] = m_ritz_val[ind[i]];
231  new_ritz_vec.col(i).noalias() = m_ritz_vec.col(ind[i]);
232  new_ritz_conv[i] = m_ritz_conv[ind[i]];
233  }
234 
235  m_ritz_val.swap(new_ritz_val);
236  m_ritz_vec.swap(new_ritz_vec);
237  m_ritz_conv.swap(new_ritz_conv);
238  }
239 
240 public:
242 
243  // If op is an lvalue
244  SymEigsBase(OpType& op, const BOpType& Bop, Index nev, Index ncv) :
245  m_op(op),
246  m_n(op.rows()),
247  m_nev(nev),
248  m_ncv(ncv > m_n ? m_n : ncv),
249  m_nmatop(0),
250  m_niter(0),
251  m_fac(ArnoldiOpType(op, Bop), m_ncv),
252  m_info(CompInfo::NotComputed)
253  {
254  if (nev < 1 || nev > m_n - 1)
255  throw std::invalid_argument("nev must satisfy 1 <= nev <= n - 1, n is the size of matrix");
256 
257  if (ncv <= nev || ncv > m_n)
258  throw std::invalid_argument("ncv must satisfy nev < ncv <= n, n is the size of matrix");
259  }
260 
261  // If op is an rvalue
262  SymEigsBase(OpType&& op, const BOpType& Bop, Index nev, Index ncv) :
263  m_op_container(create_op_container(std::move(op))),
264  m_op(m_op_container.front()),
265  m_n(m_op.rows()),
266  m_nev(nev),
267  m_ncv(ncv > m_n ? m_n : ncv),
268  m_nmatop(0),
269  m_niter(0),
270  m_fac(ArnoldiOpType(m_op, Bop), m_ncv),
271  m_info(CompInfo::NotComputed)
272  {
273  if (nev < 1 || nev > m_n - 1)
274  throw std::invalid_argument("nev must satisfy 1 <= nev <= n - 1, n is the size of matrix");
275 
276  if (ncv <= nev || ncv > m_n)
277  throw std::invalid_argument("ncv must satisfy nev < ncv <= n, n is the size of matrix");
278  }
279 
283  virtual ~SymEigsBase() {}
284 
286 
296  void init(const Scalar* init_resid)
297  {
298  // Reset all matrices/vectors to zero
299  m_ritz_val.resize(m_ncv);
300  m_ritz_vec.resize(m_ncv, m_nev);
301  m_ritz_est.resize(m_ncv);
302  m_ritz_conv.resize(m_nev);
303 
304  m_ritz_val.setZero();
305  m_ritz_vec.setZero();
306  m_ritz_est.setZero();
307  m_ritz_conv.setZero();
308 
309  m_nmatop = 0;
310  m_niter = 0;
311 
312  // Initialize the Lanczos factorization
313  MapConstVec v0(init_resid, m_n);
314  m_fac.init(v0, m_nmatop);
315  }
316 
324  void init()
325  {
326  SimpleRandom<Scalar> rng(0);
327  Vector init_resid = rng.random_vec(m_n);
328  init(init_resid.data());
329 
330  std::cout << "init_vec = " << init_resid.head(5).transpose() << std::endl
331  << std::endl;
332  }
333 
356  Index compute(SortRule selection = SortRule::LargestMagn, Index maxit = 1000,
357  Scalar tol = 1e-10, SortRule sorting = SortRule::LargestAlge)
358  {
359  // The m-step Lanczos factorization
360  m_fac.factorize_from(1, m_ncv, m_nmatop);
361  retrieve_ritzpair(selection);
362  // Restarting
363  Index i, nconv = 0, nev_adj;
364  for (i = 0; i < maxit; i++)
365  {
366  nconv = num_converged(tol);
367  if (nconv >= m_nev)
368  break;
369 
370  nev_adj = nev_adjusted(nconv);
371  restart(nev_adj, selection);
372  }
373  // Sorting results
374  sort_ritzpair(sorting);
375 
376  m_niter += i + 1;
377  m_info = (nconv >= m_nev) ? CompInfo::Successful : CompInfo::NotConverging;
378 
379  return (std::min)(m_nev, nconv);
380  }
381 
386  CompInfo info() const { return m_info; }
387 
391  Index num_iterations() const { return m_niter; }
392 
396  Index num_operations() const { return m_nmatop; }
397 
405  Vector eigenvalues() const
406  {
407  const Index nconv = m_ritz_conv.count();
408  Vector res(nconv);
409 
410  if (!nconv)
411  return res;
412 
413  Index j = 0;
414  for (Index i = 0; i < m_nev; i++)
415  {
416  if (m_ritz_conv[i])
417  {
418  res[j] = m_ritz_val[i];
419  j++;
420  }
421  }
422 
423  return res;
424  }
425 
435  virtual Matrix eigenvectors(Index nvec) const
436  {
437  const Index nconv = m_ritz_conv.count();
438  nvec = (std::min)(nvec, nconv);
439  Matrix res(m_n, nvec);
440 
441  if (!nvec)
442  return res;
443 
444  Matrix ritz_vec_conv(m_ncv, nvec);
445  Index j = 0;
446  for (Index i = 0; i < m_nev && j < nvec; i++)
447  {
448  if (m_ritz_conv[i])
449  {
450  ritz_vec_conv.col(j).noalias() = m_ritz_vec.col(i);
451  j++;
452  }
453  }
454 
455  res.noalias() = m_fac.matrix_V() * ritz_vec_conv;
456 
457  return res;
458  }
459 
463  virtual Matrix eigenvectors() const
464  {
465  return eigenvectors(m_nev);
466  }
467 };
468 
469 } // namespace Spectra
470 
471 #endif // SPECTRA_SYM_EIGS_BASE_H
Index num_operations() const
Definition: SymEigsBase.h:396
Index num_iterations() const
Definition: SymEigsBase.h:391
void init(const Scalar *init_resid)
Definition: SymEigsBase.h:296
virtual Matrix eigenvectors(Index nvec) const
Definition: SymEigsBase.h:435
CompInfo info() const
Definition: SymEigsBase.h:386
Index compute(SortRule selection=SortRule::LargestMagn, Index maxit=1000, Scalar tol=1e-10, SortRule sorting=SortRule::LargestAlge)
Definition: SymEigsBase.h:356
Vector eigenvalues() const
Definition: SymEigsBase.h:405
virtual Matrix eigenvectors() const
Definition: SymEigsBase.h:463
CompInfo
Definition: CompInfo.h:18
@ SmallestAlge
Select eigenvalues with smallest algebraic value. Only for symmetric eigen solvers.
@ Successful
Computation was successful.