Content  


Examples 5: Iterative Solvers

ViennaCL provides different iterative solvers for various classes of matrices and is not restricted to built-in types. Since the interface of ViennaCL types is compatible with uBLAS objects, the iterative solvers can directly be called with uBLAS objects. ViennaCL 1.0.x provides the following iterative solvers:

  • Conjugate Gradient (viennacl::linalg::cg_tag)
  • Stabilized Bi-Conjugate Gradient (viennacl::linalg::bicgstab_tag)
  • Generalized Minimum Residual (viennacl::linalg::gmres_tag)

An optional incomplete LU factorization with threshold can be used as preconditioner.

Iterative solvers in ViennaCL
typedef float        ScalarType;
//typedef double    ScalarType; //use this if your GPU supports double precision
 
// Set up some ublas objects:
ublas::vector<ScalarType> ublas_rhs;
ublas::vector<ScalarType> ublas_result;
ublas::compressed_matrix<ScalarType> ublas_matrix;
 
// Set up some ViennaCL objects:
viennacl::vector<ScalarType> vcl_rhs;
viennacl::vector<ScalarType> vcl_result;
viennacl::compressed_matrix<ScalarType> vcl_matrix;
 
/* Initialize and fill all objects here */
 
//
// Compute ILUT preconditioners for CPU and for GPU objects:
//
viennacl::linalg::ilut_tag ilut_conf(10, 1e-5);  //10 entries, rel. tol. 1e-5
typedef viennacl::linalg::ilut_precond< 
                    ublas::compressed_matrix<ScalarType> >     ublas_ilut_t;
//preconditioner for ublas objects:
ublas_ilut_t ublas_ilut(ublas_matrix, ilut_conf);   
 
viennacl::linalg::ilut_precond<
                    viennacl::compressed_matrix<ScalarType> >  vcl_ilut_t;
//preconditioner for ViennaCL objects:
vcl_ilut_t vcl_ilut(vcl_matrix, ilut_conf);
 
//
// Conjugate gradient solver without preconditioner:
//
ublas_result  = solve(ublas_matrix,   //using ublas objects on CPU
                      ublas_rhs,
                      viennacl::linalg::cg_tag());
vcl_result    = solve(vcl_matrix,     //using viennacl objects on GPU
                      vcl_rhs,
                      viennacl::linalg::cg_tag());
 
 
//
// Conjugate gradient solver using ILUT preconditioner
//
ublas_result  = solve(ublas_matrix,   //using ublas objects on CPU
                      ublas_rhs,
                      viennacl::linalg::cg_tag(), 
                      ublas_ilut);
vcl_result    = solve(vcl_matrix,     //using viennacl objects on GPU
                      vcl_rhs,
                      viennacl::linalg::cg_tag(),
                      vcl_ilut);
  
// for BiCGStab and GMRES, use the solver tags
// viennacl::linalg::bicgstab_tag and viennacl::linalg::gmres_tag 
// instead of viennacl::linalg::cg_tag in the calls above.