The following tutorial shows how to use the iterative solvers in ViennaCL with objects from the Eigen Library directly.
- Note
- Eigen provides its own iterative solvers in the meanwhile. Check these first.
We begin with including the necessary headers:
#include <iostream>
#ifndef NDEBUG
#define NDEBUG
#endif
#define ARMA_DONT_USE_BLAS
#define ARMA_DONT_USE_LAPACK
#include <armadillo>
#define VIENNACL_WITH_ARMADILLO 1
In the following we run the CG method, the BiCGStab method, and the GMRES method with Armadillo types directly. First, the matrices are set up, then the respective solvers are called.
Read system from file. This is a little tricky, since Armadillo does not provide a fast enough element-insertion. Therefore, we read the matrix market file to an STL-matrix and then pass the data on when creating the Armadillo sparse matrix object.
std::vector<std::map<unsigned int, ScalarType> > stl_matrix;
std::cout << "Reading matrix (this might take some time)..." << std::endl;
{
std::cout << "Error reading Matrix file. Make sure you run from the build/-folder." << std::endl;
return EXIT_FAILURE;
}
std::size_t num_nnz = 0;
for (std::size_t i=0; i<stl_matrix.size(); ++i)
num_nnz += stl_matrix[i].
size();
arma::Mat<arma::uword> arma_indices(2, num_nnz);
arma::Col<ScalarType> arma_values(num_nnz);
std::size_t index = 0;
for (std::size_t i=0; i<stl_matrix.size(); ++i)
{
for (std::map<unsigned int, ScalarType>::const_iterator it = stl_matrix[i].begin(); it != stl_matrix[i].end(); ++it)
{
arma_indices(0, index) = i;
arma_indices(1, index) = it->first;
arma_values(index) = it->second;
++index;
}
}
std::cout << "Done: reading matrix" << std::endl;
Initialize Armadillo types for iterative solvers
arma::SpMat<ScalarType> arma_matrix(arma_indices, arma_values, 65025, 65025);
arma::Col<ScalarType> arma_rhs;
arma::Col<ScalarType> arma_result;
arma::Col<ScalarType> residual;
Read the right hand side as well as the result vector from files:
{
std::cout << "Error reading RHS file" << std::endl;
return EXIT_FAILURE;
}
{
std::cout << "Error reading Result file" << std::endl;
return EXIT_FAILURE;
}
Conjugate Gradient (CG) solver:
std::cout << "----- Running CG -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
Stabilized Bi-Conjugate Gradient (BiCGStab) solver:
std::cout << "----- Running BiCGStab -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
Generalized Minimum Residual (GMRES) solver:
std::cout << "----- Running GMRES -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
That's it. Print a success message and exit.
std::cout << std::endl;
std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
std::cout << std::endl;
}
Full Example Code
#include <iostream>
#ifndef NDEBUG
#define NDEBUG
#endif
#define ARMA_DONT_USE_BLAS
#define ARMA_DONT_USE_LAPACK
#include <armadillo>
#define VIENNACL_WITH_ARMADILLO 1
{
std::vector<std::map<unsigned int, ScalarType> > stl_matrix;
std::cout << "Reading matrix (this might take some time)..." << std::endl;
{
std::cout << "Error reading Matrix file. Make sure you run from the build/-folder." << std::endl;
return EXIT_FAILURE;
}
std::size_t num_nnz = 0;
for (std::size_t i=0; i<stl_matrix.size(); ++i)
num_nnz += stl_matrix[i].
size();
arma::Mat<arma::uword> arma_indices(2, num_nnz);
arma::Col<ScalarType> arma_values(num_nnz);
std::size_t index = 0;
for (std::size_t i=0; i<stl_matrix.size(); ++i)
{
for (std::map<unsigned int, ScalarType>::const_iterator it = stl_matrix[i].begin(); it != stl_matrix[i].end(); ++it)
{
arma_indices(0, index) = i;
arma_indices(1, index) = it->first;
arma_values(index) = it->second;
++index;
}
}
std::cout << "Done: reading matrix" << std::endl;
arma::SpMat<ScalarType> arma_matrix(arma_indices, arma_values, 65025, 65025);
arma::Col<ScalarType> arma_rhs;
arma::Col<ScalarType> arma_result;
arma::Col<ScalarType> residual;
{
std::cout << "Error reading RHS file" << std::endl;
return EXIT_FAILURE;
}
{
std::cout << "Error reading Result file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "----- Running CG -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
std::cout << "----- Running BiCGStab -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
std::cout << "----- Running GMRES -----" << std::endl;
residual = arma_matrix * arma_result - arma_rhs;
std::cout << std::endl;
std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
std::cout << std::endl;
}