ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
viennacl::linalg::detail Namespace Reference

Namespace holding implementation details for linear algebra routines. Usually not of interest for a library user. More...

Namespaces

 amg
 Implementation namespace for algebraic multigrid preconditioner.
 
 fft
 
 spai
 Implementation namespace for sparse approximate inverse preconditioner.
 

Classes

class  FastMatrix
 Internal helper class representing a row-major dense matrix used for the QR method for the purpose of computing eigenvalues. More...
 
class  ilu_vector_range
 Helper range class for representing a subvector of a larger buffer. More...
 
struct  ilut_sparse_vector
 Helper struct for holding a sparse vector in linear memory. For internal use only. More...
 
struct  InputData
 In this class the input matrix is stored. More...
 
struct  op_applier
 Worker class for decomposing expression templates. More...
 
struct  op_executor
 Worker class for decomposing expression templates. More...
 
struct  ResultDataLarge
 In this class the data of the result for large matrices is stored. More...
 
struct  ResultDataSmall
 In this class the data of the result for small matrices is stored. More...
 
class  z_handler
 handles the no_precond case at minimal overhead More...
 
class  z_handler< VectorT, viennacl::linalg::no_precond >
 

Enumerations

enum  row_info_types { SPARSE_ROW_NORM_INF = 0, SPARSE_ROW_NORM_1, SPARSE_ROW_NORM_2, SPARSE_ROW_DIAGONAL }
 

Functions

template<typename NumericT >
void amg_galerkin_prod (compressed_matrix< NumericT > &A_fine, compressed_matrix< NumericT > &P, compressed_matrix< NumericT > &R, compressed_matrix< NumericT > &A_coarse)
 Sparse Galerkin product: Calculates A_coarse = trans(P)*A_fine*P = R*A_fine*P. More...
 
template<typename NumericT , typename AMGContextListT >
vcl_size_t amg_setup (std::vector< compressed_matrix< NumericT > > &list_of_A, std::vector< compressed_matrix< NumericT > > &list_of_P, std::vector< compressed_matrix< NumericT > > &list_of_R, AMGContextListT &list_of_amg_level_context, amg_tag &tag)
 Setup AMG preconditioner. More...
 
template<typename MatrixT , typename InternalT1 , typename InternalT2 >
void amg_init (MatrixT const &mat, InternalT1 &list_of_A, InternalT1 &list_of_P, InternalT1 &list_of_R, InternalT2 &list_of_amg_level_context, amg_tag &tag)
 Initialize AMG preconditioner. More...
 
template<typename InternalVectorT , typename SparseMatrixT >
void amg_setup_apply (InternalVectorT &result, InternalVectorT &result_backup, InternalVectorT &rhs, InternalVectorT &residual, SparseMatrixT const &A, vcl_size_t coarse_levels, amg_tag const &tag)
 Setup data structures for precondition phase for later use on the GPU. More...
 
template<typename NumericT , typename SparseMatrixT >
void amg_lu (viennacl::matrix< NumericT > &op, SparseMatrixT const &A, amg_tag const &tag)
 Pre-compute LU factorization for direct solve (ublas library). More...
 
template<typename MatrixT , typename NumericT >
viennacl::vector< NumericTpipelined_solve (MatrixT const &A, viennacl::vector_base< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Implementation of a pipelined stabilized Bi-conjugate gradient solver. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::compressed_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::coordinate_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::sliced_ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::hyb_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename MatrixT , typename VectorT >
VectorT solve_impl (MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)=NULL, void *monitor_data=NULL)
 Implementation of the unpreconditioned stabilized Bi-conjugate gradient solver. More...
 
template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT solve_impl (MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond, bool(*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)=NULL, void *monitor_data=NULL)
 Implementation of the preconditioned stabilized Bi-conjugate gradient solver. More...
 
template<typename NumericT , typename OtherVectorT >
void copy_vec_to_vec (viennacl::vector< NumericT > const &src, OtherVectorT &dest)
 overloaded function for copying vectors More...
 
template<typename OtherVectorT , typename NumericT >
void copy_vec_to_vec (OtherVectorT const &src, viennacl::vector< NumericT > &dest)
 
template<typename VectorT1 , typename VectorT2 >
void copy_vec_to_vec (VectorT1 const &src, VectorT2 &dest)
 
template<typename MatrixT , typename NumericT >
viennacl::vector< NumericTpipelined_solve (MatrixT const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Implementation of a pipelined conjugate gradient algorithm (no preconditioner), specialized for ViennaCL types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::compressed_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::coordinate_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::sliced_ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::hyb_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, cg_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT solve_impl (MatrixT const &matrix, VectorT const &rhs, cg_tag const &tag, PreconditionerT const &precond, bool(*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)=NULL, void *monitor_data=NULL)
 
template<typename NumericT >
void bisectSmall (const InputData< NumericT > &input, ResultDataSmall< NumericT > &result, const unsigned int mat_size, const NumericT lg, const NumericT ug, const NumericT precision)
 
template<typename NumericT >
void bisectLarge (const InputData< NumericT > &input, ResultDataLarge< NumericT > &result, const unsigned int mat_size, const NumericT lg, const NumericT ug, const NumericT precision)
 
template<typename NumericT >
void bisectLarge_OneIntervals (const InputData< NumericT > &input, ResultDataLarge< NumericT > &result, const unsigned int mat_size, const NumericT precision)
 
template<typename NumericT >
void bisectLarge_MultIntervals (const InputData< NumericT > &input, ResultDataLarge< NumericT > &result, const unsigned int mat_size, const NumericT precision)
 
template<typename NumericT >
void computeEigenvaluesLargeMatrix (InputData< NumericT > &input, ResultDataLarge< NumericT > &result, const unsigned int mat_size, const NumericT lg, const NumericT ug, const NumericT precision)
 
template<typename NumericT >
bool processResultDataLargeMatrix (ResultDataLarge< NumericT > &result, const unsigned int mat_size)
 
template<typename NumericT >
void computeEigenvaluesSmallMatrix (const InputData< NumericT > &input, ResultDataSmall< NumericT > &result, const unsigned int mat_size, const NumericT lg, const NumericT ug, const NumericT precision)
 
template<typename NumericT >
void processResultSmallMatrix (ResultDataSmall< NumericT > &result, const unsigned int mat_size)
 
template<typename NumericT >
void computeGerschgorin (std::vector< NumericT > &d, std::vector< NumericT > &s, unsigned int n, NumericT &lg, NumericT &ug)
 
template<class T >
min (const T &lhs, const T &rhs)
 Minimum. More...
 
template<class T >
max (const T &lhs, const T &rhs)
 Maximum. More...
 
float sign_f (const float &val)
 Sign of number (float) More...
 
double sign_d (const double &val)
 Sign of number (double) More...
 
unsigned int getNumBlocksLinear (const unsigned int num_threads, const unsigned int num_threads_block)
 
template<typename NumericT >
void extract_block_matrix (viennacl::compressed_matrix< NumericT > const &A, viennacl::compressed_matrix< NumericT > &diagonal_block_A, vcl_size_t start_index, vcl_size_t stop_index)
 Extracts a diagonal block from a larger system matrix. More...
 
template<typename NumericT >
void precondition (viennacl::compressed_matrix< NumericT > const &A, viennacl::compressed_matrix< NumericT > &L, viennacl::vector< NumericT > &diag_L, viennacl::compressed_matrix< NumericT > &L_trans, chow_patel_tag const &tag)
 Implementation of the parallel ICC0 factorization, Algorithm 3 in Chow-Patel paper. More...
 
template<typename NumericT >
void precondition (viennacl::compressed_matrix< NumericT > const &A, viennacl::compressed_matrix< NumericT > &L, viennacl::vector< NumericT > &diag_L, viennacl::compressed_matrix< NumericT > &U, viennacl::vector< NumericT > &diag_U, chow_patel_tag const &tag)
 Implementation of the parallel ILU0 factorization, Algorithm 2 in Chow-Patel paper. More...
 
template<typename NumericT , unsigned int AlignmentV>
void level_scheduling_setup_impl (viennacl::compressed_matrix< NumericT, AlignmentV > const &LU, viennacl::vector< NumericT > const &diagonal_LU, std::list< viennacl::backend::mem_handle > &row_index_arrays, std::list< viennacl::backend::mem_handle > &row_buffers, std::list< viennacl::backend::mem_handle > &col_buffers, std::list< viennacl::backend::mem_handle > &element_buffers, std::list< vcl_size_t > &row_elimination_num_list, bool setup_U)
 
template<typename NumericT , unsigned int AlignmentV>
void level_scheduling_setup_L (viennacl::compressed_matrix< NumericT, AlignmentV > const &LU, viennacl::vector< NumericT > const &diagonal_LU, std::list< viennacl::backend::mem_handle > &row_index_arrays, std::list< viennacl::backend::mem_handle > &row_buffers, std::list< viennacl::backend::mem_handle > &col_buffers, std::list< viennacl::backend::mem_handle > &element_buffers, std::list< vcl_size_t > &row_elimination_num_list)
 
template<typename NumericT , unsigned int AlignmentV>
void level_scheduling_setup_U (viennacl::compressed_matrix< NumericT, AlignmentV > const &LU, viennacl::vector< NumericT > const &diagonal_LU, std::list< viennacl::backend::mem_handle > &row_index_arrays, std::list< viennacl::backend::mem_handle > &row_buffers, std::list< viennacl::backend::mem_handle > &col_buffers, std::list< viennacl::backend::mem_handle > &element_buffers, std::list< vcl_size_t > &row_elimination_num_list)
 
template<typename NumericT >
void level_scheduling_substitute (viennacl::vector< NumericT > &vec, std::list< viennacl::backend::mem_handle > const &row_index_arrays, std::list< viennacl::backend::mem_handle > const &row_buffers, std::list< viennacl::backend::mem_handle > const &col_buffers, std::list< viennacl::backend::mem_handle > const &element_buffers, std::list< vcl_size_t > const &row_elimination_num_list)
 
template<typename IndexT , typename NumericT >
IndexT merge_subtract_sparse_rows (IndexT const *w_coords, NumericT const *w_elements, IndexT w_size, IndexT const *u_coords, NumericT const *u_elements, IndexT u_size, NumericT alpha, IndexT *z_coords, NumericT *z_elements)
 Subtracts a scaled sparse vector u from a sparse vector w and writes the output to z: z = w - alpha * u. More...
 
template<typename SizeT , typename NumericT >
void insert_with_value_sort (std::vector< std::pair< SizeT, NumericT > > &map, SizeT index, NumericT value)
 
template<typename NumericT , typename B >
bool op_aliasing (vector_base< NumericT > const &, B const &)
 
template<typename NumericT >
bool op_aliasing (vector_base< NumericT > const &lhs, vector_base< NumericT > const &b)
 
template<typename NumericT , typename LhsT , typename RhsT , typename OpT >
bool op_aliasing (vector_base< NumericT > const &lhs, vector_expression< const LhsT, const RhsT, OpT > const &rhs)
 
template<typename NumericT , typename B >
bool op_aliasing (matrix_base< NumericT > const &, B const &)
 
template<typename NumericT >
bool op_aliasing (matrix_base< NumericT > const &lhs, matrix_base< NumericT > const &b)
 
template<typename NumericT , typename LhsT , typename RhsT , typename OpT >
bool op_aliasing (matrix_base< NumericT > const &lhs, matrix_expression< const LhsT, const RhsT, OpT > const &rhs)
 
template<typename NumericT , typename SolverTagT >
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 \ B. More...
 
template<typename NumericT , typename SolverTagT >
void inplace_solve_vec_kernel (const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, SolverTagT)
 
template<typename MatrixT1 , typename MatrixT2 , typename SolverTagT >
void inplace_solve_lower_impl (MatrixT1 const &A, MatrixT2 &B, SolverTagT)
 
template<typename MatrixT1 , typename MatrixT2 >
void inplace_solve_impl (MatrixT1 const &A, MatrixT2 &B, viennacl::linalg::lower_tag)
 
template<typename MatrixT1 , typename MatrixT2 >
void inplace_solve_impl (MatrixT1 const &A, MatrixT2 &B, viennacl::linalg::unit_lower_tag)
 
template<typename MatrixT1 , typename MatrixT2 , typename SolverTagT >
void inplace_solve_upper_impl (MatrixT1 const &A, MatrixT2 &B, SolverTagT)
 
template<typename MatrixT1 , typename MatrixT2 >
void inplace_solve_impl (MatrixT1 const &A, MatrixT2 &B, viennacl::linalg::upper_tag)
 
template<typename MatrixT1 , typename MatrixT2 >
void inplace_solve_impl (MatrixT1 const &A, MatrixT2 &B, viennacl::linalg::unit_upper_tag)
 
template<typename MatrixT1 , typename VectorT , typename SolverTagT >
void inplace_solve_lower_vec_impl (MatrixT1 const &A, VectorT &b, SolverTagT)
 
template<typename MatrixT1 , typename VectorT >
void inplace_solve_vec_impl (MatrixT1 const &A, VectorT &B, viennacl::linalg::lower_tag)
 
template<typename MatrixT1 , typename VectorT >
void inplace_solve_vec_impl (MatrixT1 const &A, VectorT &B, viennacl::linalg::unit_lower_tag)
 
template<typename MatrixT1 , typename VectorT , typename SolverTagT >
void inplace_solve_upper_vec_impl (MatrixT1 const &A, VectorT &b, SolverTagT)
 
template<typename MatrixT1 , typename VectorT >
void inplace_solve_vec_impl (MatrixT1 const &A, VectorT &b, viennacl::linalg::upper_tag)
 
template<typename MatrixT1 , typename VectorT >
void inplace_solve_vec_impl (MatrixT1 const &A, VectorT &b, viennacl::linalg::unit_upper_tag)
 
template<typename SrcVectorT , typename DestVectorT >
void gmres_copy_helper (SrcVectorT const &src, DestVectorT &dest, vcl_size_t len, vcl_size_t start=0)
 
template<typename NumericT , typename DestVectorT >
void gmres_copy_helper (viennacl::vector< NumericT > const &src, DestVectorT &dest, vcl_size_t len, vcl_size_t start=0)
 
template<typename VectorT , typename NumericT >
void gmres_setup_householder_vector (VectorT const &input_vec, VectorT &hh_vec, NumericT &beta, NumericT &mu, vcl_size_t j)
 Computes the householder vector 'hh_vec' which rotates 'input_vec' such that all entries below the j-th entry of 'v' become zero. More...
 
template<typename VectorT , typename NumericT >
void gmres_householder_reflect (VectorT &x, VectorT const &h, NumericT beta)
 
template<typename MatrixType , typename ScalarType >
viennacl::vector< ScalarTypepipelined_solve (MatrixType const &A, viennacl::vector< ScalarType > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< ScalarType > const &, ScalarType, void *)=NULL, void *monitor_data=NULL)
 Implementation of a pipelined GMRES solver without preconditioner. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::compressed_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::coordinate_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::sliced_ell_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename NumericT >
viennacl::vector< NumericTsolve_impl (viennacl::hyb_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, gmres_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
 Overload for the pipelined CG implementation for the ViennaCL sparse matrix types. More...
 
template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT solve_impl (MatrixT const &matrix, VectorT const &rhs, gmres_tag const &tag, PreconditionerT const &precond, bool(*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)=NULL, void *monitor_data=NULL)
 Implementation of the GMRES solver. More...
 
template<typename NumericT >
void inverse_iteration (std::vector< NumericT > const &alphas, std::vector< NumericT > const &betas, NumericT &eigenvalue, std::vector< NumericT > &eigenvector)
 Inverse iteration for finding an eigenvector for an eigenvalue. More...
 
template<typename MatrixT , typename DenseMatrixT , typename NumericT >
std::vector< NumericTlanczosPRO (MatrixT const &A, vector_base< NumericT > &r, DenseMatrixT &eigenvectors_A, vcl_size_t size, lanczos_tag const &tag, bool compute_eigenvectors)
 Implementation of the Lanczos PRO algorithm (partial reorthogonalization) More...
 
template<typename MatrixT , typename DenseMatrixT , typename NumericT >
std::vector< NumericTlanczos (MatrixT const &A, vector_base< NumericT > &r, DenseMatrixT &eigenvectors_A, vcl_size_t krylov_dim, lanczos_tag const &tag, bool compute_eigenvectors)
 Implementation of the Lanczos FRO algorithm. More...
 
template<typename ScalarType >
void level_scheduling_substitute (vector< ScalarType > &vec, viennacl::backend::mem_handle const &row_index_array, viennacl::backend::mem_handle const &row_buffer, viennacl::backend::mem_handle const &col_buffer, viennacl::backend::mem_handle const &element_buffer, vcl_size_t num_rows)
 
template<typename SCALARTYPE >
SCALARTYPE pythag (SCALARTYPE a, SCALARTYPE b)
 
template<typename SCALARTYPE >
SCALARTYPE sign (SCALARTYPE val)
 
template<typename VectorType >
VectorType::value_type norm_lcl (VectorType const &x, vcl_size_t size)
 
template<typename VectorType >
void normalize (VectorType &x, vcl_size_t size)
 
template<typename VectorType >
void householder_vector (VectorType &v, vcl_size_t start)
 
template<typename SCALARTYPE >
void transpose (matrix_base< SCALARTYPE > &A)
 
template<typename T >
void cdiv (T xr, T xi, T yr, T yi, T &cdivr, T &cdivi)
 
template<typename SCALARTYPE >
void prepare_householder_vector (matrix_base< SCALARTYPE > &A, vector_base< SCALARTYPE > &D, vcl_size_t size, vcl_size_t row_start, vcl_size_t col_start, vcl_size_t start, bool is_column)
 
template<typename SCALARTYPE >
void final_iter_update_gpu (matrix_base< SCALARTYPE > &A, int n, int last_n, SCALARTYPE q, SCALARTYPE p)
 
template<typename SCALARTYPE , typename VectorType >
void update_float_QR_column_gpu (matrix_base< SCALARTYPE > &A, const VectorType &buf, viennacl::vector< SCALARTYPE > &buf_vcl, int m, int n, int last_n, bool)
 
template<typename SCALARTYPE , typename MatrixT >
void final_iter_update (MatrixT &A, int n, int last_n, SCALARTYPE q, SCALARTYPE p)
 
template<typename SCALARTYPE , typename MatrixT >
void update_float_QR_column (MatrixT &A, const std::vector< SCALARTYPE > &buf, int m, int n, int last_i, bool is_triangular)
 
template<typename SCALARTYPE , typename VectorType >
void hqr2 (viennacl::matrix< SCALARTYPE > &vcl_H, viennacl::matrix< SCALARTYPE > &V, VectorType &d, VectorType &e)
 
template<typename SCALARTYPE >
bool householder_twoside (matrix_base< SCALARTYPE > &A, matrix_base< SCALARTYPE > &Q, vector_base< SCALARTYPE > &D, vcl_size_t start)
 
template<typename SCALARTYPE >
void tridiagonal_reduction (matrix_base< SCALARTYPE > &A, matrix_base< SCALARTYPE > &Q)
 
template<typename SCALARTYPE >
void qr_method (viennacl::matrix< SCALARTYPE > &A, viennacl::matrix< SCALARTYPE > &Q, std::vector< SCALARTYPE > &D, std::vector< SCALARTYPE > &E, bool is_symmetric=true)
 
template<typename MatrixType , typename VectorType >
MatrixType::value_type setup_householder_vector_ublas (MatrixType const &A, VectorType &v, MatrixType &matrix_1x1, vcl_size_t j)
 
template<typename MatrixType , typename VectorType >
viennacl::result_of::cpu_value_type
< typename
MatrixType::value_type >::type 
setup_householder_vector_viennacl (MatrixType const &A, VectorType &v, MatrixType &matrix_1x1, vcl_size_t j)
 
template<typename MatrixType , typename VectorType , typename ScalarType >
void householder_reflect (MatrixType &A, VectorType &v, ScalarType beta, vcl_size_t j, vcl_size_t k)
 
template<typename MatrixType , typename VectorType , typename ScalarType >
void householder_reflect_ublas (MatrixType &A, VectorType &v, MatrixType &matrix_1x1, ScalarType beta, vcl_size_t j, vcl_size_t k)
 
template<typename MatrixType , typename VectorType , typename ScalarType >
void householder_reflect_viennacl (MatrixType &A, VectorType &v, MatrixType &matrix_1x1, ScalarType beta, vcl_size_t j, vcl_size_t k)
 
template<typename MatrixType , typename VectorType , typename ScalarType >
void householder_reflect (MatrixType &A, VectorType &v, ScalarType beta, vcl_size_t j)
 
template<typename MatrixType , typename VectorType >
void write_householder_to_A (MatrixType &A, VectorType const &v, vcl_size_t j)
 
template<typename MatrixType , typename VectorType >
void write_householder_to_A_ublas (MatrixType &A, VectorType const &v, vcl_size_t j)
 
template<typename MatrixType , typename VectorType >
void write_householder_to_A_viennacl (MatrixType &A, VectorType const &v, vcl_size_t j)
 
template<typename MatrixType >
std::vector< typename
MatrixType::value_type > 
inplace_qr_ublas (MatrixType &A, vcl_size_t block_size=32)
 Implementation of inplace-QR factorization for a general Boost.uBLAS compatible matrix A. More...
 
template<typename MatrixType >
std::vector< typename
viennacl::result_of::cpu_value_type
< typename
MatrixType::value_type >::type > 
inplace_qr_viennacl (MatrixType &A, vcl_size_t block_size=16)
 Implementation of a OpenCL-only QR factorization for GPUs (or multi-core CPU). DEPRECATED! Use only if you're curious and interested in playing a bit with a GPU-only implementation. More...
 
template<typename MatrixType >
std::vector< typename
viennacl::result_of::cpu_value_type
< typename
MatrixType::value_type >::type > 
inplace_qr_hybrid (MatrixType &A, vcl_size_t block_size=16)
 Implementation of a hybrid QR factorization using uBLAS on the CPU and ViennaCL for GPUs (or multi-core CPU) More...
 
template<typename SparseMatrixType , typename SCALARTYPE , unsigned int VEC_ALIGNMENT>
viennacl::enable_if
< viennacl::is_any_sparse_matrix
< SparseMatrixType >::value >
::type 
row_info (SparseMatrixType const &mat, vector< SCALARTYPE, VEC_ALIGNMENT > &vec, row_info_types info_selector)
 
template<typename SparseMatrixType , class ScalarType , typename SOLVERTAG >
viennacl::enable_if
< viennacl::is_any_sparse_matrix
< SparseMatrixType >::value >
::type 
block_inplace_solve (const matrix_expression< const SparseMatrixType, const SparseMatrixType, op_trans > &mat, viennacl::backend::mem_handle const &block_index_array, vcl_size_t num_blocks, viennacl::vector_base< ScalarType > const &mat_diagonal, viennacl::vector_base< ScalarType > &vec, SOLVERTAG tag)
 
template<typename MatrixType , typename VectorType >
void givens_prev (MatrixType &matrix, VectorType &tmp1, VectorType &tmp2, int n, int l, int k)
 
template<typename MatrixType , typename VectorType >
void change_signs (MatrixType &matrix, VectorType &signs, int n)
 
template<typename MatrixType , typename CPU_VectorType >
void svd_qr_shift (MatrixType &vcl_u, MatrixType &vcl_v, CPU_VectorType &q, CPU_VectorType &e)
 
template<typename SCALARTYPE , unsigned int ALIGNMENT>
bool householder_c (viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &A, viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &Q, viennacl::vector< SCALARTYPE, ALIGNMENT > &D, vcl_size_t row_start, vcl_size_t col_start)
 
template<typename SCALARTYPE , unsigned int ALIGNMENT>
bool householder_r (viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &A, viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &Q, viennacl::vector< SCALARTYPE, ALIGNMENT > &D, vcl_size_t row_start, vcl_size_t col_start)
 
template<typename SCALARTYPE , unsigned int ALIGNMENT>
void bidiag (viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &Ai, viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &QL, viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &QR)
 

Detailed Description

Namespace holding implementation details for linear algebra routines. Usually not of interest for a library user.

Enumeration Type Documentation

Enumerator
SPARSE_ROW_NORM_INF 
SPARSE_ROW_NORM_1 
SPARSE_ROW_NORM_2 
SPARSE_ROW_DIAGONAL 

Definition at line 837 of file forwards.h.

Function Documentation

template<typename NumericT >
void viennacl::linalg::detail::amg_galerkin_prod ( compressed_matrix< NumericT > &  A_fine,
compressed_matrix< NumericT > &  P,
compressed_matrix< NumericT > &  R,
compressed_matrix< NumericT > &  A_coarse 
)

Sparse Galerkin product: Calculates A_coarse = trans(P)*A_fine*P = R*A_fine*P.

Parameters
A_fineOperator matrix on fine grid (quadratic)
PProlongation/Interpolation matrix
RRestriction matrix
A_coarseResult matrix on coarse grid (Galerkin operator)

Definition at line 78 of file amg.hpp.

template<typename MatrixT , typename InternalT1 , typename InternalT2 >
void viennacl::linalg::detail::amg_init ( MatrixT const &  mat,
InternalT1 &  list_of_A,
InternalT1 &  list_of_P,
InternalT1 &  list_of_R,
InternalT2 &  list_of_amg_level_context,
amg_tag &  tag 
)

Initialize AMG preconditioner.

Parameters
matSystem matrix
list_of_AOperator matrices on all levels
list_of_PProlongation/Interpolation operators on all levels
list_of_RRestriction operators on all levels
list_of_amg_level_contextAuxiliary datastructures for managing the grid hierarchy (coarse nodes, etc.)
tagAMG preconditioner tag

Definition at line 169 of file amg.hpp.

template<typename NumericT , typename SparseMatrixT >
void viennacl::linalg::detail::amg_lu ( viennacl::matrix< NumericT > &  op,
SparseMatrixT const &  A,
amg_tag const &  tag 
)

Pre-compute LU factorization for direct solve (ublas library).

Speeds up precondition phase as this is computed only once overall instead of once per iteration.

Parameters
opOperator matrix for direct solve
AOperator matrix on coarsest level
tagAMG preconditioner tag

Definition at line 236 of file amg.hpp.

template<typename NumericT , typename AMGContextListT >
vcl_size_t viennacl::linalg::detail::amg_setup ( std::vector< compressed_matrix< NumericT > > &  list_of_A,
std::vector< compressed_matrix< NumericT > > &  list_of_P,
std::vector< compressed_matrix< NumericT > > &  list_of_R,
AMGContextListT &  list_of_amg_level_context,
amg_tag &  tag 
)

Setup AMG preconditioner.

Parameters
list_of_AOperator matrices on all levels
list_of_PProlongation/Interpolation operators on all levels
list_of_RRestriction operators on all levels
list_of_amg_level_contextAuxiliary datastructures for managing the grid hierarchy (coarse nodes, etc.)
tagAMG preconditioner tag

Definition at line 105 of file amg.hpp.

template<typename InternalVectorT , typename SparseMatrixT >
void viennacl::linalg::detail::amg_setup_apply ( InternalVectorT &  result,
InternalVectorT &  result_backup,
InternalVectorT &  rhs,
InternalVectorT &  residual,
SparseMatrixT const &  A,
vcl_size_t  coarse_levels,
amg_tag const &  tag 
)

Setup data structures for precondition phase for later use on the GPU.

Parameters
resultResult vector on all levels
result_backupCopy of result vector on all levels
rhsRHS vector on all levels
residualResidual vector on all levels
AOperators matrices on all levels from setup phase
coarse_levelsNumber of coarse levels for which the datastructures should be set up.
tagAMG preconditioner tag

Definition at line 199 of file amg.hpp.

template<typename SCALARTYPE , unsigned int ALIGNMENT>
void viennacl::linalg::detail::bidiag ( viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  Ai,
viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  QL,
viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  QR 
)

Definition at line 459 of file svd.hpp.

template<typename NumericT >
void viennacl::linalg::detail::bisectLarge ( const InputData< NumericT > &  input,
ResultDataLarge< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  lg,
const NumericT  ug,
const NumericT  precision 
)

Definition at line 93 of file bisect_kernel_calls.hpp.

template<typename NumericT >
void viennacl::linalg::detail::bisectLarge_MultIntervals ( const InputData< NumericT > &  input,
ResultDataLarge< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  precision 
)

Definition at line 160 of file bisect_kernel_calls.hpp.

template<typename NumericT >
void viennacl::linalg::detail::bisectLarge_OneIntervals ( const InputData< NumericT > &  input,
ResultDataLarge< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  precision 
)

Definition at line 128 of file bisect_kernel_calls.hpp.

template<typename NumericT >
void viennacl::linalg::detail::bisectSmall ( const InputData< NumericT > &  input,
ResultDataSmall< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  lg,
const NumericT  ug,
const NumericT  precision 
)

Definition at line 59 of file bisect_kernel_calls.hpp.

template<typename SparseMatrixType , class ScalarType , typename SOLVERTAG >
viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value>::type viennacl::linalg::detail::block_inplace_solve ( const matrix_expression< const SparseMatrixType, const SparseMatrixType, op_trans > &  mat,
viennacl::backend::mem_handle const &  block_index_array,
vcl_size_t  num_blocks,
viennacl::vector_base< ScalarType > const &  mat_diagonal,
viennacl::vector_base< ScalarType > &  vec,
SOLVERTAG  tag 
)

Definition at line 332 of file sparse_matrix_operations.hpp.

template<typename T >
void viennacl::linalg::detail::cdiv ( xr,
xi,
yr,
yi,
T &  cdivr,
T &  cdivi 
)

Definition at line 142 of file qr-method-common.hpp.

template<typename MatrixType , typename VectorType >
void viennacl::linalg::detail::change_signs ( MatrixType &  matrix,
VectorType &  signs,
int  n 
)

Definition at line 78 of file svd.hpp.

template<typename NumericT >
void viennacl::linalg::detail::computeEigenvaluesLargeMatrix ( InputData< NumericT > &  input,
ResultDataLarge< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  lg,
const NumericT  ug,
const NumericT  precision 
)

Run the kernels to compute the eigenvalues for large matrices

Parameters
inputhandles to input data
resulthandles to result data
mat_sizematrix size
precisiondesired precision of eigenvalues
lglower limit of Gerschgorin interval
ugupper limit of Gerschgorin interval

Definition at line 60 of file bisect_large.hpp.

template<typename NumericT >
void viennacl::linalg::detail::computeEigenvaluesSmallMatrix ( const InputData< NumericT > &  input,
ResultDataSmall< NumericT > &  result,
const unsigned int  mat_size,
const NumericT  lg,
const NumericT  ug,
const NumericT  precision 
)

Determine eigenvalues for matrices smaller than MAX_SMALL_MATRIX

Parameters
inputhandles to input data of kernel
resulthandles to result of kernel
mat_sizematrix size
lglower limit of Gerschgorin interval
ugupper limit of Gerschgorin interval
precisiondesired precision of eigenvalues

Definition at line 61 of file bisect_small.hpp.

template<typename NumericT >
void viennacl::linalg::detail::computeGerschgorin ( std::vector< NumericT > &  d,
std::vector< NumericT > &  s,
unsigned int  n,
NumericT lg,
NumericT ug 
)

Compute Gerschgorin interval for symmetric, tridiagonal matrix

Parameters
ddiagonal elements
ssuperdiagonal elements
nsize of matrix
lglower limit of Gerschgorin interval
ugupper limit of Gerschgorin interval

Definition at line 55 of file gerschgorin.hpp.

template<typename NumericT , typename OtherVectorT >
void viennacl::linalg::detail::copy_vec_to_vec ( viennacl::vector< NumericT > const &  src,
OtherVectorT &  dest 
)

overloaded function for copying vectors

Definition at line 44 of file bisect.hpp.

template<typename OtherVectorT , typename NumericT >
void viennacl::linalg::detail::copy_vec_to_vec ( OtherVectorT const &  src,
viennacl::vector< NumericT > &  dest 
)

Definition at line 50 of file bisect.hpp.

template<typename VectorT1 , typename VectorT2 >
void viennacl::linalg::detail::copy_vec_to_vec ( VectorT1 const &  src,
VectorT2 &  dest 
)

Definition at line 56 of file bisect.hpp.

template<typename NumericT >
void viennacl::linalg::detail::extract_block_matrix ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::compressed_matrix< NumericT > &  diagonal_block_A,
vcl_size_t  start_index,
vcl_size_t  stop_index 
)

Extracts a diagonal block from a larger system matrix.

Parameters
AThe full matrix
diagonal_block_AThe output matrix, to which the extracted block is written to
start_indexFirst row- and column-index of the block
stop_indexFirst row- and column-index beyond the block

Definition at line 79 of file block_ilu.hpp.

template<typename SCALARTYPE , typename MatrixT >
void viennacl::linalg::detail::final_iter_update ( MatrixT &  A,
int  n,
int  last_n,
SCALARTYPE  q,
SCALARTYPE  p 
)

Definition at line 128 of file qr-method.hpp.

template<typename SCALARTYPE >
void viennacl::linalg::detail::final_iter_update_gpu ( matrix_base< SCALARTYPE > &  A,
int  n,
int  last_n,
SCALARTYPE  q,
SCALARTYPE  p 
)

Definition at line 42 of file qr-method.hpp.

unsigned int viennacl::linalg::detail::getNumBlocksLinear ( const unsigned int  num_threads,
const unsigned int  num_threads_block 
)
inline

Get the number of blocks that are required to process num_threads with num_threads_blocks threads per block

Definition at line 96 of file util.hpp.

template<typename MatrixType , typename VectorType >
void viennacl::linalg::detail::givens_prev ( MatrixType &  matrix,
VectorType &  tmp1,
VectorType &  tmp2,
int  n,
int  l,
int  k 
)

Definition at line 48 of file svd.hpp.

template<typename SrcVectorT , typename DestVectorT >
void viennacl::linalg::detail::gmres_copy_helper ( SrcVectorT const &  src,
DestVectorT &  dest,
vcl_size_t  len,
vcl_size_t  start = 0 
)

Definition at line 107 of file gmres.hpp.

template<typename NumericT , typename DestVectorT >
void viennacl::linalg::detail::gmres_copy_helper ( viennacl::vector< NumericT > const &  src,
DestVectorT &  dest,
vcl_size_t  len,
vcl_size_t  start = 0 
)

Definition at line 114 of file gmres.hpp.

template<typename VectorT , typename NumericT >
void viennacl::linalg::detail::gmres_householder_reflect ( VectorT &  x,
VectorT const &  h,
NumericT  beta 
)

Definition at line 162 of file gmres.hpp.

template<typename VectorT , typename NumericT >
void viennacl::linalg::detail::gmres_setup_householder_vector ( VectorT const &  input_vec,
VectorT &  hh_vec,
NumericT beta,
NumericT mu,
vcl_size_t  j 
)

Computes the householder vector 'hh_vec' which rotates 'input_vec' such that all entries below the j-th entry of 'v' become zero.

Parameters
input_vecThe input vector
hh_vecThe householder vector defining the relection (I - beta * hh_vec * hh_vec^T)
betaThe coefficient beta in (I - beta * hh_vec * hh_vec^T)
muThe norm of the input vector part relevant for the reflection: norm_2(input_vec[j:size])
jIndex of the last nonzero index in 'input_vec' after applying the reflection

Definition at line 131 of file gmres.hpp.

template<typename SCALARTYPE , unsigned int ALIGNMENT>
bool viennacl::linalg::detail::householder_c ( viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  A,
viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  Q,
viennacl::vector< SCALARTYPE, ALIGNMENT > &  D,
vcl_size_t  row_start,
vcl_size_t  col_start 
)

Definition at line 328 of file svd.hpp.

template<typename SCALARTYPE , unsigned int ALIGNMENT>
bool viennacl::linalg::detail::householder_r ( viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  A,
viennacl::matrix< SCALARTYPE, row_major, ALIGNMENT > &  Q,
viennacl::vector< SCALARTYPE, ALIGNMENT > &  D,
vcl_size_t  row_start,
vcl_size_t  col_start 
)

Definition at line 415 of file svd.hpp.

template<typename MatrixType , typename VectorType , typename ScalarType >
void viennacl::linalg::detail::householder_reflect ( MatrixType &  A,
VectorType &  v,
ScalarType  beta,
vcl_size_t  j,
vcl_size_t  k 
)

Definition at line 134 of file qr.hpp.

template<typename MatrixType , typename VectorType , typename ScalarType >
void viennacl::linalg::detail::householder_reflect ( MatrixType &  A,
VectorType &  v,
ScalarType  beta,
vcl_size_t  j 
)

Definition at line 184 of file qr.hpp.

template<typename MatrixType , typename VectorType , typename ScalarType >
void viennacl::linalg::detail::householder_reflect_ublas ( MatrixType &  A,
VectorType &  v,
MatrixType &  matrix_1x1,
ScalarType  beta,
vcl_size_t  j,
vcl_size_t  k 
)

Definition at line 147 of file qr.hpp.

template<typename MatrixType , typename VectorType , typename ScalarType >
void viennacl::linalg::detail::householder_reflect_viennacl ( MatrixType &  A,
VectorType &  v,
MatrixType &  matrix_1x1,
ScalarType  beta,
vcl_size_t  j,
vcl_size_t  k 
)

Definition at line 161 of file qr.hpp.

template<typename SCALARTYPE >
bool viennacl::linalg::detail::householder_twoside ( matrix_base< SCALARTYPE > &  A,
matrix_base< SCALARTYPE > &  Q,
vector_base< SCALARTYPE > &  D,
vcl_size_t  start 
)

Definition at line 694 of file qr-method.hpp.

template<typename VectorType >
void viennacl::linalg::detail::householder_vector ( VectorType &  v,
vcl_size_t  start 
)

Definition at line 97 of file qr-method-common.hpp.

template<typename SCALARTYPE , typename VectorType >
void viennacl::linalg::detail::hqr2 ( viennacl::matrix< SCALARTYPE > &  vcl_H,
viennacl::matrix< SCALARTYPE > &  V,
VectorType &  d,
VectorType &  e 
)

Definition at line 236 of file qr-method.hpp.

template<typename MatrixType >
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixType::value_type >::type > viennacl::linalg::detail::inplace_qr_hybrid ( MatrixType &  A,
vcl_size_t  block_size = 16 
)

Implementation of a hybrid QR factorization using uBLAS on the CPU and ViennaCL for GPUs (or multi-core CPU)

Prefer the use of the convenience interface inplace_qr()

Parameters
AA dense ViennaCL matrix to be factored
block_sizeThe block size to be used. The number of columns of A must be a multiple of block_size

Definition at line 426 of file qr.hpp.

template<typename MatrixType >
std::vector<typename MatrixType::value_type> viennacl::linalg::detail::inplace_qr_ublas ( MatrixType &  A,
vcl_size_t  block_size = 32 
)

Implementation of inplace-QR factorization for a general Boost.uBLAS compatible matrix A.

Parameters
AA dense compatible to Boost.uBLAS
block_sizeThe block size to be used. The number of columns of A must be a multiple of block_size

Definition at line 228 of file qr.hpp.

template<typename MatrixType >
std::vector< typename viennacl::result_of::cpu_value_type< typename MatrixType::value_type >::type > viennacl::linalg::detail::inplace_qr_viennacl ( MatrixType &  A,
vcl_size_t  block_size = 16 
)

Implementation of a OpenCL-only QR factorization for GPUs (or multi-core CPU). DEPRECATED! Use only if you're curious and interested in playing a bit with a GPU-only implementation.

Performance is rather poor at small matrix sizes. Prefer the use of the hybrid version, which is automatically chosen using the interface function inplace_qr()

Parameters
AA dense ViennaCL matrix to be factored
block_sizeThe block size to be used. The number of columns of A must be a multiple of block_size

Definition at line 322 of file qr.hpp.

template<typename MatrixT1 , typename MatrixT2 >
void viennacl::linalg::detail::inplace_solve_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
viennacl::linalg::lower_tag   
)

Definition at line 155 of file direct_solve.hpp.

template<typename MatrixT1 , typename MatrixT2 >
void viennacl::linalg::detail::inplace_solve_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
viennacl::linalg::unit_lower_tag   
)

Definition at line 161 of file direct_solve.hpp.

template<typename MatrixT1 , typename MatrixT2 >
void viennacl::linalg::detail::inplace_solve_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
viennacl::linalg::upper_tag   
)

Definition at line 198 of file direct_solve.hpp.

template<typename MatrixT1 , typename MatrixT2 >
void viennacl::linalg::detail::inplace_solve_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
viennacl::linalg::unit_upper_tag   
)

Definition at line 204 of file direct_solve.hpp.

template<typename NumericT , typename SolverTagT >
void viennacl::linalg::detail::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 \ B.

Parameters
AThe system matrix
BThe matrix of row vectors, where the solution is directly written to

Definition at line 62 of file direct_solve.hpp.

template<typename MatrixT1 , typename MatrixT2 , typename SolverTagT >
void viennacl::linalg::detail::inplace_solve_lower_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
SolverTagT   
)

Definition at line 125 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT , typename SolverTagT >
void viennacl::linalg::detail::inplace_solve_lower_vec_impl ( MatrixT1 const &  A,
VectorT &  b,
SolverTagT   
)

Definition at line 369 of file direct_solve.hpp.

template<typename MatrixT1 , typename MatrixT2 , typename SolverTagT >
void viennacl::linalg::detail::inplace_solve_upper_impl ( MatrixT1 const &  A,
MatrixT2 &  B,
SolverTagT   
)

Definition at line 167 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT , typename SolverTagT >
void viennacl::linalg::detail::inplace_solve_upper_vec_impl ( MatrixT1 const &  A,
VectorT &  b,
SolverTagT   
)

Definition at line 407 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT >
void viennacl::linalg::detail::inplace_solve_vec_impl ( MatrixT1 const &  A,
VectorT &  B,
viennacl::linalg::lower_tag   
)

Definition at line 395 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT >
void viennacl::linalg::detail::inplace_solve_vec_impl ( MatrixT1 const &  A,
VectorT &  B,
viennacl::linalg::unit_lower_tag   
)

Definition at line 401 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT >
void viennacl::linalg::detail::inplace_solve_vec_impl ( MatrixT1 const &  A,
VectorT &  b,
viennacl::linalg::upper_tag   
)

Definition at line 433 of file direct_solve.hpp.

template<typename MatrixT1 , typename VectorT >
void viennacl::linalg::detail::inplace_solve_vec_impl ( MatrixT1 const &  A,
VectorT &  b,
viennacl::linalg::unit_upper_tag   
)

Definition at line 439 of file direct_solve.hpp.

template<typename NumericT , typename SolverTagT >
void viennacl::linalg::detail::inplace_solve_vec_kernel ( const matrix_base< NumericT > &  mat,
const vector_base< NumericT > &  vec,
SolverTagT   
)

Definition at line 94 of file direct_solve.hpp.

template<typename SizeT , typename NumericT >
void viennacl::linalg::detail::insert_with_value_sort ( std::vector< std::pair< SizeT, NumericT > > &  map,
SizeT  index,
NumericT  value 
)

Definition at line 164 of file ilut.hpp.

template<typename NumericT >
void viennacl::linalg::detail::inverse_iteration ( std::vector< NumericT > const &  alphas,
std::vector< NumericT > const &  betas,
NumericT eigenvalue,
std::vector< NumericT > &  eigenvector 
)

Inverse iteration for finding an eigenvector for an eigenvalue.

beta[0] to be ignored for consistency.

Definition at line 109 of file lanczos.hpp.

template<typename MatrixT , typename DenseMatrixT , typename NumericT >
std::vector<NumericT> viennacl::linalg::detail::lanczos ( MatrixT const &  A,
vector_base< NumericT > &  r,
DenseMatrixT &  eigenvectors_A,
vcl_size_t  krylov_dim,
lanczos_tag const &  tag,
bool  compute_eigenvectors 
)

Implementation of the Lanczos FRO algorithm.

Parameters
AThe system matrix
rRandom start vector
eigenvectors_AA dense matrix in which the eigenvectors of A will be stored. Both row- and column-major matrices are supported.
krylov_dimSize of krylov-space
tagThe Lanczos tag holding tolerances, etc.
compute_eigenvectorsBoolean flag. If true, eigenvectors are computed. Otherwise the routine returns after calculating eigenvalues.
Returns
Returns the eigenvalues (number of eigenvalues equals size of krylov-space)

Definition at line 345 of file lanczos.hpp.

template<typename MatrixT , typename DenseMatrixT , typename NumericT >
std::vector<NumericT> viennacl::linalg::detail::lanczosPRO ( MatrixT const &  A,
vector_base< NumericT > &  r,
DenseMatrixT &  eigenvectors_A,
vcl_size_t  size,
lanczos_tag const &  tag,
bool  compute_eigenvectors 
)

Implementation of the Lanczos PRO algorithm (partial reorthogonalization)

Parameters
AThe system matrix
rRandom start vector
eigenvectors_ADense matrix holding the eigenvectors of A (one eigenvector per column)
sizeSize of krylov-space
tagLanczos_tag with several options for the algorithm
compute_eigenvectorsBoolean flag. If true, eigenvectors are computed. Otherwise the routine returns after calculating eigenvalues.
Returns
Returns the eigenvalues (number of eigenvalues equals size of krylov-space)

Definition at line 158 of file lanczos.hpp.

template<typename NumericT , unsigned int AlignmentV>
void viennacl::linalg::detail::level_scheduling_setup_impl ( viennacl::compressed_matrix< NumericT, AlignmentV > const &  LU,
viennacl::vector< NumericT > const &  diagonal_LU,
std::list< viennacl::backend::mem_handle > &  row_index_arrays,
std::list< viennacl::backend::mem_handle > &  row_buffers,
std::list< viennacl::backend::mem_handle > &  col_buffers,
std::list< viennacl::backend::mem_handle > &  element_buffers,
std::list< vcl_size_t > &  row_elimination_num_list,
bool  setup_U 
)

Definition at line 51 of file common.hpp.

template<typename NumericT , unsigned int AlignmentV>
void viennacl::linalg::detail::level_scheduling_setup_L ( viennacl::compressed_matrix< NumericT, AlignmentV > const &  LU,
viennacl::vector< NumericT > const &  diagonal_LU,
std::list< viennacl::backend::mem_handle > &  row_index_arrays,
std::list< viennacl::backend::mem_handle > &  row_buffers,
std::list< viennacl::backend::mem_handle > &  col_buffers,
std::list< viennacl::backend::mem_handle > &  element_buffers,
std::list< vcl_size_t > &  row_elimination_num_list 
)

Definition at line 191 of file common.hpp.

template<typename NumericT , unsigned int AlignmentV>
void viennacl::linalg::detail::level_scheduling_setup_U ( viennacl::compressed_matrix< NumericT, AlignmentV > const &  LU,
viennacl::vector< NumericT > const &  diagonal_LU,
std::list< viennacl::backend::mem_handle > &  row_index_arrays,
std::list< viennacl::backend::mem_handle > &  row_buffers,
std::list< viennacl::backend::mem_handle > &  col_buffers,
std::list< viennacl::backend::mem_handle > &  element_buffers,
std::list< vcl_size_t > &  row_elimination_num_list 
)

Definition at line 208 of file common.hpp.

template<typename ScalarType >
void viennacl::linalg::detail::level_scheduling_substitute ( vector< ScalarType > &  vec,
viennacl::backend::mem_handle const &  row_index_array,
viennacl::backend::mem_handle const &  row_buffer,
viennacl::backend::mem_handle const &  col_buffer,
viennacl::backend::mem_handle const &  element_buffer,
vcl_size_t  num_rows 
)

Definition at line 49 of file misc_operations.hpp.

template<typename NumericT >
void viennacl::linalg::detail::level_scheduling_substitute ( viennacl::vector< NumericT > &  vec,
std::list< viennacl::backend::mem_handle > const &  row_index_arrays,
std::list< viennacl::backend::mem_handle > const &  row_buffers,
std::list< viennacl::backend::mem_handle > const &  col_buffers,
std::list< viennacl::backend::mem_handle > const &  element_buffers,
std::list< vcl_size_t > const &  row_elimination_num_list 
)

Definition at line 224 of file common.hpp.

template<class T >
T viennacl::linalg::detail::max ( const T &  lhs,
const T &  rhs 
)

Maximum.

Examples:
bandwidth-reduction.cpp.

Definition at line 59 of file util.hpp.

template<typename IndexT , typename NumericT >
IndexT viennacl::linalg::detail::merge_subtract_sparse_rows ( IndexT const *  w_coords,
NumericT const *  w_elements,
IndexT  w_size,
IndexT const *  u_coords,
NumericT const *  u_elements,
IndexT  u_size,
NumericT  alpha,
IndexT *  z_coords,
NumericT z_elements 
)

Subtracts a scaled sparse vector u from a sparse vector w and writes the output to z: z = w - alpha * u.

Sparsity pattern of u and w are usually different.

Returns
Length of new vector

Definition at line 120 of file ilut.hpp.

template<class T >
T viennacl::linalg::detail::min ( const T &  lhs,
const T &  rhs 
)

Minimum.

Definition at line 45 of file util.hpp.

template<typename VectorType >
VectorType::value_type viennacl::linalg::detail::norm_lcl ( VectorType const &  x,
vcl_size_t  size 
)

Definition at line 78 of file qr-method-common.hpp.

template<typename VectorType >
void viennacl::linalg::detail::normalize ( VectorType &  x,
vcl_size_t  size 
)

Definition at line 87 of file qr-method-common.hpp.

template<typename NumericT , typename B >
bool viennacl::linalg::detail::op_aliasing ( vector_base< NumericT > const &  ,
B const &   
)

Definition at line 36 of file op_executor.hpp.

template<typename NumericT >
bool viennacl::linalg::detail::op_aliasing ( vector_base< NumericT > const &  lhs,
vector_base< NumericT > const &  b 
)

Definition at line 42 of file op_executor.hpp.

template<typename NumericT , typename LhsT , typename RhsT , typename OpT >
bool viennacl::linalg::detail::op_aliasing ( vector_base< NumericT > const &  lhs,
vector_expression< const LhsT, const RhsT, OpT > const &  rhs 
)

Definition at line 48 of file op_executor.hpp.

template<typename NumericT , typename B >
bool viennacl::linalg::detail::op_aliasing ( matrix_base< NumericT > const &  ,
B const &   
)

Definition at line 55 of file op_executor.hpp.

template<typename NumericT >
bool viennacl::linalg::detail::op_aliasing ( matrix_base< NumericT > const &  lhs,
matrix_base< NumericT > const &  b 
)

Definition at line 61 of file op_executor.hpp.

template<typename NumericT , typename LhsT , typename RhsT , typename OpT >
bool viennacl::linalg::detail::op_aliasing ( matrix_base< NumericT > const &  lhs,
matrix_expression< const LhsT, const RhsT, OpT > const &  rhs 
)

Definition at line 67 of file op_executor.hpp.

template<typename MatrixT , typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::pipelined_solve ( MatrixT const &  A,
viennacl::vector_base< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of a pipelined stabilized Bi-conjugate gradient solver.

Definition at line 98 of file bicgstab.hpp.

template<typename MatrixT , typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::pipelined_solve ( MatrixT const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of a pipelined conjugate gradient algorithm (no preconditioner), specialized for ViennaCL types.

Pipelined version from A. T. Chronopoulos and C. W. Gear, J. Comput. Appl. Math. 25(2), 153–168 (1989)

Parameters
AThe system matrix
rhsThe load vector
tagSolver configuration tag
monitorA callback routine which is called at each GMRES restart
monitor_dataData pointer to be passed to the callback routine to pass on user-specific data
Returns
The result vector

Definition at line 129 of file cg.hpp.

template<typename MatrixType , typename ScalarType >
viennacl::vector<ScalarType> viennacl::linalg::detail::pipelined_solve ( MatrixType const &  A,
viennacl::vector< ScalarType > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< ScalarType > const &, ScalarType, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of a pipelined GMRES solver without preconditioner.

Following algorithm 2.1 proposed by Walker in "A Simpler GMRES", but uses classical Gram-Schmidt instead of modified Gram-Schmidt for better parallelization. Uses some pipelining techniques for minimizing host-device transfer

Parameters
AThe system matrix
rhsThe load vector
tagSolver configuration tag
monitorA callback routine which is called at each GMRES restart
monitor_dataData pointer to be passed to the callback routine to pass on user-specific data
Returns
The result vector

Definition at line 182 of file gmres.hpp.

template<typename NumericT >
void viennacl::linalg::detail::precondition ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::compressed_matrix< NumericT > &  L,
viennacl::vector< NumericT > &  diag_L,
viennacl::compressed_matrix< NumericT > &  L_trans,
chow_patel_tag const &  tag 
)

Implementation of the parallel ICC0 factorization, Algorithm 3 in Chow-Patel paper.

Rather than dealing with a column-major upper triangular matrix U, we use the lower-triangular matrix L such that A is approximately given by LL^T. The advantage is that L is readily available in row-major format.

Definition at line 77 of file chow_patel_ilu.hpp.

template<typename NumericT >
void viennacl::linalg::detail::precondition ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::compressed_matrix< NumericT > &  L,
viennacl::vector< NumericT > &  diag_L,
viennacl::compressed_matrix< NumericT > &  U,
viennacl::vector< NumericT > &  diag_U,
chow_patel_tag const &  tag 
)

Implementation of the parallel ILU0 factorization, Algorithm 2 in Chow-Patel paper.

Definition at line 110 of file chow_patel_ilu.hpp.

template<typename SCALARTYPE >
void viennacl::linalg::detail::prepare_householder_vector ( matrix_base< SCALARTYPE > &  A,
vector_base< SCALARTYPE > &  D,
vcl_size_t  size,
vcl_size_t  row_start,
vcl_size_t  col_start,
vcl_size_t  start,
bool  is_column 
)

Definition at line 165 of file qr-method-common.hpp.

template<typename NumericT >
bool viennacl::linalg::detail::processResultDataLargeMatrix ( ResultDataLarge< NumericT > &  result,
const unsigned int  mat_size 
)

Process the result, that is obtain result from device and do simple sanity checking

Parameters
resulthandles to result data
mat_sizematrix size

Definition at line 88 of file bisect_large.hpp.

template<typename NumericT >
void viennacl::linalg::detail::processResultSmallMatrix ( ResultDataSmall< NumericT > &  result,
const unsigned int  mat_size 
)

Process the result obtained on the device, that is transfer to host and perform basic sanity checking

Parameters
resulthandles to result data
mat_sizematrix size

Definition at line 78 of file bisect_small.hpp.

template<typename SCALARTYPE >
SCALARTYPE viennacl::linalg::detail::pythag ( SCALARTYPE  a,
SCALARTYPE  b 
)

Definition at line 65 of file qr-method-common.hpp.

template<typename SCALARTYPE >
void viennacl::linalg::detail::qr_method ( viennacl::matrix< SCALARTYPE > &  A,
viennacl::matrix< SCALARTYPE > &  Q,
std::vector< SCALARTYPE > &  D,
std::vector< SCALARTYPE > &  E,
bool  is_symmetric = true 
)

Definition at line 728 of file qr-method.hpp.

template<typename SparseMatrixType , typename SCALARTYPE , unsigned int VEC_ALIGNMENT>
viennacl::enable_if< viennacl::is_any_sparse_matrix<SparseMatrixType>::value >::type viennacl::linalg::detail::row_info ( SparseMatrixType const &  mat,
vector< SCALARTYPE, VEC_ALIGNMENT > &  vec,
row_info_types  info_selector 
)

Definition at line 50 of file sparse_matrix_operations.hpp.

template<typename MatrixType , typename VectorType >
MatrixType::value_type viennacl::linalg::detail::setup_householder_vector_ublas ( MatrixType const &  A,
VectorType &  v,
MatrixType &  matrix_1x1,
vcl_size_t  j 
)

Definition at line 53 of file qr.hpp.

template<typename MatrixType , typename VectorType >
viennacl::result_of::cpu_value_type< typename MatrixType::value_type >::type viennacl::linalg::detail::setup_householder_vector_viennacl ( MatrixType const &  A,
VectorType &  v,
MatrixType &  matrix_1x1,
vcl_size_t  j 
)

Definition at line 93 of file qr.hpp.

template<typename SCALARTYPE >
SCALARTYPE viennacl::linalg::detail::sign ( SCALARTYPE  val)

Definition at line 71 of file qr-method-common.hpp.

double viennacl::linalg::detail::sign_d ( const double &  val)
inline

Sign of number (double)

Definition at line 84 of file util.hpp.

float viennacl::linalg::detail::sign_f ( const float &  val)
inline

Sign of number (float)

Definition at line 72 of file util.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 192 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::coordinate_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 205 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 219 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 219 of file bicgstab.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::coordinate_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 232 of file bicgstab.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::sliced_ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 233 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::hyb_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
cg_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 246 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 246 of file bicgstab.hpp.

template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT viennacl::linalg::detail::solve_impl ( MatrixT const &  matrix,
VectorT const &  rhs,
cg_tag const &  tag,
PreconditionerT const &  precond,
bool(*)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Definition at line 258 of file cg.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::sliced_ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 260 of file bicgstab.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::hyb_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 273 of file bicgstab.hpp.

template<typename MatrixT , typename VectorT >
VectorT viennacl::linalg::detail::solve_impl ( MatrixT const &  matrix,
VectorT const &  rhs,
bicgstab_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of the unpreconditioned stabilized Bi-conjugate gradient solver.

Following the description in "Iterative Methods for Sparse Linear Systems" by Y. Saad

Parameters
matrixThe system matrix
rhsThe load vector
tagSolver configuration tag
monitorA callback routine which is called at each GMRES restart
monitor_dataData pointer to be passed to the callback routine to pass on user-specific data
Returns
The result vector

Definition at line 296 of file bicgstab.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::compressed_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 371 of file gmres.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::coordinate_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 384 of file gmres.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 398 of file gmres.hpp.

template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT viennacl::linalg::detail::solve_impl ( MatrixT const &  matrix,
VectorT const &  rhs,
bicgstab_tag const &  tag,
PreconditionerT const &  precond,
bool(*)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of the preconditioned stabilized Bi-conjugate gradient solver.

Following the description of the unpreconditioned case in "Iterative Methods for Sparse Linear Systems" by Y. Saad

Parameters
matrixThe system matrix
rhsThe load vector
tagSolver configuration tag
precondA preconditioner. Precondition operation is done via member function apply()
monitorA callback routine which is called at each GMRES restart
monitor_dataData pointer to be passed to the callback routine to pass on user-specific data
Returns
The result vector

Definition at line 399 of file bicgstab.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::sliced_ell_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 412 of file gmres.hpp.

template<typename NumericT >
viennacl::vector<NumericT> viennacl::linalg::detail::solve_impl ( viennacl::hyb_matrix< NumericT > const &  A,
viennacl::vector< NumericT > const &  rhs,
gmres_tag const &  tag,
viennacl::linalg::no_precond  ,
bool(*)(viennacl::vector< NumericT > const &, NumericT, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.

Definition at line 425 of file gmres.hpp.

template<typename MatrixT , typename VectorT , typename PreconditionerT >
VectorT viennacl::linalg::detail::solve_impl ( MatrixT const &  matrix,
VectorT const &  rhs,
gmres_tag const &  tag,
PreconditionerT const &  precond,
bool(*)(VectorT const &, typename viennacl::result_of::cpu_value_type< typename viennacl::result_of::value_type< VectorT >::type >::type, void *)  monitor = NULL,
void *  monitor_data = NULL 
)

Implementation of the GMRES solver.

Following the algorithm proposed by Walker in "A Simpler GMRES"

Parameters
matrixThe system matrix
rhsThe load vector
tagSolver configuration tag
precondA preconditioner. Precondition operation is done via member function apply()
monitorA callback routine which is called at each GMRES restart
monitor_dataData pointer to be passed to the callback routine to pass on user-specific data
Returns
The result vector

Definition at line 450 of file gmres.hpp.

template<typename MatrixType , typename CPU_VectorType >
void viennacl::linalg::detail::svd_qr_shift ( MatrixType &  vcl_u,
MatrixType &  vcl_v,
CPU_VectorType &  q,
CPU_VectorType &  e 
)

Definition at line 101 of file svd.hpp.

template<typename SCALARTYPE >
void viennacl::linalg::detail::transpose ( matrix_base< SCALARTYPE > &  A)

Definition at line 107 of file qr-method-common.hpp.

template<typename SCALARTYPE >
void viennacl::linalg::detail::tridiagonal_reduction ( matrix_base< SCALARTYPE > &  A,
matrix_base< SCALARTYPE > &  Q 
)

Definition at line 713 of file qr-method.hpp.

template<typename SCALARTYPE , typename MatrixT >
void viennacl::linalg::detail::update_float_QR_column ( MatrixT &  A,
const std::vector< SCALARTYPE > &  buf,
int  m,
int  n,
int  last_i,
bool  is_triangular 
)

Definition at line 145 of file qr-method.hpp.

template<typename SCALARTYPE , typename VectorType >
void viennacl::linalg::detail::update_float_QR_column_gpu ( matrix_base< SCALARTYPE > &  A,
const VectorType &  buf,
viennacl::vector< SCALARTYPE > &  buf_vcl,
int  m,
int  n,
int  last_n,
bool   
)

Definition at line 82 of file qr-method.hpp.

template<typename MatrixType , typename VectorType >
void viennacl::linalg::detail::write_householder_to_A ( MatrixType &  A,
VectorType const &  v,
vcl_size_t  j 
)

Definition at line 194 of file qr.hpp.

template<typename MatrixType , typename VectorType >
void viennacl::linalg::detail::write_householder_to_A_ublas ( MatrixType &  A,
VectorType const &  v,
vcl_size_t  j 
)

Definition at line 201 of file qr.hpp.

template<typename MatrixType , typename VectorType >
void viennacl::linalg::detail::write_householder_to_A_viennacl ( MatrixType &  A,
VectorType const &  v,
vcl_size_t  j 
)

Definition at line 211 of file qr.hpp.