ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
structured-matrices.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 
18 
19 
24 #include <iostream>
25 #include <vector>
26 #include <cmath>
27 #include <complex>
28 #include <fstream>
29 
30 //#define VIENNACL_BUILD_INFO
31 
32 //#define VIENNACL_DEBUG_ALL
33 
38 #include "viennacl/linalg/prod.hpp"
39 
40 #include "viennacl/fft.hpp"
41 
42 //
43 // A simple dense matrix class (in order to avoid an unnecessary boost dependency)
44 //
45 template<typename T>
47 {
48  public:
49  typedef std::size_t size_type;
50 
51  dense_matrix(std::size_t rows, std::size_t cols) : elements_(rows * cols), rows_(rows), cols_(cols) {}
52 
53  T & operator()(std::size_t i, std::size_t j) { return elements_[i*cols_ + j]; }
54  T const & operator()(std::size_t i, std::size_t j) const { return elements_[i*cols_ + j]; }
55 
56  std::size_t size1() const { return rows_; }
57  std::size_t size2() const { return cols_; }
58 
60  {
61  for (std::size_t i = 0; i < other.size1(); i++)
62  for (std::size_t j = 0; j < other.size2(); j++)
63  elements_[i*cols_ + j] = other.elements_[i*cols_+j];
64  return *this;
65  }
66 
67  private:
68  std::vector<T> elements_;
69  std::size_t rows_;
70  std::size_t cols_;
71 };
72 
73 template<typename T>
74 std::ostream & operator<<(std::ostream & os, dense_matrix<T> const & mat)
75 {
76  std::cout << "[" << mat.size1() << "," << mat.size2() << "](";
77  for (std::size_t i=0; i<mat.size1(); ++i)
78  {
79  std::cout << "(";
80  for (std::size_t j=0; j<mat.size2(); ++j)
81  std::cout << mat(i,j) << ",";
82  std::cout << ")";
83  }
84 
85  return os;
86 }
87 
88 
89 template<typename ScalarType>
91 {
92  ScalarType df = 0.0;
93  ScalarType d1 = 0;
94  ScalarType d2 = 0;
95 
96  for (std::size_t i = 0; i < m1.size1(); i++)
97  for (std::size_t j = 0; j < m1.size2(); j++)
98  {
99  df += (m1(i,j) - m2(i,j)) * (m1(i,j) - m2(i,j));
100  d1 += m1(i,j) * m1(i,j);
101  d2 += m2(i,j) * m2(i,j);
102  }
103 
104  if ( d1 + d2 <= 0 )
105  return 0;
106 
107  return std::sqrt(df / std::max<ScalarType>(d1, d2));
108 }
109 
110 
111 template<typename ScalarType>
112 ScalarType diff(std::vector<ScalarType>& vec, std::vector<ScalarType>& ref)
113 {
114  ScalarType df = 0.0;
115  ScalarType norm_ref = 0;
116 
117  for (std::size_t i = 0; i < vec.size(); i++)
118  {
119  df = df + pow(vec[i] - ref[i], 2);
120  norm_ref += ref[i] * ref[i];
121  }
122 
123  return std::sqrt(df / norm_ref);
124 }
125 
126 template<typename ScalarType>
127 ScalarType diff_max(std::vector<ScalarType>& vec, std::vector<ScalarType>& ref)
128 {
129  ScalarType df = 0.0;
130  ScalarType mx = 0.0;
131  ScalarType norm_max = 0;
132 
133  for (std::size_t i = 0; i < vec.size(); i++)
134  {
135  df = std::max<ScalarType>(std::fabs(vec[i] - ref[i]), df);
136  mx = std::max<ScalarType>(std::fabs(vec[i]), mx);
137 
138  if (mx > 0)
139  {
140  if (norm_max < df / mx)
141  norm_max = df / mx;
142  }
143  }
144 
145  return norm_max;
146 }
147 
148 
149 template<typename ScalarType>
151 {
152  int w = 5, h = 7;
153  std::vector<ScalarType> s_normal(2 * w * h);
154  viennacl::matrix<ScalarType> normal(w, 2 * h);
155  viennacl::matrix<ScalarType> transp(h, 2 * w);
156 
157  for (unsigned int i = 0; i < s_normal.size(); i+=2) {
158  s_normal[i] = i;
159  s_normal[i+1] = i;
160  }
161  viennacl::fast_copy(&s_normal[0], &s_normal[0] + s_normal.size(), normal);
162  std::cout << normal << std::endl;
164  std::cout << normal << std::endl;
165 }
166 
167 
168 
169 template<typename ScalarType>
171 {
172  std::size_t TOEPLITZ_SIZE = 47;
173  viennacl::toeplitz_matrix<ScalarType> vcl_toeplitz1(TOEPLITZ_SIZE, TOEPLITZ_SIZE);
174  viennacl::toeplitz_matrix<ScalarType> vcl_toeplitz2(TOEPLITZ_SIZE, TOEPLITZ_SIZE);
175 
176  viennacl::vector<ScalarType> vcl_input(TOEPLITZ_SIZE);
177  viennacl::vector<ScalarType> vcl_result(TOEPLITZ_SIZE);
178 
179  std::vector<ScalarType> input_ref(TOEPLITZ_SIZE);
180  std::vector<ScalarType> result_ref(TOEPLITZ_SIZE);
181 
182  dense_matrix<ScalarType> m1(TOEPLITZ_SIZE, TOEPLITZ_SIZE);
183  dense_matrix<ScalarType> m2(TOEPLITZ_SIZE, TOEPLITZ_SIZE);
184 
185  for (std::size_t i = 0; i < TOEPLITZ_SIZE; i++)
186  for (std::size_t j = 0; j < TOEPLITZ_SIZE; j++)
187  {
188  m1(i,j) = static_cast<ScalarType>(i) - static_cast<ScalarType>(j);
189  m2(i,j) = m1(i,j) * m1(i,j) + ScalarType(1);
190  }
191 
192  for (std::size_t i = 0; i < TOEPLITZ_SIZE; i++)
193  input_ref[i] = ScalarType(i);
194 
195  // Copy to ViennaCL
196  viennacl::copy(m1, vcl_toeplitz1);
197  viennacl::copy(m2, vcl_toeplitz2);
198  viennacl::copy(input_ref, vcl_input);
199 
200  //
201  // Matrix-Vector product:
202  //
203  vcl_result = viennacl::linalg::prod(vcl_toeplitz1, vcl_input);
204 
205  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
206  {
207  ScalarType entry = 0;
208  for (std::size_t j = 0; j < m1.size2(); j++)
209  entry += m1(i,j) * input_ref[j];
210 
211  result_ref[i] = entry;
212  }
213 
214  viennacl::copy(vcl_result, input_ref);
215  std::cout << "Matrix-Vector Product: " << diff_max(input_ref, result_ref);
216  if (diff_max(input_ref, result_ref) < epsilon)
217  std::cout << " [OK]" << std::endl;
218  else
219  {
220  for (std::size_t i=0; i<input_ref.size(); ++i)
221  std::cout << "Should: " << result_ref[i] << ", is: " << input_ref[i] << std::endl;
222  std::cout << " [FAILED]" << std::endl;
223  return EXIT_FAILURE;
224  }
225 
226 
227  //
228  // Matrix addition:
229  //
230  vcl_toeplitz1 += vcl_toeplitz2;
231 
232  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
233  for (std::size_t j = 0; j < m1.size2(); j++)
234  m1(i,j) += m2(i,j);
235 
236  viennacl::copy(vcl_toeplitz1, m2);
237  std::cout << "Matrix Addition: " << diff(m1, m2);
238  if (diff(m1, m2) < epsilon)
239  std::cout << " [OK]" << std::endl;
240  else
241  {
242  std::cout << " [FAILED]" << std::endl;
243  return EXIT_FAILURE;
244  }
245 
246  //
247  // Per-Element access:
248  //
249  vcl_toeplitz1(2,4) = 42;
250 
251  for (std::size_t i=0; i<m1.size1(); ++i) //reference calculation
252  {
253  if (i + 2 < m1.size2())
254  m1(i, i+2) = 42;
255  }
256 
257  viennacl::copy(vcl_toeplitz1, m2);
258  std::cout << "Element manipulation: " << diff(m1, m2);
259  if (diff(m1, m2) < epsilon)
260  std::cout << " [OK]" << std::endl;
261  else
262  {
263  std::cout << " [FAILED]" << std::endl;
264  return EXIT_FAILURE;
265  }
266 
267  return EXIT_SUCCESS;
268 }
269 
270 template<typename ScalarType>
272 {
273  std::size_t CIRCULANT_SIZE = 53;
274  viennacl::circulant_matrix<ScalarType> vcl_circulant1(CIRCULANT_SIZE, CIRCULANT_SIZE);
275  viennacl::circulant_matrix<ScalarType> vcl_circulant2(CIRCULANT_SIZE, CIRCULANT_SIZE);
276 
277  viennacl::vector<ScalarType> vcl_input(CIRCULANT_SIZE);
278  viennacl::vector<ScalarType> vcl_result(CIRCULANT_SIZE);
279 
280  std::vector<ScalarType> input_ref(CIRCULANT_SIZE);
281  std::vector<ScalarType> result_ref(CIRCULANT_SIZE);
282 
283  dense_matrix<ScalarType> m1(vcl_circulant1.size1(), vcl_circulant1.size2());
284  dense_matrix<ScalarType> m2(vcl_circulant1.size1(), vcl_circulant1.size2());
285 
286  for (std::size_t i = 0; i < m1.size1(); i++)
287  for (std::size_t j = 0; j < m1.size2(); j++)
288  {
289  m1(i,j) = static_cast<ScalarType>((i - j + m1.size1()) % m1.size1());
290  m2(i,j) = m1(i,j) * m1(i,j) + ScalarType(1);
291  }
292 
293  for (std::size_t i = 0; i < input_ref.size(); i++)
294  input_ref[i] = ScalarType(i);
295 
296  // Copy to ViennaCL
297  viennacl::copy(m1, vcl_circulant1);
298  viennacl::copy(m2, vcl_circulant2);
299  viennacl::copy(input_ref, vcl_input);
300 
301  //
302  // Matrix-Vector product:
303  //
304  vcl_result = viennacl::linalg::prod(vcl_circulant1, vcl_input);
305 
306  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
307  {
308  ScalarType entry = 0;
309  for (std::size_t j = 0; j < m1.size2(); j++)
310  entry += m1(i,j) * input_ref[j];
311 
312  result_ref[i] = entry;
313  }
314 
315  viennacl::copy(vcl_result, input_ref);
316  std::cout << "Matrix-Vector Product: " << diff_max(input_ref, result_ref);
317  if (diff_max(input_ref, result_ref) < epsilon)
318  std::cout << " [OK]" << std::endl;
319  else
320  {
321  for (std::size_t i=0; i<input_ref.size(); ++i)
322  std::cout << "Should: " << result_ref[i] << ", is: " << input_ref[i] << std::endl;
323  std::cout << " [FAILED]" << std::endl;
324  return EXIT_FAILURE;
325  }
326 
327 
328  //
329  // Matrix addition:
330  //
331  vcl_circulant1 += vcl_circulant2;
332 
333  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
334  for (std::size_t j = 0; j < m1.size2(); j++)
335  m1(i,j) += m2(i,j);
336 
337  viennacl::copy(vcl_circulant1, m2);
338  std::cout << "Matrix Addition: " << diff(m1, m2);
339  if (diff(m1, m2) < epsilon)
340  std::cout << " [OK]" << std::endl;
341  else
342  {
343  std::cout << " [FAILED]" << std::endl;
344  return EXIT_FAILURE;
345  }
346 
347  //
348  // Per-Element access:
349  //
350  vcl_circulant1(4,2) = 42;
351 
352  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
353  for (std::size_t j = 0; j < m1.size2(); j++)
354  {
355  if ((i - j + m1.size1()) % m1.size1() == 2)
356  m1(i, j) = 42;
357  }
358 
359  viennacl::copy(vcl_circulant1, m2);
360  std::cout << "Element manipulation: " << diff(m1, m2);
361  if (diff(m1, m2) < epsilon)
362  std::cout << " [OK]" << std::endl;
363  else
364  {
365  std::cout << " [FAILED]" << std::endl;
366  return EXIT_FAILURE;
367  }
368 
369  return EXIT_SUCCESS;
370 }
371 
372 template<typename ScalarType>
374 {
375  std::size_t VANDERMONDE_SIZE = 61;
376 
377  viennacl::vandermonde_matrix<ScalarType> vcl_vandermonde1(VANDERMONDE_SIZE, VANDERMONDE_SIZE);
378  viennacl::vandermonde_matrix<ScalarType> vcl_vandermonde2(VANDERMONDE_SIZE, VANDERMONDE_SIZE);
379 
380  viennacl::vector<ScalarType> vcl_input(VANDERMONDE_SIZE);
381  viennacl::vector<ScalarType> vcl_result(VANDERMONDE_SIZE);
382 
383  std::vector<ScalarType> input_ref(VANDERMONDE_SIZE);
384  std::vector<ScalarType> result_ref(VANDERMONDE_SIZE);
385 
386  dense_matrix<ScalarType> m1(vcl_vandermonde1.size1(), vcl_vandermonde1.size2());
387  dense_matrix<ScalarType> m2(m1.size1(), m1.size2());
388 
389  for (std::size_t i = 0; i < m1.size1(); i++)
390  for (std::size_t j = 0; j < m1.size2(); j++)
391  {
392  m1(i,j) = std::pow(ScalarType(1.0) + ScalarType(i)/ScalarType(1000.0), ScalarType(j));
393  m2(i,j) = std::pow(ScalarType(1.0) - ScalarType(i)/ScalarType(2000.0), ScalarType(j));
394  }
395 
396  for (std::size_t i = 0; i < input_ref.size(); i++)
397  input_ref[i] = ScalarType(i);
398 
399  // Copy to ViennaCL
400  viennacl::copy(m1, vcl_vandermonde1);
401  viennacl::copy(m2, vcl_vandermonde2);
402  viennacl::copy(input_ref, vcl_input);
403 
404  //
405  // Matrix-Vector product:
406  //
407  vcl_result = viennacl::linalg::prod(vcl_vandermonde1, vcl_input);
408 
409  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
410  {
411  ScalarType entry = 0;
412  for (std::size_t j = 0; j < m1.size2(); j++)
413  entry += m1(i,j) * input_ref[j];
414 
415  result_ref[i] = entry;
416  }
417 
418  viennacl::copy(vcl_result, input_ref);
419  std::cout << "Matrix-Vector Product: " << diff_max(input_ref, result_ref);
420  if (diff_max(input_ref, result_ref) < epsilon)
421  std::cout << " [OK]" << std::endl;
422  else
423  {
424  for (std::size_t i=0; i<input_ref.size(); ++i)
425  std::cout << "Should: " << result_ref[i] << ", is: " << input_ref[i] << std::endl;
426  std::cout << " [FAILED]" << std::endl;
427  return EXIT_FAILURE;
428  }
429 
430 
431  //
432  // Note: Matrix addition does not make sense for a Vandermonde matrix
433  //
434 
435 
436  //
437  // Per-Element access:
438  //
439  vcl_vandermonde1(4) = static_cast<ScalarType>(1.0001);
440 
441  for (std::size_t j = 0; j < m1.size2(); j++)
442  {
443  m1(4, j) = std::pow(ScalarType(1.0001), ScalarType(j));
444  }
445 
446  viennacl::copy(vcl_vandermonde1, m2);
447  std::cout << "Element manipulation: " << diff(m1, m2);
448  if (diff(m1, m2) < epsilon)
449  std::cout << " [OK]" << std::endl;
450  else
451  {
452  std::cout << " [FAILED]" << std::endl;
453  return EXIT_FAILURE;
454  }
455 
456  return EXIT_SUCCESS;
457 }
458 
459 template<typename ScalarType>
461 {
462  std::size_t HANKEL_SIZE = 7;
463  viennacl::hankel_matrix<ScalarType> vcl_hankel1(HANKEL_SIZE, HANKEL_SIZE);
464  viennacl::hankel_matrix<ScalarType> vcl_hankel2(HANKEL_SIZE, HANKEL_SIZE);
465 
466  viennacl::vector<ScalarType> vcl_input(HANKEL_SIZE);
467  viennacl::vector<ScalarType> vcl_result(HANKEL_SIZE);
468 
469  std::vector<ScalarType> input_ref(HANKEL_SIZE);
470  std::vector<ScalarType> result_ref(HANKEL_SIZE);
471 
472  dense_matrix<ScalarType> m1(vcl_hankel1.size1(), vcl_hankel1.size2());
473  dense_matrix<ScalarType> m2(m1.size1(), m1.size2());
474 
475  for (std::size_t i = 0; i < m1.size1(); i++)
476  for (std::size_t j = 0; j < m1.size2(); j++)
477  {
478  m1(i,j) = static_cast<ScalarType>((i + j) % (2 * m1.size1()));
479  m2(i,j) = m1(i,j) * m1(i,j) + ScalarType(1);
480  }
481 
482  for (std::size_t i = 0; i < input_ref.size(); i++)
483  input_ref[i] = ScalarType(i);
484 
485  // Copy to ViennaCL
486  viennacl::copy(m1, vcl_hankel1);
487  viennacl::copy(m2, vcl_hankel2);
488  viennacl::copy(input_ref, vcl_input);
489 
490  //
491  // Matrix-Vector product:
492  //
493  vcl_result = viennacl::linalg::prod(vcl_hankel1, vcl_input);
494 
495  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
496  {
497  ScalarType entry = 0;
498  for (std::size_t j = 0; j < m1.size2(); j++)
499  entry += m1(i,j) * input_ref[j];
500 
501  result_ref[i] = entry;
502  }
503 
504  viennacl::copy(vcl_result, input_ref);
505  std::cout << "Matrix-Vector Product: " << diff_max(input_ref, result_ref);
506  if (diff_max(input_ref, result_ref) < epsilon)
507  std::cout << " [OK]" << std::endl;
508  else
509  {
510  for (std::size_t i=0; i<input_ref.size(); ++i)
511  std::cout << "Should: " << result_ref[i] << ", is: " << input_ref[i] << std::endl;
512  std::cout << " [FAILED]" << std::endl;
513  return EXIT_FAILURE;
514  }
515 
516 
517  //
518  // Matrix addition:
519  //
520  vcl_hankel1 += vcl_hankel2;
521 
522  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
523  for (std::size_t j = 0; j < m1.size2(); j++)
524  m1(i,j) += m2(i,j);
525 
526  viennacl::copy(vcl_hankel1, m2);
527  std::cout << "Matrix Addition: " << diff(m1, m2);
528  if (diff(m1, m2) < epsilon)
529  std::cout << " [OK]" << std::endl;
530  else
531  {
532  std::cout << " [FAILED]" << std::endl;
533  return EXIT_FAILURE;
534  }
535 
536  //
537  // Per-Element access:
538  //
539  vcl_hankel1(4,2) = 42;
540 
541  for (std::size_t i = 0; i < m1.size1(); i++) //reference calculation
542  for (std::size_t j = 0; j < m1.size2(); j++)
543  {
544  if ((i + j) % (2*m1.size1()) == 6)
545  m1(i, j) = 42;
546  }
547 
548  viennacl::copy(vcl_hankel1, m2);
549  std::cout << "Element manipulation: " << diff(m1, m2);
550  if (diff(m1, m2) < epsilon)
551  std::cout << " [OK]" << std::endl;
552  else
553  {
554  std::cout << " [FAILED]" << std::endl;
555  return EXIT_FAILURE;
556  }
557 
558  return EXIT_SUCCESS;
559 }
560 
561 int main()
562 {
563  std::cout << std::endl;
564  std::cout << "----------------------------------------------" << std::endl;
565  std::cout << "----------------------------------------------" << std::endl;
566  std::cout << "## Test :: Structured Matrices" << std::endl;
567  std::cout << "----------------------------------------------" << std::endl;
568  std::cout << "----------------------------------------------" << std::endl;
569  std::cout << std::endl;
570 
571  double eps = 1e-3;
572 
573  std::cout << "# Testing setup:" << std::endl;
574  std::cout << " eps: " << eps << std::endl;
575  std::cout << " numeric: float" << std::endl;
576  std::cout << std::endl;
577  std::cout << " -- Vandermonde matrix -- " << std::endl;
578  if (vandermonde_test<float>(static_cast<float>(eps)) == EXIT_FAILURE)
579  return EXIT_FAILURE;
580 
581  std::cout << " -- Circulant matrix -- " << std::endl;
582  if (circulant_test<float>(static_cast<float>(eps)) == EXIT_FAILURE)
583  return EXIT_FAILURE;
584 
585  std::cout << " -- Toeplitz matrix -- " << std::endl;
586  if (toeplitz_test<float>(static_cast<float>(eps)) == EXIT_FAILURE)
587  return EXIT_FAILURE;
588 
589  std::cout << " -- Hankel matrix -- " << std::endl;
590  if (hankel_test<float>(static_cast<float>(eps)) == EXIT_FAILURE)
591  return EXIT_FAILURE;
592 
593 
594  std::cout << std::endl;
595 
596  if ( viennacl::ocl::current_device().double_support() )
597  {
598  eps = 1e-10;
599 
600  std::cout << std::endl;
601  std::cout << "# Testing setup:" << std::endl;
602  std::cout << " eps: " << eps << std::endl;
603  std::cout << " numeric: double" << std::endl;
604  std::cout << std::endl;
605 
606  std::cout << " -- Vandermonde matrix -- " << std::endl;
607  if (vandermonde_test<double>(eps) == EXIT_FAILURE)
608  return EXIT_FAILURE;
609 
610  std::cout << " -- Circulant matrix -- " << std::endl;
611  if (circulant_test<double>(eps) == EXIT_FAILURE)
612  return EXIT_FAILURE;
613 
614  std::cout << " -- Toeplitz matrix -- " << std::endl;
615  if (toeplitz_test<double>(eps) == EXIT_FAILURE)
616  return EXIT_FAILURE;
617 
618  std::cout << " -- Hankel matrix -- " << std::endl;
619  if (hankel_test<double>(eps) == EXIT_FAILURE)
620  return EXIT_FAILURE;
621  }
622 
623  std::cout << std::endl;
624  std::cout << "------- Test completed --------" << std::endl;
625  std::cout << std::endl;
626 
627 
628  return EXIT_SUCCESS;
629 }
dense_matrix(std::size_t rows, std::size_t cols)
std::size_t size2() const
int vandermonde_test(ScalarType epsilon)
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
ScalarType diff(dense_matrix< ScalarType > const &m1, dense_matrix< ScalarType > const &m2)
int hankel_test(ScalarType epsilon)
vcl_size_t size2() const
Returns the number of columns of the matrix.
vcl_size_t size2() const
Returns the number of columns of the matrix.
A Hankel matrix class.
Definition: forwards.h:412
A dense matrix class.
Definition: forwards.h:375
void transpose_test()
int toeplitz_test(ScalarType epsilon)
vcl_size_t size1() const
Returns the number of rows of the matrix.
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
Definition: backend.hpp:351
std::size_t size1() const
A Vandermonde matrix class.
Definition: forwards.h:418
dense_matrix & operator+=(dense_matrix const &other)
int circulant_test(ScalarType epsilon)
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
Implementation of the hankel_matrix class for efficient manipulation of Hankel matrices. Experimental.
Implementation of the circulant_matrix class for efficient manipulation of circulant matrices...
vcl_size_t size1() const
Returns the number of rows of the matrix.
ScalarType diff_max(std::vector< ScalarType > &vec, std::vector< ScalarType > &ref)
vcl_size_t size1() const
Returns the number of rows of the matrix.
T & operator()(std::size_t i, std::size_t j)
Implementation of the vandermonde_matrix class for efficient manipulation of Vandermonde matrices...
void transpose(viennacl::matrix< NumericT, viennacl::row_major, AlignmentV > &input)
Inplace_transpose matrix.
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) ...
float ScalarType
Definition: fft_1d.cpp:42
All routines related to the Fast Fourier Transform. Experimental.
A Toeplitz matrix class.
Definition: forwards.h:415
A Circulant matrix class.
viennacl::matrix< float > m1
vcl_size_t size2() const
Returns the number of columns of the matrix.
std::size_t size_type
int main()
T const & operator()(std::size_t i, std::size_t j) const
Implementation of the toeplitz_matrix class for efficient manipulation of Toeplitz matrices...
void fast_copy(const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_begin, const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_end, CPU_ITERATOR cpu_begin)