SparseSymShiftSolve.h
Go to the documentation of this file.
1 // Copyright (C) 2016-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_SYM_SHIFT_SOLVE_H
8 #define SPARSE_SYM_SHIFT_SOLVE_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseLU>
13 #include <stdexcept>
14 
15 namespace Spectra {
16 
24 template <typename Scalar, int Uplo = Eigen::Lower, int Flags = 0, typename StorageIndex = int>
26 {
27 private:
34 
35  ConstGenericSparseMatrix m_mat;
36  const int m_n;
38 
39 public:
47  SparseSymShiftSolve(ConstGenericSparseMatrix& mat) :
48  m_mat(mat), m_n(mat.rows())
49  {
50  if (mat.rows() != mat.cols())
51  throw std::invalid_argument("SparseSymShiftSolve: matrix must be square");
52  }
53 
57  Index rows() const { return m_n; }
61  Index cols() const { return m_n; }
62 
67  {
68  SparseMatrix mat = m_mat.template selfadjointView<Uplo>();
69  SparseMatrix identity(m_n, m_n);
70  identity.setIdentity();
71  mat = mat - sigma * identity;
72  m_solver.isSymmetric(true);
73  m_solver.compute(mat);
74  if (m_solver.info() != Eigen::Success)
75  throw std::invalid_argument("SparseSymShiftSolve: factorization failed with the given shift");
76  }
77 
84  // y_out = inv(A - sigma * I) * x_in
85  void perform_op(const Scalar* x_in, Scalar* y_out) const
86  {
87  MapConstVec x(x_in, m_n);
88  MapVec y(y_out, m_n);
89  y.noalias() = m_solver.solve(x);
90  }
91 };
92 
93 } // namespace Spectra
94 
95 #endif // SPARSE_SYM_SHIFT_SOLVE_H
SCALAR Scalar
Definition: bench_gemm.cpp:33
Scalar * y
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
Eigen::Map< const Vector > MapConstVec
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
static const double sigma
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
const Solve< SparseLU< _MatrixType, _OrderingType >, Rhs > solve(const MatrixBase< Rhs > &b) const
ConstGenericSparseMatrix m_mat
const Eigen::Ref< const SparseMatrix > ConstGenericSparseMatrix
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseLU.h:202
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:17
SparseSymShiftSolve(ConstGenericSparseMatrix &mat)
void compute(const MatrixType &matrix)
Definition: SparseLU.h:124
void isSymmetric(bool sym)
Definition: SparseLU.h:135
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Eigen::SparseMatrix< Scalar, Flags, StorageIndex > SparseMatrix
void perform_op(const Scalar *x_in, Scalar *y_out) const
Eigen::SparseLU< SparseMatrix > m_solver
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:39