LOBPCGSolver.h
Go to the documentation of this file.
1 // Written by Anna Araslanova
2 // Modified by Yixuan Qiu
3 // License: MIT
4 
5 #ifndef LOBPCG_SOLVER
6 #define LOBPCG_SOLVER
7 
8 #include <functional>
9 #include <map>
10 
11 #include <Eigen/Core>
12 #include <Eigen/SparseCore>
13 #include <Eigen/Eigenvalues>
14 #include <Eigen/SVD>
15 #include <Eigen/SparseCholesky>
16 
17 #include "../SymGEigsSolver.h"
18 
19 namespace Spectra {
20 
24 
81 
82 template <typename Scalar = long double>
84 {
85 private:
88 
89  typedef std::complex<Scalar> Complex;
92 
95 
96  const int m_n; // dimension of matrix A
97  const int m_nev; // number of eigenvalues requested
101 
102 public:
106  int m_info;
107 
108 private:
109  // B-orthonormalize matrix M
111  SparseMatrix& true_BM, bool has_true_BM = false)
112  {
113  SparseMatrix BM;
114 
115  if (has_true_BM == false)
116  {
117  if (flag_with_B)
118  {
119  BM = B * M;
120  }
121  else
122  {
123  BM = M;
124  }
125  }
126  else
127  {
128  BM = true_BM;
129  }
130 
131  Eigen::SimplicialLDLT<SparseMatrix> chol_MBM(M.transpose() * BM);
132 
133  if (chol_MBM.info() != SUCCESSFUL)
134  {
135  // LDLT decomposition fail
136  m_info = chol_MBM.info();
137  return chol_MBM.info();
138  }
139 
140  SparseComplexMatrix Upper_MBM = chol_MBM.matrixU().template cast<Complex>();
141  ComplexVector D_MBM_vec = chol_MBM.vectorD().template cast<Complex>();
142 
143  D_MBM_vec = D_MBM_vec.cwiseSqrt();
144 
145  for (int i = 0; i < D_MBM_vec.rows(); i++)
146  {
147  D_MBM_vec(i) = Complex(1.0, 0.0) / D_MBM_vec(i);
148  }
149 
150  SparseComplexMatrix D_MBM_mat(D_MBM_vec.asDiagonal());
151 
152  SparseComplexMatrix U_inv(Upper_MBM.rows(), Upper_MBM.cols());
153  U_inv.setIdentity();
154  Upper_MBM.template triangularView<Eigen::Upper>().solveInPlace(U_inv);
155 
156  SparseComplexMatrix right_product = U_inv * D_MBM_mat;
157  M = M * right_product.real();
158  if (flag_with_B)
159  {
160  true_BM = B * M;
161  }
162  else
163  {
164  true_BM = M;
165  }
166 
167  return SUCCESSFUL;
168  }
169 
171  SparseMatrix& B)
172  {
173  SparseMatrix BY;
174  if (flag_with_B)
175  {
176  BY = B * Y;
177  }
178  else
179  {
180  BY = Y;
181  }
182 
183  SparseMatrix YBY = Y.transpose() * BY;
184  SparseMatrix BYX = BY.transpose() * X;
185 
186  SparseMatrix YBY_XYX = (Matrix(YBY).bdcSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(Matrix(BYX))).sparseView();
187  X = X - Y * YBY_XYX;
188  }
189 
190  /*
191  return
192  'AB
193  CD'
194  */
196  Matrix C, Matrix D)
197  {
198  Matrix result(A.rows() + C.rows(), A.cols() + B.cols());
199  result.topLeftCorner(A.rows(), A.cols()) = A;
200  result.topRightCorner(B.rows(), B.cols()) = B;
201  result.bottomLeftCorner(C.rows(), C.cols()) = C;
202  result.bottomRightCorner(D.rows(), D.cols()) = D;
203  return result;
204  }
205 
207  Matrix D, Matrix E, Matrix F,
208  Matrix G, Matrix H, Matrix I)
209  {
210  Matrix result(A.rows() + D.rows() + G.rows(), A.cols() + B.cols() + C.cols());
211  result.block(0, 0, A.rows(), A.cols()) = A;
212  result.block(0, A.cols(), B.rows(), B.cols()) = B;
213  result.block(0, A.cols() + B.cols(), C.rows(), C.cols()) = C;
214  result.block(A.rows(), 0, D.rows(), D.cols()) = D;
215  result.block(A.rows(), A.cols(), E.rows(), E.cols()) = E;
216  result.block(A.rows(), A.cols() + B.cols(), F.rows(), F.cols()) = F;
217  result.block(A.rows() + D.rows(), 0, G.rows(), G.cols()) = G;
218  result.block(A.rows() + D.rows(), A.cols(), H.rows(), H.cols()) = H;
219  result.block(A.rows() + D.rows(), A.cols() + B.cols(), I.rows(), I.cols()) = I;
220 
221  return result;
222  }
223 
224  void sort_epairs(Vector& evalues, Matrix& evectors, int SelectionRule)
225  {
226  std::function<bool(Scalar, Scalar)> cmp;
227  if (SelectionRule == SMALLEST_ALGE)
228  cmp = std::less<Scalar>{};
229  else
230  cmp = std::greater<Scalar>{};
231 
232  std::map<Scalar, Vector, decltype(cmp)> epairs(cmp);
233  for (int i = 0; i < m_evectors.cols(); ++i)
234  epairs.insert(std::make_pair(evalues(i), evectors.col(i)));
235 
236  int i = 0;
237  for (auto& epair : epairs)
238  {
239  evectors.col(i) = epair.second;
240  evalues(i) = epair.first;
241  i++;
242  }
243  }
244 
245  void removeColumns(SparseMatrix& matrix, std::vector<int>& colToRemove)
246  {
247  // remove columns through matrix multiplication
248  SparseMatrix new_matrix(matrix.cols(), matrix.cols() - int(colToRemove.size()));
249  int iCol = 0;
250  std::vector<Eigen::Triplet<Scalar>> tripletList;
251  tripletList.reserve(matrix.cols() - int(colToRemove.size()));
252 
253  for (int iRow = 0; iRow < matrix.cols(); iRow++)
254  {
255  if (std::find(colToRemove.begin(), colToRemove.end(), iRow) == colToRemove.end())
256  {
257  tripletList.push_back(Eigen::Triplet<Scalar>(iRow, iCol, 1));
258  iCol++;
259  }
260  }
261 
262  new_matrix.setFromTriplets(tripletList.begin(), tripletList.end());
263  matrix = matrix * new_matrix;
264  }
265 
266  int checkConvergence_getBlocksize(SparseMatrix& m_residuals, Scalar tolerance_L2, std::vector<int>& columnsToDelete)
267  {
268  // square roots from sum of squares by column
269  int BlockSize = m_nev;
270  Scalar sum, buffer;
271 
272  for (int iCol = 0; iCol < m_nev; iCol++)
273  {
274  sum = 0;
275  for (int iRow = 0; iRow < m_n; iRow++)
276  {
277  buffer = m_residuals.coeff(iRow, iCol);
278  sum += buffer * buffer;
279  }
280 
281  if (sqrt(sum) < tolerance_L2)
282  {
283  BlockSize--;
284  columnsToDelete.push_back(iCol);
285  }
286  }
287  return BlockSize;
288  }
289 
290 public:
292  m_n(A.rows()),
293  m_nev(X.cols()),
295  flag_with_constraints(false),
296  flag_with_B(false),
298  A(A),
299  X(X)
300  {
301  if (A.rows() != X.rows() || A.rows() != A.cols())
302  throw std::invalid_argument("Wrong size");
303 
304  //if (m_n < 5* m_nev)
305  // throw std::invalid_argument("The problem size is small compared to the block size. Use standard eigensolver");
306  }
307 
309  {
310  m_Y = Y;
311  flag_with_constraints = true;
312  }
313 
314  void setB(const SparseMatrix& B)
315  {
316  if (B.rows() != A.rows() || B.cols() != A.cols())
317  throw std::invalid_argument("Wrong size");
318  m_B = B;
319  flag_with_B = true;
320  }
321 
322  void setPreconditioner(const SparseMatrix& preconditioner)
323  {
324  m_preconditioner = preconditioner;
326  }
327 
328  void compute(int maxit = 10, Scalar tol_div_n = 1e-7)
329  {
330  Scalar tolerance_L2 = tol_div_n * m_n;
331  int BlockSize;
332  int max_iter = std::min(m_n, maxit);
333 
334  SparseMatrix directions, AX, AR, BX, AD, ADD, DD, BDD, BD, XAD, RAD, DAD, XBD, RBD, BR, sparse_eVecX, sparse_eVecR, sparse_eVecD, inverse_matrix;
335  Matrix XAR, RAR, XBR, gramA, gramB, eVecX, eVecR, eVecD;
336  std::vector<int> columnsToDelete;
337 
339  {
340  // Apply the constraints Y to X
342  }
343 
344  // Make initial vectors orthonormal
345  // implicit BX declaration
346  if (orthogonalizeInPlace(X, m_B, BX) != SUCCESSFUL)
347  {
348  max_iter = 0;
349  }
350 
351  AX = A * X;
352  // Solve the following NxN eigenvalue problem for all N eigenvalues and -vectors:
353  // first approximation via a dense problem
354  Eigen::EigenSolver<Matrix> eigs(Matrix(X.transpose() * AX));
355 
356  if (eigs.info() != SUCCESSFUL)
357  {
358  m_info = eigs.info();
359  max_iter = 0;
360  }
361  else
362  {
363  m_evalues = eigs.eigenvalues().real();
364  m_evectors = eigs.eigenvectors().real();
366  sparse_eVecX = m_evectors.sparseView();
367 
368  X = X * sparse_eVecX;
369  AX = AX * sparse_eVecX;
370  BX = BX * sparse_eVecX;
371  }
372 
373  for (int iter_num = 0; iter_num < max_iter; iter_num++)
374  {
376  for (int i = 0; i < m_nev; i++)
377  {
378  m_residuals.col(i) = AX.col(i) - m_evalues(i) * BX.col(i);
379  }
380  BlockSize = checkConvergence_getBlocksize(m_residuals, tolerance_L2, columnsToDelete);
381 
382  if (BlockSize == 0)
383  {
384  m_info = SUCCESSFUL;
385  break;
386  }
387 
388  // substitution of the original active mask
389  if (columnsToDelete.size() > 0)
390  {
391  removeColumns(m_residuals, columnsToDelete);
392  if (iter_num > 0)
393  {
394  removeColumns(directions, columnsToDelete);
395  removeColumns(AD, columnsToDelete);
396  removeColumns(BD, columnsToDelete);
397  }
398  columnsToDelete.clear(); // for next iteration
399  }
400 
402  {
403  // Apply the preconditioner to the residuals
405  }
406 
408  {
409  // Apply the constraints Y to residuals
411  }
412 
414  {
415  break;
416  }
417  AR = A * m_residuals;
418 
419  // Orthonormalize conjugate directions
420  if (iter_num > 0)
421  {
422  if (orthogonalizeInPlace(directions, m_B, BD, true) != SUCCESSFUL)
423  {
424  break;
425  }
426  AD = A * directions;
427  }
428 
429  // Perform the Rayleigh Ritz Procedure
430  XAR = Matrix(X.transpose() * AR);
431  RAR = Matrix(m_residuals.transpose() * AR);
432  XBR = Matrix(X.transpose() * BR);
433 
434  if (iter_num > 0)
435  {
436  XAD = X.transpose() * AD;
437  RAD = m_residuals.transpose() * AD;
438  DAD = directions.transpose() * AD;
439  XBD = X.transpose() * BD;
440  RBD = m_residuals.transpose() * BD;
441 
442  gramA = stack_9_matricies(m_evalues.asDiagonal(), XAR, XAD, XAR.transpose(), RAR, RAD, XAD.transpose(), RAD.transpose(), DAD.transpose());
443  gramB = stack_9_matricies(Matrix::Identity(m_nev, m_nev), XBR, XBD, XBR.transpose(), Matrix::Identity(BlockSize, BlockSize), RBD, XBD.transpose(), RBD.transpose(), Matrix::Identity(BlockSize, BlockSize));
444  }
445  else
446  {
447  gramA = stack_4_matricies(m_evalues.asDiagonal(), XAR, XAR.transpose(), RAR);
448  gramB = stack_4_matricies(Matrix::Identity(m_nev, m_nev), XBR, XBR.transpose(), Matrix::Identity(BlockSize, BlockSize));
449  }
450 
451  //calculate the lowest/largest m eigenpairs; Solve the generalized eigenvalue problem.
452  DenseSymMatProd<Scalar> Aop(gramA);
453  DenseCholesky<Scalar> Bop(gramB);
454 
457  geigs(&Aop, &Bop, m_nev, std::min(10, int(gramA.rows()) - 1));
458 
459  geigs.init();
460  int nconv = geigs.compute();
461 
462  //Mat evecs;
463  if (geigs.info() == SUCCESSFUL)
464  {
465  m_evalues = geigs.eigenvalues();
466  m_evectors = geigs.eigenvectors();
468  }
469  else
470  {
471  // Problem With General EgenVec
472  m_info = geigs.info();
473  break;
474  }
475 
476  // Compute Ritz vectors
477  if (iter_num > 0)
478  {
479  eVecX = m_evectors.block(0, 0, m_nev, m_nev);
480  eVecR = m_evectors.block(m_nev, 0, BlockSize, m_nev);
481  eVecD = m_evectors.block(m_nev + BlockSize, 0, BlockSize, m_nev);
482 
483  sparse_eVecX = eVecX.sparseView();
484  sparse_eVecR = eVecR.sparseView();
485  sparse_eVecD = eVecD.sparseView();
486 
487  DD = m_residuals * sparse_eVecR; // new conjugate directions
488  ADD = AR * sparse_eVecR;
489  BDD = BR * sparse_eVecR;
490 
491  DD = DD + directions * sparse_eVecD;
492  ADD = ADD + AD * sparse_eVecD;
493  BDD = BDD + BD * sparse_eVecD;
494  }
495  else
496  {
497  eVecX = m_evectors.block(0, 0, m_nev, m_nev);
498  eVecR = m_evectors.block(m_nev, 0, BlockSize, m_nev);
499 
500  sparse_eVecX = eVecX.sparseView();
501  sparse_eVecR = eVecR.sparseView();
502 
503  DD = m_residuals * sparse_eVecR;
504  ADD = AR * sparse_eVecR;
505  BDD = BR * sparse_eVecR;
506  }
507 
508  X = X * sparse_eVecX + DD;
509  AX = AX * sparse_eVecX + ADD;
510  BX = BX * sparse_eVecX + BDD;
511 
512  directions = DD;
513  AD = ADD;
514  BD = BDD;
515 
516  } // iteration loop
517 
518  // calculate last residuals
520  for (int i = 0; i < m_nev; i++)
521  {
522  m_residuals.col(i) = AX.col(i) - m_evalues(i) * BX.col(i);
523  }
524  BlockSize = checkConvergence_getBlocksize(m_residuals, tolerance_L2, columnsToDelete);
525 
526  if (BlockSize == 0)
527  {
528  m_info = SUCCESSFUL;
529  }
530  } // compute
531 
533  {
534  return m_evalues;
535  }
536 
538  {
539  return m_evectors;
540  }
541 
543  {
544  return Matrix(m_residuals);
545  }
546 
547  int info() { return m_info; }
548 };
549 
550 } // namespace Spectra
551 
552 #endif // LOBPCG_SOLVER
Eigen::SparseMatrix::resize
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:626
Spectra::LOBPCGSolver::ComplexMatrix
Eigen::Matrix< Complex, Eigen::Dynamic, Eigen::Dynamic > ComplexMatrix
Definition: LOBPCGSolver.h:90
Eigen::SparseMatrix::cols
Index cols() const
Definition: SparseMatrix.h:140
Spectra::LOBPCGSolver::m_evectors
Matrix m_evectors
Definition: LOBPCGSolver.h:104
H
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate set rrange[*:*] noreverse nowriteback set trange[*:*] noreverse nowriteback set urange[*:*] noreverse nowriteback set vrange[*:*] noreverse nowriteback set xlabel matrix size set x2label set timefmt d m y n H
Definition: gnuplot_common_settings.hh:74
Spectra::LOBPCGSolver::m_B
SparseMatrix m_B
Definition: LOBPCGSolver.h:99
Spectra::LOBPCGSolver::flag_with_constraints
bool flag_with_constraints
Definition: LOBPCGSolver.h:100
Eigen::SparseMatrix< Scalar >
D
MatrixXcd D
Definition: EigenSolver_EigenSolver_MatrixType.cpp:14
Spectra::LOBPCGSolver::m_info
int m_info
Definition: LOBPCGSolver.h:106
Spectra::LOBPCGSolver::LOBPCGSolver
LOBPCGSolver(const SparseMatrix &A, const SparseMatrix X)
Definition: LOBPCGSolver.h:291
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Eigen::SimplicialLDLT::matrixU
const MatrixU matrixU() const
Definition: SimplicialCholesky.h:463
Spectra::LOBPCGSolver::removeColumns
void removeColumns(SparseMatrix &matrix, std::vector< int > &colToRemove)
Definition: LOBPCGSolver.h:245
Eigen::SimplicialLDLT
A direct sparse LDLT Cholesky factorizations without square root.
Definition: SimplicialCholesky.h:275
gtsam::Y
GaussianFactorGraphValuePair Y
Definition: HybridGaussianProductFactor.cpp:29
Spectra::LOBPCGSolver::ComplexVector
Eigen::Matrix< Complex, Eigen::Dynamic, 1 > ComplexVector
Definition: LOBPCGSolver.h:91
Spectra::LOBPCGSolver::m_residuals
SparseMatrix m_residuals
Definition: LOBPCGSolver.h:103
Spectra::LOBPCGSolver::setPreconditioner
void setPreconditioner(const SparseMatrix &preconditioner)
Definition: LOBPCGSolver.h:322
B
Definition: test_numpy_dtypes.cpp:299
Spectra::LOBPCGSolver::Complex
std::complex< Scalar > Complex
Definition: LOBPCGSolver.h:89
Spectra::DenseSymMatProd
Definition: DenseSymMatProd.h:22
Spectra::LOBPCGSolver::orthogonalizeInPlace
int orthogonalizeInPlace(SparseMatrix &M, SparseMatrix &B, SparseMatrix &true_BM, bool has_true_BM=false)
Definition: LOBPCGSolver.h:110
buffer
Definition: pytypes.h:2270
Spectra::LOBPCGSolver::eigenvectors
Matrix eigenvectors()
Definition: LOBPCGSolver.h:537
Eigen::SimplicialLDLT::vectorD
const VectorType vectorD() const
Definition: SimplicialCholesky.h:452
Spectra::LOBPCGSolver::X
SparseMatrix X
Definition: LOBPCGSolver.h:98
result
Values result
Definition: OdometryOptimize.cpp:8
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
Spectra::SymGEigsSolver
Definition: SymGEigsSolver.h:49
Spectra::LOBPCGSolver::A
SparseMatrix A
Definition: LOBPCGSolver.h:98
Spectra::LOBPCGSolver::m_n
const int m_n
Definition: LOBPCGSolver.h:96
Spectra::LOBPCGSolver::stack_9_matricies
Matrix stack_9_matricies(Matrix A, Matrix B, Matrix C, Matrix D, Matrix E, Matrix F, Matrix G, Matrix H, Matrix I)
Definition: LOBPCGSolver.h:206
AD
static double AD[8]
Definition: airy.c:83
Spectra::LOBPCGSolver::flag_with_preconditioner
bool flag_with_preconditioner
Definition: LOBPCGSolver.h:100
Spectra::LOBPCGSolver::setB
void setB(const SparseMatrix &B)
Definition: LOBPCGSolver.h:314
A
Definition: test_numpy_dtypes.cpp:298
BR
BearingRangeFactor< Pose, Point2 > BR
Definition: SolverComparer.cpp:72
Spectra::LOBPCGSolver::Matrix
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Definition: LOBPCGSolver.h:86
I
#define I
Definition: main.h:112
Eigen::ComputeThinU
@ ComputeThinU
Definition: Constants.h:395
Eigen::SparseMatrix::setIdentity
void setIdentity()
Definition: SparseMatrix.h:749
Spectra::NOT_COMPUTED
@ NOT_COMPUTED
Definition: CompInfo.h:21
Eigen::PlainObjectBase::rows
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:143
Spectra::DenseCholesky
Definition: DenseCholesky.h:26
BD
static double BD[10]
Definition: dawsn.c:93
Spectra::LOBPCGSolver::compute
void compute(int maxit=10, Scalar tol_div_n=1e-7)
Definition: LOBPCGSolver.h:328
Eigen::Triplet
A small structure to hold a non zero as a triplet (i,j,value).
Definition: SparseUtil.h:162
Eigen::EigenSolver::eigenvectors
EigenvectorsType eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: EigenSolver.h:345
gtsam::symbol_shorthand::F
Key F(std::uint64_t j)
Definition: inference/Symbol.h:153
Eigen::SparseMatrix::coeff
Scalar coeff(Index row, Index col) const
Definition: SparseMatrix.h:190
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
Spectra::LOBPCGSolver::sort_epairs
void sort_epairs(Vector &evalues, Matrix &evectors, int SelectionRule)
Definition: LOBPCGSolver.h:224
E
DiscreteKey E(5, 2)
Spectra::LOBPCGSolver::m_nev
const int m_nev
Definition: LOBPCGSolver.h:97
Spectra::LOBPCGSolver::residuals
Matrix residuals()
Definition: LOBPCGSolver.h:542
Spectra::LOBPCGSolver::Vector
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Definition: LOBPCGSolver.h:87
Eigen::ComputeThinV
@ ComputeThinV
Definition: Constants.h:399
Spectra::GEIGS_CHOLESKY
@ GEIGS_CHOLESKY
Using Cholesky decomposition to solve generalized eigenvalues.
Definition: GEigsMode.h:19
Spectra::LOBPCGSolver::SparseComplexMatrix
Eigen::SparseMatrix< Complex > SparseComplexMatrix
Definition: LOBPCGSolver.h:94
C
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:50
Spectra::LOBPCGSolver::m_evalues
Vector m_evalues
Definition: LOBPCGSolver.h:105
Eigen::PlainObjectBase::cols
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:145
Spectra::LOBPCGSolver::flag_with_B
bool flag_with_B
Definition: LOBPCGSolver.h:100
Spectra::LOBPCGSolver::m_Y
SparseMatrix m_Y
Definition: LOBPCGSolver.h:99
Spectra::LOBPCGSolver::eigenvalues
Vector eigenvalues()
Definition: LOBPCGSolver.h:532
G
JacobiRotation< float > G
Definition: Jacobi_makeGivens.cpp:2
Eigen::SparseMatrix::setFromTriplets
void setFromTriplets(const InputIterators &begin, const InputIterators &end)
Definition: SparseMatrix.h:1108
Spectra
Definition: LOBPCGSolver.h:19
Spectra::LOBPCGSolver::SparseMatrix
Eigen::SparseMatrix< Scalar > SparseMatrix
Definition: LOBPCGSolver.h:93
Eigen::EigenSolver
Computes eigenvalues and eigenvectors of general matrices.
Definition: EigenSolver.h:64
Spectra::LOBPCGSolver::checkConvergence_getBlocksize
int checkConvergence_getBlocksize(SparseMatrix &m_residuals, Scalar tolerance_L2, std::vector< int > &columnsToDelete)
Definition: LOBPCGSolver.h:266
min
#define min(a, b)
Definition: datatypes.h:19
Spectra::LOBPCGSolver::stack_4_matricies
Matrix stack_4_matricies(Matrix A, Matrix B, Matrix C, Matrix D)
Definition: LOBPCGSolver.h:195
Spectra::SUCCESSFUL
@ SUCCESSFUL
Computation was successful.
Definition: CompInfo.h:19
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic >
Spectra::SMALLEST_ALGE
@ SMALLEST_ALGE
Select eigenvalues with smallest algebraic value. Only for symmetric eigen solvers.
Definition: SelectionRule.h:51
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
Eigen::SparseMatrix::rows
Index rows() const
Definition: SparseMatrix.h:138
Eigen::EigenSolver::info
ComputationInfo info() const
Definition: EigenSolver.h:281
Spectra::LOBPCGSolver::m_preconditioner
SparseMatrix m_preconditioner
Definition: LOBPCGSolver.h:99
Spectra::LOBPCGSolver::applyConstraintsInPlace
void applyConstraintsInPlace(SparseMatrix &X, SparseMatrix &Y, SparseMatrix &B)
Definition: LOBPCGSolver.h:170
ceres::sqrt
Jet< T, N > sqrt(const Jet< T, N > &f)
Definition: jet.h:418
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Spectra::LOBPCGSolver::info
int info()
Definition: LOBPCGSolver.h:547
Spectra::LOBPCGSolver::setConstraints
void setConstraints(const SparseMatrix &Y)
Definition: LOBPCGSolver.h:308
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Spectra::LOBPCGSolver
Definition: LOBPCGSolver.h:83
M
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:51
Eigen::EigenSolver::eigenvalues
const EigenvalueType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: EigenSolver.h:244


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:02:43