1 #ifndef VIENNACL_LINALG_MATRIX_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_MATRIX_OPERATIONS_HPP_
40 #ifdef VIENNACL_WITH_OPENCL
44 #ifdef VIENNACL_WITH_CUDA
53 template<
typename DestNumericT,
typename SrcNumericT>
64 #ifdef VIENNACL_WITH_OPENCL
69 #ifdef VIENNACL_WITH_CUDA
82 typename SizeT,
typename DistanceT>
91 #ifdef VIENNACL_WITH_OPENCL
96 #ifdef VIENNACL_WITH_CUDA
110 typename ScalarType1>
119 #ifdef VIENNACL_WITH_OPENCL
124 #ifdef VIENNACL_WITH_CUDA
138 typename ScalarType1,
typename ScalarType2>
147 mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
148 mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
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);
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);
173 typename ScalarType1,
typename ScalarType2>
182 mat2, alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
183 mat3, beta, len_beta, reciprocal_beta, flip_sign_beta);
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);
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);
207 template<
typename NumericT>
215 #ifdef VIENNACL_WITH_OPENCL
220 #ifdef VIENNACL_WITH_CUDA
233 template<
typename NumericT>
241 #ifdef VIENNACL_WITH_OPENCL
246 #ifdef VIENNACL_WITH_CUDA
260 template<
typename NumericT>
268 #ifdef VIENNACL_WITH_OPENCL
273 #ifdef VIENNACL_WITH_CUDA
286 template<
typename NumericT>
294 #ifdef VIENNACL_WITH_OPENCL
299 #ifdef VIENNACL_WITH_CUDA
311 template<
typename NumericT>
319 #ifdef VIENNACL_WITH_OPENCL
324 #ifdef VIENNACL_WITH_CUDA
336 template<
typename NumericT>
344 #ifdef VIENNACL_WITH_OPENCL
349 #ifdef VIENNACL_WITH_CUDA
437 template<
typename NumericT>
450 #ifdef VIENNACL_WITH_OPENCL
455 #ifdef VIENNACL_WITH_CUDA
478 template<
typename NumericT>
491 #ifdef VIENNACL_WITH_OPENCL
496 #ifdef VIENNACL_WITH_CUDA
518 template<
typename NumericT,
typename ScalarType >
535 #ifdef VIENNACL_WITH_OPENCL
540 #ifdef VIENNACL_WITH_CUDA
559 template<
typename NumericT,
typename ScalarType >
577 #ifdef VIENNACL_WITH_OPENCL
582 #ifdef VIENNACL_WITH_CUDA
602 template<
typename NumericT,
typename ScalarType >
618 #ifdef VIENNACL_WITH_OPENCL
623 #ifdef VIENNACL_WITH_CUDA
642 template<
typename NumericT,
typename ScalarType >
658 #ifdef VIENNACL_WITH_OPENCL
663 #ifdef VIENNACL_WITH_CUDA
678 template<
typename NumericT>
685 template<
typename NumericT>
701 template<
typename T,
typename OP>
713 #ifdef VIENNACL_WITH_OPENCL
718 #ifdef VIENNACL_WITH_CUDA
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)\
736 return viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_binary<op_##OPNAME> >(A, B);\
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)\
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);\
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)\
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);\
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)\
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);\
785 #undef VIENNACL_GENERATE_BINARY_OP_OVERLOADS
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) \
794 return viennacl::matrix_expression<const matrix_base<T>, const matrix_base<T>, op_element_unary<op_##funcname> >(A, A); \
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) \
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); \
825 #undef VIENNACL_MAKE_UNARY_ELEMENT_OP
838 template<
typename NumericT>
858 template<
typename NumericT,
typename S1>
860 S1
const & alpha,
vcl_size_t len_alpha,
bool reciprocal_alpha,
bool flip_sign_alpha,
868 alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
871 #ifdef VIENNACL_WITH_OPENCL
874 alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
878 #ifdef VIENNACL_WITH_CUDA
881 alpha, len_alpha, reciprocal_alpha, flip_sign_alpha,
900 template <
typename NumericT,
typename VectorType>
911 #ifdef VIENNACL_WITH_OPENCL
917 #ifdef VIENNACL_WITH_CUDA
941 template <
typename SCALARTYPE>
954 #ifdef VIENNACL_WITH_OPENCL
960 #ifdef VIENNACL_WITH_CUDA
980 template <
typename NumericT>
990 #ifdef VIENNACL_WITH_OPENCL
996 #ifdef VIENNACL_WITH_CUDA
1017 template <
typename NumericT>
1026 #ifdef VIENNACL_WITH_OPENCL
1032 #ifdef VIENNACL_WITH_CUDA
1052 template <
typename NumericT>
1062 #ifdef VIENNACL_WITH_OPENCL
1068 #ifdef VIENNACL_WITH_CUDA
1091 template<
typename NumericT>
1104 #ifdef VIENNACL_WITH_OPENCL
1110 #ifdef VIENNACL_WITH_CUDA
1139 template<
typename NumericT>
1157 template<
typename NumericT>
1180 template<
typename NumericT>
1198 template<
typename NumericT>
1207 result = v1 - result;
1221 template<
typename NumericT>
1228 assert(
viennacl::traits::size2(proxy.lhs()) == v1.
size() && bool(
"Size check failed in v1 += trans(A) * v2: size2(A) != size(v1)"));
1242 template<
typename NumericT>
1249 assert(
viennacl::traits::size2(proxy.lhs()) == v1.
size() && bool(
"Size check failed in v1 += trans(A) * v2: size2(A) != size(v1)"));
1264 template<
typename NumericT>
1284 template<
typename NumericT>
1295 result = v1 - result;
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 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_assign(matrix_base< NumericT > &mat, NumericT s, bool clear=false)
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 matrix_diag_from_vector(const vector_base< NumericT > &vec, int k, matrix_base< NumericT > &mat)
This class represents a single scalar value on the GPU and behaves mostly like a built-in scalar type...
Implementations of dense matrix related operations, including matrix-vector products, using OpenCL.
void matrix_diag_to_vector(const matrix_base< NumericT > &mat, int k, vector_base< NumericT > &vec)
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 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)
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)
void prod_impl(const matrix_base< NumericT > &mat, bool trans_A, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
void matrix_column(const matrix_base< NumericT > &mat, unsigned int j, vector_base< NumericT > &vec)
Exception class in case of memory errors.
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)
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.)
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...
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.
void clear(VectorType &vec)
Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...
This file provides the forward declarations for the main types used within ViennaCL.
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.
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...
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...
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 norm_2_impl(vector_base< T > const &vec, scalar< T > &result)
Computes the l^2-norm of a vector - dispatcher interface.
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)
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
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)
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)
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
size_type stride1() const
Returns the number of rows.
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 am(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarT1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha)
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)
size_type size2() const
Returns the number of columns.
handle_type & handle()
Returns the OpenCL handle, non-const-version.
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)
void element_op(matrix_base< T > &A, matrix_expression< const matrix_base< T >, const matrix_base< T >, op_element_binary< OP > > const &proxy)
Implementation of binary element-wise operations A = OP(B,C)
size_type size1() const
Returns the number of rows.
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. ...
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.
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 ...
Represents a vector consisting of scalars 's' only, i.e. v[i] = s for all i. To be used as an initial...
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)
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...
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.
size_type start2() const
Returns the number of columns.
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 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)
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
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.
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 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 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.
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...
void matrix_row(matrix_base< NumericT > const &mat, unsigned int i, vector_base< NumericT > &vec)