Eigen  3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
 
Loading...
Searching...
No Matches
Matrix-free solvers

Iterative solvers such as ConjugateGradient and BiCGSTAB can be used in a matrix free context. To this end, user must provide a wrapper class inheriting EigenBase<> and implementing the following methods:

Eigen::internal::traits<> must also be specialized for the wrapper type.

Here is a complete example wrapping an Eigen::SparseMatrix:

#include <iostream>
#include <Eigen/Core>
#include <Eigen/Dense>
#include <Eigen/IterativeLinearSolvers>
#include <unsupported/Eigen/IterativeSolvers>
class MatrixReplacement;
namespace Eigen {
namespace internal {
// MatrixReplacement looks-like a SparseMatrix, so let's inherit its traits:
template <>
struct traits<MatrixReplacement> : public Eigen::internal::traits<Eigen::SparseMatrix<double> > {};
} // namespace internal
} // namespace Eigen
// Example of a matrix-free wrapper from a user type to Eigen's compatible type
// For the sake of simplicity, this example simply wrap a Eigen::SparseMatrix.
class MatrixReplacement : public Eigen::EigenBase<MatrixReplacement> {
public:
// Required typedefs, constants, and method:
typedef double Scalar;
typedef double RealScalar;
typedef int StorageIndex;
enum { ColsAtCompileTime = Eigen::Dynamic, MaxColsAtCompileTime = Eigen::Dynamic, IsRowMajor = false };
Index rows() const { return mp_mat->rows(); }
Index cols() const { return mp_mat->cols(); }
template <typename Rhs>
Eigen::Product<MatrixReplacement, Rhs, Eigen::AliasFreeProduct> operator*(const Eigen::MatrixBase<Rhs>& x) const {
return Eigen::Product<MatrixReplacement, Rhs, Eigen::AliasFreeProduct>(*this, x.derived());
}
// Custom API:
MatrixReplacement() : mp_mat(0) {}
void attachMyMatrix(const SparseMatrix<double>& mat) { mp_mat = &mat; }
const SparseMatrix<double> my_matrix() const { return *mp_mat; }
private:
const SparseMatrix<double>* mp_mat;
};
// Implementation of MatrixReplacement * Eigen::DenseVector though a specialization of internal::generic_product_impl:
namespace Eigen {
namespace internal {
template <typename Rhs>
struct generic_product_impl<MatrixReplacement, Rhs, SparseShape, DenseShape,
GemvProduct> // GEMV stands for matrix-vector
: generic_product_impl_base<MatrixReplacement, Rhs, generic_product_impl<MatrixReplacement, Rhs> > {
typedef typename Product<MatrixReplacement, Rhs>::Scalar Scalar;
template <typename Dest>
static void scaleAndAddTo(Dest& dst, const MatrixReplacement& lhs, const Rhs& rhs, const Scalar& alpha) {
// This method should implement "dst += alpha * lhs * rhs" inplace,
// however, for iterative solvers, alpha is always equal to 1, so let's not bother about it.
eigen_assert(alpha == Scalar(1) && "scaling is not implemented");
EIGEN_ONLY_USED_FOR_DEBUG(alpha);
// Here we could simply call dst.noalias() += lhs.my_matrix() * rhs,
// but let's do something fancier (and less efficient):
for (Index i = 0; i < lhs.cols(); ++i) dst += rhs(i) * lhs.my_matrix().col(i);
}
};
} // namespace internal
} // namespace Eigen
int main() {
int n = 10;
Eigen::SparseMatrix<double> S = Eigen::MatrixXd::Random(n, n).sparseView(0.5, 1);
S = S.transpose() * S;
MatrixReplacement A;
A.attachMyMatrix(S);
Eigen::VectorXd b(n), x;
b.setRandom();
// Solve Ax = b using various iterative solver with matrix-free version:
{
cg.compute(A);
x = cg.solve(b);
std::cout << "CG: #iterations: " << cg.iterations() << ", estimated error: " << cg.error() << std::endl;
}
{
bicg.compute(A);
x = bicg.solve(b);
std::cout << "BiCGSTAB: #iterations: " << bicg.iterations() << ", estimated error: " << bicg.error() << std::endl;
}
{
Eigen::GMRES<MatrixReplacement, Eigen::IdentityPreconditioner> gmres;
gmres.compute(A);
x = gmres.solve(b);
std::cout << "GMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
}
{
Eigen::DGMRES<MatrixReplacement, Eigen::IdentityPreconditioner> gmres;
gmres.compute(A);
x = gmres.solve(b);
std::cout << "DGMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
}
{
Eigen::MINRES<MatrixReplacement, Eigen::Lower | Eigen::Upper, Eigen::IdentityPreconditioner> minres;
minres.compute(A);
x = minres.solve(b);
std::cout << "MINRES: #iterations: " << minres.iterations() << ", estimated error: " << minres.error()
<< std::endl;
}
}
A bi conjugate gradient stabilized solver for sparse square problems.
Definition BiCGSTAB.h:153
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition ConjugateGradient.h:152
RealScalar error() const
Definition IterativeSolverBase.h:270
Derived & compute(const EigenBase< MatrixDerived > &A)
Definition IterativeSolverBase.h:210
Index iterations() const
Definition IterativeSolverBase.h:262
A versatible sparse matrix representation.
Definition SparseUtil.h:47
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:84
Matrix< double, Dynamic, 1 > VectorXd
DynamicĂ—1 vector of type double.
Definition Matrix.h:479
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const Product< MatrixDerived, PermutationDerived, AliasFreeProduct > operator*(const MatrixBase< MatrixDerived > &matrix, const PermutationBase< PermutationDerived > &permutation)
Definition PermutationMatrix.h:471
const int Dynamic
Definition Constants.h:25
Definition ForwardDeclarations.h:57
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition EigenBase.h:61
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition EigenBase.h:59

Output:

CG:       #iterations: 19, estimated error: 1.0666e-18
BiCGSTAB: #iterations: 20, estimated error: 1.56382e-24
GMRES:    #iterations: 10, estimated error: 0
DGMRES:   #iterations: 20, estimated error: 5.63173e-29
MINRES:   #iterations: 19, estimated error: 5.04856e-18