Implementation namespace for sparse approximate inverse preconditioner. More...
Classes | |
class | block_matrix |
Represents contigious matrices on GPU. More... | |
class | block_vector |
Represents a contiguous vector on the GPU to represent a concatentation of small vectors. More... | |
struct | CompareSecond |
Helper functor for comparing std::pair<> based on the second member. More... | |
class | fspai_tag |
A tag for FSPAI. Experimental. More... | |
class | spai_tag |
A tag for SPAI. More... | |
class | sparse_vector |
Represents a sparse vector based on std::map<unsigned int, NumericT> More... | |
Functions | |
template<typename MatrixT , typename NumericT > | |
void | sym_sparse_matrix_to_stl (MatrixT const &A, std::vector< std::map< unsigned int, NumericT > > &STL_A) |
template<typename MatrixT > | |
void | generateJ (MatrixT const &A, std::vector< std::vector< vcl_size_t > > &J) |
template<typename NumericT , typename MatrixT , typename VectorT > | |
void | fill_blocks (std::vector< std::map< unsigned int, NumericT > > &A, std::vector< MatrixT > &blocks, std::vector< std::vector< vcl_size_t > > const &J, std::vector< VectorT > &Y) |
template<typename MatrixT > | |
void | cholesky_decompose (MatrixT &A) |
template<typename MatrixT , typename VectorT > | |
void | cholesky_solve (MatrixT const &L, VectorT &b) |
template<typename MatrixT , typename VectorT > | |
void | computeL (MatrixT const &A, MatrixT &L, MatrixT &L_trans, std::vector< VectorT > &Y, std::vector< std::vector< vcl_size_t > > &J) |
template<typename MatrixT > | |
void | computeFSPAI (MatrixT const &A, MatrixT const &PatternA, MatrixT &L, MatrixT &L_trans, fspai_tag) |
template<typename T , typename InputIteratorT > | |
void | Print (std::ostream &ostr, InputIteratorT it_begin, InputIteratorT it_end) |
template<typename VectorT , typename MatrixT > | |
void | write_to_block (VectorT &con_A_I_J, unsigned int start_ind, std::vector< unsigned int > const &I, std::vector< unsigned int > const &J, MatrixT &m) |
template<typename VectorT > | |
void | print_continious_matrix (VectorT &con_A_I_J, std::vector< cl_uint > &blocks_ind, std::vector< std::vector< unsigned int > > const &g_I, std::vector< std::vector< unsigned int > > const &g_J) |
template<typename VectorT > | |
void | print_continious_vector (VectorT &con_v, std::vector< cl_uint > &block_ind, std::vector< std::vector< unsigned int > > const &g_J) |
void | compute_blocks_size (std::vector< std::vector< unsigned int > > const &g_I, std::vector< std::vector< unsigned int > > const &g_J, unsigned int &sz, std::vector< cl_uint > &blocks_ind, std::vector< cl_uint > &matrix_dims) |
**************************************** BLOCK FUNCTIONS ************************************// More... | |
template<typename SizeT > | |
void | get_size (std::vector< std::vector< SizeT > > const &inds, SizeT &size) |
Computes size of particular container of index set. More... | |
template<typename SizeT > | |
void | init_start_inds (std::vector< std::vector< SizeT > > const &inds, std::vector< cl_uint > &start_inds) |
Initializes start indices of particular index set. More... | |
template<typename MatrixT , typename NumericT > | |
void | dot_prod (MatrixT const &A, unsigned int beg_ind, NumericT &res) |
Dot prod of particular column of martix A with it's self starting at a certain index beg_ind. More... | |
template<typename MatrixT , typename VectorT , typename NumericT > | |
void | custom_inner_prod (MatrixT const &A, VectorT const &v, unsigned int col_ind, unsigned int start_ind, NumericT &res) |
Dot prod of particular matrix column with arbitrary vector: A(:, col_ind) More... | |
template<typename MatrixT , typename VectorT > | |
void | copy_vector (MatrixT const &A, VectorT &v, unsigned int beg_ind) |
Copying part of matrix column. More... | |
template<typename MatrixT , typename VectorT , typename NumericT > | |
void | householder_vector (MatrixT const &A, unsigned int j, VectorT &v, NumericT &b) |
Computation of Householder vector, householder reflection c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.210. More... | |
template<typename MatrixT , typename VectorT , typename NumericT > | |
void | apply_householder_reflection (MatrixT &A, unsigned int iter_cnt, VectorT &v, NumericT b) |
Inplace application of Householder vector to a matrix A. More... | |
template<typename MatrixT , typename VectorT > | |
void | store_householder_vector (MatrixT &A, unsigned int ind, VectorT &v) |
Storage of vector v in column(A, ind), starting from ind-1 index of a column. More... | |
template<typename MatrixT , typename VectorT > | |
void | single_qr (MatrixT &R, VectorT &b_v) |
Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224. More... | |
template<typename SizeT > | |
void | get_max_block_size (std::vector< std::vector< SizeT > > const &inds, SizeT &max_size) |
Getting max size of rows/columns from container of index set. More... | |
template<typename MatrixT , typename VectorT , typename NumericT > | |
void | custom_dot_prod (MatrixT const &A, VectorT const &v, unsigned int ind, NumericT &res) |
Dot_prod(column(A, ind), v) starting from index ind+1. More... | |
template<typename MatrixT , typename VectorT > | |
void | apply_q_trans_vec (MatrixT const &R, VectorT const &b_v, VectorT &y) |
Recovery Q from matrix R and vector of betas b_v. More... | |
template<typename MatrixT , typename VectorT > | |
void | apply_q_trans_mat (MatrixT const &R, VectorT const &b_v, MatrixT &A) |
Multiplication of Q'*A, where Q is in implicit for lower part of R and vector of betas - b_v. More... | |
template<typename NumericT > | |
void | block_qr (std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, block_matrix &g_A_I_J_vcl, block_vector &g_bv_vcl, std::vector< cl_uint > &g_is_update, viennacl::context ctx) |
Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224 performed on GPU. More... | |
template<typename MatrixT > | |
void | make_rotation_matrix (MatrixT &mat, vcl_size_t new_size, vcl_size_t off_diagonal_distance=4) |
template<typename MatrixT > | |
double | determinant (boost::numeric::ublas::matrix_expression< MatrixT > const &mat_r) |
template<typename MatrixT > | |
void | composeNewR (MatrixT const &A, MatrixT const &R_n, MatrixT &R) |
Composition of new matrix R, that is going to be used in Least Square problem solving. More... | |
template<typename VectorT > | |
void | composeNewVector (VectorT const &v_n, VectorT &v) |
Composition of new vector of coefficients beta from QR factorizations(necessary for Q recovery) More... | |
template<typename SparseVectorT , typename NumericT > | |
void | sparse_norm_2 (SparseVectorT const &v, NumericT &norm) |
Computation of Euclidean norm for sparse vector. More... | |
template<typename SparseVectorT , typename NumericT > | |
void | sparse_inner_prod (SparseVectorT const &v1, SparseVectorT const &v2, NumericT &res_v) |
Dot product of two sparse vectors. More... | |
template<typename SparseVectorT , typename NumericT > | |
bool | buildAugmentedIndexSet (std::vector< SparseVectorT > const &A_v_c, SparseVectorT const &res, std::vector< unsigned int > &J, std::vector< unsigned int > &J_u, spai_tag const &tag) |
Building a new set of column indices J_u, cf. Kallischko dissertation p.31. More... | |
template<typename SparseVectorT > | |
void | buildNewRowSet (std::vector< SparseVectorT > const &A_v_c, std::vector< unsigned int > const &I, std::vector< unsigned int > const &J_n, std::vector< unsigned int > &I_n) |
Building a new indices to current set of row indices I_n, cf. Kallischko dissertation p.32. More... | |
template<typename MatrixT > | |
void | QRBlockComposition (MatrixT const &A_I_J, MatrixT const &A_I_J_u, MatrixT &A_I_u_J_u) |
Composition of new block for QR factorization cf. Kallischko dissertation p.82, figure 4.7. More... | |
template<typename SparseMatrixT , typename SparseVectorT , typename DenseMatrixT , typename VectorT > | |
void | block_update (SparseMatrixT const &A, std::vector< SparseVectorT > const &A_v_c, std::vector< SparseVectorT > &g_res, std::vector< bool > &g_is_update, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, std::vector< VectorT > &g_b_v, std::vector< DenseMatrixT > &g_A_I_J, spai_tag const &tag) |
CPU-based dynamic update for SPAI preconditioner. More... | |
template<typename NumericT > | |
void | block_q_multiplication (std::vector< std::vector< unsigned int > > const &g_J_u, std::vector< std::vector< unsigned int > > const &g_I, block_matrix &g_A_I_J_vcl, block_vector &g_bv_vcl, block_matrix &g_A_I_J_u_vcl, std::vector< cl_uint > &g_is_update, viennacl::context ctx) |
Performs multiplication Q'*A(I, \tilde J) on GPU. More... | |
template<typename SizeT > | |
void | assemble_qr_row_inds (std::vector< std::vector< SizeT > > const &g_I, std::vector< std::vector< SizeT > > const &g_J, std::vector< std::vector< SizeT > > const &g_I_u, std::vector< std::vector< SizeT > > &g_I_q) |
Assembly of container of index row sets: I_q, row indices for new "QR block". More... | |
template<typename NumericT > | |
void | assemble_qr_block (std::vector< std::vector< unsigned int > > const &g_J, std::vector< std::vector< unsigned int > > const &g_I, std::vector< std::vector< unsigned int > > const &g_J_u, std::vector< std::vector< unsigned int > > const &g_I_u, std::vector< std::vector< unsigned int > > &g_I_q, block_matrix &g_A_I_J_u_vcl, viennacl::ocl::handle< cl_mem > &matrix_dimensions, block_matrix &g_A_I_u_J_u_vcl, std::vector< cl_uint > &g_is_update, bool is_empty_block, viennacl::context ctx) |
Performs assembly for new QR block. More... | |
template<typename NumericT > | |
void | assemble_r (std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, block_matrix &g_A_I_J_vcl, block_matrix &g_A_I_J_u_vcl, block_matrix &g_A_I_u_J_u_vcl, block_vector &g_bv_vcl, block_vector &g_bv_vcl_u, std::vector< cl_uint > &g_is_update, viennacl::context ctx) |
Performs assembly for new R matrix on GPU. More... | |
template<typename NumericT , unsigned int AlignmentV, typename SparseVectorT > | |
void | block_update (viennacl::compressed_matrix< NumericT, AlignmentV > const &A, std::vector< SparseVectorT > const &A_v_c, std::vector< cl_uint > &g_is_update, std::vector< SparseVectorT > &g_res, std::vector< std::vector< unsigned int > > &g_J, std::vector< std::vector< unsigned int > > &g_I, block_matrix &g_A_I_J_vcl, block_vector &g_bv_vcl, spai_tag const &tag) |
GPU-based block update. More... | |
template<typename SizeT > | |
bool | isInIndexSet (std::vector< SizeT > const &J, SizeT ind) |
Determines if element ind is in set {J}. More... | |
template<typename VectorT , typename SparseVectorT > | |
void | fanOutVector (VectorT const &m_in, std::vector< unsigned int > const &J, SparseVectorT &m) |
Projects solution of LS problem onto original column m. More... | |
template<typename MatrixT , typename VectorT > | |
void | backwardSolve (MatrixT const &R, VectorT const &y, VectorT &x) |
Solution of linear:R*x=y system by backward substitution. More... | |
template<typename VectorT , typename NumericT > | |
void | projectI (std::vector< unsigned int > const &I, VectorT &y, unsigned int ind) |
Perform projection of set I on the unit-vector. More... | |
template<typename SparseVectorT > | |
void | buildColumnIndexSet (SparseVectorT const &v, std::vector< unsigned int > &J) |
Builds index set of projected columns for current column of preconditioner. More... | |
template<typename SparseMatrixT > | |
void | initPreconditioner (SparseMatrixT const &A, SparseMatrixT &M) |
Initialize preconditioner with sparcity pattern = p(A) More... | |
template<typename SparseVectorT > | |
void | projectRows (std::vector< SparseVectorT > const &A_v_c, std::vector< unsigned int > const &J, std::vector< unsigned int > &I) |
Row projection for matrix A(:,J) -> A(I,J), building index set of non-zero rows. More... | |
template<typename SparseVectorT > | |
void | print_sparse_vector (SparseVectorT const &v) |
template<typename DenseMatrixT > | |
void | print_matrix (DenseMatrixT &m) |
template<typename SparseVectorT , typename NumericT > | |
void | add_sparse_vectors (SparseVectorT const &v, NumericT b, SparseVectorT &res_v) |
Add two sparse vectors res_v = b*v. More... | |
template<typename SparseVectorT , typename NumericT > | |
void | compute_spai_residual (std::vector< SparseVectorT > const &A_v_c, SparseVectorT const &v, unsigned int ind, SparseVectorT &res) |
Computation of residual res = A*v - e. More... | |
template<typename SparseVectorT > | |
void | build_index_set (std::vector< SparseVectorT > const &A_v_c, SparseVectorT const &v, std::vector< unsigned int > &J, std::vector< unsigned int > &I) |
Setting up index set of columns and rows for certain column. More... | |
template<typename SparseMatrixT , typename DenseMatrixT > | |
void | initProjectSubMatrix (SparseMatrixT const &A_in, std::vector< unsigned int > const &J, std::vector< unsigned int > &I, DenseMatrixT &A_out) |
Initializes a dense matrix from a sparse one. More... | |
template<typename SparseMatrixT , typename DenseMatrixT , typename SparseVectorT , typename VectorT > | |
void | block_set_up (SparseMatrixT const &A, std::vector< SparseVectorT > const &A_v_c, std::vector< SparseVectorT > const &M_v, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, std::vector< DenseMatrixT > &g_A_I_J, std::vector< VectorT > &g_b_v) |
Setting up blocks and QR factorizing them on CPU. More... | |
template<typename SparseVectorT > | |
void | index_set_up (std::vector< SparseVectorT > const &A_v_c, std::vector< SparseVectorT > const &M_v, std::vector< std::vector< unsigned int > > &g_J, std::vector< std::vector< unsigned int > > &g_I) |
Setting up index set of columns and rows for all columns. More... | |
template<typename NumericT , unsigned int AlignmentV, typename SparseVectorT > | |
void | block_set_up (viennacl::compressed_matrix< NumericT, AlignmentV > const &A, std::vector< SparseVectorT > const &A_v_c, std::vector< SparseVectorT > const &M_v, std::vector< cl_uint > g_is_update, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, block_matrix &g_A_I_J, block_vector &g_bv) |
Setting up blocks and QR factorizing them on GPU. More... | |
template<typename NumericT , typename SparseVectorT > | |
void | custom_fan_out (std::vector< NumericT > const &m_in, unsigned int start_m_ind, std::vector< unsigned int > const &J, SparseVectorT &m) |
Elicitation of sparse vector m for particular column from m_in - contigious vector for all columns. More... | |
template<typename SparseVectorT , typename NumericT > | |
void | least_square_solve (std::vector< SparseVectorT > &A_v_c, std::vector< SparseVectorT > &M_v, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, block_matrix &g_A_I_J_vcl, block_vector &g_bv_vcl, std::vector< SparseVectorT > &g_res, std::vector< cl_uint > &g_is_update, const spai_tag &tag, viennacl::context ctx) |
Solution of Least square problem on GPU. More... | |
template<typename SparseVectorT , typename DenseMatrixT , typename VectorT > | |
void | least_square_solve (std::vector< SparseVectorT > const &A_v_c, std::vector< DenseMatrixT > &g_R, std::vector< VectorT > &g_b_v, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, std::vector< SparseVectorT > &g_res, std::vector< bool > &g_is_update, std::vector< SparseVectorT > &M_v, spai_tag const &tag) |
Solution of Least square problem on CPU. More... | |
template<typename VectorType > | |
bool | is_all_update (VectorType ¶llel_is_update) |
template<typename SparseMatrixT , typename SparseVectorT > | |
void | vectorize_column_matrix (SparseMatrixT const &M_in, std::vector< SparseVectorT > &M_v) |
Solution of Least square problem on CPU. More... | |
template<typename SparseMatrixT , typename SparseVectorT > | |
void | vectorize_row_matrix (SparseMatrixT const &M_in, std::vector< SparseVectorT > &M_v) |
template<typename SizeT > | |
void | write_set_to_array (std::vector< std::vector< SizeT > > const &ind_set, std::vector< cl_uint > &a) |
template<typename NumericT , unsigned int AlignmentV> | |
void | block_assembly (viennacl::compressed_matrix< NumericT, AlignmentV > const &A, std::vector< std::vector< unsigned int > > const &g_J, std::vector< std::vector< unsigned int > > const &g_I, block_matrix &g_A_I_J_vcl, std::vector< cl_uint > &g_is_update, bool &is_empty_block) |
Assembly of blocks on GPU by a gived set of row indices: g_I and column indices: g_J. More... | |
template<typename SparseMatrixT , typename SparseVectorT > | |
void | insert_sparse_columns (std::vector< SparseVectorT > const &M_v, SparseMatrixT &M, bool is_right) |
Insertion of vectorized matrix column into original sparse matrix. More... | |
template<typename MatrixT > | |
void | sparse_transpose (MatrixT const &A_in, MatrixT &A) |
Transposition of sparse matrix. More... | |
template<typename MatrixT > | |
void | computeSPAI (MatrixT const &A, MatrixT &M, spai_tag &tag) |
Construction of SPAI preconditioner on CPU. More... | |
template<typename NumericT , unsigned int AlignmentV> | |
void | computeSPAI (viennacl::compressed_matrix< NumericT, AlignmentV > const &A, boost::numeric::ublas::compressed_matrix< NumericT > const &cpu_A, boost::numeric::ublas::compressed_matrix< NumericT > &cpu_M, viennacl::compressed_matrix< NumericT, AlignmentV > &M, spai_tag const &tag) |
Construction of SPAI preconditioner on GPU. More... | |
Implementation namespace for sparse approximate inverse preconditioner.
void viennacl::linalg::detail::spai::add_sparse_vectors | ( | SparseVectorT const & | v, |
NumericT | b, | ||
SparseVectorT & | res_v | ||
) |
void viennacl::linalg::detail::spai::apply_householder_reflection | ( | MatrixT & | A, |
unsigned int | iter_cnt, | ||
VectorT & | v, | ||
NumericT | b | ||
) |
void viennacl::linalg::detail::spai::apply_q_trans_mat | ( | MatrixT const & | R, |
VectorT const & | b_v, | ||
MatrixT & | A | ||
) |
void viennacl::linalg::detail::spai::apply_q_trans_vec | ( | MatrixT const & | R, |
VectorT const & | b_v, | ||
VectorT & | y | ||
) |
void viennacl::linalg::detail::spai::assemble_qr_block | ( | std::vector< std::vector< unsigned int > > const & | g_J, |
std::vector< std::vector< unsigned int > > const & | g_I, | ||
std::vector< std::vector< unsigned int > > const & | g_J_u, | ||
std::vector< std::vector< unsigned int > > const & | g_I_u, | ||
std::vector< std::vector< unsigned int > > & | g_I_q, | ||
block_matrix & | g_A_I_J_u_vcl, | ||
viennacl::ocl::handle< cl_mem > & | matrix_dimensions, | ||
block_matrix & | g_A_I_u_J_u_vcl, | ||
std::vector< cl_uint > & | g_is_update, | ||
bool | is_empty_block, | ||
viennacl::context | ctx | ||
) |
Performs assembly for new QR block.
g_J | container of column indices |
g_I | container of row indices |
g_J_u | container of new column indices |
g_I_u | container of new row indices |
g_I_q | container of row indices for new QR blocks |
g_A_I_J_u_vcl | blocks of Q'*A(I, \tilde J) |
matrix_dimensions | array with matrix dimensions for all blocks |
g_A_I_u_J_u_vcl | blocks A(\tilde I, \tilde J) |
g_is_update | container with update indicators |
is_empty_block | indicator if all previous blocks A(\tilde I, \tilde J) - are empty, in case if they are empty kernel with smaller number of arguments is used |
ctx | Optional context in which the matrix is created (one out of multiple OpenCL contexts, CUDA, host) |
Definition at line 439 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::assemble_qr_row_inds | ( | std::vector< std::vector< SizeT > > const & | g_I, |
std::vector< std::vector< SizeT > > const & | g_J, | ||
std::vector< std::vector< SizeT > > const & | g_I_u, | ||
std::vector< std::vector< SizeT > > & | g_I_q | ||
) |
Assembly of container of index row sets: I_q, row indices for new "QR block".
g_I | container of row indices |
g_J | container of column indices |
g_I_u | container of new row indices |
g_I_q | container of row indices for new QR blocks |
Definition at line 406 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::assemble_r | ( | std::vector< std::vector< unsigned int > > & | g_I, |
std::vector< std::vector< unsigned int > > & | g_J, | ||
block_matrix & | g_A_I_J_vcl, | ||
block_matrix & | g_A_I_J_u_vcl, | ||
block_matrix & | g_A_I_u_J_u_vcl, | ||
block_vector & | g_bv_vcl, | ||
block_vector & | g_bv_vcl_u, | ||
std::vector< cl_uint > & | g_is_update, | ||
viennacl::context | ctx | ||
) |
Performs assembly for new R matrix on GPU.
g_I | container of row indices |
g_J | container of column indices |
g_A_I_J_vcl | container of block matrices from previous update |
g_A_I_J_u_vcl | container of block matrices Q'*A(I, \tilde J) |
g_A_I_u_J_u_vcl | container of block matrices QR factored on current iteration |
g_bv_vcl | block of beta vectors from previous iteration |
g_bv_vcl_u | block of updated beta vectors got after recent QR factorization |
g_is_update | container with identificators that shows which block should be modified |
ctx | Optional context in which the auxiliary data is created (one out of multiple OpenCL contexts, CUDA, host) |
Definition at line 533 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::backwardSolve | ( | MatrixT const & | R, |
VectorT const & | y, | ||
VectorT & | x | ||
) |
Solution of linear:R*x=y system by backward substitution.
R | uppertriangular matrix |
y | right handside vector |
x | solution vector |
Definition at line 102 of file spai-static.hpp.
void viennacl::linalg::detail::spai::block_assembly | ( | viennacl::compressed_matrix< NumericT, AlignmentV > const & | A, |
std::vector< std::vector< unsigned int > > const & | g_J, | ||
std::vector< std::vector< unsigned int > > const & | g_I, | ||
block_matrix & | g_A_I_J_vcl, | ||
std::vector< cl_uint > & | g_is_update, | ||
bool & | is_empty_block | ||
) |
Assembly of blocks on GPU by a gived set of row indices: g_I and column indices: g_J.
A | intial sparse matrix |
g_J | container of column index set |
g_I | container of row index set |
g_A_I_J_vcl | contigious blocks A(I, J) using GPU memory |
g_is_update | container with indicators which blocks are active |
is_empty_block | parameter that indicates if no block were assembled |
void viennacl::linalg::detail::spai::block_q_multiplication | ( | std::vector< std::vector< unsigned int > > const & | g_J_u, |
std::vector< std::vector< unsigned int > > const & | g_I, | ||
block_matrix & | g_A_I_J_vcl, | ||
block_vector & | g_bv_vcl, | ||
block_matrix & | g_A_I_J_u_vcl, | ||
std::vector< cl_uint > & | g_is_update, | ||
viennacl::context | ctx | ||
) |
Performs multiplication Q'*A(I, \tilde J) on GPU.
g_J_u | container of sets of new column indices |
g_I | container of row indices |
g_A_I_J_vcl | block matrix composed from previous blocks, they are blocks of R |
g_bv_vcl | block of beta vectors |
g_A_I_J_u_vcl | block of matrices A(I, \tilde J) |
g_is_update | indicators, that show if a certain block should be processed |
ctx | Optional context in which the auxiliary data is created (one out of multiple OpenCL contexts, CUDA, host) |
Definition at line 361 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::block_qr | ( | std::vector< std::vector< unsigned int > > & | g_I, |
std::vector< std::vector< unsigned int > > & | g_J, | ||
block_matrix & | g_A_I_J_vcl, | ||
block_vector & | g_bv_vcl, | ||
std::vector< cl_uint > & | g_is_update, | ||
viennacl::context | ctx | ||
) |
Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224 performed on GPU.
g_I | container of row indices |
g_J | container of column indices |
g_A_I_J_vcl | contigious matrices, GPU memory is used |
g_bv_vcl | contigiuos vectors beta, GPU memory is used |
g_is_update | container of indicators that show active blocks |
ctx | Optional context in which the auxiliary data is created (one out of multiple OpenCL contexts, CUDA, host) |
void viennacl::linalg::detail::spai::block_set_up | ( | SparseMatrixT const & | A, |
std::vector< SparseVectorT > const & | A_v_c, | ||
std::vector< SparseVectorT > const & | M_v, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
std::vector< DenseMatrixT > & | g_A_I_J, | ||
std::vector< VectorT > & | g_b_v | ||
) |
Setting up blocks and QR factorizing them on CPU.
A | initial sparse matrix |
A_v_c | column major vectorized initial sparse matrix |
M_v | initialized preconditioner |
g_I | container of row indices |
g_J | container of column indices |
g_A_I_J | container of dense matrices -> R matrices after QR factorization |
g_b_v | container of vectors beta, necessary for Q recovery |
void viennacl::linalg::detail::spai::block_set_up | ( | viennacl::compressed_matrix< NumericT, AlignmentV > const & | A, |
std::vector< SparseVectorT > const & | A_v_c, | ||
std::vector< SparseVectorT > const & | M_v, | ||
std::vector< cl_uint > | g_is_update, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
block_matrix & | g_A_I_J, | ||
block_vector & | g_bv | ||
) |
Setting up blocks and QR factorizing them on GPU.
A | initial sparse matrix |
A_v_c | column major vectorized initial sparse matrix |
M_v | initialized preconditioner |
g_is_update | container that indicates which blocks are active |
g_I | container of row indices |
g_J | container of column indices |
g_A_I_J | container of dense matrices -> R matrices after QR factorization |
g_bv | container of vectors beta, necessary for Q recovery |
void viennacl::linalg::detail::spai::block_update | ( | SparseMatrixT const & | A, |
std::vector< SparseVectorT > const & | A_v_c, | ||
std::vector< SparseVectorT > & | g_res, | ||
std::vector< bool > & | g_is_update, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
std::vector< VectorT > & | g_b_v, | ||
std::vector< DenseMatrixT > & | g_A_I_J, | ||
spai_tag const & | tag | ||
) |
CPU-based dynamic update for SPAI preconditioner.
A | initial sparse matrix |
A_v_c | vectorized column-wise initial matrix |
g_res | container of residuals for all columns |
g_is_update | container with identificators that shows which block should be modified |
g_I | container of row index sets for all columns |
g_J | container of column index sets for all columns |
g_b_v | container of vectors of beta for Q recovery(cf. Golub Van Loan "Matrix Computations", 3rd edition p.211) |
g_A_I_J | container of block matrices from previous update |
tag | SPAI configuration tag |
Definition at line 291 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::block_update | ( | viennacl::compressed_matrix< NumericT, AlignmentV > const & | A, |
std::vector< SparseVectorT > const & | A_v_c, | ||
std::vector< cl_uint > & | g_is_update, | ||
std::vector< SparseVectorT > & | g_res, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
block_matrix & | g_A_I_J_vcl, | ||
block_vector & | g_bv_vcl, | ||
spai_tag const & | tag | ||
) |
GPU-based block update.
A | sparse matrix |
A_v_c | vectorized column-wise initial matrix |
g_is_update | container with identificators that shows which block should be modified |
g_res | container of residuals for all columns |
g_J | container of column index sets for all columns |
g_I | container of row index sets for all columns |
g_A_I_J_vcl | container of block matrices from previous update |
g_bv_vcl | block of beta vectors from previous iteration |
tag | SPAI configuration tag |
Definition at line 625 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::build_index_set | ( | std::vector< SparseVectorT > const & | A_v_c, |
SparseVectorT const & | v, | ||
std::vector< unsigned int > & | J, | ||
std::vector< unsigned int > & | I | ||
) |
bool viennacl::linalg::detail::spai::buildAugmentedIndexSet | ( | std::vector< SparseVectorT > const & | A_v_c, |
SparseVectorT const & | res, | ||
std::vector< unsigned int > & | J, | ||
std::vector< unsigned int > & | J_u, | ||
spai_tag const & | tag | ||
) |
Building a new set of column indices J_u, cf. Kallischko dissertation p.31.
A_v_c | vectorized column-wise initial matrix |
res | residual vector |
J | set of column indices |
J_u | set of new column indices |
tag | SPAI tag with parameters |
Definition at line 188 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::buildColumnIndexSet | ( | SparseVectorT const & | v, |
std::vector< unsigned int > & | J | ||
) |
Builds index set of projected columns for current column of preconditioner.
v | current column of preconditioner |
J | output - index set of non-zero columns |
Definition at line 140 of file spai-static.hpp.
void viennacl::linalg::detail::spai::buildNewRowSet | ( | std::vector< SparseVectorT > const & | A_v_c, |
std::vector< unsigned int > const & | I, | ||
std::vector< unsigned int > const & | J_n, | ||
std::vector< unsigned int > & | I_n | ||
) |
Building a new indices to current set of row indices I_n, cf. Kallischko dissertation p.32.
A_v_c | vectorized column-wise initial matrix |
I | set of previous determined row indices |
J_n | set of new column indices |
I_n | set of new indices |
Definition at line 228 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::cholesky_decompose | ( | MatrixT & | A | ) |
void viennacl::linalg::detail::spai::cholesky_solve | ( | MatrixT const & | L, |
VectorT & | b | ||
) |
void viennacl::linalg::detail::spai::composeNewR | ( | MatrixT const & | A, |
MatrixT const & | R_n, | ||
MatrixT & | R | ||
) |
Composition of new matrix R, that is going to be used in Least Square problem solving.
A | matrix Q'*A(I, \tilde J), where \tilde J - set of new column indices |
R_n | matrix A_Iu_J_u after QR factorization |
R | previously composed matrix R |
Definition at line 94 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::composeNewVector | ( | VectorT const & | v_n, |
VectorT & | v | ||
) |
Composition of new vector of coefficients beta from QR factorizations(necessary for Q recovery)
v_n | new vector from last QR factorization |
v | composition of previous vectors from QR factorizations |
Definition at line 124 of file spai-dynamic.hpp.
|
inline |
**************************************** BLOCK FUNCTIONS ************************************//
Computes size of elements, start indices and matrix dimensions for a certain block
g_I | container of row indices |
g_J | container of column indices |
sz | general size for all elements in a certain block |
blocks_ind | start indices in a certain |
matrix_dims | matrix dimensions for each block |
void viennacl::linalg::detail::spai::compute_spai_residual | ( | std::vector< SparseVectorT > const & | A_v_c, |
SparseVectorT const & | v, | ||
unsigned int | ind, | ||
SparseVectorT & | res | ||
) |
void viennacl::linalg::detail::spai::computeFSPAI | ( | MatrixT const & | A, |
MatrixT const & | PatternA, | ||
MatrixT & | L, | ||
MatrixT & | L_trans, | ||
fspai_tag | |||
) |
void viennacl::linalg::detail::spai::computeL | ( | MatrixT const & | A, |
MatrixT & | L, | ||
MatrixT & | L_trans, | ||
std::vector< VectorT > & | Y, | ||
std::vector< std::vector< vcl_size_t > > & | J | ||
) |
void viennacl::linalg::detail::spai::computeSPAI | ( | MatrixT const & | A, |
MatrixT & | M, | ||
spai_tag & | tag | ||
) |
void viennacl::linalg::detail::spai::computeSPAI | ( | viennacl::compressed_matrix< NumericT, AlignmentV > const & | A, |
boost::numeric::ublas::compressed_matrix< NumericT > const & | cpu_A, | ||
boost::numeric::ublas::compressed_matrix< NumericT > & | cpu_M, | ||
viennacl::compressed_matrix< NumericT, AlignmentV > & | M, | ||
spai_tag const & | tag | ||
) |
void viennacl::linalg::detail::spai::copy_vector | ( | MatrixT const & | A, |
VectorT & | v, | ||
unsigned int | beg_ind | ||
) |
void viennacl::linalg::detail::spai::custom_dot_prod | ( | MatrixT const & | A, |
VectorT const & | v, | ||
unsigned int | ind, | ||
NumericT & | res | ||
) |
void viennacl::linalg::detail::spai::custom_fan_out | ( | std::vector< NumericT > const & | m_in, |
unsigned int | start_m_ind, | ||
std::vector< unsigned int > const & | J, | ||
SparseVectorT & | m | ||
) |
void viennacl::linalg::detail::spai::custom_inner_prod | ( | MatrixT const & | A, |
VectorT const & | v, | ||
unsigned int | col_ind, | ||
unsigned int | start_ind, | ||
NumericT & | res | ||
) |
double viennacl::linalg::detail::spai::determinant | ( | boost::numeric::ublas::matrix_expression< MatrixT > const & | mat_r | ) |
Definition at line 84 of file small_matrix.hpp.
void viennacl::linalg::detail::spai::dot_prod | ( | MatrixT const & | A, |
unsigned int | beg_ind, | ||
NumericT & | res | ||
) |
void viennacl::linalg::detail::spai::fanOutVector | ( | VectorT const & | m_in, |
std::vector< unsigned int > const & | J, | ||
SparseVectorT & | m | ||
) |
Projects solution of LS problem onto original column m.
m_in | solution of LS |
J | set of non-zero columns |
m | original column of M |
Definition at line 88 of file spai-static.hpp.
void viennacl::linalg::detail::spai::fill_blocks | ( | std::vector< std::map< unsigned int, NumericT > > & | A, |
std::vector< MatrixT > & | blocks, | ||
std::vector< std::vector< vcl_size_t > > const & | J, | ||
std::vector< VectorT > & | Y | ||
) |
void viennacl::linalg::detail::spai::generateJ | ( | MatrixT const & | A, |
std::vector< std::vector< vcl_size_t > > & | J | ||
) |
void viennacl::linalg::detail::spai::get_max_block_size | ( | std::vector< std::vector< SizeT > > const & | inds, |
SizeT & | max_size | ||
) |
void viennacl::linalg::detail::spai::get_size | ( | std::vector< std::vector< SizeT > > const & | inds, |
SizeT & | size | ||
) |
void viennacl::linalg::detail::spai::householder_vector | ( | MatrixT const & | A, |
unsigned int | j, | ||
VectorT & | v, | ||
NumericT & | b | ||
) |
void viennacl::linalg::detail::spai::index_set_up | ( | std::vector< SparseVectorT > const & | A_v_c, |
std::vector< SparseVectorT > const & | M_v, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
std::vector< std::vector< unsigned int > > & | g_I | ||
) |
void viennacl::linalg::detail::spai::init_start_inds | ( | std::vector< std::vector< SizeT > > const & | inds, |
std::vector< cl_uint > & | start_inds | ||
) |
void viennacl::linalg::detail::spai::initPreconditioner | ( | SparseMatrixT const & | A, |
SparseMatrixT & | M | ||
) |
Initialize preconditioner with sparcity pattern = p(A)
A | input matrix |
M | output matrix - initialized preconditioner |
Definition at line 154 of file spai-static.hpp.
void viennacl::linalg::detail::spai::initProjectSubMatrix | ( | SparseMatrixT const & | A_in, |
std::vector< unsigned int > const & | J, | ||
std::vector< unsigned int > & | I, | ||
DenseMatrixT & | A_out | ||
) |
void viennacl::linalg::detail::spai::insert_sparse_columns | ( | std::vector< SparseVectorT > const & | M_v, |
SparseMatrixT & | M, | ||
bool | is_right | ||
) |
bool viennacl::linalg::detail::spai::is_all_update | ( | VectorType & | parallel_is_update | ) |
bool viennacl::linalg::detail::spai::isInIndexSet | ( | std::vector< SizeT > const & | J, |
SizeT | ind | ||
) |
Determines if element ind is in set {J}.
J | current set |
ind | current element |
Definition at line 72 of file spai-static.hpp.
void viennacl::linalg::detail::spai::least_square_solve | ( | std::vector< SparseVectorT > & | A_v_c, |
std::vector< SparseVectorT > & | M_v, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
block_matrix & | g_A_I_J_vcl, | ||
block_vector & | g_bv_vcl, | ||
std::vector< SparseVectorT > & | g_res, | ||
std::vector< cl_uint > & | g_is_update, | ||
const spai_tag & | tag, | ||
viennacl::context | ctx | ||
) |
Solution of Least square problem on GPU.
A_v_c | column-major vectorized initial sparse matrix |
M_v | column-major vectorized sparse preconditioner matrix |
g_I | container of row set indices |
g_J | container of column set indices |
g_A_I_J_vcl | contigious matrix that consists of blocks A(I_k, J_k) |
g_bv_vcl | contigious vector that consists of betas, necessary for Q recovery |
g_res | container of residuals |
g_is_update | container with indicators which blocks are active |
tag | spai tag |
ctx | Optional context in which the auxiliary data is created (one out of multiple OpenCL contexts, CUDA, host) |
void viennacl::linalg::detail::spai::least_square_solve | ( | std::vector< SparseVectorT > const & | A_v_c, |
std::vector< DenseMatrixT > & | g_R, | ||
std::vector< VectorT > & | g_b_v, | ||
std::vector< std::vector< unsigned int > > & | g_I, | ||
std::vector< std::vector< unsigned int > > & | g_J, | ||
std::vector< SparseVectorT > & | g_res, | ||
std::vector< bool > & | g_is_update, | ||
std::vector< SparseVectorT > & | M_v, | ||
spai_tag const & | tag | ||
) |
Solution of Least square problem on CPU.
A_v_c | column-major vectorized initial sparse matrix |
g_R | blocks for least square solution |
g_b_v | vectors beta, necessary for Q recovery |
g_I | container of row index set for all columns of matrix M |
g_J | container of column index set for all columns of matrix M |
g_res | container of residuals |
g_is_update | container with indicators which blocks are active |
M_v | column-major vectorized sparse matrix, final preconditioner |
tag | spai tag |
void viennacl::linalg::detail::spai::make_rotation_matrix | ( | MatrixT & | mat, |
vcl_size_t | new_size, | ||
vcl_size_t | off_diagonal_distance = 4 |
||
) |
Definition at line 61 of file small_matrix.hpp.
void viennacl::linalg::detail::spai::Print | ( | std::ostream & | ostr, |
InputIteratorT | it_begin, | ||
InputIteratorT | it_end | ||
) |
void viennacl::linalg::detail::spai::print_continious_matrix | ( | VectorT & | con_A_I_J, |
std::vector< cl_uint > & | blocks_ind, | ||
std::vector< std::vector< unsigned int > > const & | g_I, | ||
std::vector< std::vector< unsigned int > > const & | g_J | ||
) |
void viennacl::linalg::detail::spai::print_continious_vector | ( | VectorT & | con_v, |
std::vector< cl_uint > & | block_ind, | ||
std::vector< std::vector< unsigned int > > const & | g_J | ||
) |
void viennacl::linalg::detail::spai::print_matrix | ( | DenseMatrixT & | m | ) |
void viennacl::linalg::detail::spai::print_sparse_vector | ( | SparseVectorT const & | v | ) |
void viennacl::linalg::detail::spai::projectI | ( | std::vector< unsigned int > const & | I, |
VectorT & | y, | ||
unsigned int | ind | ||
) |
Perform projection of set I on the unit-vector.
I | set of non-zero rows |
y | result vector |
ind | index of unit vector |
Definition at line 122 of file spai-static.hpp.
void viennacl::linalg::detail::spai::projectRows | ( | std::vector< SparseVectorT > const & | A_v_c, |
std::vector< unsigned int > const & | J, | ||
std::vector< unsigned int > & | I | ||
) |
Row projection for matrix A(:,J) -> A(I,J), building index set of non-zero rows.
A_v_c | input matrix |
J | set of non-zero rows |
I | output matrix |
Definition at line 171 of file spai-static.hpp.
void viennacl::linalg::detail::spai::QRBlockComposition | ( | MatrixT const & | A_I_J, |
MatrixT const & | A_I_J_u, | ||
MatrixT & | A_I_u_J_u | ||
) |
Composition of new block for QR factorization cf. Kallischko dissertation p.82, figure 4.7.
A_I_J | previously composed block |
A_I_J_u | matrix Q'*A(I, \tilde J), where \tilde J - set of new column indices |
A_I_u_J_u | is composition of lower part A(I, \tilde J) and A(\tilde I, \tilde J) - new block for QR decomposition |
Definition at line 250 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::single_qr | ( | MatrixT & | R, |
VectorT & | b_v | ||
) |
void viennacl::linalg::detail::spai::sparse_inner_prod | ( | SparseVectorT const & | v1, |
SparseVectorT const & | v2, | ||
NumericT & | res_v | ||
) |
Dot product of two sparse vectors.
v1 | initial sparse vector |
v2 | initial sparse vector |
res_v | scalar that represents dot product result |
Definition at line 157 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::sparse_norm_2 | ( | SparseVectorT const & | v, |
NumericT & | norm | ||
) |
Computation of Euclidean norm for sparse vector.
v | initial sparse vector |
norm | scalar that represents Euclidean norm |
Definition at line 141 of file spai-dynamic.hpp.
void viennacl::linalg::detail::spai::sparse_transpose | ( | MatrixT const & | A_in, |
MatrixT & | A | ||
) |
void viennacl::linalg::detail::spai::store_householder_vector | ( | MatrixT & | A, |
unsigned int | ind, | ||
VectorT & | v | ||
) |
void viennacl::linalg::detail::spai::sym_sparse_matrix_to_stl | ( | MatrixT const & | A, |
std::vector< std::map< unsigned int, NumericT > > & | STL_A | ||
) |
void viennacl::linalg::detail::spai::vectorize_column_matrix | ( | SparseMatrixT const & | M_in, |
std::vector< SparseVectorT > & | M_v | ||
) |
void viennacl::linalg::detail::spai::vectorize_row_matrix | ( | SparseMatrixT const & | M_in, |
std::vector< SparseVectorT > & | M_v | ||
) |
void viennacl::linalg::detail::spai::write_set_to_array | ( | std::vector< std::vector< SizeT > > const & | ind_set, |
std::vector< cl_uint > & | a | ||
) |