Implementation of a dynamic SPAI. Provides the routines for automatic pattern updates Experimental in 1.2.x. More...
#include <utility>
#include <iostream>
#include <fstream>
#include <string>
#include <algorithm>
#include <vector>
#include <math.h>
#include <map>
#include "boost/numeric/ublas/vector.hpp"
#include "boost/numeric/ublas/matrix.hpp"
#include "boost/numeric/ublas/matrix_proxy.hpp"
#include "boost/numeric/ublas/vector_proxy.hpp"
#include "boost/numeric/ublas/storage.hpp"
#include "boost/numeric/ublas/io.hpp"
#include "boost/numeric/ublas/lu.hpp"
#include "boost/numeric/ublas/triangular.hpp"
#include "boost/numeric/ublas/matrix_expression.hpp"
#include "viennacl/linalg/prod.hpp"
#include "viennacl/matrix.hpp"
#include "viennacl/compressed_matrix.hpp"
#include "viennacl/linalg/compressed_matrix_operations.hpp"
#include "viennacl/linalg/matrix_operations.hpp"
#include "viennacl/scalar.hpp"
#include "viennacl/linalg/cg.hpp"
#include "viennacl/linalg/inner_prod.hpp"
#include "viennacl/linalg/ilu.hpp"
#include "viennacl/ocl/backend.hpp"
#include "viennacl/linalg/detail/spai/block_matrix.hpp"
#include "viennacl/linalg/detail/spai/block_vector.hpp"
#include "viennacl/linalg/detail/spai/qr.hpp"
#include "viennacl/linalg/detail/spai/spai_tag.hpp"
#include "viennacl/linalg/kernels/spai_source.h"
#include "viennacl/linalg/kernels/spai_kernels.h"
Go to the source code of this file.
Data Structures | |
struct | CompareSecond |
Namespaces | |
namespace | viennacl |
namespace | viennacl::linalg |
namespace | viennacl::linalg::detail |
namespace | viennacl::linalg::detail::spai |
Typedefs | |
typedef std::pair< unsigned int, double > | PairT |
Functions | |
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 Dense matrix from sparse one. | |
bool | isInIndexSet (const std::vector< unsigned int > &J, const unsigned int &ind) |
Determines if element ind is in set {J}. | |
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, const unsigned int cur_iter) |
Performs multiplication Q'*A(I, \tilde J) on GPU. | |
void | assemble_qr_row_inds (const std::vector< std::vector< unsigned int > > &g_I, const std::vector< std::vector< unsigned int > > g_J, const std::vector< std::vector< unsigned int > > &g_I_u, std::vector< std::vector< unsigned int > > &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, const unsigned int cur_iter) |
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, const unsigned int cur_iter) |
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, const unsigned int cur_iter) |
GPU-based block update. |
Implementation of a dynamic SPAI. Provides the routines for automatic pattern updates Experimental in 1.2.x.
SPAI code contributed by Nikolay Lukash