1 #ifndef VIENNACL_LINALG_AMG_HPP_
2 #define VIENNACL_LINALG_AMG_HPP_
27 #include <boost/numeric/ublas/matrix.hpp>
28 #include <boost/numeric/ublas/lu.hpp>
29 #include <boost/numeric/ublas/operation.hpp>
30 #include <boost/numeric/ublas/vector_proxy.hpp>
31 #include <boost/numeric/ublas/matrix_proxy.hpp>
32 #include <boost/numeric/ublas/vector.hpp>
33 #include <boost/numeric/ublas/triangular.hpp>
47 #ifdef VIENNACL_WITH_OPENMP
53 #define VIENNACL_AMG_COARSE_LIMIT 50
54 #define VIENNACL_AMG_MAX_LEVELS 100
71 template<
typename InternalT1,
typename InternalT2>
72 void amg_setup(InternalT1 & A, InternalT1 & P, InternalT2 & pointvector,
amg_tag & tag)
74 typedef typename InternalT2::value_type PointVectorType;
76 unsigned int i, iterations, c_points, f_points;
86 slicing.
init(iterations);
88 for (i=0; i<iterations; ++i)
91 pointvector[i] = PointVectorType(static_cast<unsigned int>(A[i].
size1()));
92 pointvector[i].init_points();
98 c_points = pointvector[i].get_cpoints();
99 f_points = pointvector[i].get_fpoints();
101 #if defined (VIENNACL_AMG_DEBUG) //or defined(VIENNACL_AMG_DEBUGBENCH)
102 std::cout <<
"Level " << i <<
": ";
103 std::cout <<
"No of C points = " << c_points <<
", ";
104 std::cout <<
"No of F points = " << f_points << std::endl;
108 if (c_points == 0 || f_points == 0)
120 pointvector[i].delete_points();
122 #ifdef VIENNACL_AMG_DEBUG
123 std::cout <<
"Coarse Grid Operator Matrix:" << std::endl;
145 template<
typename MatrixT,
typename InternalT1,
typename InternalT2>
146 void amg_init(MatrixT
const & mat, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector,
amg_tag & tag)
149 typedef typename InternalT1::value_type SparseMatrixType;
165 SparseMatrixType A0(mat);
166 A.insert_element(0, A0);
178 template<
typename InternalT1,
typename InternalT2>
191 A[i].resize(A_setup[i].
size1(),A_setup[i].
size2(),
false);
196 P[i].resize(P_setup[i].
size1(),P_setup[i].
size2(),
false);
201 R[i].resize(P_setup[i].
size2(),P_setup[i].
size1(),
false);
202 P_setup[i].set_trans(
true);
204 P_setup[i].set_trans(
false);
218 template<
typename InternalT1,
typename InternalT2>
221 typedef typename InternalT2::value_type::value_type NumericType;
239 for (
vcl_size_t j=0; j<P_setup[i].get_internal_pointer()->size(); ++j)
240 nonzeros += (*P_setup[i].get_internal_pointer())[j].
size();
249 P_setup[i].set_trans(
true);
251 P_setup[i].set_trans(
false);
263 template<
typename InternalVectorT,
typename SparseMatrixT>
264 void amg_setup_apply(InternalVectorT & result, InternalVectorT & rhs, InternalVectorT & residual, SparseMatrixT
const & A,
amg_tag const & tag)
266 typedef typename InternalVectorT::value_type VectorType;
274 result[level] = VectorType(A[level].
size1());
275 result[level].clear();
276 rhs[level] = VectorType(A[level].
size1());
281 residual[level] = VectorType(A[level].
size1());
282 residual[level].clear();
296 template<
typename InternalVectorT,
typename SparseMatrixT>
299 typedef typename InternalVectorT::value_type VectorType;
307 result[level] = VectorType(A[level].
size1(), ctx);
308 rhs[level] = VectorType(A[level].
size1(), ctx);
312 residual[level] = VectorType(A[level].
size1(), ctx);
324 template<
typename NumericT,
typename SparseMatrixT>
325 void amg_lu(boost::numeric::ublas::compressed_matrix<NumericT> & op, boost::numeric::ublas::permutation_matrix<> & permutation, SparseMatrixT
const & A)
327 typedef typename SparseMatrixT::const_iterator1 ConstRowIterator;
328 typedef typename SparseMatrixT::const_iterator2 ConstColIterator;
331 op.resize(A.size1(),A.size2(),
false);
332 for (ConstRowIterator row_iter = A.begin1(); row_iter != A.end1(); ++row_iter)
333 for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
334 op (col_iter.index1(), col_iter.index2()) = *col_iter;
337 permutation = boost::numeric::ublas::permutation_matrix<> (op.size1());
343 template<
typename MatrixT>
346 typedef typename MatrixT::value_type NumericType;
347 typedef boost::numeric::ublas::vector<NumericType> VectorType;
356 boost::numeric::ublas::vector<SparseMatrixType> A_setup_;
357 boost::numeric::ublas::vector<SparseMatrixType> P_setup_;
358 boost::numeric::ublas::vector<MatrixT> A_;
359 boost::numeric::ublas::vector<MatrixT> P_;
360 boost::numeric::ublas::vector<MatrixT> R_;
361 boost::numeric::ublas::vector<PointVectorType> pointvector_;
363 mutable boost::numeric::ublas::compressed_matrix<NumericType> op_;
364 mutable boost::numeric::ublas::permutation_matrix<> permutation_;
366 mutable boost::numeric::ublas::vector<VectorType> result_;
367 mutable boost::numeric::ublas::vector<VectorType> rhs_;
368 mutable boost::numeric::ublas::vector<VectorType> residual_;
370 mutable bool done_init_apply_;
385 amg_init (mat, A_setup_, P_setup_, pointvector_, tag_);
387 done_init_apply_ =
false;
395 amg_setup(A_setup_, P_setup_, pointvector_, tag_);
399 done_init_apply_ =
false;
413 done_init_apply_ =
true;
421 template<
typename VectorT>
425 unsigned int nonzero=0, systemmat_nonzero=0, level_coefficients=0;
429 level_coefficients = 0;
430 for (
InternalRowIterator row_iter = A_setup_[level].begin1(); row_iter != A_setup_[level].end1(); ++row_iter)
437 level_coefficients++;
440 avgstencil[level] =
static_cast<NumericType
>(level_coefficients)/static_cast<NumericType>(A_setup_[level].
size1());
442 return static_cast<NumericType
>(nonzero) / static_cast<NumericType>(systemmat_nonzero);
449 template<
typename VectorT>
453 if (!done_init_apply_)
462 result_[level].clear();
467 #ifdef VIENNACL_AMG_DEBUG
468 std::cout <<
"After presmooth:" << std::endl;
475 #ifdef VIENNACL_AMG_DEBUG
476 std::cout <<
"Residual:" << std::endl;
483 #ifdef VIENNACL_AMG_DEBUG
484 std::cout <<
"Restricted Residual: " << std::endl;
490 result_[level] = rhs_[level];
493 #ifdef VIENNACL_AMG_DEBUG
494 std::cout <<
"After direct solve: " << std::endl;
500 #ifdef VIENNACL_AMG_DEBUG
501 std::cout <<
"Coarse Error: " << std::endl;
508 #ifdef VIENNACL_AMG_DEBUG
509 std::cout <<
"Corrected Result: " << std::endl;
516 #ifdef VIENNACL_AMG_DEBUG
517 std::cout <<
"After postsmooth: " << std::endl;
530 template<
typename VectorT>
531 void smooth_jacobi(
int level,
int const iterations, VectorT & x, VectorT
const & rhs_smooth)
const
533 VectorT old_result(x.size());
536 for (
int i=0; i<iterations; ++i)
540 #ifdef VIENNACL_WITH_OPENMP
541 #pragma omp parallel for
543 for (index=0; index < static_cast<long>(A_setup_[level].size1()); ++index)
548 NumericType
diag = 1;
551 if (col_iter.index1() == col_iter.index2())
554 sum += *col_iter * old_result[col_iter.index2()];
556 x[index]=
static_cast<NumericType
>(tag_.
get_jacobiweight()) * (rhs_smooth[index] - sum) / diag + (1-
static_cast<NumericType
>(tag_.
get_jacobiweight())) * old_result[index];
568 template<
typename NumericT,
unsigned int AlignmentV>
581 boost::numeric::ublas::vector<SparseMatrixType> A_setup_;
582 boost::numeric::ublas::vector<SparseMatrixType> P_setup_;
583 boost::numeric::ublas::vector<MatrixType> A_;
584 boost::numeric::ublas::vector<MatrixType> P_;
585 boost::numeric::ublas::vector<MatrixType> R_;
586 boost::numeric::ublas::vector<PointVectorType> pointvector_;
588 mutable boost::numeric::ublas::compressed_matrix<NumericT> op_;
589 mutable boost::numeric::ublas::permutation_matrix<> permutation_;
591 mutable boost::numeric::ublas::vector<VectorType> result_;
592 mutable boost::numeric::ublas::vector<VectorType> rhs_;
593 mutable boost::numeric::ublas::vector<VectorType> residual_;
597 mutable bool done_init_apply_;
615 std::vector<std::map<unsigned int, NumericT> > mat2 = std::vector<std::map<unsigned int, NumericT> >(mat.
size1());
619 amg_init (mat2, A_setup_, P_setup_, pointvector_, tag_);
621 done_init_apply_ =
false;
629 amg_setup(A_setup_, P_setup_, pointvector_, tag_);
633 done_init_apply_ =
false;
647 done_init_apply_ =
true;
655 template<
typename VectorT>
659 unsigned int nonzero=0, systemmat_nonzero=0, level_coefficients=0;
663 level_coefficients = 0;
664 for (
InternalRowIterator row_iter = A_setup_[level].begin1(); row_iter != A_setup_[level].end1(); ++row_iter)
671 level_coefficients++;
674 avgstencil[level] = level_coefficients/
static_cast<double>(A_[level].size1());
676 return nonzero/
static_cast<double>(systemmat_nonzero);
683 template<
typename VectorT>
686 if (!done_init_apply_)
695 result_[level].clear();
700 #ifdef VIENNACL_AMG_DEBUG
701 std::cout <<
"After presmooth: " << std::endl;
708 residual_[level] = rhs_[level] - residual_[level];
710 #ifdef VIENNACL_AMG_DEBUG
711 std::cout <<
"Residual: " << std::endl;
719 #ifdef VIENNACL_AMG_DEBUG
720 std::cout <<
"Restricted Residual: " << std::endl;
727 result_[level] = rhs_[level];
728 boost::numeric::ublas::vector<NumericT> result_cpu(result_[level].
size());
734 #ifdef VIENNACL_AMG_DEBUG
735 std::cout <<
"After direct solve: " << std::endl;
739 for (
int level2 = static_cast<int>(tag_.
get_coarselevels()-1); level2 >= 0; level2--)
743 #ifdef VIENNACL_AMG_DEBUG
744 std::cout <<
"Coarse Error: " << std::endl;
751 #ifdef VIENNACL_AMG_DEBUG
752 std::cout <<
"Corrected Result: " << std::endl;
759 #ifdef VIENNACL_AMG_DEBUG
760 std::cout <<
"After postsmooth: " << std::endl;
773 template<
typename VectorT>
782 for (
unsigned int i=0; i<iterations; ++i)
789 viennacl::traits::opencl_handle(old_result),
790 viennacl::traits::opencl_handle(x),
791 viennacl::traits::opencl_handle(rhs_smooth),
792 static_cast<cl_uint
>(rhs_smooth.size())));
void apply(VectorT &vec) const
Precondition Operation.
Debug functionality for AMG. To be removed.
unsigned int get_coarselevels() const
double get_jacobiweight() const
void amg_lu(boost::numeric::ublas::compressed_matrix< NumericT > &op, boost::numeric::ublas::permutation_matrix<> &permutation, SparseMatrixT const &A)
Pre-compute LU factorization for direct solve (ublas library).
const vcl_size_t & size1() const
Returns the number of rows.
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Represents an OpenCL kernel within ViennaCL.
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
void lu_substitute(matrix< NumericT, F1, AlignmentV1 > const &A, matrix< NumericT, F2, AlignmentV2 > &B)
LU substitution for the system LU = rhs.
void init_apply() const
Prepare data structures for preconditioning: Build data structures for precondition phase...
#define VIENNACL_AMG_COARSE_RS3
NumericT calc_complexity(VectorT &avgstencil)
Returns complexity measures.
Manages an OpenCL context and provides the respective convenience functions for creating buffers...
This file provides the forward declarations for the main types used within ViennaCL.
unsigned int get_postsmooth() const
#define VIENNACL_AMG_COARSE_LIMIT
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
void amg_setup(InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
Setup AMG preconditioner.
statement sum(scalar< NumericT > const *s, vector_base< NumericT > const *x)
#define VIENNACL_AMG_COARSE_RS0
void smooth_jacobi(int level, int const iterations, VectorT &x, VectorT const &rhs_smooth) const
(Weighted) Jacobi Smoother (CPU version)
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
void init_apply() const
Prepare data structures for preconditioning: Build data structures for precondition phase...
#define VIENNACL_AMG_MAX_LEVELS
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
Main kernel class for generating OpenCL kernels for compressed_matrix.
viennacl::ocl::kernel & get_kernel(std::string const &program_name, std::string const &kernel_name)
Convenience function for retrieving the kernel of a program directly from the context.
Implementations of several variants of the AMG coarsening procedure (setup phase). Experimental.
void init(unsigned int levels, unsigned int threads=0)
void amg_interpol(unsigned int level, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
Calls the right function to build interpolation matrix.
unsigned int get_presmooth() const
AMG preconditioner class, can be supplied to solve()-routines.
void amg_galerkin_prod(SparseMatrixT &A, SparseMatrixT &P, SparseMatrixT &RES)
Sparse Galerkin product: Calculates RES = trans(P)*A*P.
void prod(const MatrixT1 &A, bool transposed_A, const MatrixT2 &B, bool transposed_B, MatrixT3 &C, ScalarT alpha, ScalarT beta)
vector_expression< const matrix_base< NumericT >, const int, op_matrix_diag > diag(const matrix_base< NumericT > &A, int k=0)
void setup()
Start setup phase for this class and copy data structures.
void printvector(VectorT const &)
A class for the sparse matrix type. Uses vector of maps as data structure for higher performance and ...
Implementations of dense direct solvers are found here.
void amg_init(MatrixT const &mat, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
Initialize AMG preconditioner.
void copy_impl(const CPUMatrixT &cpu_matrix, compressed_compressed_matrix< NumericT > &gpu_matrix, vcl_size_t nonzero_rows, vcl_size_t nonzeros)
amg_precond(compressed_matrix< NumericT, AlignmentV > const &mat, amg_tag const &tag)
The constructor. Builds data structures.
void smooth_jacobi(vcl_size_t level, unsigned int iterations, VectorT &x, VectorT const &rhs_smooth) const
Jacobi Smoother (GPU version)
NumericType calc_complexity(VectorT &avgstencil)
Returns complexity measures.
void amg_coarse(unsigned int level, InternalT1 &A, InternalT2 &pointvector, InternalT3 &slicing, amg_tag &tag)
Calls the right coarsening procedure.
unsigned int get_coarse() const
void clear()
Resets all entries to zero. Does not change the size of the vector.
void enqueue(KernelType &k, viennacl::ocl::command_queue const &queue)
Enqueues a kernel in the provided queue.
static void init(viennacl::ocl::context &ctx)
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
void amg_transform_cpu(InternalT1 &A, InternalT1 &P, InternalT1 &R, InternalT2 &A_setup, InternalT2 &P_setup, amg_tag &tag)
Save operators after setup phase for CPU computation.
void amg_transform_gpu(InternalT1 &A, InternalT1 &P, InternalT1 &R, InternalT2 &A_setup, InternalT2 &P_setup, amg_tag &tag, viennacl::context ctx)
Save operators after setup phase for GPU computation.
void setup()
Start setup phase for this class and copy data structures.
A sparse square matrix in compressed sparse rows format.
Helper classes and functions for the AMG preconditioner. Experimental.
amg_precond(MatrixT const &mat, amg_tag const &tag)
The constructor. Saves system matrix, tag and builds data structures for setup.
void printmatrix(MatrixT &, int)
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
detail::amg::amg_tag amg_tag
void apply(VectorT &vec) const
Precondition Operation.
void set_coarselevels(unsigned int coarselevels)
A class for the matrix slicing for parallel coarsening schemes (RS0/RS3).
void lu_factorize(matrix< NumericT, viennacl::row_major > &A)
LU factorization of a row-major dense matrix.
A class for the AMG points. Holds pointers of type amg_point in a vector that can be accessed using [...
void amg_setup_apply(InternalVectorT &result, InternalVectorT &rhs, InternalVectorT &residual, SparseMatrixT const &A, amg_tag const &tag)
Setup data structures for precondition phase.
Implementations of several variants of the AMG interpolation operators (setup phase). Experimental.
void switch_memory_context(T &obj, viennacl::context new_ctx)
Generic convenience routine for migrating data of an object to a new memory domain.