ViennaCL - The Vienna Computing Library  1.6.2
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
eigen-with-viennacl.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2014, 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 "Random.hpp"
45 #include "vector-io.hpp"
46 #include "../benchmarks/benchmark-utils.hpp"
47 
53 //dense matrix:
54 template<typename T>
55 struct Eigen_dense_matrix
56 {
57  typedef typename T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;
58 };
59 
60 template<>
61 struct Eigen_dense_matrix<float>
62 {
63  typedef Eigen::MatrixXf type;
64 };
65 
66 template<>
67 struct Eigen_dense_matrix<double>
68 {
69  typedef Eigen::MatrixXd type;
70 };
71 
72 
73 //sparse matrix
74 template<typename T>
75 struct Eigen_vector
76 {
77  typedef typename T::ERROR_NO_EIGEN_TYPE_AVAILABLE error_type;
78 };
79 
80 template<>
81 struct Eigen_vector<float>
82 {
83  typedef Eigen::VectorXf type;
84 };
85 
86 template<>
87 struct Eigen_vector<double>
88 {
89  typedef Eigen::VectorXd type;
90 };
91 
92 
93 
104 template<typename ScalarType>
105 void run_tutorial()
106 {
111  typedef typename Eigen_dense_matrix<ScalarType>::type EigenMatrix;
112  typedef typename Eigen_vector<ScalarType>::type EigenVector;
113 
117  EigenMatrix eigen_densemat(6, 5);
118  EigenMatrix eigen_densemat2(6, 5);
119  eigen_densemat(0,0) = 2.0; eigen_densemat(0,1) = -1.0;
120  eigen_densemat(1,0) = -1.0; eigen_densemat(1,1) = 2.0; eigen_densemat(1,2) = -1.0;
121  eigen_densemat(2,1) = -1.0; eigen_densemat(2,2) = -1.0; eigen_densemat(2,3) = -1.0;
122  eigen_densemat(3,2) = -1.0; eigen_densemat(3,3) = 2.0; eigen_densemat(3,4) = -1.0;
123  eigen_densemat(5,4) = -1.0; eigen_densemat(4,4) = -1.0;
124 
128  Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat(6, 5);
129  Eigen::SparseMatrix<ScalarType, Eigen::RowMajor> eigen_sparsemat2(6, 5);
130  eigen_sparsemat.reserve(5*2);
131  eigen_sparsemat.insert(0,0) = 2.0; eigen_sparsemat.insert(0,1) = -1.0;
132  eigen_sparsemat.insert(1,1) = 2.0; eigen_sparsemat.insert(1,2) = -1.0;
133  eigen_sparsemat.insert(2,2) = -1.0; eigen_sparsemat.insert(2,3) = -1.0;
134  eigen_sparsemat.insert(3,3) = 2.0; eigen_sparsemat.insert(3,4) = -1.0;
135  eigen_sparsemat.insert(5,4) = -1.0;
136  //eigen_sparsemat.endFill();
137 
141  EigenVector eigen_rhs(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 
167  viennacl::copy(eigen_densemat, vcl_densemat);
168  viennacl::copy(eigen_sparsemat, vcl_sparsemat);
169  std::cout << "VCL sparsematrix dimensions: " << vcl_sparsemat.size1() << ", " << vcl_sparsemat.size2() << std::endl;
170 
171  // For completeness: Copy matrices from ViennaCL back to Eigen:
172  viennacl::copy(vcl_densemat, eigen_densemat2);
173  viennacl::copy(vcl_sparsemat, eigen_sparsemat2);
174 
175 
179  eigen_result = eigen_densemat * eigen_rhs;
180  vcl_result = viennacl::linalg::prod(vcl_densemat, vcl_rhs);
181  viennacl::copy(vcl_result, eigen_temp);
182  std::cout << "Difference for dense matrix-vector product: " << (eigen_result - eigen_temp).norm() << std::endl;
183  std::cout << "Difference for dense matrix-vector product (Eigen->ViennaCL->Eigen): "
184  << (eigen_densemat2 * eigen_rhs - eigen_temp).norm() << std::endl;
185 
189  eigen_result = eigen_sparsemat * eigen_rhs;
190  vcl_result = viennacl::linalg::prod(vcl_sparsemat, vcl_rhs);
191  viennacl::copy(vcl_result, eigen_temp);
192  std::cout << "Difference for sparse matrix-vector product: " << (eigen_result - eigen_temp).norm() << std::endl;
193  std::cout << "Difference for sparse matrix-vector product (Eigen->ViennaCL->Eigen): "
194  << (eigen_sparsemat2 * eigen_rhs - eigen_temp).norm() << std::endl;
195 }
196 
197 
201 int main(int, char *[])
202 {
203  std::cout << "----------------------------------------------" << std::endl;
204  std::cout << "## Single precision" << std::endl;
205  std::cout << "----------------------------------------------" << std::endl;
206  run_tutorial<float>();
207 
208 #ifdef VIENNACL_HAVE_OPENCL
210 #endif
211  {
212  std::cout << "----------------------------------------------" << std::endl;
213  std::cout << "## Double precision" << std::endl;
214  std::cout << "----------------------------------------------" << std::endl;
215  run_tutorial<double>();
216  }
217 
221  std::cout << std::endl;
222  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
223  std::cout << std::endl;
224 
225 }
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:160
A dense matrix class.
Definition: forwards.h:374
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:91
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.