ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
iterative.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2016, 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 //
28 // include necessary system headers
29 //
30 #include <iostream>
31 
32 //
33 // Necessary to obtain a suitable performance in ublas
34 #ifndef NDEBUG
35  #define BOOST_UBLAS_NDEBUG
36 #endif
37 
38 //
39 // ublas includes
40 //
41 #include <boost/numeric/ublas/io.hpp>
42 #include <boost/numeric/ublas/triangular.hpp>
43 #include <boost/numeric/ublas/matrix_sparse.hpp>
44 #include <boost/numeric/ublas/matrix.hpp>
45 #include <boost/numeric/ublas/matrix_proxy.hpp>
46 #include <boost/numeric/ublas/operation.hpp>
47 #include <boost/numeric/ublas/operation_sparse.hpp>
48 #include <boost/numeric/ublas/io.hpp>
49 #include <boost/numeric/ublas/lu.hpp>
50 
51 // Must be set if you want to use ViennaCL algorithms on ublas objects
52 #define VIENNACL_WITH_UBLAS 1
53 
54 
55 //
56 // ViennaCL includes
57 //
58 #include "viennacl/scalar.hpp"
59 #include "viennacl/vector.hpp"
62 #include "viennacl/linalg/prod.hpp"
63 #include "viennacl/linalg/ilu.hpp"
65 #include "viennacl/linalg/cg.hpp"
69 
70 
71 // Some helper functions for this tutorial:
72 #include "vector-io.hpp"
73 
77 using namespace boost::numeric;
78 
79 
85 int main()
86 {
87  typedef float ScalarType;
88 
92  ublas::vector<ScalarType> rhs;
93  ublas::vector<ScalarType> ref_result;
94  ublas::vector<ScalarType> result;
95  ublas::compressed_matrix<ScalarType> ublas_matrix;
96 
100  if (!viennacl::io::read_matrix_market_file(ublas_matrix, "../examples/testdata/mat65k.mtx"))
101  {
102  std::cout << "Error reading Matrix file" << std::endl;
103  return EXIT_FAILURE;
104  }
105 
106  if (!readVectorFromFile("../examples/testdata/rhs65025.txt", rhs))
107  {
108  std::cout << "Error reading RHS file" << std::endl;
109  return EXIT_FAILURE;
110  }
111 
112  if (!readVectorFromFile("../examples/testdata/result65025.txt", ref_result))
113  {
114  std::cout << "Error reading Result file" << std::endl;
115  return EXIT_FAILURE;
116  }
117 
121  std::size_t vcl_size = rhs.size();
122  viennacl::compressed_matrix<ScalarType> vcl_compressed_matrix;
123  viennacl::coordinate_matrix<ScalarType> vcl_coordinate_matrix;
124  viennacl::vector<ScalarType> vcl_rhs(vcl_size);
125  viennacl::vector<ScalarType> vcl_result(vcl_size);
126  viennacl::vector<ScalarType> vcl_ref_result(vcl_size);
127 
128  viennacl::copy(rhs.begin(), rhs.end(), vcl_rhs.begin());
129  viennacl::copy(ref_result.begin(), ref_result.end(), vcl_ref_result.begin());
130 
131 
135  viennacl::copy(ublas_matrix, vcl_compressed_matrix);
136 
143  std::vector< std::map< unsigned int, ScalarType> > stl_matrix(rhs.size());
144  for (ublas::compressed_matrix<ScalarType>::iterator1 iter1 = ublas_matrix.begin1();
145  iter1 != ublas_matrix.end1();
146  ++iter1)
147  {
148  for (ublas::compressed_matrix<ScalarType>::iterator2 iter2 = iter1.begin();
149  iter2 != iter1.end();
150  ++iter2)
151  stl_matrix[iter2.index1()][static_cast<unsigned int>(iter2.index2())] = *iter2;
152  }
153  std::vector<ScalarType> stl_rhs(rhs.size()), stl_result(result.size());
154  std::copy(rhs.begin(), rhs.end(), stl_rhs.begin());
155  viennacl::copy(stl_matrix, vcl_coordinate_matrix);
156  viennacl::copy(vcl_coordinate_matrix, stl_matrix);
157 
161  std::cout << "Setting up preconditioners for uBLAS-matrix..." << std::endl;
165  viennacl::linalg::ilu0_tag> ublas_block_ilu0(ublas_matrix, viennacl::linalg::ilu0_tag());
166 
167  std::cout << "Setting up preconditioners for ViennaCL-matrix..." << std::endl;
171  viennacl::linalg::ilu0_tag> vcl_block_ilu0(vcl_compressed_matrix, viennacl::linalg::ilu0_tag());
172 
178 
182  std::cout << "----- CG Method -----" << std::endl;
183 
189  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag());
190  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag(1e-6, 20), ublas_ilut);
191  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag(1e-6, 20), ublas_jacobi);
192 
193 
198  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::cg_tag());
199  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::cg_tag(1e-6, 20), vcl_ilut);
200  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::cg_tag(1e-6, 20), vcl_jacobi);
201 
206  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::cg_tag());
207  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::cg_tag(1e-6, 20), vcl_ilut);
208  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::cg_tag(1e-6, 20), vcl_jacobi);
209 
210 
214  std::cout << "----- BiCGStab Method -----" << std::endl;
215 
219  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag()); //without preconditioner
220  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), ublas_ilut); //with preconditioner
221  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), ublas_jacobi); //with preconditioner
222 
227  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::bicgstab_tag()); //without preconditioner
228  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), vcl_ilut); //with preconditioner
229  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), vcl_jacobi); //with preconditioner
230 
235  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::bicgstab_tag());
236  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), vcl_ilut);
237  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), vcl_jacobi);
238 
239 
243  std::cout << "----- GMRES Method -----" << std::endl;
244 
248  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag()); //without preconditioner
249  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag(1e-6, 20), ublas_ilut);//with preconditioner
250  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag(1e-6, 20), ublas_jacobi);//with preconditioner
251 
256  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::gmres_tag()); //without preconditioner
257  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::gmres_tag(1e-6, 20), vcl_ilut);//with preconditioner
258  vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::gmres_tag(1e-6, 20), vcl_jacobi);//with preconditioner
259 
264  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::gmres_tag());
265  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::gmres_tag(1e-6, 20), vcl_ilut);
266  stl_result = viennacl::linalg::solve(stl_matrix, stl_rhs, viennacl::linalg::gmres_tag(1e-6, 20), vcl_jacobi);
267 
268 
272  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
273 
274  return EXIT_SUCCESS;
275 }
276 
A reader and writer for the matrix market format is implemented here.
ILU0 preconditioner class, can be supplied to solve()-routines.
Definition: ilu0.hpp:146
Jacobi preconditioner class, can be supplied to solve()-routines. Generic version for non-ViennaCL ma...
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
int main()
Definition: bisect.cpp:91
A tag for incomplete LU factorization with static pattern (ILU0)
Definition: ilu0.hpp:58
The stabilized bi-conjugate gradient method is implemented here.
A tag for a jacobi preconditioner.
A block ILU preconditioner class, can be supplied to solve()-routines.
Definition: block_ilu.hpp:132
Implementation of a simple Jacobi preconditioner.
A tag for the solver GMRES. Used for supplying solver parameters and for dispatching the solve() func...
Definition: gmres.hpp:49
Implementation of the coordinate_matrix class.
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
Definition: bicgstab.hpp:496
Implementations of the generalized minimum residual method are in this file.
Implementations of incomplete factorization preconditioners. Convenience header file.
A tag for incomplete LU factorization with threshold (ILUT)
Definition: ilut.hpp:45
Implementation of the compressed_matrix class.
bool readVectorFromFile(const std::string &filename, VectorType &vec)
Definition: vector-io.hpp:104
ILUT preconditioner class, can be supplied to solve()-routines.
Definition: ilut.hpp:352
The conjugate gradient method is implemented here.
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
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
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
Definition: cg.hpp:48
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)
Implementation of the ViennaCL scalar class.
A sparse square matrix, where entries are stored as triplets (i,j, val), where i and j are the row an...