This tutorial shows the use of algebraic multigrid (AMG) preconditioners.
- Warning
- AMG is currently only experimentally available with the OpenCL backend and depends on Boost.uBLAS
We start with some rather general includes and preprocessor variables:
#ifndef NDEBUG //without NDEBUG the performance of sparse ublas matrices is poor.
#define NDEBUG
#endif
#include <boost/numeric/ublas/matrix_sparse.hpp>
#include <boost/numeric/ublas/operation_sparse.hpp>
#define VIENNACL_WITH_UBLAS 1
Import the AMG functionality:
Some more includes:
#include <iostream>
#include <vector>
#include <ctime>
Part 1: Worker routines
Run the Solver
Runs the provided solver specified in the solver
object with the provided preconditioner precond
template<typename MatrixType, typename VectorType, typename SolverTag, typename PrecondTag>
void run_solver(MatrixType
const & matrix, VectorType
const & rhs, VectorType
const & ref_result, SolverTag
const & solver, PrecondTag
const & precond)
{
VectorType result(rhs);
VectorType residual(rhs);
std::cout << " > Iterations: " << solver.iters() << std::endl;
result -= ref_result;
}
Compare AMG preconditioner for uBLAS and ViennaCL types
The AMG implementations in ViennaCL can be used with uBLAS types as well as ViennaCL types. This function compares the two in terms of execution time.
template<typename ScalarType>
boost::numeric::ublas::vector<ScalarType> & ,
boost::numeric::ublas::vector<ScalarType> & ,
boost::numeric::ublas::compressed_matrix<ScalarType> & ublas_matrix,
std::string info,
{
boost::numeric::ublas::vector<ScalarType> avgstencil;
std::cout << "-- CG with AMG preconditioner, " << info << " --" << std::endl;
std::cout << " * Setup phase (ublas types)..." << std::endl;
ublas_amg.
tag().set_coarselevels(coarselevels);
std::cout <<
" * Operator complexity: " << ublas_amg.
calc_complexity(avgstencil) << std::endl;
std::cout << " * Setup phase (ViennaCL types)..." << std::endl;
vcl_amg.
tag().set_coarselevels(coarselevels);
std::cout << " * CG solver (ublas types)..." << std::endl;
std::cout << " * CG solver (ViennaCL types)..." << std::endl;
run_solver(vcl_compressed_matrix, vcl_vec, vcl_result, cg_solver, vcl_amg);
}
Part 2: Run Solvers with AMG Preconditioners
In this
Print some device info at the beginning. If there is more than one OpenCL device available, use the second device.
std::cout << std::endl;
std::cout << "----------------------------------------------" << std::endl;
std::cout << " Device Info" << std::endl;
std::cout << "----------------------------------------------" << std::endl;
#ifdef VIENNACL_WITH_OPENCL
std::vector<viennacl::ocl::device>
const & devices = pf.
devices();
if (devices.size() > 1)
else
#else
#endif
Set up the matrices and vectors for the iterative solvers (cf. iterative.cpp)
boost::numeric::ublas::vector<ScalarType> ublas_vec, ublas_result;
boost::numeric::ublas::compressed_matrix<ScalarType> ublas_matrix;
{
std::cout << "Error reading Matrix file" << std::endl;
return EXIT_FAILURE;
}
{
std::cout << "Error reading RHS file" << std::endl;
return 0;
}
{
std::cout << "Error reading Result file" << std::endl;
return 0;
}
Instantiate a tag for the conjugate gradient solver, the AMG preconditioner tag, and create an AMG preconditioner object:
Run solver without preconditioner. This serves as a baseline for comparison. Note that iterative solvers without preconditioner on GPUs can be very efficient because they map well to the massively parallel hardware.
std::cout << "-- CG solver (CPU, no preconditioner) --" << std::endl;
std::cout << "-- CG solver (GPU, no preconditioner) --" << std::endl;
Generate the setup for an AMG preconditioner of Ruge-Stueben type with direct interpolation (RS+DIRECT) and run the solver:
0.25,
0.2,
0.67,
3,
3,
0);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, DIRECT INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner of Ruge-Stueben type with classic interpolation (RS+CLASSIC) and run the solver:
run_amg ( cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, CLASSIC INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner of Ruge-Stueben type with only one pass and direct interpolation (ONEPASS+DIRECT)
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "ONEPASS COARSENING, DIRECT INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner of parallel Ruge-Stueben type with direct interpolation (RS0+DIRECT)
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS0 COARSENING, DIRECT INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner of parallel Ruge-Stueben type (with communication across domains) and use direct interpolation (RS3+DIRECT)
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS3 COARSENING, DIRECT INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner which as aggregation-based (AG)
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, AG INTERPOLATION", amg_tag);
Generate the setup for an AMG preconditioner with smoothed aggregation (SA)
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, SA INTERPOLATION",amg_tag);
That's it.
std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
return EXIT_SUCCESS;
}
Full Example Code
#ifndef NDEBUG //without NDEBUG the performance of sparse ublas matrices is poor.
#define NDEBUG
#endif
#include <boost/numeric/ublas/matrix_sparse.hpp>
#include <boost/numeric/ublas/operation_sparse.hpp>
#define VIENNACL_WITH_UBLAS 1
#include <iostream>
#include <vector>
#include <ctime>
template<typename MatrixType, typename VectorType, typename SolverTag, typename PrecondTag>
void run_solver(MatrixType
const & matrix, VectorType
const & rhs, VectorType
const & ref_result, SolverTag
const & solver, PrecondTag
const & precond)
{
VectorType result(rhs);
VectorType residual(rhs);
std::cout << " > Iterations: " << solver.iters() << std::endl;
result -= ref_result;
}
template<typename ScalarType>
boost::numeric::ublas::vector<ScalarType> & ,
boost::numeric::ublas::vector<ScalarType> & ,
boost::numeric::ublas::compressed_matrix<ScalarType> & ublas_matrix,
std::string info,
{
boost::numeric::ublas::vector<ScalarType> avgstencil;
std::cout << "-- CG with AMG preconditioner, " << info << " --" << std::endl;
std::cout << " * Setup phase (ublas types)..." << std::endl;
ublas_amg.
tag().set_coarselevels(coarselevels);
std::cout <<
" * Operator complexity: " << ublas_amg.
calc_complexity(avgstencil) << std::endl;
std::cout << " * Setup phase (ViennaCL types)..." << std::endl;
vcl_amg.
tag().set_coarselevels(coarselevels);
std::cout << " * CG solver (ublas types)..." << std::endl;
std::cout << " * CG solver (ViennaCL types)..." << std::endl;
run_solver(vcl_compressed_matrix, vcl_vec, vcl_result, cg_solver, vcl_amg);
}
{
std::cout << std::endl;
std::cout << "----------------------------------------------" << std::endl;
std::cout << " Device Info" << std::endl;
std::cout << "----------------------------------------------" << std::endl;
#ifdef VIENNACL_WITH_OPENCL
std::vector<viennacl::ocl::device>
const & devices = pf.
devices();
if (devices.size() > 1)
else
#else
#endif
boost::numeric::ublas::vector<ScalarType> ublas_vec, ublas_result;
boost::numeric::ublas::compressed_matrix<ScalarType> ublas_matrix;
{
std::cout << "Error reading Matrix file" << std::endl;
return EXIT_FAILURE;
}
{
std::cout << "Error reading RHS file" << std::endl;
return 0;
}
{
std::cout << "Error reading Result file" << std::endl;
return 0;
}
std::cout << "-- CG solver (CPU, no preconditioner) --" << std::endl;
std::cout << "-- CG solver (GPU, no preconditioner) --" << std::endl;
0.25,
0.2,
0.67,
3,
3,
0);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, DIRECT INTERPOLATION", amg_tag);
run_amg ( cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, CLASSIC INTERPOLATION", amg_tag);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "ONEPASS COARSENING, DIRECT INTERPOLATION", amg_tag);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS0 COARSENING, DIRECT INTERPOLATION", amg_tag);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS3 COARSENING, DIRECT INTERPOLATION", amg_tag);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, AG INTERPOLATION", amg_tag);
run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, SA INTERPOLATION",amg_tag);
std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
return EXIT_SUCCESS;
}