1 #ifndef VIENNACL_LINALG_DIRECT_SOLVE_HPP_
2 #define VIENNACL_LINALG_DIRECT_SOLVE_HPP_
34 #ifdef VIENNACL_WITH_OPENCL
38 #ifdef VIENNACL_WITH_CUDA
42 #define VIENNACL_DIRECT_SOLVE_BLOCKSIZE 128
61 template<
typename NumericT,
typename SolverTagT>
71 #ifdef VIENNACL_WITH_OPENCL
76 #ifdef VIENNACL_WITH_CUDA
93 template<
typename NumericT,
typename SolverTagT>
98 assert( (mat.
size1() == vec.
size()) &&
bool(
"Size check failed in inplace_solve(): size1(A) != size(b)"));
99 assert( (mat.
size2() == vec.
size()) &&
bool(
"Size check failed in inplace_solve(): size2(A) != size(b)"));
106 #ifdef VIENNACL_WITH_OPENCL
111 #ifdef VIENNACL_WITH_CUDA
124 template<
typename MatrixT1,
typename MatrixT2,
typename SolverTagT>
130 if (A.size1() <= blockSize)
134 for (
vcl_size_t i = 0; i < A.size1(); i = i + blockSize)
137 vcl_size_t Apos2 = std::min<vcl_size_t>(A.size1(), i + blockSize);
142 if (Apos2 < A.size1())
148 NumericType(-1.0), NumericType(1.0));
154 template<
typename MatrixT1,
typename MatrixT2>
160 template<
typename MatrixT1,
typename MatrixT2>
166 template<
typename MatrixT1,
typename MatrixT2,
typename SolverTagT>
172 if (static_cast<int>(A.size1()) <= blockSize)
176 for (
int i = static_cast<int>(A.size1()); i > 0; i = i - blockSize)
191 NumericType(-1.0), NumericType(1.0));
197 template<
typename MatrixT1,
typename MatrixT2>
203 template<
typename MatrixT1,
typename MatrixT2>
216 template<
typename NumericT,
typename SolverTagT>
229 template<
typename NumericT,
typename SolverTagT>
237 proxy_B.lhs().size2(), proxy_B.lhs().start2(), proxy_B.lhs().stride2(), proxy_B.lhs().internal_size2(),
238 proxy_B.lhs().size1(), proxy_B.lhs().start1(), proxy_B.lhs().stride1(), proxy_B.lhs().internal_size1(),
239 !proxy_B.lhs().row_major());
250 template<
typename NumericT,
typename SolverTagT>
258 proxy_A.lhs().size2(), proxy_A.lhs().start2(), proxy_A.lhs().stride2(), proxy_A.lhs().internal_size2(),
259 proxy_A.lhs().size1(), proxy_A.lhs().start1(), proxy_A.lhs().stride1(), proxy_A.lhs().internal_size1(),
260 !proxy_A.lhs().row_major());
270 template<
typename NumericT,
typename SolverTagT>
278 proxy_A.lhs().size2(), proxy_A.lhs().start2(), proxy_A.lhs().stride2(), proxy_A.lhs().internal_size2(),
279 proxy_A.lhs().size1(), proxy_A.lhs().start1(), proxy_A.lhs().stride1(), proxy_A.lhs().internal_size1(),
280 !proxy_A.lhs().row_major());
283 proxy_B.lhs().size2(), proxy_B.lhs().start2(), proxy_B.lhs().stride2(), proxy_B.lhs().internal_size2(),
284 proxy_B.lhs().size1(), proxy_B.lhs().start1(), proxy_B.lhs().stride1(), proxy_B.lhs().internal_size1(),
285 !proxy_B.lhs().row_major());
300 template<
typename NumericT,
typename SolverTagT>
317 template<
typename NumericT,
typename SolverTagT>
334 template<
typename NumericT,
typename SolverTagT>
351 template<
typename NumericT,
typename SolverTagT>
368 template<
typename MatrixT1,
typename VectorT,
typename SolverTagT>
372 if (A.size1() <= blockSize)
377 for (
vcl_size_t i = 0; i < A.size1(); i = i + blockSize)
380 vcl_size_t Apos2 = std::min<vcl_size_t>(A.size1(), i + blockSize);
384 if (Apos2 < A.size1())
394 template<
typename MatrixT1,
typename VectorT>
400 template<
typename MatrixT1,
typename VectorT>
406 template<
typename MatrixT1,
typename VectorT,
typename SolverTagT>
410 if (static_cast<int>(A.size1()) <= blockSize)
415 for (
int i = static_cast<int>(A.size1()); i > 0; i = i - blockSize)
432 template<
typename MatrixT1,
typename VectorT>
438 template<
typename MatrixT1,
typename VectorT>
452 template<
typename NumericT,
typename SolverTagT>
455 SolverTagT
const & tag)
467 template<
typename NumericT,
typename SolverTagT>
470 SolverTagT
const & tag)
476 proxy.lhs().size2(), proxy.lhs().start2(), proxy.lhs().stride2(), proxy.lhs().internal_size2(),
477 proxy.lhs().size1(), proxy.lhs().start1(), proxy.lhs().stride1(), proxy.lhs().internal_size1(),
478 !proxy.lhs().row_major());
491 template<
typename NumericT>
510 template<
typename NumericT>
529 template<
typename NumericT>
548 template<
typename NumericT>
565 template<
typename NumericT,
typename SolverTagT>
568 SolverTagT
const & tag)
viennacl::tools::shared_ptr< char > handle_type
Implementations of dense direct triangular solvers are found here.
void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT tag)
Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...
void inplace_solve(const matrix_base< NumericT > &A, matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...
void inplace_solve_upper_impl(MatrixT1 const &A, MatrixT2 &B, SolverTagT)
Exception class in case of memory errors.
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Implementation of the dense matrix class.
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
A tag class representing a lower triangular matrix.
void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...
Expression template class for representing a tree of expressions which ultimately result in a matrix...
This file provides the forward declarations for the main types used within ViennaCL.
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
void inplace_solve_upper_vec_impl(MatrixT1 const &A, VectorT &b, SolverTagT)
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
void inplace_solve_vec_impl(MatrixT1 const &A, VectorT &B, viennacl::linalg::lower_tag)
A tag class representing an upper triangular matrix.
void inplace_solve_kernel(const matrix_base< NumericT > &A, const matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for dense triangular systems using a single kernel launch. Matlab notation: A \...
matrix_range< MatrixType > project(MatrixType const &A, viennacl::range const &r1, viennacl::range const &r2)
size_type size2() const
Returns the number of columns.
void inplace_solve_impl(MatrixT1 const &A, MatrixT2 &B, viennacl::linalg::lower_tag)
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
size_type size1() const
Returns the number of rows.
Proxy classes for vectors.
#define VIENNACL_DIRECT_SOLVE_BLOCKSIZE
void inplace_solve(matrix_base< NumericT > const &A, matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for dense triangular systems. Matlab notation: A \ B.
Proxy classes for matrices.
void inplace_solve_lower_impl(MatrixT1 const &A, MatrixT2 &B, SolverTagT)
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
void inplace_solve_lower_vec_impl(MatrixT1 const &A, VectorT &b, SolverTagT)
size_type size() const
Returns the length of the vector (cf. std::vector)
Implementations of dense direct solvers using CUDA are found here.
A range class that refers to an interval [start, stop), where 'start' is included, and 'stop' is excluded.
A tag class representing a lower triangular matrix with unit diagonal.
A tag class representing transposed matrices.
Class for representing non-strided submatrices of a bigger matrix A.
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 inplace_solve_vec_kernel(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, SolverTagT)
Implementations of dense direct solvers are found here.
Simple enable-if variant that uses the SFINAE pattern.
A tag class representing an upper triangular matrix with unit diagonal.
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...