SparseRegularInverse.h
Go to the documentation of this file.
1 // Copyright (C) 2017-2025 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_SPARSE_REGULAR_INVERSE_H
8 #define SPECTRA_SPARSE_REGULAR_INVERSE_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/SparseCore>
12 #include <Eigen/IterativeLinearSolvers>
13 #include <stdexcept>
14 
15 namespace Spectra {
16 
36 template <typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor, typename StorageIndex = int>
38 {
39 public:
43  using Scalar = Scalar_;
44 
45 private:
52 
54  const Index m_n;
56  mutable CompInfo m_info;
57 
58 public:
66  template <typename Derived>
68  m_mat(mat), m_n(mat.rows())
69  {
70  static_assert(
71  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(SparseMatrix::IsRowMajor),
72  "SparseRegularInverse: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
73 
74  if (mat.rows() != mat.cols())
75  throw std::invalid_argument("SparseRegularInverse: matrix must be square");
76 
77  m_cg.compute(mat);
78  m_info = (m_cg.info() == Eigen::Success) ?
81  }
82 
86  Index rows() const { return m_n; }
90  Index cols() const { return m_n; }
91 
96  CompInfo info() const { return m_info; }
97 
104  // y_out = inv(B) * x_in
105  void solve(const Scalar* x_in, Scalar* y_out) const
106  {
107  MapConstVec x(x_in, m_n);
108  MapVec y(y_out, m_n);
109  y.noalias() = m_cg.solve(x);
110 
111  m_info = (m_cg.info() == Eigen::Success) ?
115  throw std::runtime_error("SparseRegularInverse: CG solver does not converge");
116  }
117 
124  // y_out = B * x_in
125  void perform_op(const Scalar* x_in, Scalar* y_out) const
126  {
127  MapConstVec x(x_in, m_n);
128  MapVec y(y_out, m_n);
129  y.noalias() = m_mat.template selfadjointView<Uplo>() * x;
130  }
131 };
132 
133 } // namespace Spectra
134 
135 #endif // SPECTRA_SPARSE_REGULAR_INVERSE_H
Eigen::SparseMatrix
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
Spectra::SparseRegularInverse
Definition: SparseRegularInverse.h:37
Spectra::SparseRegularInverse::m_info
CompInfo m_info
Definition: SparseRegularInverse.h:56
Spectra::CompInfo::Successful
@ Successful
Computation was successful.
Spectra::SparseRegularInverse::SparseRegularInverse
SparseRegularInverse(const Eigen::SparseMatrixBase< Derived > &mat)
Definition: SparseRegularInverse.h:67
Spectra::SparseRegularInverse::perform_op
void perform_op(const Scalar *x_in, Scalar *y_out) const
Definition: SparseRegularInverse.h:125
x
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 x
Definition: gnuplot_common_settings.hh:12
Eigen::Success
@ Success
Definition: Constants.h:442
mat
MatrixXf mat
Definition: Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
Spectra::SparseRegularInverse< double >::Scalar
double Scalar
Definition: SparseRegularInverse.h:43
Spectra::SparseRegularInverse::m_n
const Index m_n
Definition: SparseRegularInverse.h:54
Spectra::CompInfo::NotConverging
@ NotConverging
Spectra::CompInfo
CompInfo
Definition: CompInfo.h:17
Eigen::ConjugateGradient
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition: ConjugateGradient.h:97
Spectra::SparseRegularInverse::solve
void solve(const Scalar *x_in, Scalar *y_out) const
Definition: SparseRegularInverse.h:105
Spectra::SparseRegularInverse::m_cg
Eigen::ConjugateGradient< SparseMatrix > m_cg
Definition: SparseRegularInverse.h:55
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
y
Scalar * y
Definition: level1_cplx_impl.h:124
Spectra::CompInfo::NumericalIssue
@ NumericalIssue
Spectra::SparseRegularInverse::rows
Index rows() const
Definition: SparseRegularInverse.h:86
Eigen::Ref< const SparseMatrix >
Spectra::SparseRegularInverse::info
CompInfo info() const
Definition: SparseRegularInverse.h:96
Spectra
Definition: LOBPCGSolver.h:19
Eigen::SparseMatrixBase
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:301
Spectra::SparseRegularInverse::m_mat
ConstGenericSparseMatrix m_mat
Definition: SparseRegularInverse.h:53
Spectra::SparseRegularInverse::cols
Index cols() const
Definition: SparseRegularInverse.h:90
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 >
Spectra::SparseRegularInverse< double >::Index
Eigen::Index Index
Definition: SparseRegularInverse.h:46
Eigen::SparseCompressedBase< SparseMatrix< _Scalar, _Options, _StorageIndex > >::IsRowMajor
@ IsRowMajor
Definition: SparseMatrixBase.h:100
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74


gtsam
Author(s):
autogenerated on Thu Apr 10 2025 03:03:49