SparseRegularInverse.h
Go to the documentation of this file.
1 // Copyright (C) 2017-2019 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 SPARSE_REGULAR_INVERSE_H
8 #define 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 
28 template <typename Scalar, int Uplo = Eigen::Lower, int Flags = 0, typename StorageIndex = int>
30 {
31 private:
38 
39  ConstGenericSparseMatrix m_mat;
40  const int m_n;
42 
43 public:
51  SparseRegularInverse(ConstGenericSparseMatrix& mat) :
52  m_mat(mat), m_n(mat.rows())
53  {
54  if (mat.rows() != mat.cols())
55  throw std::invalid_argument("SparseRegularInverse: matrix must be square");
56 
57  m_cg.compute(mat);
58  }
59 
63  Index rows() const { return m_n; }
67  Index cols() const { return m_n; }
68 
75  // y_out = inv(B) * x_in
76  void solve(const Scalar* x_in, Scalar* y_out) const
77  {
78  MapConstVec x(x_in, m_n);
79  MapVec y(y_out, m_n);
80  y.noalias() = m_cg.solve(x);
81  }
82 
89  // y_out = B * x_in
90  void mat_prod(const Scalar* x_in, Scalar* y_out) const
91  {
92  MapConstVec x(x_in, m_n);
93  MapVec y(y_out, m_n);
94  y.noalias() = m_mat.template selfadjointView<Uplo>() * x;
95  }
96 };
97 
98 } // namespace Spectra
99 
100 #endif // SPARSE_REGULAR_INVERSE_H
SCALAR Scalar
Definition: bench_gemm.cpp:33
Scalar * y
SparseRegularInverse(ConstGenericSparseMatrix &mat)
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
const Solve< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs > solve(const MatrixBase< Rhs > &b) const
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
void solve(const Scalar *x_in, Scalar *y_out) const
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute(const EigenBase< MatrixDerived > &A)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
void mat_prod(const Scalar *x_in, Scalar *y_out) const
Eigen::Map< const Vector > MapConstVec
Eigen::SparseMatrix< Scalar, Flags, StorageIndex > SparseMatrix
Eigen::ConjugateGradient< SparseMatrix > m_cg
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
const Eigen::Ref< const SparseMatrix > ConstGenericSparseMatrix
ConstGenericSparseMatrix m_mat
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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:44:37