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-mtl4.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 // include necessary system headers
28 #include <iostream>
29 
30 // MTL4 headers
31 #include <boost/numeric/mtl/mtl.hpp>
32 #include <boost/numeric/itl/itl.hpp>
33 
34 // Must be set prior to any ViennaCL includes if you want to use ViennaCL algorithms on Eigen objects
35 #define VIENNACL_WITH_MTL4 1
36 
37 // ViennaCL includes
38 #include "viennacl/linalg/ilu.hpp"
39 #include "viennacl/linalg/cg.hpp"
43 
44 
45 // Some helper functions for this tutorial:
46 #include "vector-io.hpp"
47 
52 int main(int, char *[])
53 {
54  typedef double ScalarType;
55 
56  mtl::compressed2D<ScalarType> mtl4_matrix;
57  mtl4_matrix.change_dim(65025, 65025);
58  set_to_zero(mtl4_matrix);
59 
60  mtl::dense_vector<ScalarType> mtl4_rhs(65025, 1.0);
61  mtl::dense_vector<ScalarType> mtl4_result(65025, 0.0);
62  mtl::dense_vector<ScalarType> mtl4_residual(65025, 0.0);
63 
68  mtl::io::matrix_market_istream("../examples/testdata/mat65k.mtx") >> mtl4_matrix;
69 
73  std::cout << "----- Running CG -----" << std::endl;
74  mtl4_result = viennacl::linalg::solve(mtl4_matrix, mtl4_rhs, viennacl::linalg::cg_tag());
75 
76  mtl4_residual = mtl4_matrix * mtl4_result - mtl4_rhs;
77  std::cout << "Relative residual: " << viennacl::linalg::norm_2(mtl4_residual) / viennacl::linalg::norm_2(mtl4_rhs) << std::endl;
78 
82  std::cout << "----- Running BiCGStab -----" << std::endl;
83  mtl4_result = viennacl::linalg::solve(mtl4_matrix, mtl4_rhs, viennacl::linalg::bicgstab_tag());
84 
85  mtl4_residual = mtl4_matrix * mtl4_result - mtl4_rhs;
86  std::cout << "Relative residual: " << viennacl::linalg::norm_2(mtl4_residual) / viennacl::linalg::norm_2(mtl4_rhs) << std::endl;
87 
91  std::cout << "----- Running GMRES -----" << std::endl;
92  mtl4_result = viennacl::linalg::solve(mtl4_matrix, mtl4_rhs, viennacl::linalg::gmres_tag());
93 
94  mtl4_residual = mtl4_matrix * mtl4_result - mtl4_rhs;
95  std::cout << "Relative residual: " << viennacl::linalg::norm_2(mtl4_residual) / viennacl::linalg::norm_2(mtl4_rhs) << std::endl;
96 
100  std::cout << std::endl;
101  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
102  std::cout << std::endl;
103 }
104 
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.
int main()
Definition: bisect.cpp:91
The stabilized bi-conjugate gradient method is implemented here.
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.
The conjugate gradient method is implemented here.
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