ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
matrix_operations.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_MATRIX_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_MATRIX_OPERATIONS_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2015, 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 
25 #include "viennacl/forwards.h"
26 #include "viennacl/scalar.hpp"
27 #include "viennacl/vector.hpp"
29 #include "viennacl/tools/tools.hpp"
33 #include "viennacl/traits/size.hpp"
37 #include "viennacl/vector.hpp"
39 
40 #ifdef VIENNACL_WITH_OPENCL
42 #endif
43 
44 #ifdef VIENNACL_WITH_CUDA
46 #endif
47 
48 namespace viennacl
49 {
50  namespace linalg
51  {
52 
53  template<typename DestNumericT, typename SrcNumericT>
55  {
56  assert(viennacl::traits::size1(dest) == viennacl::traits::size1(src) && bool("Incompatible matrix sizes in m1 = m2 (convert): size1(m1) != size1(m2)"));
57  assert(viennacl::traits::size2(dest) == viennacl::traits::size2(src) && bool("Incompatible matrix sizes in m1 = m2 (convert): size2(m1) != size2(m2)"));
58 
59  switch (viennacl::traits::handle(dest).get_active_handle_id())
60  {
63  break;
64 #ifdef VIENNACL_WITH_OPENCL
67  break;
68 #endif
69 #ifdef VIENNACL_WITH_CUDA
72  break;
73 #endif
75  throw memory_exception("not initialised!");
76  default:
77  throw memory_exception("not implemented");
78  }
79  }
80 
81  template<typename NumericT,
82  typename SizeT, typename DistanceT>
84  matrix_base<NumericT> & temp_trans)
85  {
86  switch (viennacl::traits::handle(proxy).get_active_handle_id())
87  {
89  viennacl::linalg::host_based::trans(proxy, temp_trans);
90  break;
91 #ifdef VIENNACL_WITH_OPENCL
93  viennacl::linalg::opencl::trans(proxy,temp_trans);
94  break;
95 #endif
96 #ifdef VIENNACL_WITH_CUDA
98  viennacl::linalg::cuda::trans(proxy,temp_trans);
99  break;
100 #endif
102  throw memory_exception("not initialised!");
103  default:
104  throw memory_exception("not implemented");
105  }
106  }
107 
108 
109  template<typename NumericT,
110  typename ScalarType1>
112  matrix_base<NumericT> const & mat2, ScalarType1 const & alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha)
113  {
114  switch (viennacl::traits::handle(mat1).get_active_handle_id())
115  {
117  viennacl::linalg::host_based::am(mat1, mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha);
118  break;
119 #ifdef VIENNACL_WITH_OPENCL
121  viennacl::linalg::opencl::am(mat1, mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha);
122  break;
123 #endif
124 #ifdef VIENNACL_WITH_CUDA
126  viennacl::linalg::cuda::am(mat1, mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha);
127  break;
128 #endif
130  throw memory_exception("not initialised!");
131  default:
132  throw memory_exception("not implemented");
133  }
134  }
135 
136 
137  template<typename NumericT,
138  typename ScalarType1, typename ScalarType2>
140  matrix_base<NumericT> const & mat2, ScalarType1 const & alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha,
141  matrix_base<NumericT> const & mat3, ScalarType2 const & beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
142  {
143  switch (viennacl::traits::handle(mat1).get_active_handle_id())
144  {
147  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
148  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
149  break;
150 #ifdef VIENNACL_WITH_OPENCL
153  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
154  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
155  break;
156 #endif
157 #ifdef VIENNACL_WITH_CUDA
160  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
161  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
162  break;
163 #endif
165  throw memory_exception("not initialised!");
166  default:
167  throw memory_exception("not implemented");
168  }
169  }
170 
171 
172  template<typename NumericT,
173  typename ScalarType1, typename ScalarType2>
175  matrix_base<NumericT> const & mat2, ScalarType1 const & alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha,
176  matrix_base<NumericT> const & mat3, ScalarType2 const & beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
177  {
178  switch (viennacl::traits::handle(mat1).get_active_handle_id())
179  {
182  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
183  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
184  break;
185 #ifdef VIENNACL_WITH_OPENCL
188  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
189  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
190  break;
191 #endif
192 #ifdef VIENNACL_WITH_CUDA
195  mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
196  mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
197  break;
198 #endif
200  throw memory_exception("not initialised!");
201  default:
202  throw memory_exception("not implemented");
203  }
204  }
205 
206 
207  template<typename NumericT>
208  void matrix_assign(matrix_base<NumericT> & mat, NumericT s, bool clear = false)
209  {
210  switch (viennacl::traits::handle(mat).get_active_handle_id())
211  {
214  break;
215 #ifdef VIENNACL_WITH_OPENCL
218  break;
219 #endif
220 #ifdef VIENNACL_WITH_CUDA
223  break;
224 #endif
226  throw memory_exception("not initialised!");
227  default:
228  throw memory_exception("not implemented");
229  }
230  }
231 
232 
233  template<typename NumericT>
235  {
236  switch (viennacl::traits::handle(mat).get_active_handle_id())
237  {
240  break;
241 #ifdef VIENNACL_WITH_OPENCL
244  break;
245 #endif
246 #ifdef VIENNACL_WITH_CUDA
249  break;
250 #endif
252  throw memory_exception("not initialised!");
253  default:
254  throw memory_exception("not implemented");
255  }
256  }
257 
258 
260  template<typename NumericT>
262  {
263  switch (viennacl::traits::handle(v).get_active_handle_id())
264  {
267  break;
268 #ifdef VIENNACL_WITH_OPENCL
271  break;
272 #endif
273 #ifdef VIENNACL_WITH_CUDA
276  break;
277 #endif
279  throw memory_exception("not initialised!");
280  default:
281  throw memory_exception("not implemented");
282  }
283  }
284 
286  template<typename NumericT>
288  {
289  switch (viennacl::traits::handle(A).get_active_handle_id())
290  {
293  break;
294 #ifdef VIENNACL_WITH_OPENCL
297  break;
298 #endif
299 #ifdef VIENNACL_WITH_CUDA
302  break;
303 #endif
305  throw memory_exception("not initialised!");
306  default:
307  throw memory_exception("not implemented");
308  }
309  }
310 
311  template<typename NumericT>
312  void matrix_row(const matrix_base<NumericT> & A, unsigned int i, vector_base<NumericT> & v)
313  {
314  switch (viennacl::traits::handle(A).get_active_handle_id())
315  {
318  break;
319 #ifdef VIENNACL_WITH_OPENCL
322  break;
323 #endif
324 #ifdef VIENNACL_WITH_CUDA
327  break;
328 #endif
330  throw memory_exception("not initialised!");
331  default:
332  throw memory_exception("not implemented");
333  }
334  }
335 
336  template<typename NumericT>
337  void matrix_column(const matrix_base<NumericT> & A, unsigned int j, vector_base<NumericT> & v)
338  {
339  switch (viennacl::traits::handle(A).get_active_handle_id())
340  {
343  break;
344 #ifdef VIENNACL_WITH_OPENCL
347  break;
348 #endif
349 #ifdef VIENNACL_WITH_CUDA
352  break;
353 #endif
355  throw memory_exception("not initialised!");
356  default:
357  throw memory_exception("not implemented");
358  }
359  }
360 
368  template<typename T>
370  scalar<T> & result)
371  {
372  typedef typename matrix_base<T>::handle_type HandleType;
373 
374  if ((A.start1() > 0) || (A.start2() > 0) || (A.stride1() > 1) || (A.stride2() > 1)) {
375  if (A.row_major()) {
377  viennacl::vector_base<T> temp(const_cast<HandleType &>(temp_A.handle()), temp_A.internal_size(), 0, 1);
378  norm_2_impl(temp, result);
379  } else {
381  viennacl::vector_base<T> temp(const_cast<HandleType &>(temp_A.handle()), temp_A.internal_size(), 0, 1);
382  norm_2_impl(temp, result);
383  }
384  } else {
385  viennacl::vector_base<T> temp(const_cast<HandleType &>(A.handle()), A.internal_size(), 0, 1);
386  norm_2_impl(temp, result);
387  }
388 
389  }
390 
398  template<typename T>
400  T & result)
401  {
402  typedef typename matrix_base<T>::handle_type HandleType;
403 
404  if ((A.start1() > 0) || (A.start2() > 0) || (A.stride1() > 1) || (A.stride2() > 1)) {
405  if (A.row_major()) {
407  viennacl::vector_base<T> temp(const_cast<HandleType &>(temp_A.handle()), temp_A.internal_size(), 0, 1);
408  norm_2_cpu(temp, result);
409  } else {
411  viennacl::vector_base<T> temp(const_cast<HandleType &>(temp_A.handle()), temp_A.internal_size(), 0, 1);
412  norm_2_cpu(temp, result);
413  }
414  } else {
415  viennacl::vector_base<T> temp(const_cast<HandleType &>(A.handle()), A.internal_size(), 0, 1);
416  norm_2_cpu(temp, result);
417  }
418 
419  }
420 
421  //
423  //
424 
425 
426 
427  // A * x
428 
437  template<typename NumericT>
439  const vector_base<NumericT> & vec,
440  vector_base<NumericT> & result)
441  {
442  assert( (viennacl::traits::size1(mat) == viennacl::traits::size(result)) && bool("Size check failed at v1 = prod(A, v2): size1(A) != size(v1)"));
443  assert( (viennacl::traits::size2(mat) == viennacl::traits::size(vec)) && bool("Size check failed at v1 = prod(A, v2): size2(A) != size(v2)"));
444 
445  switch (viennacl::traits::handle(mat).get_active_handle_id())
446  {
448  viennacl::linalg::host_based::prod_impl(mat, false, vec, result);
449  break;
450 #ifdef VIENNACL_WITH_OPENCL
452  viennacl::linalg::opencl::prod_impl(mat, false, vec, result);
453  break;
454 #endif
455 #ifdef VIENNACL_WITH_CUDA
457  viennacl::linalg::cuda::prod_impl(mat, false, vec, result);
458  break;
459 #endif
461  throw memory_exception("not initialised!");
462  default:
463  throw memory_exception("not implemented");
464  }
465  }
466 
467 
468  // trans(A) * x
469 
478  template<typename NumericT>
480  const vector_base<NumericT> & vec,
481  vector_base<NumericT> & result)
482  {
483  assert( (viennacl::traits::size1(mat_trans.lhs()) == viennacl::traits::size(vec)) && bool("Size check failed at v1 = trans(A) * v2: size1(A) != size(v2)"));
484  assert( (viennacl::traits::size2(mat_trans.lhs()) == viennacl::traits::size(result)) && bool("Size check failed at v1 = trans(A) * v2: size2(A) != size(v1)"));
485 
486  switch (viennacl::traits::handle(mat_trans.lhs()).get_active_handle_id())
487  {
489  viennacl::linalg::host_based::prod_impl(mat_trans.lhs(), true, vec, result);
490  break;
491 #ifdef VIENNACL_WITH_OPENCL
493  viennacl::linalg::opencl::prod_impl(mat_trans.lhs(), true, vec, result);
494  break;
495 #endif
496 #ifdef VIENNACL_WITH_CUDA
498  viennacl::linalg::cuda::prod_impl(mat_trans.lhs(), true, vec, result);
499  break;
500 #endif
502  throw memory_exception("not initialised!");
503  default:
504  throw memory_exception("not implemented");
505  }
506  }
507 
508 
509  //
511  //
512 
518  template<typename NumericT, typename ScalarType >
520  const matrix_base<NumericT> & B,
522  ScalarType alpha,
523  ScalarType beta)
524  {
525  assert( (viennacl::traits::size1(A) == viennacl::traits::size1(C)) && bool("Size check failed at C = prod(A, B): size1(A) != size1(C)"));
526  assert( (viennacl::traits::size2(A) == viennacl::traits::size1(B)) && bool("Size check failed at C = prod(A, B): size2(A) != size1(B)"));
527  assert( (viennacl::traits::size2(B) == viennacl::traits::size2(C)) && bool("Size check failed at C = prod(A, B): size2(B) != size2(C)"));
528 
529 
530  switch (viennacl::traits::handle(A).get_active_handle_id())
531  {
533  viennacl::linalg::host_based::prod_impl(A, false, B, false, C, alpha, beta);
534  break;
535 #ifdef VIENNACL_WITH_OPENCL
537  viennacl::linalg::opencl::prod_impl(A, false, B, false, C, alpha, beta);
538  break;
539 #endif
540 #ifdef VIENNACL_WITH_CUDA
542  viennacl::linalg::cuda::prod_impl(A, false, B, false, C, alpha, beta);
543  break;
544 #endif
546  throw memory_exception("not initialised!");
547  default:
548  throw memory_exception("not implemented");
549  }
550  }
551 
552 
553 
559  template<typename NumericT, typename ScalarType >
561  const matrix_base<NumericT>,
562  op_trans> & A,
563  const matrix_base<NumericT> & B,
565  ScalarType alpha,
566  ScalarType beta)
567  {
568  assert(viennacl::traits::size2(A.lhs()) == viennacl::traits::size1(C) && bool("Size check failed at C = prod(trans(A), B): size2(A) != size1(C)"));
569  assert(viennacl::traits::size1(A.lhs()) == viennacl::traits::size1(B) && bool("Size check failed at C = prod(trans(A), B): size1(A) != size1(B)"));
570  assert(viennacl::traits::size2(B) == viennacl::traits::size2(C) && bool("Size check failed at C = prod(trans(A), B): size2(B) != size2(C)"));
571 
572  switch (viennacl::traits::handle(A.lhs()).get_active_handle_id())
573  {
575  viennacl::linalg::host_based::prod_impl(A.lhs(), true, B, false, C, alpha, beta);
576  break;
577 #ifdef VIENNACL_WITH_OPENCL
579  viennacl::linalg::opencl::prod_impl(A.lhs(), true, B, false, C, alpha, beta);
580  break;
581 #endif
582 #ifdef VIENNACL_WITH_CUDA
584  viennacl::linalg::cuda::prod_impl(A.lhs(), true, B, false, C, alpha, beta);
585  break;
586 #endif
588  throw memory_exception("not initialised!");
589  default:
590  throw memory_exception("not implemented");
591  }
592  }
593 
594 
595 
596 
602  template<typename NumericT, typename ScalarType >
606  ScalarType alpha,
607  ScalarType beta)
608  {
609  assert(viennacl::traits::size1(A) == viennacl::traits::size1(C) && bool("Size check failed at C = prod(A, trans(B)): size1(A) != size1(C)"));
610  assert(viennacl::traits::size2(A) == viennacl::traits::size2(B.lhs()) && bool("Size check failed at C = prod(A, trans(B)): size2(A) != size2(B)"));
611  assert(viennacl::traits::size1(B.lhs()) == viennacl::traits::size2(C) && bool("Size check failed at C = prod(A, trans(B)): size1(B) != size2(C)"));
612 
614  {
616  viennacl::linalg::host_based::prod_impl(A, false, B.lhs(), true, C, alpha, beta);
617  break;
618 #ifdef VIENNACL_WITH_OPENCL
620  viennacl::linalg::opencl::prod_impl(A, false, B.lhs(), true, C, alpha, beta);
621  break;
622 #endif
623 #ifdef VIENNACL_WITH_CUDA
625  viennacl::linalg::cuda::prod_impl(A, false, B.lhs(), true, C, alpha, beta);
626  break;
627 #endif
629  throw memory_exception("not initialised!");
630  default:
631  throw memory_exception("not implemented");
632  }
633  }
634 
635 
636 
642  template<typename NumericT, typename ScalarType >
646  ScalarType alpha,
647  ScalarType beta)
648  {
649  assert(viennacl::traits::size2(A.lhs()) == viennacl::traits::size1(C) && bool("Size check failed at C = prod(trans(A), trans(B)): size2(A) != size1(C)"));
650  assert(viennacl::traits::size1(A.lhs()) == viennacl::traits::size2(B.lhs()) && bool("Size check failed at C = prod(trans(A), trans(B)): size1(A) != size2(B)"));
651  assert(viennacl::traits::size1(B.lhs()) == viennacl::traits::size2(C) && bool("Size check failed at C = prod(trans(A), trans(B)): size1(B) != size2(C)"));
652 
653  switch (viennacl::traits::handle(A.lhs()).get_active_handle_id())
654  {
656  viennacl::linalg::host_based::prod_impl(A.lhs(), true, B.lhs(), true, C, alpha, beta);
657  break;
658 #ifdef VIENNACL_WITH_OPENCL
660  viennacl::linalg::opencl::prod_impl(A.lhs(), true, B.lhs(), true, C, alpha, beta);
661  break;
662 #endif
663 #ifdef VIENNACL_WITH_CUDA
665  viennacl::linalg::cuda::prod_impl(A.lhs(), true, B.lhs(), true, C, alpha, beta);
666  break;
667 #endif
669  throw memory_exception("not initialised!");
670  default:
671  throw memory_exception("not implemented");
672  }
673  }
674 
675 
677 
678  template<typename NumericT>
680  {
682  viennacl::linalg::prod_impl(A, all_ones, result);
683  }
684 
685  template<typename NumericT>
687  {
690  }
691 
693 
694 
695 
701  template<typename T, typename OP>
703  matrix_expression<const matrix_base<T>, const matrix_base<T>, OP> const & proxy)
704  {
705  assert( (viennacl::traits::size1(A) == viennacl::traits::size1(proxy)) && bool("Size check failed at A = element_op(B): size1(A) != size1(B)"));
706  assert( (viennacl::traits::size2(A) == viennacl::traits::size2(proxy)) && bool("Size check failed at A = element_op(B): size2(A) != size2(B)"));
707 
708  switch (viennacl::traits::handle(A).get_active_handle_id())
709  {
712  break;
713 #ifdef VIENNACL_WITH_OPENCL
716  break;
717 #endif
718 #ifdef VIENNACL_WITH_CUDA
721  break;
722 #endif
724  throw memory_exception("not initialised!");
725  default:
726  throw memory_exception("not implemented");
727  }
728  }
729 
730 
731 #define VIENNACL_MAKE_BINARY_OP(OPNAME)\
732  template<typename T>\
733  viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_binary<op_##OPNAME> >\
734  element_##OPNAME(matrix_base<T> const & A, matrix_base<T> const & B)\
735  {\
736  return viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_binary<op_##OPNAME> >(A, B);\
737  }\
738 \
739  template<typename M1, typename M2, typename OP, typename T>\
740  viennacl::matrix_expression<const matrix_expression<const M1, const M2, OP>,\
741  const matrix_base<T>,\
742  op_element_binary<op_##OPNAME> >\
743  element_##OPNAME(matrix_expression<const M1, const M2, OP> const & proxy, matrix_base<T> const & B)\
744  {\
745  return viennacl::matrix_expression<const matrix_expression<const M1, const M2, OP>,\
746  const matrix_base<T>,\
747  op_element_binary<op_##OPNAME> >(proxy, B);\
748  }\
749 \
750  template<typename T, typename M2, typename M3, typename OP>\
751  viennacl::matrix_expression<const matrix_base<T>,\
752  const matrix_expression<const M2, const M3, OP>,\
753  op_element_binary<op_##OPNAME> >\
754  element_##OPNAME(matrix_base<T> const & A, matrix_expression<const M2, const M3, OP> const & proxy)\
755  {\
756  return viennacl::matrix_expression<const matrix_base<T>,\
757  const matrix_expression<const M2, const M3, OP>,\
758  op_element_binary<op_##OPNAME> >(A, proxy);\
759  }\
760 \
761  template<typename M1, typename M2, typename OP1,\
762  typename M3, typename M4, typename OP2>\
763  viennacl::matrix_expression<const matrix_expression<const M1, const M2, OP1>,\
764  const matrix_expression<const M3, const M4, OP2>,\
765  op_element_binary<op_##OPNAME> >\
766  element_##OPNAME(matrix_expression<const M1, const M2, OP1> const & proxy1,\
767  matrix_expression<const M3, const M4, OP2> const & proxy2)\
768  {\
769  return viennacl::matrix_expression<const matrix_expression<const M1, const M2, OP1>,\
770  const matrix_expression<const M3, const M4, OP2>,\
771  op_element_binary<op_##OPNAME> >(proxy1, proxy2);\
772  }
773 
777 
780  VIENNACL_MAKE_BINARY_OP(greater)
784 
785 #undef VIENNACL_GENERATE_BINARY_OP_OVERLOADS
786 
787 
788 
789 #define VIENNACL_MAKE_UNARY_ELEMENT_OP(funcname) \
790  template<typename T> \
791  viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_unary<op_##funcname> > \
792  element_##funcname(matrix_base<T> const & A) \
793  { \
794  return viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_unary<op_##funcname> >(A, A); \
795  } \
796  template<typename LHS, typename RHS, typename OP> \
797  viennacl::matrix_expression<const matrix_expression<const LHS, const RHS, OP>, \
798  const matrix_expression<const LHS, const RHS, OP>, \
799  op_element_unary<op_##funcname> > \
800  element_##funcname(matrix_expression<const LHS, const RHS, OP> const & proxy) \
801  { \
802  return viennacl::matrix_expression<const matrix_expression<const LHS, const RHS, OP>, \
803  const matrix_expression<const LHS, const RHS, OP>, \
804  op_element_unary<op_##funcname> >(proxy, proxy); \
805  } \
806 
824 
825 #undef VIENNACL_MAKE_UNARY_ELEMENT_OP
826 
827 
828  //
830  //
831 
832 
838  template<typename NumericT>
839  viennacl::matrix_expression<const vector_base<NumericT>, const vector_base<NumericT>, op_prod>
840  outer_prod(const vector_base<NumericT> & vec1, const vector_base<NumericT> & vec2)
841  {
842  return viennacl::matrix_expression< const vector_base<NumericT>, const vector_base<NumericT>, op_prod>(vec1, vec2);
843  }
844 
845 
858  template<typename NumericT, typename S1>
860  S1 const & alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha,
861  const vector_base<NumericT> & vec1,
862  const vector_base<NumericT> & vec2)
863  {
864  switch (viennacl::traits::handle(mat1).get_active_handle_id())
865  {
868  alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
869  vec1, vec2);
870  break;
871 #ifdef VIENNACL_WITH_OPENCL
874  alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
875  vec1, vec2);
876  break;
877 #endif
878 #ifdef VIENNACL_WITH_CUDA
881  alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
882  vec1, vec2);
883  break;
884 #endif
886  throw memory_exception("not initialised!");
887  default:
888  throw memory_exception("not implemented");
889  }
890  }
891 
900  template <typename NumericT, typename VectorType>
902  VectorType & dh,
903  VectorType & sh
904  )
905  {
906  switch (viennacl::traits::handle(A).get_active_handle_id())
907  {
910  break;
911 #ifdef VIENNACL_WITH_OPENCL
914  break;
915 #endif
916 
917 #ifdef VIENNACL_WITH_CUDA
920  break;
921 #endif
922 
924  throw memory_exception("not initialised!");
925  default:
926  throw memory_exception("not implemented");
927  }
928 
929 
930  }
941  template <typename SCALARTYPE>
944  vcl_size_t row_start,
945  vcl_size_t col_start,
946  bool copy_col
947  )
948  {
949  switch (viennacl::traits::handle(A).get_active_handle_id())
950  {
952  viennacl::linalg::host_based::copy_vec(A, V, row_start, col_start, copy_col);
953  break;
954 #ifdef VIENNACL_WITH_OPENCL
956  viennacl::linalg::opencl::copy_vec(A, V, row_start, col_start, copy_col);
957  break;
958 #endif
959 
960 #ifdef VIENNACL_WITH_CUDA
962  viennacl::linalg::cuda::copy_vec(A, V, row_start, col_start, copy_col);
963  break;
964 #endif
965 
967  throw memory_exception("not initialised!");
968  default:
969  throw memory_exception("not implemented");
970  }
971 
972  }
973 
980  template <typename NumericT>
982  vector_base<NumericT> & D,
984  {
985  switch (viennacl::traits::handle(A).get_active_handle_id())
986  {
989  break;
990 #ifdef VIENNACL_WITH_OPENCL
993  break;
994 #endif
995 
996 #ifdef VIENNACL_WITH_CUDA
999  break;
1000 #endif
1001 
1003  throw memory_exception("not initialised!");
1004  default:
1005  throw memory_exception("not implemented");
1006  }
1007  }
1008 
1009 
1017  template <typename NumericT>
1019  vector_base<NumericT> & D)
1020  {
1021  switch (viennacl::traits::handle(A).get_active_handle_id())
1022  {
1023  case viennacl::MAIN_MEMORY:
1025  break;
1026 #ifdef VIENNACL_WITH_OPENCL
1029  break;
1030 #endif
1031 
1032 #ifdef VIENNACL_WITH_CUDA
1033  case viennacl::CUDA_MEMORY:
1035  break;
1036 #endif
1037 
1039  throw memory_exception("not initialised!");
1040  default:
1041  throw memory_exception("not implemented");
1042  }
1043  }
1044 
1052  template <typename NumericT>
1054  vector_base<NumericT> & D,
1055  vcl_size_t A_size1)
1056  {
1057  switch (viennacl::traits::handle(Q).get_active_handle_id())
1058  {
1059  case viennacl::MAIN_MEMORY:
1061  break;
1062 #ifdef VIENNACL_WITH_OPENCL
1065  break;
1066 #endif
1067 
1068 #ifdef VIENNACL_WITH_CUDA
1069  case viennacl::CUDA_MEMORY:
1071  break;
1072 #endif
1073 
1075  throw memory_exception("not initialised!");
1076  default:
1077  throw memory_exception("not implemented");
1078  }
1079  }
1080 
1081 
1091  template<typename NumericT>
1093  vector_base<NumericT> & tmp1,
1094  vector_base<NumericT> & tmp2,
1095  int l,
1096  int m
1097  )
1098  {
1099  switch (viennacl::traits::handle(Q).get_active_handle_id())
1100  {
1101  case viennacl::MAIN_MEMORY:
1102  viennacl::linalg::host_based::givens_next(Q, tmp1, tmp2, l, m);
1103  break;
1104 #ifdef VIENNACL_WITH_OPENCL
1106  viennacl::linalg::opencl::givens_next(Q, tmp1, tmp2, l, m);
1107  break;
1108 #endif
1109 
1110 #ifdef VIENNACL_WITH_CUDA
1111  case viennacl::CUDA_MEMORY:
1112  viennacl::linalg::cuda::givens_next(Q, tmp1, tmp2, l, m);
1113  break;
1114 #endif
1115 
1117  throw memory_exception("not initialised!");
1118  default:
1119  throw memory_exception("not implemented");
1120  }
1121  }
1122 
1123  } //namespace linalg
1124 
1125 
1126 
1127 
1128  //
1130  //
1131 
1132 
1133  //v += A * x
1139  template<typename NumericT>
1140  vector<NumericT>
1143  {
1144  assert(viennacl::traits::size1(proxy.lhs()) == v1.size() && bool("Size check failed for v1 += A * v2: size1(A) != size(v1)"));
1145 
1146  vector<NumericT> result(viennacl::traits::size1(proxy.lhs()));
1147  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1148  v1 += result;
1149  return v1;
1150  }
1151 
1157  template<typename NumericT>
1158  vector<NumericT>
1161  {
1162  assert(viennacl::traits::size1(proxy.lhs()) == v1.size() && bool("Size check failed for v1 -= A * v2: size1(A) != size(v1)"));
1163 
1164  vector<NumericT> result(viennacl::traits::size1(proxy.lhs()));
1165  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1166  v1 -= result;
1167  return v1;
1168  }
1169 
1170 
1171 
1172 
1173 
1174  //free functions:
1180  template<typename NumericT>
1184  {
1185  assert(viennacl::traits::size1(proxy.lhs()) == viennacl::traits::size(v1) && bool("Size check failed for v1 + A * v2: size1(A) != size(v1)"));
1186 
1188  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1189  result += v1;
1190  return result;
1191  }
1192 
1198  template<typename NumericT>
1202  {
1203  assert(viennacl::traits::size1(proxy.lhs()) == viennacl::traits::size(v1) && bool("Size check failed for v1 - A * v2: size1(A) != size(v1)"));
1204 
1206  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1207  result = v1 - result;
1208  return result;
1209  }
1210 
1211 
1213 
1214 
1215  //v += A^T * x
1221  template<typename NumericT>
1222  vector<NumericT>
1225  const vector_base<NumericT>,
1226  op_prod> & proxy)
1227  {
1228  assert(viennacl::traits::size2(proxy.lhs()) == v1.size() && bool("Size check failed in v1 += trans(A) * v2: size2(A) != size(v1)"));
1229 
1230  vector<NumericT> result(viennacl::traits::size2(proxy.lhs()));
1231  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1232  v1 += result;
1233  return v1;
1234  }
1235 
1236  //v -= A^T * x
1242  template<typename NumericT>
1243  vector<NumericT>
1246  const vector_base<NumericT>,
1247  op_prod> & proxy)
1248  {
1249  assert(viennacl::traits::size2(proxy.lhs()) == v1.size() && bool("Size check failed in v1 += trans(A) * v2: size2(A) != size(v1)"));
1250 
1251  vector<NumericT> result(viennacl::traits::size2(proxy.lhs()));
1252  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1253  v1 -= result;
1254  return v1;
1255  }
1256 
1257 
1258  //free functions:
1264  template<typename NumericT>
1265  vector<NumericT>
1268  const vector_base<NumericT>,
1269  op_prod> & proxy)
1270  {
1271  assert(viennacl::traits::size2(proxy.lhs()) == viennacl::traits::size(v1) && bool("Size check failed in v1 + trans(A) * v2: size2(A) != size(v1)"));
1272 
1274  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1275  result += v1;
1276  return result;
1277  }
1278 
1284  template<typename NumericT>
1285  vector<NumericT>
1288  const vector_base<NumericT>,
1289  op_prod> & proxy)
1290  {
1291  assert(viennacl::traits::size2(proxy.lhs()) == viennacl::traits::size(v1) && bool("Size check failed in v1 - trans(A) * v2: size2(A) != size(v1)"));
1292 
1294  viennacl::linalg::prod_impl(proxy.lhs(), proxy.rhs(), result);
1295  result = v1 - result;
1296  return result;
1297  }
1298 
1299 
1300 } //namespace viennacl
1301 
1302 
1303 #endif
void norm_frobenius_cpu(matrix_base< T > const &vec, T &result)
Computes the Frobenius norm of a vector with final reduction on the CPU.
void row_sum_impl(const matrix_base< NumericT > &A, vector_base< NumericT > &result)
Implementations of dense matrix related operations, including matrix-vector products, using a plain single-threaded or OpenMP-enabled execution on CPU.
void matrix_diag_from_vector(const vector_base< NumericT > &vec, int k, matrix_base< NumericT > &A)
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
void matrix_diag_to_vector(const matrix_base< NumericT > &A, int k, vector_base< NumericT > &v)
Dispatcher interface for v = diag(A, k)
void house_update_QL(matrix_base< NumericT > &Q, vector_base< NumericT > &D, vcl_size_t A_size1)
This function updates the matrix Q, which is needed for the computation of the eigenvectors.
void am(matrix_base< NumericT > &A, matrix_base< NumericT > const &B, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha)
void matrix_column(const matrix_base< NumericT > &A, unsigned int j, vector_base< NumericT > &vec)
This class represents a single scalar value on the GPU and behaves mostly like a built-in scalar type...
Definition: forwards.h:227
Implementations of dense matrix related operations, including matrix-vector products, using OpenCL.
void convert(matrix_base< DestNumericT > &mat1, matrix_base< SrcNumericT > const &mat2)
void house_update_A_right(matrix_base< NumericT > &A, vector_base< NumericT > &D)
This function applies a householder transformation to a matrix: A <- A * P with a householder reflect...
vector< NumericT > operator-=(vector_base< NumericT > &v1, const viennacl::vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, viennacl::op_prod > &proxy)
Implementation of the operation v1 -= A * v2, where A is a matrix.
void matrix_diag_to_vector(const matrix_base< NumericT > &mat, int k, vector_base< NumericT > &vec)
size_type internal_size() const
Returns the total amount of allocated memory in multiples of sizeof(NumericT)
Definition: matrix_def.hpp:242
Exception class in case of memory errors.
Definition: forwards.h:572
void copy_vec(matrix_base< SCALARTYPE > &A, vector_base< SCALARTYPE > &V, vcl_size_t row_start, vcl_size_t col_start, bool copy_col)
This function copies a row or a column from a matrix to a vector.
Generic size and resize functionality for different vector and matrix types.
void trans(matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > const &proxy, matrix_base< NumericT > &temp_trans)
Extracts the underlying OpenCL start index handle from a vector, a matrix, an expression etc...
void matrix_assign(matrix_base< NumericT > &mat, NumericT s, bool clear=false)
Various little tools used here and there in ViennaCL.
void house_update_A_right(matrix_base< NumericT > &A, vector_base< NumericT > &D)
This function applies a householder transformation to a matrix: A <- A * P with a householder reflect...
vector< NumericT > operator+=(vector_base< NumericT > &v1, const viennacl::vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, viennacl::op_prod > &proxy)
Implementation of the operation v1 += A * v2, where A is a matrix.
void matrix_assign(matrix_base< NumericT > &mat, NumericT s, bool clear=false)
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
void norm_2_cpu(vector_base< T > const &vec, T &result)
Computes the l^2-norm of a vector with final reduction on the CPU - dispatcher interface.
void house_update_QL(matrix_base< NumericT > &Q, vector_base< NumericT > &D, vcl_size_t A_size1)
This function updates the matrix Q, which is needed for the computation of the eigenvectors.
void am(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha)
Expression template class for representing a tree of expressions which ultimately result in a matrix...
Definition: forwards.h:341
void ambm(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarType2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
size_type stride2() const
Returns the number of columns.
Definition: matrix_def.hpp:234
void clear(VectorType &vec)
Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...
Definition: clear.hpp:43
This file provides the forward declarations for the main types used within ViennaCL.
A dense matrix class.
Definition: forwards.h:375
void ambm(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarT2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
Determines row and column increments for matrices and matrix proxies.
void bidiag_pack(matrix_base< NumericT > &A, viennacl::vector< NumericT > &dh, viennacl::vector< NumericT > &sh)
void prod_impl(const matrix_base< NumericT > &mat, bool mat_transpose, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
An expression template class that represents a binary operation that yields a vector.
Definition: forwards.h:239
void element_op(matrix_base< T > &A, matrix_expression< const matrix_base< T >, const matrix_base< T >, OP > const &proxy)
Implementation of the element-wise operation A = B .* C and A = B ./ C for matrices (using MATLAB syn...
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
void scaled_rank_1_update(matrix_base< NumericT > &A, ScalarT1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
The implementation of the operation mat += alpha * vec1 * vec2^T, i.e. a scaled rank 1 update...
viennacl::vector< NumericT > operator-(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)
Implementation of the operation 'result = v1 - A * v2', where A is a matrix.
void house_update_A_right(matrix_base< NumericT > &A, vector_base< NumericT > &D)
This function applies a householder transformation to a matrix: A <- A * P with a householder reflect...
float NumericT
Definition: bisect.cpp:40
void copy_vec(matrix_base< NumericT > &A, vector_base< NumericT > &V, vcl_size_t row_start, vcl_size_t col_start, bool copy_col)
void house_update_A_right(matrix_base< NumericT > &A, vector_base< NumericT > &D)
void element_op(matrix_base< NumericT > &A, matrix_expression< const matrix_base< NumericT >, const matrix_base< NumericT >, op_element_binary< OpT > > const &proxy)
Implementation of binary element-wise operations A = OP(B,C)
void norm_2_impl(vector_base< T > const &vec, scalar< T > &result)
Computes the l^2-norm of a vector - dispatcher interface.
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:34
void ambm(matrix_base< NumericT > &A, matrix_base< NumericT > const &B, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &C, ScalarT2 const &beta, vcl_size_t, bool reciprocal_beta, bool flip_sign_beta)
void norm_frobenius_impl(matrix_base< T > const &vec, scalar< T > &result)
Computes the Frobenius norm of a matrix - dispatcher interface.
void house_update_QL(matrix_base< NumericT > &Q, vector_base< NumericT > &D, vcl_size_t A_size1)
void convert(matrix_base< DestNumericT > &mat1, matrix_base< SrcNumericT > const &mat2)
viennacl::vector< float > v1
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:235
void convert(matrix_base< DestNumericT > &dest, matrix_base< SrcNumericT > const &src)
void matrix_row(const matrix_base< NumericT > &A, unsigned int i, vector_base< NumericT > &v)
void copy_vec(matrix_base< NumericT > &A, vector_base< S1 > &V, vcl_size_t row_start, vcl_size_t col_start, bool copy_col)
This function copies a row or a column from a matrix to a vector.
void prod_impl(const matrix_base< NumericT > &mat, bool trans, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
void ambm_m(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarT2 const &beta, vcl_size_t, bool reciprocal_beta, bool flip_sign_beta)
void am(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha)
void matrix_diag_to_vector(matrix_base< NumericT > const &mat, int k, vector_base< NumericT > &vec)
#define VIENNACL_MAKE_UNARY_ELEMENT_OP(funcname)
void matrix_column(const matrix_base< NumericT > &mat, unsigned int j, vector_base< NumericT > &vec)
void bidiag_pack(matrix_base< NumericT > &A, VectorType &dh, VectorType &sh)
This function stores the diagonal and the superdiagonal of a matrix in two vectors.
void house_update_A_left(matrix_base< NumericT > &A, vector_base< NumericT > &D, vcl_size_t start)
This function applies a householder transformation to a matrix. A <- P * A with a householder reflect...
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
void element_op(matrix_base< NumericT > &A, matrix_expression< const matrix_base< NumericT >, const matrix_base< NumericT >, op_element_binary< OpT > > const &proxy)
Implementation of the element-wise operations A = B .* C and A = B ./ C (using MATLAB syntax) ...
result_of::size_type< T >::type start(T const &obj)
Definition: start.hpp:44
void house_update_QL(matrix_base< NumericT > &Q, vector_base< NumericT > &D, vcl_size_t A_size1)
This function updates the matrix Q, which is needed for the computation of the eigenvectors.
void scaled_rank_1_update(matrix_base< NumericT > &mat1, ScalarT const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
The implementation of the operation mat += alpha * vec1 * vec2^T, i.e. a scaled rank 1 update...
void ambm_m(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarT2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
size_type stride1() const
Returns the number of rows.
Definition: matrix_def.hpp:232
void house_update_A_left(matrix_base< NumericT > &A, vector_base< NumericT > &D, vcl_size_t start)
This function applies a householder transformation to a matrix. A <- P * A with a householder reflect...
#define VIENNACL_MAKE_BINARY_OP(OPNAME)
void scaled_rank_1_update(matrix_base< NumericT > &mat1, S1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
The implementation of the operation mat += alpha * vec1 * vec2^T, i.e. a scaled rank 1 update...
void matrix_diag_from_vector(const vector_base< NumericT > &vec, int k, matrix_base< NumericT > &mat)
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
void ambm(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarT2 const &beta, vcl_size_t, bool reciprocal_beta, bool flip_sign_beta)
std::size_t vcl_size_t
Definition: forwards.h:75
size_type size2() const
Returns the number of columns.
Definition: matrix_def.hpp:226
handle_type & handle()
Returns the OpenCL handle, non-const-version.
Definition: matrix_def.hpp:244
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
void house_update_A_left(matrix_base< NumericT > &A, vector_base< NumericT > &D, vcl_size_t start)
void matrix_column(const matrix_base< NumericT > &A, unsigned int j, vector_base< NumericT > &v)
size_type size1() const
Returns the number of rows.
Definition: matrix_def.hpp:224
Proxy classes for vectors.
void scaled_rank_1_update(matrix_base< NumericT > &mat1, ScalarT const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha, const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
The implementation of the operation mat += alpha * vec1 * vec2^T, i.e. a scaled rank 1 update...
All the predicates used within ViennaCL. Checks for expressions to be vectors, etc.
void matrix_row(const matrix_base< NumericT > &mat, unsigned int i, vector_base< NumericT > &vec)
void convert(matrix_base< DestNumericT > &dest, matrix_base< SrcNumericT > const &src)
viennacl::vector< NumericT > operator+(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)
Implementation of the operation 'result = v1 + A * v2', where A is a matrix.
void am(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha)
void matrix_column(const matrix_base< NumericT > &mat, unsigned int j, vector_base< NumericT > &vec)
viennacl::matrix_expression< const vector_base< NumericT >, const vector_base< NumericT >, op_prod > outer_prod(const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
Returns a proxy class for the operation mat += vec1 * vec2^T, i.e. a rank 1 update.
A tag class representing matrix-vector products and element-wise multiplications. ...
Definition: forwards.h:94
void matrix_diag_to_vector(const matrix_base< NumericT > &A, int k, vector_base< NumericT > &vec)
Implementations of dense matrix related operations, including matrix-vector products, using CUDA.
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
void element_op(matrix_base< NumericT, SizeT > &A, matrix_expression< const matrix_base< NumericT, SizeT >, const matrix_base< NumericT, SizeT >, op_element_binary< OpT > > const &proxy)
void matrix_row(matrix_base< NumericT > const &mat, unsigned int i, vector_base< NumericT > &vec)
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
bool row_major() const
Definition: matrix_def.hpp:248
Represents a vector consisting of scalars 's' only, i.e. v[i] = s for all i. To be used as an initial...
Definition: vector_def.hpp:87
void matrix_assign(matrix_base< NumericT > &A, NumericT s, bool up_to_internal_size=false)
void givens_next(matrix_base< NumericT > &Q, vector_base< NumericT > &tmp1, vector_base< NumericT > &tmp2, int l, int m)
This function updates the matrix Q. It is part of the tql2 algorithm.
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
size_type size() const
Returns the length of the vector (cf. std::vector)
Definition: vector_def.hpp:118
float ScalarType
Definition: fft_1d.cpp:42
void bidiag_pack(matrix_base< NumericT > &A, VectorType &dh, VectorType &sh)
This function stores the diagonal and the superdiagonal of a matrix in two vectors.
void bidiag_pack(matrix_base< NumericT > &A, VectorType &dh, VectorType &sh)
Main abstraction class for multiple memory domains. Represents a buffer in either main RAM...
Definition: mem_handle.hpp:89
void givens_next(matrix_base< NumericT > &matrix, vector_base< NumericT > &tmp1, vector_base< NumericT > &tmp2, int l, int m)
A tag class representing transposed matrices.
Definition: forwards.h:220
size_type start2() const
Returns the number of columns.
Definition: matrix_def.hpp:230
void house_update_A_left(matrix_base< NumericT > &A, vector_base< NumericT > &D, vcl_size_t start)
This function applies a householder transformation to a matrix. A <- P * A with a householder reflect...
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
void matrix_row(const matrix_base< NumericT > &A, unsigned int i, vector_base< NumericT > &vec)
void givens_next(matrix_base< NumericT > &Q, vector_base< NumericT > &tmp1, vector_base< NumericT > &tmp2, int l, int m)
This function updates the matrix Q. It is part of the tql2 algorithm.
Extracts the underlying OpenCL handle from a vector, a matrix, an expression etc. ...
void prod_impl(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
Definition: handle.hpp:41
void matrix_diag_from_vector(const vector_base< NumericT > &v, int k, matrix_base< NumericT > &A)
Dispatcher interface for A = diag(v, k)
void matrix_diag_from_vector(const vector_base< NumericT > &vec, int k, matrix_base< NumericT > &mat)
void ambm_m(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarType2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
Implementation of the ViennaCL scalar class.
A collection of compile time type deductions.
void matrix_assign(matrix_base< NumericT > &mat, NumericT s, bool clear=false)
void prod_impl(const matrix_base< NumericT > &A, bool trans_A, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
void ambm_m(matrix_base< NumericT > &A, matrix_base< NumericT > const &B, ScalarT1 const &alpha, vcl_size_t, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &C, ScalarT2 const &beta, vcl_size_t, bool reciprocal_beta, bool flip_sign_beta)
void givens_next(matrix_base< NumericT > &Q, vector_base< NumericT > &tmp1, vector_base< NumericT > &tmp2, int l, int m)
This function updates the matrix Q. It is part of the tql2 algorithm.
void matrix_diagonal_assign(matrix_base< NumericT > &A, NumericT s)
void copy_vec(matrix_base< NumericT > &A, vector_base< NumericT > &V, vcl_size_t row_start, vcl_size_t col_start, bool copy_col)
This function copies a row or a column from a matrix to a vector.
Simple enable-if variant that uses the SFINAE pattern.
size_type start1() const
Returns the number of rows.
Definition: matrix_def.hpp:228
void column_sum_impl(const matrix_base< NumericT > &A, vector_base< NumericT > &result)
memory_types get_active_handle_id() const
Returns an ID for the currently active memory buffer. Other memory buffers might contain old or no da...
Definition: mem_handle.hpp:118