ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
eigen-with-viennacl.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 // System headers
26 #include <iostream>
27 
28 // Eigen headers
29 #include <Eigen/Core>
30 #include <Eigen/Sparse>
31 
32 // IMPORTANT: Must be set prior to any ViennaCL includes if you want to use ViennaCL algorithms on Eigen objects
33 #define VIENNACL_WITH_EIGEN 1
34 
35 
36 // ViennaCL includes
37 #include "viennacl/vector.hpp"
38 #include "viennacl/matrix.hpp"
40 #include "viennacl/linalg/prod.hpp"
41 
42 
43 // Helper functions for this tutorial:
44 #include "vector-io.hpp"
45 
51 //dense matrix:
52 template<typename T>
53 struct Eigen_dense_matrix
54 {
55  typedef typename T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;
56 };
57 
58 template<>
59 struct Eigen_dense_matrix<float>
60 {
61  typedef Eigen::MatrixXf type;
62 };
63 
64 template<>
65 struct Eigen_dense_matrix<double>
66 {
67  typedef Eigen::MatrixXd type;
68 };
69 
70 
71 //sparse matrix
72 template<typename T>
73 struct Eigen_vector
74 {
75  typedef typename T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;
76 };
77 
78 template<>
79 struct Eigen_vector<float>
80 {
81  typedef Eigen::VectorXf type;
82 };
83 
84 template<>
85 struct Eigen_vector<double>
86 {
87  typedef Eigen::VectorXd type;
88 };
89 
90 
91 
102 template<typename ScalarType>
103 void run_tutorial()
104 {
109  typedef typename Eigen_dense_matrix<ScalarType>::type EigenMatrix;
110  typedef typename Eigen_vector<ScalarType>::type EigenVector;
111 
115  EigenMatrix eigen_densemat(6, 5);
116  EigenMatrix eigen_densemat2(6, 5);
117  eigen_densemat(0,0) = 2.0; eigen_densemat(0,1) = -1.0;
118  eigen_densemat(1,0) = -1.0; eigen_densemat(1,1) = 2.0; eigen_densemat(1,2) = -1.0;
119  eigen_densemat(2,1) = -1.0; eigen_densemat(2,2) = -1.0; eigen_densemat(2,3) = -1.0;
120  eigen_densemat(3,2) = -1.0; eigen_densemat(3,3) = 2.0; eigen_densemat(3,4) = -1.0;
121  eigen_densemat(5,4) = -1.0; eigen_densemat(4,4) = -1.0;
122  Eigen::Map<EigenMatrix> eigen_densemat_map(eigen_densemat.data(), 6, 5); // same as eigen_densemat, but emulating user-provided buffer
123 
127  Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat(6, 5);
128  Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat2(6, 5);
129  eigen_sparsemat.reserve(5*2);
130  eigen_sparsemat.insert(0,0) = 2.0; eigen_sparsemat.insert(0,1) = -1.0;
131  eigen_sparsemat.insert(1,1) = 2.0; eigen_sparsemat.insert(1,2) = -1.0;
132  eigen_sparsemat.insert(2,2) = -1.0; eigen_sparsemat.insert(2,3) = -1.0;
133  eigen_sparsemat.insert(3,3) = 2.0; eigen_sparsemat.insert(3,4) = -1.0;
134  eigen_sparsemat.insert(5,4) = -1.0;
135  //eigen_sparsemat.endFill();
136 
140  EigenVector eigen_rhs(5);
141  Eigen::Map<EigenVector> eigen_rhs_map(eigen_rhs.data(), 5);
142  EigenVector eigen_result(6);
143  EigenVector eigen_temp(6);
144 
145  eigen_rhs(0) = 10.0;
146  eigen_rhs(1) = 11.0;
147  eigen_rhs(2) = 12.0;
148  eigen_rhs(3) = 13.0;
149  eigen_rhs(4) = 14.0;
150 
151 
155  viennacl::vector<ScalarType> vcl_rhs(5);
156  viennacl::vector<ScalarType> vcl_result(6);
157  viennacl::matrix<ScalarType> vcl_densemat(6, 5);
158  viennacl::compressed_matrix<ScalarType> vcl_sparsemat(6, 5);
159 
160 
164  viennacl::copy(&(eigen_rhs[0]), &(eigen_rhs[0]) + 5, vcl_rhs.begin()); // Method 1: via iterator interface (cf. std::copy())
165  viennacl::copy(eigen_rhs, vcl_rhs); // Method 2: via built-in wrappers (convenience layer)
166  viennacl::copy(eigen_rhs_map, vcl_rhs); // Same as method 2, but for a mapped vector
167 
168  viennacl::copy(eigen_densemat, vcl_densemat);
169  viennacl::copy(eigen_densemat_map, vcl_densemat); //same as above, using mapped matrix
170  viennacl::copy(eigen_sparsemat, vcl_sparsemat);
171  std::cout << "VCL sparsematrix dimensions: " << vcl_sparsemat.size1() << ", " << vcl_sparsemat.size2() << std::endl;
172 
173  // For completeness: Copy matrices from ViennaCL back to Eigen:
174  viennacl::copy(vcl_densemat, eigen_densemat2);
175  viennacl::copy(vcl_sparsemat, eigen_sparsemat2);
176 
177 
181  eigen_result = eigen_densemat * eigen_rhs;
182  vcl_result = viennacl::linalg::prod(vcl_densemat, vcl_rhs);
183  viennacl::copy(vcl_result, eigen_temp);
184  std::cout << "Difference for dense matrix-vector product: " << (eigen_result - eigen_temp).norm() << std::endl;
185  std::cout << "Difference for dense matrix-vector product (Eigen->ViennaCL->Eigen): "
186  << (eigen_densemat2 * eigen_rhs - eigen_temp).norm() << std::endl;
187 
191  eigen_result = eigen_sparsemat * eigen_rhs;
192  vcl_result = viennacl::linalg::prod(vcl_sparsemat, vcl_rhs);
193  viennacl::copy(vcl_result, eigen_temp);
194  std::cout << "Difference for sparse matrix-vector product: " << (eigen_result - eigen_temp).norm() << std::endl;
195  std::cout << "Difference for sparse matrix-vector product (Eigen->ViennaCL->Eigen): "
196  << (eigen_sparsemat2 * eigen_rhs - eigen_temp).norm() << std::endl;
197 }
198 
199 
203 int main(int, char *[])
204 {
205  std::cout << "----------------------------------------------" << std::endl;
206  std::cout << "## Single precision" << std::endl;
207  std::cout << "----------------------------------------------" << std::endl;
208  run_tutorial<float>();
209 
210 #ifdef VIENNACL_HAVE_OPENCL
212 #endif
213  {
214  std::cout << "----------------------------------------------" << std::endl;
215  std::cout << "## Double precision" << std::endl;
216  std::cout << "----------------------------------------------" << std::endl;
217  run_tutorial<double>();
218  }
219 
223  std::cout << std::endl;
224  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
225  std::cout << std::endl;
226 
227 }
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
A dense matrix class.
Definition: forwards.h:375
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
Definition: backend.hpp:351
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
Implementation of the compressed_matrix class.
bool double_support() const
ViennaCL convenience function: Returns true if the device supports double precision.
Definition: device.hpp:956
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) ...
A sparse square matrix in compressed sparse rows format.