ViennaCL - The Vienna Computing Library  1.4.2
Data Structures | Typedefs | Functions
viennacl::linalg::detail::spai Namespace Reference

Implementation namespace for sparse approximate inverse preconditioner. More...

Data Structures

class  block_matrix
 Represents a contigious matrices on GPU. More...
class  block_vector
 Represents a contigious vector on GPU. More...
class  fspai_tag
 A tag for FSPAI. Experimental. Contains values for the algorithm. Must be passed to spai_precond constructor. More...
struct  CompareSecond
class  spai_tag
 A tag for SPAI Contains values for the algorithm. Must be passed to spai_precond constructor. More...
class  sparse_vector
 Represents sparse vector based on std::map<unsigned int, ScalarType> More...

Typedefs

typedef std::pair< unsigned
int, double > 
PairT

Functions

template<typename MatrixType , typename ScalarType >
void sym_sparse_matrix_to_stl (MatrixType const &A, std::vector< std::map< unsigned int, ScalarType > > &STL_A)
template<typename MatrixType >
void generateJ (MatrixType const &A, std::vector< std::vector< std::size_t > > &J)
template<typename ScalarType , typename MatrixType , typename VectorType >
void fill_blocks (std::vector< std::map< unsigned int, ScalarType > > &A, std::vector< MatrixType > &blocks, std::vector< std::vector< std::size_t > > const &J, std::vector< VectorType > &Y)
template<typename MatrixType >
void cholesky_decompose (MatrixType &A)
template<typename MatrixType , typename VectorType >
void cholesky_solve (MatrixType const &L, VectorType &b)
template<typename MatrixType , typename VectorType1 >
void computeL (MatrixType const &A, MatrixType &L, MatrixType &L_trans, std::vector< VectorType1 > &Y, std::vector< std::vector< std::size_t > > &J)
template<typename MatrixType >
void computeFSPAI (MatrixType const &A, MatrixType const &PatternA, MatrixType &L, MatrixType &L_trans, fspai_tag)
template<typename T , typename InputIterator >
void Print (std::ostream &ostr, InputIterator it_begin, InputIterator it_end)
template<typename VectorType , typename MatrixType >
void write_to_block (VectorType &con_A_I_J, unsigned int start_ind, const std::vector< unsigned int > &I, const std::vector< unsigned int > &J, MatrixType &m)
template<typename VectorType >
void print_continious_matrix (VectorType &con_A_I_J, std::vector< cl_uint > &blocks_ind, const std::vector< std::vector< unsigned int > > &g_I, const std::vector< std::vector< unsigned int > > &g_J)
template<typename VectorType >
void print_continious_vector (VectorType &con_v, std::vector< cl_uint > &block_ind, const std::vector< std::vector< unsigned int > > &g_J)
void compute_blocks_size (const std::vector< std::vector< unsigned int > > &g_I, const std::vector< std::vector< unsigned int > > &g_J, unsigned int &sz, std::vector< cl_uint > &blocks_ind, std::vector< cl_uint > &matrix_dims)
 **************************************** BLOCK FUNCTIONS ************************************//
template<typename SizeType >
void get_size (const std::vector< std::vector< SizeType > > &inds, SizeType &size)
 Computes size of particular container of index set.
template<typename SizeType >
void init_start_inds (const std::vector< std::vector< SizeType > > &inds, std::vector< cl_uint > &start_inds)
 Initializes start indices of particular index set.
template<typename MatrixType , typename ScalarType >
void dot_prod (const MatrixType &A, unsigned int beg_ind, ScalarType &res)
 Dot prod of particular column of martix A with it's self starting at a certain index beg_ind.
template<typename MatrixType , typename VectorType , typename ScalarType >
void custom_inner_prod (const MatrixType &A, const VectorType &v, unsigned int col_ind, unsigned int start_ind, ScalarType &res)
 Dot prod of particular matrix column with arbitrary vector: A(:, col_ind)
template<typename MatrixType , typename VectorType >
void copy_vector (const MatrixType &A, VectorType &v, const unsigned int beg_ind)
 Copying part of matrix column.
template<typename MatrixType , typename VectorType , typename ScalarType >
void householder_vector (const MatrixType &A, unsigned int j, VectorType &v, ScalarType &b)
 Coputation of Householder vector, householder reflection c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.210.
template<typename MatrixType , typename VectorType , typename ScalarType >
void apply_householder_reflection (MatrixType &A, unsigned int iter_cnt, VectorType &v, ScalarType b)
 Inplace application of Householder vector to a matrix A.
template<typename MatrixType , typename VectorType >
void store_householder_vector (MatrixType &A, unsigned int ind, VectorType &v)
 Storage of vector v in column(A, ind), starting from ind-1 index of a column.
template<typename MatrixType , typename VectorType >
void single_qr (MatrixType &R, VectorType &b_v)
 Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224.
template<typename SizeType >
void get_max_block_size (const std::vector< std::vector< SizeType > > &inds, SizeType max_size)
 Getting max size of rows/columns from container of index set.
template<typename MatrixType , typename VectorType , typename ScalarType >
void custom_dot_prod (const MatrixType &A, const VectorType &v, unsigned int ind, ScalarType &res)
 Dot_prod(column(A, ind), v) starting from index ind+1.
template<typename MatrixType , typename VectorType >
void apply_q_trans_vec (const MatrixType &R, const VectorType &b_v, VectorType &y)
 Recovery Q from matrix R and vector of betas b_v.
template<typename MatrixType , typename VectorType >
void apply_q_trans_mat (const MatrixType &R, const VectorType &b_v, MatrixType &A)
 Multiplication of Q'*A, where Q is in implicit for lower part of R and vector of betas - b_v.
template<typename ScalarType >
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)
 Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224 performed on GPU.
template<typename MatrixType >
void make_rotation_matrix (MatrixType &mat, std::size_t new_size, std::size_t off_diagonal_distance=4)
template<typename MatrixType >
double determinant (boost::numeric::ublas::matrix_expression< MatrixType > const &mat_r)
template<typename MatrixType >
void composeNewR (const MatrixType &A, const MatrixType &R_n, MatrixType &R)
 Composition of new matrix R, that is going to be used in Least Square problem solving.
template<typename VectorType >
void composeNewVector (const VectorType &v_n, VectorType &v)
 Composition of new vector of coefficients beta from QR factorizations(necessary for Q recovery)
template<typename SparseVectorType , typename ScalarType >
void sparse_norm_2 (const SparseVectorType &v, ScalarType &norm)
 Computation of Euclidean norm for sparse vector.
template<typename SparseVectorType , typename ScalarType >
void sparse_inner_prod (const SparseVectorType &v1, const SparseVectorType &v2, ScalarType &res_v)
 Dot product of two sparse vectors.
template<typename SparseVectorType , typename ScalarType >
bool buildAugmentedIndexSet (const std::vector< SparseVectorType > &A_v_c, const SparseVectorType &res, std::vector< unsigned int > &J, std::vector< unsigned int > &J_u, const spai_tag &tag)
 Building a new set of column indices J_u, cf. Kallischko dissertation p.31.
template<typename SparseVectorType >
void buildNewRowSet (const std::vector< SparseVectorType > &A_v_c, const std::vector< unsigned int > &I, const std::vector< unsigned int > &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.
template<typename MatrixType >
void QRBlockComposition (const MatrixType &A_I_J, const MatrixType &A_I_J_u, MatrixType &A_I_u_J_u)
 Composition of new block for QR factorization cf. Kallischko dissertation p.82, figure 4.7.
template<typename SparseMatrixType , typename SparseVectorType , typename DenseMatrixType , typename VectorType >
void block_update (const SparseMatrixType &A, const std::vector< SparseVectorType > &A_v_c, std::vector< SparseVectorType > &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< VectorType > &g_b_v, std::vector< DenseMatrixType > &g_A_I_J, spai_tag const &tag)
 CPU-based dynamic update for SPAI preconditioner.
template<typename ScalarType >
void block_q_multiplication (const std::vector< std::vector< unsigned int > > &g_J_u, const std::vector< std::vector< unsigned int > > &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)
 Performs multiplication Q'*A(I, \tilde J) on GPU.
template<typename SizeType >
void assemble_qr_row_inds (const std::vector< std::vector< SizeType > > &g_I, const std::vector< std::vector< SizeType > > g_J, const std::vector< std::vector< SizeType > > &g_I_u, std::vector< std::vector< SizeType > > &g_I_q)
 Assembly of container of index row sets: I_q, row indices for new "QR block".
template<typename ScalarType >
void assemble_qr_block (const std::vector< std::vector< unsigned int > > &g_J, const std::vector< std::vector< unsigned int > > &g_I, const std::vector< std::vector< unsigned int > > &g_J_u, const std::vector< std::vector< unsigned int > > &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, const bool is_empty_block)
 Performs assembly for new QR block.
template<typename ScalarType >
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)
 Performs assembly for new R matrix on GPU.
template<typename ScalarType , unsigned int MAT_ALIGNMENT, typename SparseVectorType >
void block_update (const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &A, const std::vector< SparseVectorType > &A_v_c, std::vector< cl_uint > &g_is_update, std::vector< SparseVectorType > &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.
template<typename SizeType >
bool isInIndexSet (const std::vector< SizeType > &J, SizeType ind)
 Determines if element ind is in set {J}.
template<typename VectorType , typename SparseVectorType >
void fanOutVector (const VectorType &m_in, const std::vector< unsigned int > &J, SparseVectorType &m)
 Projects solution of LS problem onto original column m.
template<typename MatrixType , typename VectorType >
void backwardSolve (const MatrixType &R, const VectorType &y, VectorType &x)
 Solution of linear:R*x=y system by backward substitution.
template<typename VectorType , typename ScalarType >
void projectI (const std::vector< unsigned int > &I, VectorType &y, unsigned int ind)
 Perform projection of set I on the unit-vector.
template<typename SparseVectorType >
void buildColumnIndexSet (const SparseVectorType &v, std::vector< unsigned int > &J)
 Builds index set of projected columns for current column of preconditioner.
template<typename SparseMatrixType >
void initPreconditioner (const SparseMatrixType &A, SparseMatrixType &M)
 Initialize preconditioner with sparcity pattern = p(A)
template<typename SparseVectorType >
void projectRows (const std::vector< SparseVectorType > &A_v_c, const std::vector< unsigned int > &J, std::vector< unsigned int > &I)
 Row projection for matrix A(:,J) -> A(I,J), building index set of non-zero rows.
template<typename SparseVectorType >
void print_sparse_vector (const SparseVectorType &v)
template<typename DenseMatrixType >
void print_matrix (DenseMatrixType &m)
template<typename SparseVectorType , typename ScalarType >
void add_sparse_vectors (const SparseVectorType &v, const ScalarType b, SparseVectorType &res_v)
 Add two sparse vectors res_v = b*v.
template<typename SparseVectorType , typename ScalarType >
void compute_spai_residual (const std::vector< SparseVectorType > &A_v_c, const SparseVectorType &v, const unsigned int ind, SparseVectorType &res)
 Computation of residual res = A*v - e.
template<typename SparseVectorType >
void build_index_set (const std::vector< SparseVectorType > &A_v_c, const SparseVectorType &v, std::vector< unsigned int > &J, std::vector< unsigned int > &I)
 Setting up index set of columns and rows for certain column.
template<typename SparseMatrixType , typename DenseMatrixType >
void initProjectSubMatrix (const SparseMatrixType &A_in, const std::vector< unsigned int > &J, std::vector< unsigned int > &I, DenseMatrixType &A_out)
 Initializes a dense matrix from a sparse one.
template<typename SparseMatrixType , typename DenseMatrixType , typename SparseVectorType , typename VectorType >
void block_set_up (const SparseMatrixType &A, const std::vector< SparseVectorType > &A_v_c, const std::vector< SparseVectorType > &M_v, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, std::vector< DenseMatrixType > &g_A_I_J, std::vector< VectorType > &g_b_v)
 Setting up blocks and QR factorizing them on CPU.
template<typename SparseVectorType >
void index_set_up (const std::vector< SparseVectorType > &A_v_c, const std::vector< SparseVectorType > &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.
template<typename ScalarType , unsigned int MAT_ALIGNMENT, typename SparseVectorType >
void block_set_up (const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &A, const std::vector< SparseVectorType > &A_v_c, const std::vector< SparseVectorType > &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.
template<typename ScalarType , typename SparseVectorType >
void custom_fan_out (const std::vector< ScalarType > &m_in, unsigned int start_m_ind, const std::vector< unsigned int > &J, SparseVectorType &m)
 Elicitation of sparse vector m for particular column from m_in - contigious vector for all columns.
template<typename SparseVectorType , typename ScalarType >
void least_square_solve (std::vector< SparseVectorType > &A_v_c, std::vector< SparseVectorType > &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< SparseVectorType > &g_res, std::vector< cl_uint > &g_is_update, const spai_tag &tag)
 Solution of Least square problem on GPU.
template<typename SparseVectorType , typename DenseMatrixType , typename VectorType >
void least_square_solve (const std::vector< SparseVectorType > &A_v_c, std::vector< DenseMatrixType > &g_R, std::vector< VectorType > &g_b_v, std::vector< std::vector< unsigned int > > &g_I, std::vector< std::vector< unsigned int > > &g_J, std::vector< SparseVectorType > &g_res, std::vector< bool > &g_is_update, std::vector< SparseVectorType > &M_v, const spai_tag &tag)
 Solution of Least square problem on CPU.
template<typename VectorType >
bool is_all_update (VectorType &parallel_is_update)
template<typename SparseMatrixType , typename SparseVectorType >
void vectorize_column_matrix (const SparseMatrixType &M_in, std::vector< SparseVectorType > &M_v)
 Solution of Least square problem on CPU.
template<typename SparseMatrixType , typename SparseVectorType >
void vectorize_row_matrix (const SparseMatrixType &M_in, std::vector< SparseVectorType > &M_v)
template<typename SizeType >
void write_set_to_array (const std::vector< std::vector< SizeType > > &ind_set, std::vector< cl_uint > &a)
template<typename ScalarType , unsigned int MAT_ALIGNMENT>
void block_assembly (const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &A, const std::vector< std::vector< unsigned int > > &g_J, const std::vector< std::vector< unsigned int > > &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.
template<typename SparseMatrixType , typename SparseVectorType >
void insert_sparse_columns (const std::vector< SparseVectorType > &M_v, SparseMatrixType &M, bool is_right)
 Insertion of vectorized matrix column into original sparse matrix.
template<typename MatrixType >
void sparse_transpose (const MatrixType &A_in, MatrixType &A)
 Transposition of sparse matrix.
template<typename MatrixType >
void computeSPAI (const MatrixType &A, MatrixType &M, spai_tag &tag)
 Construction of SPAI preconditioner on CPU.
template<typename ScalarType , unsigned int MAT_ALIGNMENT>
void computeSPAI (const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &A, const boost::numeric::ublas::compressed_matrix< ScalarType > &cpu_A, boost::numeric::ublas::compressed_matrix< ScalarType > &cpu_M, viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &M, const spai_tag &tag)
 Construction of SPAI preconditioner on GPU.

Detailed Description

Implementation namespace for sparse approximate inverse preconditioner.


Typedef Documentation

typedef std::pair<unsigned int, double> PairT

Function Documentation

void viennacl::linalg::detail::spai::add_sparse_vectors ( const SparseVectorType &  v,
const ScalarType  b,
SparseVectorType &  res_v 
)

Add two sparse vectors res_v = b*v.

Parameters:
vinitial sparse vector
bscalar
res_voutput vector
void viennacl::linalg::detail::spai::apply_householder_reflection ( MatrixType &  A,
unsigned int  iter_cnt,
VectorType &  v,
ScalarType  b 
)

Inplace application of Householder vector to a matrix A.

Parameters:
Ainit matrix
iter_cntcurrent iteration
vHouseholder vector
bbeta
void viennacl::linalg::detail::spai::apply_q_trans_mat ( const MatrixType &  R,
const VectorType &  b_v,
MatrixType &  A 
)

Multiplication of Q'*A, where Q is in implicit for lower part of R and vector of betas - b_v.

Parameters:
Rinput matrix
b_vvector of betas
Aoutput matrix
void viennacl::linalg::detail::spai::apply_q_trans_vec ( const MatrixType &  R,
const VectorType &  b_v,
VectorType &  y 
)

Recovery Q from matrix R and vector of betas b_v.

Parameters:
Rinput matrix
b_vvector of betas
youtput vector
void viennacl::linalg::detail::spai::assemble_qr_block ( const std::vector< std::vector< unsigned int > > &  g_J,
const std::vector< std::vector< unsigned int > > &  g_I,
const std::vector< std::vector< unsigned int > > &  g_J_u,
const std::vector< std::vector< unsigned int > > &  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,
const bool  is_empty_block 
)

Performs assembly for new QR block.

Parameters:
g_Jcontainer of column indices
g_Icontainer of row indices
g_J_ucontainer of new column indices
g_I_ucontainer of new row indices
g_I_qcontainer of row indices for new QR blocks
g_A_I_J_u_vclblocks of Q'*A(I, \tilde J)
matrix_dimensionsarray with matrix dimensions for all blocks
g_A_I_u_J_u_vclblocks A(\tilde I, \tilde J)
g_is_updatecontainer with update indicators
is_empty_blockindicator 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
void viennacl::linalg::detail::spai::assemble_qr_row_inds ( const std::vector< std::vector< SizeType > > &  g_I,
const std::vector< std::vector< SizeType > >  g_J,
const std::vector< std::vector< SizeType > > &  g_I_u,
std::vector< std::vector< SizeType > > &  g_I_q 
)

Assembly of container of index row sets: I_q, row indices for new "QR block".

Parameters:
g_Icontainer of row indices
g_Jcontainer of column indices
g_I_ucontainer of new row indices
g_I_qcontainer of row indices for new QR blocks
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 
)

Performs assembly for new R matrix on GPU.

Parameters:
g_Icontainer of row indices
g_Jcontainer of column indices
g_A_I_J_vclcontainer of block matrices from previous update
g_A_I_J_u_vclcontainer of block matrices Q'*A(I, \tilde J)
g_A_I_u_J_u_vclcontainer of block matrices QR factored on current iteration
g_bv_vclblock of beta vectors from previous iteration
g_bv_vcl_ublock of updated beta vectors got after recent QR factorization
g_is_updatecontainer with identificators that shows which block should be modified
void viennacl::linalg::detail::spai::backwardSolve ( const MatrixType &  R,
const VectorType &  y,
VectorType &  x 
)

Solution of linear:R*x=y system by backward substitution.

Parameters:
Ruppertriangular matrix
yright handside vector
xsolution vector
void viennacl::linalg::detail::spai::block_assembly ( const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &  A,
const std::vector< std::vector< unsigned int > > &  g_J,
const std::vector< std::vector< unsigned int > > &  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.

Parameters:
Aintial sparse matrix
g_Jcontainer of column index set
g_Icontainer of row index set
g_A_I_J_vclcontigious blocks A(I, J) using GPU memory
g_is_updatecontainer with indicators which blocks are active
is_empty_blockparameter that indicates if no block were assembled
void viennacl::linalg::detail::spai::block_q_multiplication ( const std::vector< std::vector< unsigned int > > &  g_J_u,
const std::vector< std::vector< unsigned int > > &  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 
)

Performs multiplication Q'*A(I, \tilde J) on GPU.

Parameters:
g_J_ucontainer of sets of new column indices
g_Icontainer of row indices
g_A_I_J_vclblock matrix composed from previous blocks, they are blocks of R
g_bv_vclblock of beta vectors
g_A_I_J_u_vclblock of matrices A(I, \tilde J)
g_is_updateindicators, that show if a certain block should be processed
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 
)

Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224 performed on GPU.

Parameters:
g_Icontainer of row indices
g_Jcontainer of column indices
g_A_I_J_vclcontigious matrices, GPU memory is used
g_bv_vclcontigios vectors beta, GPU memory is used
g_is_updatecontainer of indicators that show active blocks
void viennacl::linalg::detail::spai::block_set_up ( const SparseMatrixType &  A,
const std::vector< SparseVectorType > &  A_v_c,
const std::vector< SparseVectorType > &  M_v,
std::vector< std::vector< unsigned int > > &  g_I,
std::vector< std::vector< unsigned int > > &  g_J,
std::vector< DenseMatrixType > &  g_A_I_J,
std::vector< VectorType > &  g_b_v 
)

Setting up blocks and QR factorizing them on CPU.

Parameters:
Ainitial sparse matrix
A_v_ccolumn major vectorized initial sparse matrix
M_vinitialized preconditioner
g_Icontainer of row indices
g_Jcontainer of column indices
g_A_I_Jcontainer of dense matrices -> R matrices after QR factorization
g_b_vcontainer of vectors beta, necessary for Q recovery
void viennacl::linalg::detail::spai::block_set_up ( const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &  A,
const std::vector< SparseVectorType > &  A_v_c,
const std::vector< SparseVectorType > &  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.

Parameters:
Ainitial sparse matrix
A_v_ccolumn major vectorized initial sparse matrix
M_vinitialized preconditioner
g_is_updatecontainer that indicates which blocks are active
g_Icontainer of row indices
g_Jcontainer of column indices
g_A_I_Jcontainer of dense matrices -> R matrices after QR factorization
g_bvcontainer of vectors beta, necessary for Q recovery
void viennacl::linalg::detail::spai::block_update ( const SparseMatrixType &  A,
const std::vector< SparseVectorType > &  A_v_c,
std::vector< SparseVectorType > &  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< VectorType > &  g_b_v,
std::vector< DenseMatrixType > &  g_A_I_J,
spai_tag const &  tag 
)

CPU-based dynamic update for SPAI preconditioner.

Parameters:
Ainitial sparse matrix
A_v_cvectorized column-wise initial matrix
g_rescontainer of residuals for all columns
g_is_updatecontainer with identificators that shows which block should be modified
g_Icontainer of row index sets for all columns
g_Jcontainer of column index sets for all columns
g_b_vcontainer of vectors of beta for Q recovery(cf. Golub Van Loan "Matrix Computations", 3rd edition p.211)
g_A_I_Jcontainer of block matrices from previous update
tagSPAI configuration tag
void viennacl::linalg::detail::spai::block_update ( const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &  A,
const std::vector< SparseVectorType > &  A_v_c,
std::vector< cl_uint > &  g_is_update,
std::vector< SparseVectorType > &  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.

Parameters:
Asparse matrix
A_v_cvectorized column-wise initial matrix
g_is_updatecontainer with identificators that shows which block should be modified
g_rescontainer of residuals for all columns
g_Jcontainer of column index sets for all columns
g_Icontainer of row index sets for all columns
g_A_I_J_vclcontainer of block matrices from previous update
g_bv_vclblock of beta vectors from previous iteration
tagSPAI configuration tag
void viennacl::linalg::detail::spai::build_index_set ( const std::vector< SparseVectorType > &  A_v_c,
const SparseVectorType &  v,
std::vector< unsigned int > &  J,
std::vector< unsigned int > &  I 
)

Setting up index set of columns and rows for certain column.

Parameters:
A_v_ccolumn major vectorized initial sparse matrix
vcurrent column of preconditioner matrix
Jset of column indices
Iset of row indices
bool viennacl::linalg::detail::spai::buildAugmentedIndexSet ( const std::vector< SparseVectorType > &  A_v_c,
const SparseVectorType &  res,
std::vector< unsigned int > &  J,
std::vector< unsigned int > &  J_u,
const spai_tag &  tag 
)

Building a new set of column indices J_u, cf. Kallischko dissertation p.31.

Parameters:
A_v_cvectorized column-wise initial matrix
resresidual vector
Jset of column indices
J_uset of new column indices
tagSPAI tag with parameters
void viennacl::linalg::detail::spai::buildColumnIndexSet ( const SparseVectorType &  v,
std::vector< unsigned int > &  J 
)

Builds index set of projected columns for current column of preconditioner.

Parameters:
vcurrent column of preconditioner
Joutput - index set of non-zero columns
void viennacl::linalg::detail::spai::buildNewRowSet ( const std::vector< SparseVectorType > &  A_v_c,
const std::vector< unsigned int > &  I,
const std::vector< unsigned int > &  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.

Parameters:
A_v_cvectorized column-wise initial matrix
Iset of previous determined row indices
J_nset of new column indices
I_nset of new indices
void viennacl::linalg::detail::spai::cholesky_solve ( MatrixType const &  L,
VectorType &  b 
)
void viennacl::linalg::detail::spai::composeNewR ( const MatrixType &  A,
const MatrixType &  R_n,
MatrixType &  R 
)

Composition of new matrix R, that is going to be used in Least Square problem solving.

Parameters:
Amatrix Q'*A(I, \tilde J), where \tilde J - set of new column indices
R_nmatrix A_Iu_J_u after QR factorization
Rpreviously composed matrix R
void viennacl::linalg::detail::spai::composeNewVector ( const VectorType &  v_n,
VectorType &  v 
)

Composition of new vector of coefficients beta from QR factorizations(necessary for Q recovery)

Parameters:
v_nnew vector from last QR factorization
vcomposition of previous vectors from QR factorizations
void viennacl::linalg::detail::spai::compute_blocks_size ( const std::vector< std::vector< unsigned int > > &  g_I,
const std::vector< std::vector< unsigned int > > &  g_J,
unsigned int &  sz,
std::vector< cl_uint > &  blocks_ind,
std::vector< cl_uint > &  matrix_dims 
) [inline]

**************************************** BLOCK FUNCTIONS ************************************//

Computes size of elements, start indices and matrix dimensions for a certain block

Parameters:
g_Icontainer of row indices
g_Jcontainer of column indices
szgeneral size for all elements in a certain block
blocks_indstart indices in a certain
matrix_dimsmatrix dimensions for each block
void viennacl::linalg::detail::spai::compute_spai_residual ( const std::vector< SparseVectorType > &  A_v_c,
const SparseVectorType &  v,
const unsigned int  ind,
SparseVectorType &  res 
)

Computation of residual res = A*v - e.

Parameters:
A_v_ccolumn major vectorized input sparse matrix
vsparse vector, in this case new column of preconditioner matrix
indindex for current column
resresidual
void viennacl::linalg::detail::spai::computeFSPAI ( MatrixType const &  A,
MatrixType const &  PatternA,
MatrixType &  L,
MatrixType &  L_trans,
fspai_tag   
)
void viennacl::linalg::detail::spai::computeL ( MatrixType const &  A,
MatrixType &  L,
MatrixType &  L_trans,
std::vector< VectorType1 > &  Y,
std::vector< std::vector< std::size_t > > &  J 
)
void viennacl::linalg::detail::spai::computeSPAI ( const MatrixType &  A,
MatrixType &  M,
spai_tag &  tag 
)

Construction of SPAI preconditioner on CPU.

Parameters:
Ainitial sparse matrix
Moutput preconditioner
tagspai tag
void viennacl::linalg::detail::spai::computeSPAI ( const viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &  A,
const boost::numeric::ublas::compressed_matrix< ScalarType > &  cpu_A,
boost::numeric::ublas::compressed_matrix< ScalarType > &  cpu_M,
viennacl::compressed_matrix< ScalarType, MAT_ALIGNMENT > &  M,
const spai_tag &  tag 
)

Construction of SPAI preconditioner on GPU.

Parameters:
Ainitial sparse matrix
cpu_Acopy of initial matrix on CPU
cpu_Moutput preconditioner on CPU
Moutput preconditioner
tagSPAI tag class with parameters
void viennacl::linalg::detail::spai::copy_vector ( const MatrixType &  A,
VectorType &  v,
const unsigned int  beg_ind 
)

Copying part of matrix column.

Parameters:
Ainit matrix
voutput vector
beg_indstart index for copying
void viennacl::linalg::detail::spai::custom_dot_prod ( const MatrixType &  A,
const VectorType &  v,
unsigned int  ind,
ScalarType &  res 
)

Dot_prod(column(A, ind), v) starting from index ind+1.

Parameters:
Ainput matrix
vinput vector
indindex
resresult value
void viennacl::linalg::detail::spai::custom_fan_out ( const std::vector< ScalarType > &  m_in,
unsigned int  start_m_ind,
const std::vector< unsigned int > &  J,
SparseVectorType &  m 
)

Elicitation of sparse vector m for particular column from m_in - contigious vector for all columns.

Parameters:
m_incontigious sparse vector for all columns
start_m_indstart index of particular vector
Jcolumn index set
msparse vector for particular column
void viennacl::linalg::detail::spai::custom_inner_prod ( const MatrixType &  A,
const VectorType &  v,
unsigned int  col_ind,
unsigned int  start_ind,
ScalarType &  res 
)

Dot prod of particular matrix column with arbitrary vector: A(:, col_ind)

Parameters:
Ainit matrix
vinput vector
col_indstarting column index
start_indstarting index inside column
resresult of dot product
double viennacl::linalg::detail::spai::determinant ( boost::numeric::ublas::matrix_expression< MatrixType > const &  mat_r)
void viennacl::linalg::detail::spai::dot_prod ( const MatrixType &  A,
unsigned int  beg_ind,
ScalarType &  res 
)

Dot prod of particular column of martix A with it's self starting at a certain index beg_ind.

Parameters:
Ainit matrix
beg_indstarting index
resresult of dot product
void viennacl::linalg::detail::spai::fanOutVector ( const VectorType &  m_in,
const std::vector< unsigned int > &  J,
SparseVectorType &  m 
)

Projects solution of LS problem onto original column m.

Parameters:
m_insolution of LS
Jset of non-zero columns
moriginal column of M
void viennacl::linalg::detail::spai::fill_blocks ( std::vector< std::map< unsigned int, ScalarType > > &  A,
std::vector< MatrixType > &  blocks,
std::vector< std::vector< std::size_t > > const &  J,
std::vector< VectorType > &  Y 
)
void viennacl::linalg::detail::spai::generateJ ( MatrixType const &  A,
std::vector< std::vector< std::size_t > > &  J 
)
void viennacl::linalg::detail::spai::get_max_block_size ( const std::vector< std::vector< SizeType > > &  inds,
SizeType  max_size 
)

Getting max size of rows/columns from container of index set.

Parameters:
indscontainer of index set
max_sizemax size that corresponds to that container
void viennacl::linalg::detail::spai::get_size ( const std::vector< std::vector< SizeType > > &  inds,
SizeType &  size 
)

Computes size of particular container of index set.

Parameters:
indscontainer of index sets
sizeoutput size
void viennacl::linalg::detail::spai::householder_vector ( const MatrixType &  A,
unsigned int  j,
VectorType &  v,
ScalarType &  b 
)

Coputation of Householder vector, householder reflection c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.210.

Parameters:
Ainit matrix
jstart index for computations
voutput Householder vector
bbeta
void viennacl::linalg::detail::spai::index_set_up ( const std::vector< SparseVectorType > &  A_v_c,
const std::vector< SparseVectorType > &  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.

Parameters:
A_v_ccolumn major vectorized initial sparse matrix
M_vinitialized preconditioner
g_Jcontainer of column indices
g_Icontainer of row indices
void viennacl::linalg::detail::spai::init_start_inds ( const std::vector< std::vector< SizeType > > &  inds,
std::vector< cl_uint > &  start_inds 
)

Initializes start indices of particular index set.

Parameters:
indscontainer of index sets
start_indsoutput index set
void viennacl::linalg::detail::spai::initPreconditioner ( const SparseMatrixType &  A,
SparseMatrixType &  M 
)

Initialize preconditioner with sparcity pattern = p(A)

Parameters:
Ainput matrix
Moutput matrix - initialized preconditioner
void viennacl::linalg::detail::spai::initProjectSubMatrix ( const SparseMatrixType &  A_in,
const std::vector< unsigned int > &  J,
std::vector< unsigned int > &  I,
DenseMatrixType &  A_out 
)

Initializes a dense matrix from a sparse one.

Parameters:
A_inRiginal sparse matrix
JSet of column indices
ISet of row indices
A_outdense matrix output
void viennacl::linalg::detail::spai::insert_sparse_columns ( const std::vector< SparseVectorType > &  M_v,
SparseMatrixType &  M,
bool  is_right 
)

Insertion of vectorized matrix column into original sparse matrix.

Parameters:
M_vcolumn-major vectorized matrix
Moriginal sparse matrix
is_rightindicates if matrix should be transposed in the output
bool viennacl::linalg::detail::spai::is_all_update ( VectorType &  parallel_is_update)
bool viennacl::linalg::detail::spai::isInIndexSet ( const std::vector< SizeType > &  J,
SizeType  ind 
)

Determines if element ind is in set {J}.

Parameters:
Jcurrent set
indcurrent element
void viennacl::linalg::detail::spai::least_square_solve ( std::vector< SparseVectorType > &  A_v_c,
std::vector< SparseVectorType > &  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< SparseVectorType > &  g_res,
std::vector< cl_uint > &  g_is_update,
const spai_tag &  tag 
)

Solution of Least square problem on GPU.

Parameters:
A_v_ccolumn-major vectorized initial sparse matrix
M_vcolumn-major vectorized sparse preconditioner matrix
g_Icontainer of row set indices
g_Jcontainer of column set indices
g_A_I_J_vclcontigious matrix that consists of blocks A(I_k, J_k)
g_bv_vclcontigious vector that consists of betas, necessary for Q recovery
g_rescontainer of residuals
g_is_updatecontainer with indicators which blocks are active
tagspai tag
void viennacl::linalg::detail::spai::least_square_solve ( const std::vector< SparseVectorType > &  A_v_c,
std::vector< DenseMatrixType > &  g_R,
std::vector< VectorType > &  g_b_v,
std::vector< std::vector< unsigned int > > &  g_I,
std::vector< std::vector< unsigned int > > &  g_J,
std::vector< SparseVectorType > &  g_res,
std::vector< bool > &  g_is_update,
std::vector< SparseVectorType > &  M_v,
const spai_tag &  tag 
)

Solution of Least square problem on CPU.

Parameters:
A_v_ccolumn-major vectorized initial sparse matrix
g_Rblocks for least square solution
g_b_vvectors beta, necessary for Q recovery
g_Icontainer of row index set for all columns of matrix M
g_Jcontainer of column index set for all columns of matrix M
g_rescontainer of residuals
g_is_updatecontainer with indicators which blocks are active
M_vcolumn-major vectorized sparse matrix, final preconditioner
tagspai tag
void viennacl::linalg::detail::spai::make_rotation_matrix ( MatrixType &  mat,
std::size_t  new_size,
std::size_t  off_diagonal_distance = 4 
)
void viennacl::linalg::detail::spai::Print ( std::ostream &  ostr,
InputIterator  it_begin,
InputIterator  it_end 
)
void viennacl::linalg::detail::spai::print_continious_matrix ( VectorType &  con_A_I_J,
std::vector< cl_uint > &  blocks_ind,
const std::vector< std::vector< unsigned int > > &  g_I,
const std::vector< std::vector< unsigned int > > &  g_J 
)
void viennacl::linalg::detail::spai::print_continious_vector ( VectorType &  con_v,
std::vector< cl_uint > &  block_ind,
const std::vector< std::vector< unsigned int > > &  g_J 
)
void viennacl::linalg::detail::spai::print_matrix ( DenseMatrixType &  m)
void viennacl::linalg::detail::spai::print_sparse_vector ( const SparseVectorType &  v)
void viennacl::linalg::detail::spai::projectI ( const std::vector< unsigned int > &  I,
VectorType &  y,
unsigned int  ind 
)

Perform projection of set I on the unit-vector.

Parameters:
Iset of non-zero rows
yresult vector
indindex of unit vector
void viennacl::linalg::detail::spai::projectRows ( const std::vector< SparseVectorType > &  A_v_c,
const std::vector< unsigned int > &  J,
std::vector< unsigned int > &  I 
)

Row projection for matrix A(:,J) -> A(I,J), building index set of non-zero rows.

Parameters:
A_v_cinput matrix
Jset of non-zero rows
Ioutput matrix
void viennacl::linalg::detail::spai::QRBlockComposition ( const MatrixType &  A_I_J,
const MatrixType &  A_I_J_u,
MatrixType &  A_I_u_J_u 
)

Composition of new block for QR factorization cf. Kallischko dissertation p.82, figure 4.7.

Parameters:
A_I_Jpreviously composed block
A_I_J_umatrix Q'*A(I, \tilde J), where \tilde J - set of new column indices
A_I_u_J_uis composition of lower part A(I, \tilde J) and A(\tilde I, \tilde J) - new block for QR decomposition
void viennacl::linalg::detail::spai::single_qr ( MatrixType &  R,
VectorType &  b_v 
)

Inplace QR factorization via Householder reflections c.f. Gene H. Golub, Charles F. Van Loan "Matrix Computations" 3rd edition p.224.

Parameters:
Rinput matrix
b_vvector of betas
void viennacl::linalg::detail::spai::sparse_inner_prod ( const SparseVectorType &  v1,
const SparseVectorType &  v2,
ScalarType &  res_v 
)

Dot product of two sparse vectors.

Parameters:
v1initial sparse vector
v2initial sparse vector
res_vscalar that represents dot product result
void viennacl::linalg::detail::spai::sparse_norm_2 ( const SparseVectorType &  v,
ScalarType &  norm 
)

Computation of Euclidean norm for sparse vector.

Parameters:
vinitial sparse vector
normscalar that represents Euclidean norm
void viennacl::linalg::detail::spai::sparse_transpose ( const MatrixType &  A_in,
MatrixType &  A 
)

Transposition of sparse matrix.

Parameters:
A_inintial sparse matrix
Aoutput transposed matrix
void viennacl::linalg::detail::spai::store_householder_vector ( MatrixType &  A,
unsigned int  ind,
VectorType &  v 
)

Storage of vector v in column(A, ind), starting from ind-1 index of a column.

Parameters:
Ainit matrix
indindex of a column
vvector that should be stored
void viennacl::linalg::detail::spai::sym_sparse_matrix_to_stl ( MatrixType const &  A,
std::vector< std::map< unsigned int, ScalarType > > &  STL_A 
)
void viennacl::linalg::detail::spai::vectorize_column_matrix ( const SparseMatrixType &  M_in,
std::vector< SparseVectorType > &  M_v 
)

Solution of Least square problem on CPU.

Parameters:
M_ininput sparse, boost::numeric::ublas::compressed_matrix
M_varray of sparse vectors
void viennacl::linalg::detail::spai::vectorize_row_matrix ( const SparseMatrixType &  M_in,
std::vector< SparseVectorType > &  M_v 
)
void viennacl::linalg::detail::spai::write_set_to_array ( const std::vector< std::vector< SizeType > > &  ind_set,
std::vector< cl_uint > &  a 
)
void viennacl::linalg::detail::spai::write_to_block ( VectorType &  con_A_I_J,
unsigned int  start_ind,
const std::vector< unsigned int > &  I,
const std::vector< unsigned int > &  J,
MatrixType &  m 
)