ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
spai.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2015, Institute for Microelectronics,
3  Institute for Analysis and Scientific Computing,
4  TU Wien.
5  Portions of this software are copyright by UChicago Argonne, LLC.
6 
7  -----------------
8  ViennaCL - The Vienna Computing Library
9  -----------------
10 
11  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
12 
13  (A list of authors and contributors can be found in the PDF manual)
14 
15  License: MIT (X11), see file LICENSE in the base directory
16 ============================================================================= */
17 
27 // enable Boost.uBLAS support
28 #define VIENNACL_WITH_UBLAS
29 
30 #ifndef NDEBUG
31  #define BOOST_UBLAS_NDEBUG
32 #endif
33 
34 // System headers:
35 #include <utility>
36 #include <iostream>
37 #include <fstream>
38 #include <string>
39 #include <cmath>
40 #include <algorithm>
41 #include <stdio.h>
42 
43 // ViennaCL headers:
44 #include "viennacl/scalar.hpp"
45 #include "viennacl/matrix.hpp"
47 #include "viennacl/linalg/cg.hpp"
49 #include "viennacl/linalg/prod.hpp"
51 #include "viennacl/linalg/ilu.hpp"
54 #include "viennacl/linalg/spai.hpp"
55 
56 // Boost headers:
57 #include "boost/numeric/ublas/vector.hpp"
58 #include "boost/numeric/ublas/matrix.hpp"
59 #include "boost/numeric/ublas/io.hpp"
60 
61 // auxiliary functionality:
62 #include "vector-io.hpp"
63 
67 template<typename MatrixT, typename VectorT, typename SolverTagT, typename PreconditionerT>
68 void run_solver(MatrixT const & A, VectorT const & b, SolverTagT const & solver_tag, PreconditionerT const & precond)
69 {
70  VectorT result = viennacl::linalg::solve(A, b, solver_tag, precond);
71  std::cout << " * Solver iterations: " << solver_tag.iters() << std::endl;
72  VectorT residual = viennacl::linalg::prod(A, result);
73  residual -= b;
74  std::cout << " * Rel. Residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(b) << std::endl;
75 }
76 
87 int main (int, const char **)
88 {
89  typedef float ScalarType;
90  typedef boost::numeric::ublas::compressed_matrix<ScalarType> MatrixType;
91  typedef boost::numeric::ublas::vector<ScalarType> VectorType;
92  typedef viennacl::compressed_matrix<ScalarType> GPUMatrixType;
93  typedef viennacl::vector<ScalarType> GPUVectorType;
94 
98 #ifdef VIENNACL_WITH_OPENCL
99  // Optional: Customize OpenCL backend
101  std::vector<viennacl::ocl::device> const & devices = pf.devices();
102 
103  // Optional: Set first device to first context:
104  viennacl::ocl::setup_context(0, devices[0]);
105 
106  // Optional: Set second device for second context (use the same device for the second context if only one device available):
107  if (devices.size() > 1)
108  viennacl::ocl::setup_context(1, devices[1]);
109  else
110  viennacl::ocl::setup_context(1, devices[0]);
111 
112  std::cout << viennacl::ocl::current_device().info() << std::endl;
114 #else
115  viennacl::context ctx;
116 #endif
117 
121  MatrixType M;
122 
123  if (!viennacl::io::read_matrix_market_file(M, "../examples/testdata/mat65k.mtx"))
124  {
125  std::cerr<<"ERROR: Could not read matrix file " << std::endl;
126  exit(EXIT_FAILURE);
127  }
128 
129  std::cout << "Size of matrix: " << M.size1() << std::endl;
130  std::cout << "Avg. Entries per row: " << double(M.nnz()) / static_cast<double>(M.size1()) << std::endl;
131 
135  VectorType rhs(M.size2());
136  for (std::size_t i=0; i<rhs.size(); ++i)
137  rhs(i) = ScalarType(1);
138 
142  GPUMatrixType gpu_M(M.size1(), M.size2(), ctx);
143  GPUVectorType gpu_rhs(M.size1(), ctx);
144  viennacl::copy(M, gpu_M);
145  viennacl::copy(rhs, gpu_rhs);
146 
153  viennacl::linalg::bicgstab_tag solver_tag(1e-10, 50); //for simplicity and reasonably short execution times we use only 50 iterations here
154 
158  std::cout << "--- Reference 1: Pure BiCGStab on CPU ---" << std::endl;
159  VectorType result = viennacl::linalg::solve(M, rhs, solver_tag);
160  std::cout << " * Solver iterations: " << solver_tag.iters() << std::endl;
161  VectorType residual = viennacl::linalg::prod(M, result) - rhs;
162  std::cout << " * Rel. Residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(rhs) << std::endl;
163 
164  std::cout << "--- Reference 2: Pure BiCGStab on GPU ---" << std::endl;
165  GPUVectorType gpu_result = viennacl::linalg::solve(gpu_M, gpu_rhs, solver_tag);
166  std::cout << " * Solver iterations: " << solver_tag.iters() << std::endl;
167  GPUVectorType gpu_residual = viennacl::linalg::prod(gpu_M, gpu_result);
168  gpu_residual -= gpu_rhs;
169  std::cout << " * Rel. Residual: " << viennacl::linalg::norm_2(gpu_residual) / viennacl::linalg::norm_2(gpu_rhs) << std::endl;
170 
171 
175  std::cout << "--- Reference 2: BiCGStab with ILUT on CPU ---" << std::endl;
176  std::cout << " * Preconditioner setup..." << std::endl;
178  std::cout << " * Iterative solver run..." << std::endl;
179  run_solver(M, rhs, solver_tag, ilut);
180 
181 
185  std::cout << "--- Test 1: CPU-based SPAI ---" << std::endl;
186  std::cout << " * Preconditioner setup..." << std::endl;
188  std::cout << " * Iterative solver run..." << std::endl;
189  run_solver(M, rhs, solver_tag, spai_cpu);
190 
194  std::cout << "--- Test 2: CPU-based FSPAI ---" << std::endl;
195  std::cout << " * Preconditioner setup..." << std::endl;
197  std::cout << " * Iterative solver run..." << std::endl;
198  run_solver(M, rhs, solver_tag, fspai_cpu);
199 
203  std::cout << "--- Test 3: GPU-based SPAI ---" << std::endl;
204  std::cout << " * Preconditioner setup..." << std::endl;
206  std::cout << " * Iterative solver run..." << std::endl;
207  run_solver(gpu_M, gpu_rhs, solver_tag, spai_gpu);
208 
212  std::cout << "--- Test 4: GPU-based FSPAI ---" << std::endl;
213  std::cout << " * Preconditioner setup..." << std::endl;
215  std::cout << " * Iterative solver run..." << std::endl;
216  run_solver(gpu_M, gpu_rhs, solver_tag, fspai_gpu);
217 
221  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
222 
223  return EXIT_SUCCESS;
224 }
225 
T norm_2(std::vector< T, A > const &v1)
Definition: norm_2.hpp:96
A reader and writer for the matrix market format is implemented here.
std::vector< platform > get_platforms()
Definition: platform.hpp:124
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
Wrapper class for an OpenCL platform.
Definition: platform.hpp:45
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Implementation of the dense matrix class.
int main()
Definition: bisect.cpp:91
std::vector< device > devices(cl_device_type dtype=CL_DEVICE_TYPE_DEFAULT)
Returns the available devices of the supplied device type.
Definition: platform.hpp:91
The stabilized bi-conjugate gradient method is implemented here.
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
Definition: backend.hpp:351
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations.
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
Definition: bicgstab.hpp:496
std::string info(vcl_size_t indent=0, char indent_char= ' ') const
Returns an info string with a few properties of the device. Use full_info() to get all details...
Definition: device.hpp:995
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Definition: context.hpp:39
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
A tag for FSPAI. Experimental.
Definition: fspai.hpp:71
Implementations of incomplete factorization preconditioners. Convenience header file.
Implementation of the Factored SParse Approximate Inverse Algorithm for a generic, uBLAS-compatible matrix type.
Definition: spai.hpp:189
A tag for incomplete LU factorization with threshold (ILUT)
Definition: ilut.hpp:45
Implementation of the compressed_matrix class.
Implementation of the SParse Approximate Inverse Algorithm for a generic, uBLAS-compatible matrix typ...
Definition: spai.hpp:75
void run_solver(MatrixType const &matrix, VectorType const &rhs, VectorType const &ref_result, SolverTag const &solver, PrecondTag const &precond, long ops)
Definition: solver.cpp:101
ILUT preconditioner class, can be supplied to solve()-routines.
Definition: ilut.hpp:352
The conjugate gradient method is implemented here.
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
float ScalarType
Definition: fft_1d.cpp:42
Main include file for the sparse approximate inverse preconditioner family (SPAI and FSPAI)...
A sparse square matrix in compressed sparse rows format.
A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...
Definition: bicgstab.hpp:47
long read_matrix_market_file(MatrixT &mat, const char *file, long index_base=1)
Reads a sparse matrix from a file (MatrixMarket format)
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
Definition: backend.hpp:225
Implementation of the ViennaCL scalar class.
void setup_context(long i, std::vector< cl_device_id > const &devices)
Convenience function for setting devices for a context.
Definition: backend.hpp:231