Getting started

The easiest way to solve a problem with AMGCL is to use the amgcl::make_solver class. It has two template parameters: the first one specifies a preconditioner to use, and the second chooses an iterative solver. The class constructor takes the system matrix in one of supported formats and parameters for the chosen algorithms and for the backend.

Solving Poisson’s equation

Let us consider a simple example of Poisson’s equation in a unit square. Here is how the problem may be solved with AMGCL. We will use BiCGStab solver preconditioned with smoothed aggregation multigrid with SPAI(0) for relaxation (smoothing). First, we include the necessary headers. Each of those brings in the corresponding component of the method:

#include <amgcl/make_solver.hpp>
#include <amgcl/solver/bicgstab.hpp>
#include <amgcl/amg.hpp>
#include <amgcl/coarsening/smoothed_aggregation.hpp>
#include <amgcl/relaxation/spai0.hpp>
#include <amgcl/adapter/crs_tuple.hpp>

Next, we assemble sparse matrix for the Poisson’s equation on a uniform 1000x1000 grid. See Assembling matrix for Poisson’s equation for the source code of the poisson() function:

std::vector<int>    ptr, col;
std::vector<double> val, rhs;
int n = poisson(1000, ptr, col, val, rhs);

For this example, we select the builtin backend with double precision numbers as value type:

typedef amgcl::backend::builtin<double> Backend;

Now we can construct the solver for our system matrix. We use the convenient adapter for boost tuples here and just tie together the matrix size and its CRS components:

typedef amgcl::make_solver<
    // Use AMG as preconditioner:
    amgcl::amg<
        Backend,
        amgcl::coarsening::smoothed_aggregation,
        amgcl::relaxation::spai0
        >,
    // And BiCGStab as iterative solver:
    amgcl::solver::bicgstab<Backend>
    > Solver;

Solver solve( boost::tie(n, ptr, col, val) );

Once the solver is constructed, we can apply it to the right-hand side to obtain the solution. This may be repeated multiple times for different right-hand sides. Here we start with a zero initial approximation. The solver returns a boost tuple with number of iterations and norm of the achieved residual:

std::vector<double> x(n, 0.0);
int    iters;
double error;
boost::tie(iters, error) = solve(rhs, x);

That’s it! Vector x contains the solution of our problem now.

Input formats

We used STL vectors to store the matrix components in the above axample. This may seem too restrictive if you want to use AMGCL with your own types. But the crs_tuple adapter will take anything that the Boost.Range library recognizes as a random access range. For example, you can wrap raw pointers to your data into a boost::iterator_range:

Solver solve( boost::make_tuple(
    n,
    boost::make_iterator_range(ptr.data(), ptr.data() + ptr.size()),
    boost::make_iterator_range(col.data(), col.data() + col.size()),
    boost::make_iterator_range(val.data(), val.data() + val.size())
    ) );

Same applies to the right-hand side and the solution vectors. And if that is still not general enough, you can provide your own adapter for your matrix type. See Matrix adapters for further information on this.

Setting parameters

Any component in AMGCL defines its own parameters by declaring a param subtype. When a class wraps several subclasses, it includes parameters of its children into its own param. For example, parameters for the amgcl::make_solver<Precond, Solver> are declared as

struct params {
    typename Precond::params precond;
    typename Solver::params solver;
};

Knowing that, we can easily set the parameters for individual components. For example, we can set the desired tolerance for the iterative solver in the above example like this:

Solver::params prm;
prm.solver.tol = 1e-3;
Solver solve( boost::tie(n, ptr, col, val), prm );

Parameters may also be initialized with a boost::property_tree::ptree. This is especially convenient when Runtime interface is used, and the exact structure of the parameters is not known at compile time:

boost::property_tree::ptree prm;
prm.put("solver.tol", 1e-3);
Solver solve( boost::tie(n, ptr, col, val), prm );

The make_solver class

template <class Precond, class IterativeSolver>
class amgcl::make_solver

Convenience class that bundles together a preconditioner and an iterative solver.

Public Functions

template <class Matrix>
make_solver(const Matrix &A, const params &prm = params (), const backend_params &bprm = backend_params())

Sets up the preconditioner and creates the iterative solver.

template <class Matrix, class Vec1, class Vec2>
boost::tuple<size_t, scalar_type> operator()(Matrix const &A, Vec1 const &rhs, Vec2 &&x) const

Computes the solution for the given system matrix A and the right-hand side rhs. Returns the number of iterations made and the achieved residual as a boost::tuple. The solution vector x provides initial approximation in input and holds the computed solution on output.

The system matrix may differ from the matrix used during initialization. This may be used for the solution of non-stationary problems with slowly changing coefficients. There is a strong chance that a preconditioner built for a time step will act as a reasonably good preconditioner for several subsequent time steps [DeSh12].

template <class Vec1, class Vec2>
boost::tuple<size_t, scalar_type> operator()(Vec1 const &rhs, Vec2 &&x) const

Computes the solution for the given right-hand side rhs. Returns the number of iterations made and the achieved residual as a boost::tuple. The solution vector x provides initial approximation in input and holds the computed solution on output.

template <class Vec1, class Vec2>
void apply(const Vec1 &rhs, Vec2 &&x) const

Acts as a preconditioner. That is, applies the solver to the right-hand side rhs to get the solution x with zero initial approximation. Iterative methods usually use estimated residual for exit condition. For some problems the value of the estimated residual can get too far from the true residual due to round-off errors. Nesting iterative solvers in this way may allow to shave the last bits off the error. The method should not be used directly but rather allows nesting make_solver classes as in the following example:

typedef amgcl::make_solver<
  amgcl::make_solver<
    amgcl::amg<
      Backend, amgcl::coarsening::smoothed_aggregation, amgcl::relaxation::spai0
      >,
    amgcl::solver::cg<Backend>
    >,
  amgcl::solver::cg<Backend>
  > NestedSolver;

const Precond &precond() const

Returns reference to the constructed preconditioner.

const IterativeSolver &solver() const

Returns reference to the constructed iterative solver.

Precond::matrix const &system_matrix() const

Returns the system matrix in the backend format.

void get_params(boost::property_tree::ptree &p) const

Stores the parameters used during construction into the property tree p.

size_t size() const

Returns the size of the system matrix.

struct params

Combined parameters of the bundled preconditioner and the iterative solver.

Public Members

template<>
Precond::params precond

Preconditioner parameters.

template<>
IterativeSolver::params solver

Iterative solver parameters.