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>
179 void amg_transform_cpu(InternalT1 & A, InternalT1 & P, InternalT1 & R, InternalT2 & A_setup, InternalT2 & P_setup, amg_tag & tag)
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>
244 P_setup[i].set_trans(
true);
246 P_setup[i].set_trans(
false);
258 template<
typename InternalVectorT,
typename SparseMatrixT>
259 void amg_setup_apply(InternalVectorT & result, InternalVectorT & rhs, InternalVectorT & residual, SparseMatrixT
const & A, amg_tag
const & tag)
261 typedef typename InternalVectorT::value_type VectorType;
269 result[level] = VectorType(A[level].
size1());
270 result[level].clear();
271 rhs[level] = VectorType(A[level].
size1());
276 residual[level] = VectorType(A[level].
size1());
277 residual[level].clear();
291 template<
typename InternalVectorT,
typename SparseMatrixT>
294 typedef typename InternalVectorT::value_type VectorType;
302 result[level] = VectorType(A[level].
size1(), ctx);
303 rhs[level] = VectorType(A[level].
size1(), ctx);
307 residual[level] = VectorType(A[level].
size1(), ctx);
319 template<
typename NumericT,
typename SparseMatrixT>
320 void amg_lu(boost::numeric::ublas::compressed_matrix<NumericT> & op, boost::numeric::ublas::permutation_matrix<> & permutation, SparseMatrixT
const & A)
322 typedef typename SparseMatrixT::const_iterator1 ConstRowIterator;
323 typedef typename SparseMatrixT::const_iterator2 ConstColIterator;
326 op.resize(A.size1(),A.size2(),
false);
327 for (ConstRowIterator row_iter = A.begin1(); row_iter != A.end1(); ++row_iter)
328 for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
329 op (col_iter.index1(), col_iter.index2()) = *col_iter;
332 permutation = boost::numeric::ublas::permutation_matrix<> (op.size1());
338 template<
typename MatrixT>
341 typedef typename MatrixT::value_type NumericType;
342 typedef boost::numeric::ublas::vector<NumericType> VectorType;
351 boost::numeric::ublas::vector<SparseMatrixType> A_setup_;
352 boost::numeric::ublas::vector<SparseMatrixType> P_setup_;
353 boost::numeric::ublas::vector<MatrixT> A_;
354 boost::numeric::ublas::vector<MatrixT> P_;
355 boost::numeric::ublas::vector<MatrixT> R_;
356 boost::numeric::ublas::vector<PointVectorType> pointvector_;
358 mutable boost::numeric::ublas::compressed_matrix<NumericType> op_;
359 mutable boost::numeric::ublas::permutation_matrix<> permutation_;
361 mutable boost::numeric::ublas::vector<VectorType> result_;
362 mutable boost::numeric::ublas::vector<VectorType> rhs_;
363 mutable boost::numeric::ublas::vector<VectorType> residual_;
365 mutable bool done_init_apply_;
380 amg_init (mat, A_setup_, P_setup_, pointvector_, tag_);
382 done_init_apply_ =
false;
390 amg_setup(A_setup_, P_setup_, pointvector_, tag_);
394 done_init_apply_ =
false;
408 done_init_apply_ =
true;
416 template<
typename VectorT>
420 unsigned int nonzero=0, systemmat_nonzero=0, level_coefficients=0;
424 level_coefficients = 0;
425 for (InternalRowIterator row_iter = A_setup_[level].begin1(); row_iter != A_setup_[level].end1(); ++row_iter)
427 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
432 level_coefficients++;
435 avgstencil[level] =
static_cast<NumericType
>(level_coefficients)/static_cast<NumericType>(A_setup_[level].
size1());
437 return static_cast<NumericType
>(nonzero) / static_cast<NumericType>(systemmat_nonzero);
444 template<
typename VectorT>
448 if (!done_init_apply_)
457 result_[level].clear();
462 #ifdef VIENNACL_AMG_DEBUG
463 std::cout <<
"After presmooth:" << std::endl;
470 #ifdef VIENNACL_AMG_DEBUG
471 std::cout <<
"Residual:" << std::endl;
478 #ifdef VIENNACL_AMG_DEBUG
479 std::cout <<
"Restricted Residual: " << std::endl;
485 result_[level] = rhs_[level];
488 #ifdef VIENNACL_AMG_DEBUG
489 std::cout <<
"After direct solve: " << std::endl;
495 #ifdef VIENNACL_AMG_DEBUG
496 std::cout <<
"Coarse Error: " << std::endl;
503 #ifdef VIENNACL_AMG_DEBUG
504 std::cout <<
"Corrected Result: " << std::endl;
511 #ifdef VIENNACL_AMG_DEBUG
512 std::cout <<
"After postsmooth: " << std::endl;
525 template<
typename VectorT>
526 void smooth_jacobi(
int level,
int const iterations, VectorT & x, VectorT
const & rhs_smooth)
const
528 VectorT old_result(x.size());
532 for (
int i=0; i<iterations; ++i)
536 #ifdef VIENNACL_WITH_OPENMP
537 #pragma omp parallel for private (sum,diag) shared (rhs_smooth, x)
539 for (index=0; index < static_cast<long>(A_setup_[level].size1()); ++index)
541 InternalConstRowIterator row_iter = A_setup_[level].begin1();
545 for (InternalConstColIterator col_iter = row_iter.
begin(); col_iter != row_iter.
end(); ++col_iter)
547 if (col_iter.index1() == col_iter.index2())
550 sum += *col_iter * old_result[col_iter.index2()];
557 amg_tag &
tag() {
return tag_; }
564 template<
typename NumericT,
unsigned int AlignmentV>
577 boost::numeric::ublas::vector<SparseMatrixType> A_setup_;
578 boost::numeric::ublas::vector<SparseMatrixType> P_setup_;
579 boost::numeric::ublas::vector<MatrixType> A_;
580 boost::numeric::ublas::vector<MatrixType> P_;
581 boost::numeric::ublas::vector<MatrixType> R_;
582 boost::numeric::ublas::vector<PointVectorType> pointvector_;
584 mutable boost::numeric::ublas::compressed_matrix<NumericT> op_;
585 mutable boost::numeric::ublas::permutation_matrix<> permutation_;
587 mutable boost::numeric::ublas::vector<VectorType> result_;
588 mutable boost::numeric::ublas::vector<VectorType> rhs_;
589 mutable boost::numeric::ublas::vector<VectorType> residual_;
593 mutable bool done_init_apply_;
611 std::vector<std::map<unsigned int, NumericT> > mat2 = std::vector<std::map<unsigned int, NumericT> >(mat.
size1());
615 amg_init (mat2, A_setup_, P_setup_, pointvector_, tag_);
617 done_init_apply_ =
false;
625 amg_setup(A_setup_, P_setup_, pointvector_, tag_);
629 done_init_apply_ =
false;
643 done_init_apply_ =
true;
651 template<
typename VectorT>
655 unsigned int nonzero=0, systemmat_nonzero=0, level_coefficients=0;
659 level_coefficients = 0;
660 for (InternalRowIterator row_iter = A_setup_[level].begin1(); row_iter != A_setup_[level].end1(); ++row_iter)
662 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
667 level_coefficients++;
670 avgstencil[level] = level_coefficients/
static_cast<double>(A_[level].size1());
672 return nonzero/
static_cast<double>(systemmat_nonzero);
679 template<
typename VectorT>
682 if (!done_init_apply_)
691 result_[level].clear();
696 #ifdef VIENNACL_AMG_DEBUG
697 std::cout <<
"After presmooth: " << std::endl;
704 residual_[level] = rhs_[level] - residual_[level];
706 #ifdef VIENNACL_AMG_DEBUG
707 std::cout <<
"Residual: " << std::endl;
715 #ifdef VIENNACL_AMG_DEBUG
716 std::cout <<
"Restricted Residual: " << std::endl;
723 result_[level] = rhs_[level];
724 boost::numeric::ublas::vector<NumericT> result_cpu(result_[level].
size());
730 #ifdef VIENNACL_AMG_DEBUG
731 std::cout <<
"After direct solve: " << std::endl;
735 for (
int level2 = static_cast<int>(tag_.
get_coarselevels()-1); level2 >= 0; level2--)
739 #ifdef VIENNACL_AMG_DEBUG
740 std::cout <<
"Coarse Error: " << std::endl;
747 #ifdef VIENNACL_AMG_DEBUG
748 std::cout <<
"Corrected Result: " << std::endl;
755 #ifdef VIENNACL_AMG_DEBUG
756 std::cout <<
"After postsmooth: " << std::endl;
769 template<
typename VectorT>
772 VectorType old_result = x;
778 for (
unsigned int i=0; i<iterations; ++i)
785 viennacl::traits::opencl_handle(old_result),
786 viennacl::traits::opencl_handle(x),
787 viennacl::traits::opencl_handle(rhs_smooth),
788 static_cast<cl_uint
>(rhs_smooth.size())));
793 amg_tag &
tag() {
return tag_; }
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...
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
#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.
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.