ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
iterative-ublas.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 
25 //
26 // include necessary system headers
27 //
28 #include <iostream>
29 
30 //
31 // Necessary to obtain a suitable performance in ublas
32 #ifndef NDEBUG
33  #define BOOST_UBLAS_NDEBUG
34 #endif
35 
36 
37 //
38 // ublas includes
39 //
40 #include <boost/numeric/ublas/io.hpp>
41 #include <boost/numeric/ublas/triangular.hpp>
42 #include <boost/numeric/ublas/matrix_sparse.hpp>
43 #include <boost/numeric/ublas/matrix.hpp>
44 #include <boost/numeric/ublas/matrix_proxy.hpp>
45 #include <boost/numeric/ublas/operation.hpp>
46 #include <boost/numeric/ublas/operation_sparse.hpp>
47 #include <boost/numeric/ublas/io.hpp>
48 #include <boost/numeric/ublas/lu.hpp>
49 
50 // Must be set if you want to use ViennaCL algorithms on ublas objects
51 #define VIENNACL_WITH_UBLAS 1
52 
53 //
54 // ViennaCL includes
55 //
56 #include "viennacl/linalg/ilu.hpp"
57 #include "viennacl/linalg/cg.hpp"
61 
62 // Some helper functions for this tutorial:
63 #include "vector-io.hpp"
64 
68 using namespace boost::numeric;
69 
73 int main()
74 {
75  typedef float ScalarType;
76 
77  //
78  // Set up some ublas objects
79  //
80  ublas::vector<ScalarType> rhs;
81  ublas::vector<ScalarType> ref_result;
82  ublas::vector<ScalarType> result;
83  ublas::compressed_matrix<ScalarType> ublas_matrix;
84 
88  if (!viennacl::io::read_matrix_market_file(ublas_matrix, "../examples/testdata/mat65k.mtx"))
89  {
90  std::cout << "Error reading Matrix file" << std::endl;
91  return EXIT_FAILURE;
92  }
93 
97  if (!readVectorFromFile("../examples/testdata/rhs65025.txt", rhs))
98  {
99  std::cout << "Error reading RHS file" << std::endl;
100  return EXIT_FAILURE;
101  }
102  //std::cout << "done reading rhs" << std::endl;
103 
104  if (!readVectorFromFile("../examples/testdata/result65025.txt", ref_result))
105  {
106  std::cout << "Error reading Result file" << std::endl;
107  return EXIT_FAILURE;
108  }
109  //std::cout << "done reading result" << std::endl;
110 
111 
118  viennacl::linalg::ilu0_tag> ublas_block_ilu0(ublas_matrix, viennacl::linalg::ilu0_tag());
119 
123  std::cout << "----- CG Test -----" << std::endl;
124 
125  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag());
126  std::cout << "Residual norm: " << norm_2(prod(ublas_matrix, result) - rhs) << std::endl;
127  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag(1e-6, 20), ublas_ilut);
128  std::cout << "Residual norm: " << norm_2(prod(ublas_matrix, result) - rhs) << std::endl;
129  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag(1e-6, 20), ublas_ilu0);
130  std::cout << "Residual norm: " << norm_2(prod(ublas_matrix, result) - rhs) << std::endl;
131  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::cg_tag(1e-6, 20), ublas_block_ilu0);
132  std::cout << "Residual norm: " << norm_2(prod(ublas_matrix, result) - rhs) << std::endl;
133 
137  std::cout << "----- BiCGStab Test -----" << std::endl;
138 
139  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag()); //without preconditioner
140  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), ublas_ilut); //with preconditioner
141  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::bicgstab_tag(1e-6, 20), ublas_ilu0); //with preconditioner
142 
146  std::cout << "----- GMRES Test -----" << std::endl;
147 
148  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag()); //without preconditioner
149  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag(1e-6, 20), ublas_ilut);//with preconditioner
150  result = viennacl::linalg::solve(ublas_matrix, rhs, viennacl::linalg::gmres_tag(1e-6, 20), ublas_ilu0);//with preconditioner
151 
155  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
156 
157  return EXIT_SUCCESS;
158 }
159 
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.
ILU0 preconditioner class, can be supplied to solve()-routines.
Definition: ilu0.hpp:146
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 block ILU preconditioner class, can be supplied to solve()-routines.
Definition: block_ilu.hpp:132
A tag for the solver GMRES. Used for supplying solver parameters and for dispatching the solve() func...
Definition: gmres.hpp:49
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
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.
void prod(std::vector< std::map< IndexT, NumericT > > const &stl_A, std::vector< std::map< IndexT, NumericT > > const &stl_B, std::vector< std::map< IndexT, NumericT > > &stl_C)
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 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)