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-eigen.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 // System headers
28 #include <iostream>
29 
30 #ifndef NDEBUG
31  #define NDEBUG
32 #endif
33 
34 
35 // Eigen headers
36 #include <Eigen/Core>
37 #include <Eigen/Sparse>
38 
39 // Must be set prior to any ViennaCL includes if you want to use ViennaCL algorithms on Eigen objects
40 #define VIENNACL_WITH_EIGEN 1
41 
42 // ViennaCL headers
43 #include "viennacl/linalg/ilu.hpp"
44 #include "viennacl/linalg/cg.hpp"
48 
49 
50 // Some helper functions for this tutorial:
51 #include "vector-io.hpp"
52 
57 int main(int, char *[])
58 {
59  typedef float ScalarType;
60 
61  Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_matrix(65025, 65025);
62  Eigen::VectorXf eigen_rhs;
63  Eigen::VectorXf eigen_result;
64  Eigen::VectorXf ref_result;
65  Eigen::VectorXf residual;
66 
70  std::cout << "Reading matrix (this might take some time)..." << std::endl;
71  eigen_matrix.reserve(65025 * 7);
72  if (!viennacl::io::read_matrix_market_file(eigen_matrix, "../examples/testdata/mat65k.mtx"))
73  {
74  std::cout << "Error reading Matrix file. Make sure you run from the build/-folder." << std::endl;
75  return EXIT_FAILURE;
76  }
77  //eigen_matrix.endFill();
78  std::cout << "Done: reading matrix" << std::endl;
79 
80  if (!readVectorFromFile("../examples/testdata/rhs65025.txt", eigen_rhs))
81  {
82  std::cout << "Error reading RHS file" << std::endl;
83  return EXIT_FAILURE;
84  }
85 
86  if (!readVectorFromFile("../examples/testdata/result65025.txt", ref_result))
87  {
88  std::cout << "Error reading Result file" << std::endl;
89  return EXIT_FAILURE;
90  }
91 
95  std::cout << "----- Running CG -----" << std::endl;
96  eigen_result = viennacl::linalg::solve(eigen_matrix, eigen_rhs, viennacl::linalg::cg_tag());
97 
98  residual = eigen_matrix * eigen_result - eigen_rhs;
99  std::cout << "Relative residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(eigen_rhs) << std::endl;
100 
104  std::cout << "----- Running BiCGStab -----" << std::endl;
105  eigen_result = viennacl::linalg::solve(eigen_matrix, eigen_rhs, viennacl::linalg::bicgstab_tag());
106 
107  residual = eigen_matrix * eigen_result - eigen_rhs;
108  std::cout << "Relative residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(eigen_rhs) << std::endl;
109 
113  std::cout << "----- Running GMRES -----" << std::endl;
114  eigen_result = viennacl::linalg::solve(eigen_matrix, eigen_rhs, viennacl::linalg::gmres_tag());
115 
116  residual = eigen_matrix * eigen_result - eigen_rhs;
117  std::cout << "Relative residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(eigen_rhs) << std::endl;
118 
122  std::cout << std::endl;
123  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
124  std::cout << std::endl;
125 }
126 
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.
bool readVectorFromFile(const std::string &filename, VectorType &vec)
Definition: vector-io.hpp:104
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
long read_matrix_market_file(MatrixT &mat, const char *file, long index_base=1)
Reads a sparse matrix from a file (MatrixMarket format)