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
vandermonde_matrix.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_VANDERMONDE_MATRIX_HPP
2 #define VIENNACL_VANDERMONDE_MATRIX_HPP
3 
4 /* =========================================================================
5  Copyright (c) 2010-2016, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
21 #include <cmath>
22 
27 #include "viennacl/forwards.h"
28 #include "viennacl/vector.hpp"
29 #include "viennacl/ocl/backend.hpp"
30 
31 #include "viennacl/fft.hpp"
32 
34 
35 namespace viennacl
36 {
37 
43 template<class NumericT, unsigned int AlignmentV>
44 class vandermonde_matrix
45 {
46 public:
49 
54  explicit vandermonde_matrix() {}
55 
62  explicit vandermonde_matrix(vcl_size_t rows, vcl_size_t cols) : elements_(rows)
63  {
64  assert(rows == cols && bool("Vandermonde matrix must be square in this release!"));
65  (void)cols; // avoid 'unused parameter' warning in optimized builds
66  }
67 
74  void resize(vcl_size_t sz, bool preserve = true)
75  {
76  elements_.resize(sz, preserve);
77  }
78 
83  handle_type const & handle() const { return elements_.handle(); }
84 
90  viennacl::vector<NumericT, AlignmentV> const & elements() const { return elements_; }
91 
95  vcl_size_t size1() const { return elements_.size(); }
96 
100  vcl_size_t size2() const { return elements_.size(); }
101 
107  vcl_size_t internal_size() const { return elements_.internal_size(); }
108 
116  {
117  return elements_[row_index];
118  }
119 
127  NumericT operator()(vcl_size_t row_index, vcl_size_t col_index) const
128  {
129  assert(row_index < size1() && col_index < size2() && bool("Invalid access"));
130 
131  return pow(elements_[row_index], static_cast<int>(col_index));
132  }
133 
134 private:
136  vandermonde_matrix & operator=(vandermonde_matrix const & t);
137 
139 };
140 
147 template<typename NumericT, unsigned int AlignmentV>
148 void copy(std::vector<NumericT>& cpu_vec, vandermonde_matrix<NumericT, AlignmentV>& gpu_mat)
149 {
150  assert(cpu_vec.size() == gpu_mat.size1() && bool("Size mismatch"));
151  copy(cpu_vec, gpu_mat.elements());
152 }
153 
160 template<typename NumericT, unsigned int AlignmentV>
161 void copy(vandermonde_matrix<NumericT, AlignmentV>& gpu_mat, std::vector<NumericT>& cpu_vec)
162 {
163  assert(cpu_vec.size() == gpu_mat.size1() && bool("Size mismatch"));
164  copy(gpu_mat.elements(), cpu_vec);
165 }
166 
173 template<typename NumericT, unsigned int AlignmentV, typename MatrixT>
174 void copy(vandermonde_matrix<NumericT, AlignmentV>& vander_src, MatrixT& com_dst)
175 {
176  assert(vander_src.size1() == viennacl::traits::size1(com_dst) && bool("Size mismatch"));
177  assert(vander_src.size2() == viennacl::traits::size2(com_dst) && bool("Size mismatch"));
178 
179  vcl_size_t size = vander_src.size1();
180  std::vector<NumericT> tmp(size);
181  copy(vander_src, tmp);
182 
183  for (vcl_size_t i = 0; i < size; i++)
184  for (vcl_size_t j = 0; j < size; j++)
185  com_dst(i, j) = std::pow(tmp[i], static_cast<int>(j));
186 
187 }
188 
195 template<typename NumericT, unsigned int AlignmentV, typename MatrixT>
196 void copy(MatrixT& com_src, vandermonde_matrix<NumericT, AlignmentV>& vander_dst)
197 {
198  assert( (vander_dst.size1() == 0 || vander_dst.size1() == viennacl::traits::size1(com_src)) && bool("Size mismatch"));
199  assert( (vander_dst.size2() == 0 || vander_dst.size2() == viennacl::traits::size2(com_src)) && bool("Size mismatch"));
200 
201  vcl_size_t size = vander_dst.size1();
202  std::vector<NumericT> tmp(size);
203 
204  for (vcl_size_t i = 0; i < size; i++)
205  tmp[i] = com_src(i, 1);
206 
207  copy(tmp, vander_dst);
208 }
209 
210 /*template<typename NumericT, unsigned int AlignmentV, unsigned int VECTOR_AlignmentV>
211  void prod_impl(vandermonde_matrix<NumericT, AlignmentV>& mat,
212  vector<NumericT, VECTOR_AlignmentV>& vec,
213  vector<NumericT, VECTOR_AlignmentV>& result) {
214  assert(mat.size1() == vec.size());
215 
216  fft::vandermonde_prod<NumericT>(mat.handle(), vec.handle(), result.handle(), mat.size1());
217  } */
218 
224 template<class NumericT, unsigned int AlignmentV>
225 std::ostream & operator<<(std::ostream& s, vandermonde_matrix<NumericT, AlignmentV>& gpu_matrix)
226 {
227  vcl_size_t size = gpu_matrix.size1();
228  std::vector<NumericT> tmp(size);
229  copy(gpu_matrix, tmp);
230  s << "[" << size << "," << size << "](\n";
231 
232  for (vcl_size_t i = 0; i < size; i++)
233  {
234  s << "(";
235  for (vcl_size_t j = 0; j < size; j++)
236  {
237  s << pow(tmp[i], static_cast<NumericT>(j));
238  if (j < (size - 1))
239  s << ",";
240  }
241  s << ")";
242  }
243  s << ")";
244  return s;
245 }
246 
247 
248 //
249 // Specify available operations:
250 //
251 
254 namespace linalg
255 {
256 namespace detail
257 {
258  // x = A * y
259  template<typename T, unsigned int A>
260  struct op_executor<vector_base<T>, op_assign, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> >
261  {
262  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
263  {
264  // check for the special case x = A * x
265  if (viennacl::traits::handle(lhs) == viennacl::traits::handle(rhs.rhs()))
266  {
267  viennacl::vector<T> temp(lhs);
268  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
269  lhs = temp;
270  }
271  else
272  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), lhs);
273  }
274  };
275 
276  template<typename T, unsigned int A>
277  struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> >
278  {
279  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
280  {
281  viennacl::vector<T> temp(lhs);
282  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
283  lhs += temp;
284  }
285  };
286 
287  template<typename T, unsigned int A>
288  struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> >
289  {
290  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
291  {
292  viennacl::vector<T> temp(lhs);
293  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
294  lhs -= temp;
295  }
296  };
297 
298 
299  // x = A * vec_op
300  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
301  struct op_executor<vector_base<T>, op_assign, vector_expression<const vandermonde_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> >
302  {
303  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
304  {
305  viennacl::vector<T> temp(rhs.rhs());
306  viennacl::linalg::prod_impl(rhs.lhs(), temp, lhs);
307  }
308  };
309 
310  // x = A * vec_op
311  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
312  struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const vandermonde_matrix<T, A>, vector_expression<const LHS, const RHS, OP>, op_prod> >
313  {
314  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
315  {
316  viennacl::vector<T> temp(rhs.rhs());
317  viennacl::vector<T> temp_result(lhs);
318  viennacl::linalg::prod_impl(rhs.lhs(), temp, temp_result);
319  lhs += temp_result;
320  }
321  };
322 
323  // x = A * vec_op
324  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
325  struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const vandermonde_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> >
326  {
327  static void apply(vector_base<T> & lhs, vector_expression<const vandermonde_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
328  {
329  viennacl::vector<T> temp(rhs.rhs());
330  viennacl::vector<T> temp_result(lhs);
331  viennacl::linalg::prod_impl(rhs.lhs(), temp, temp_result);
332  lhs -= temp_result;
333  }
334  };
335 
336 } // namespace detail
337 } // namespace linalg
338 
340 }
341 
342 #endif // VIENNACL_VANDERMONDE_MATRIX_HPP
This class represents a single scalar value on the GPU and behaves mostly like a built-in scalar type...
Definition: forwards.h:227
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:163
viennacl::backend::mem_handle handle_type
scalar< typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT< NumericT >::ResultType > value_type
This file provides the forward declarations for the main types used within ViennaCL.
vandermonde_matrix(vcl_size_t rows, vcl_size_t cols)
Creates the matrix with the given size.
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
Definition: size.hpp:201
handle_type const & handle() const
Returns the OpenCL handle.
float NumericT
Definition: bisect.cpp:40
A Vandermonde matrix class.
Definition: forwards.h:418
Implementations of operations using vandermonde_matrix. Experimental.
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:239
void resize(size_type new_size, bool preserve=true)
Resizes the allocated memory for the vector. Pads the memory to be a multiple of 'AlignmentV'.
Definition: vector.hpp:1046
NumericT operator()(vcl_size_t row_index, vcl_size_t col_index) const
Read access to a element of the matrix.
std::size_t vcl_size_t
Definition: forwards.h:75
entry_proxy< NumericT > operator()(vcl_size_t row_index)
Read-write access to a base element of the matrix.
Implementations of the OpenCL backend, where all contexts are stored in.
vcl_size_t internal_size() const
Returns the internal size of matrix representtion. Usually required for launching OpenCL kernels only...
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) ...
size_type size() const
Returns the length of the vector (cf. std::vector)
Definition: vector_def.hpp:118
All routines related to the Fast Fourier Transform. Experimental.
vcl_size_t size2() const
Returns the number of columns of the matrix.
Main abstraction class for multiple memory domains. Represents a buffer in either main RAM...
Definition: mem_handle.hpp:89
viennacl::vector< NumericT, AlignmentV > & elements()
Returns an internal viennacl::vector, which represents a Vandermonde matrix elements.
viennacl::vector< NumericT, AlignmentV > const & elements() const
vandermonde_matrix()
The default constructor. Does not allocate any memory.
void resize(vcl_size_t sz, bool preserve=true)
Resizes the matrix. Existing entries can be preserved.
void prod_impl(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
size_type internal_size() const
Returns the internal length of the vector, which is given by size() plus the extra memory due to padd...
Definition: vector_def.hpp:120
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
Definition: handle.hpp:41
vcl_size_t size1() const
Returns the number of rows of the matrix.
A proxy class for a single element of a vector or matrix. This proxy should not be noticed by end-use...
Definition: forwards.h:233
const handle_type & handle() const
Returns the memory handle.
Definition: vector_def.hpp:128