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 <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> >
{};
}
}
// 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>
}
// 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");
// 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);
}
};
}
}
int main()
{
int n = 10;
S = S.transpose()*S;
MatrixReplacement A;
A.attachMyMatrix(S);
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;
}
}
Array< int, 3, 1 > b
int n
MatrixXcf A
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:914
#define eigen_assert(x)
Definition: Macros.h:902
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:161
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
static const RandomReturnType Random()
Definition: Random.h:114
Derived & compute(const EigenBase< MatrixDerived > &A)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
Expression of the product of two arbitrary matrices or vectors.
Definition: Product.h:77
TransposeReturnType transpose()
A versatible sparse matrix representation.
Definition: SparseMatrix.h:125
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
int main(int, char **)
Definition: class_Block.cpp:18
: InteropHeaders
Definition: Core:139
@ GemvProduct
Definition: Constants.h:504
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)
const int Dynamic
Definition: Constants.h:24

Output:

CG:       #iterations: 20, estimated error: 1.46043e-14
BiCGSTAB: #iterations: 17, estimated error: 9.10876e-17
GMRES:    #iterations: 10, estimated error: 0
DGMRES:   #iterations: 20, estimated error: 2.83433e-26
MINRES:   #iterations: 20, estimated error: 3.81447e-15