ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
amg.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_AMG_HPP_
2 #define VIENNACL_LINALG_AMG_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2015, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
27 #include <vector>
28 #include <cmath>
29 #include "viennacl/forwards.h"
30 #include "viennacl/tools/tools.hpp"
31 #include "viennacl/linalg/prod.hpp"
34 
38 #include "viennacl/tools/timer.hpp"
40 #include "viennacl/linalg/lu.hpp"
41 
42 #include <map>
43 
44 #ifdef VIENNACL_WITH_OPENMP
45  #include <omp.h>
46 #endif
47 
48 #define VIENNACL_AMG_MAX_LEVELS 20
49 
50 namespace viennacl
51 {
52 namespace linalg
53 {
54 
55 class amg_coarse_problem_too_large_exception : public std::runtime_error
56 {
57 public:
58  amg_coarse_problem_too_large_exception(std::string const & msg, vcl_size_t num_points) : std::runtime_error(msg), c_points_(num_points) {}
59 
61  vcl_size_t coarse_points() const { return c_points_; }
62 
63 private:
64  vcl_size_t c_points_;
65 };
66 
67 
68 namespace detail
69 {
77  template<typename NumericT>
81  compressed_matrix<NumericT> & A_coarse)
82  {
83 
85 
86  // transpose P in memory (no known way of efficiently multiplying P^T * B for CSR-matrices P and B):
88 
89  // compute Galerkin product using a temporary for the result of A_fine * P
90  A_fine_times_P = viennacl::linalg::prod(A_fine, P);
91  A_coarse = viennacl::linalg::prod(R, A_fine_times_P);
92 
93  }
94 
95 
104  template<typename NumericT, typename AMGContextListT>
106  std::vector<compressed_matrix<NumericT> > & list_of_P,
107  std::vector<compressed_matrix<NumericT> > & list_of_R,
108  AMGContextListT & list_of_amg_level_context,
109  amg_tag & tag)
110  {
111  // Set number of iterations. If automatic coarse grid construction is chosen (0), then set a maximum size and stop during the process.
112  vcl_size_t iterations = tag.get_coarse_levels();
113  if (iterations == 0)
114  iterations = VIENNACL_AMG_MAX_LEVELS;
115 
116  for (vcl_size_t i=0; i<iterations; ++i)
117  {
118  list_of_amg_level_context[i].switch_context(tag.get_setup_context());
119  list_of_amg_level_context[i].resize(list_of_A[i].size1(), list_of_A[i].nnz());
120 
121  // Construct C and F points on coarse level (i is fine level, i+1 coarse level).
122  detail::amg::amg_coarse(list_of_A[i], list_of_amg_level_context[i], tag);
123 
124  // Calculate number of C and F points on level i.
125  unsigned int c_points = list_of_amg_level_context[i].num_coarse_;
126  unsigned int f_points = static_cast<unsigned int>(list_of_A[i].size1()) - c_points;
127 
128  if (f_points == 0 && c_points > tag.get_coarsening_cutoff())
129  {
130  std::stringstream ss;
131  ss << "No further coarsening possible (" << c_points << " coarse points). Consider changing the strong connection threshold or increasing the coarsening cutoff." << std::endl;
132  throw amg_coarse_problem_too_large_exception(ss.str(), c_points);
133  }
134 
135  // Stop routine when the maximal coarse level is found (no C or F point). Coarsest level is level i.
136  if (c_points == 0 || f_points == 0)
137  break;
138 
139  // Construct interpolation matrix for level i.
140  detail::amg::amg_interpol(list_of_A[i], list_of_P[i], list_of_amg_level_context[i], tag);
141 
142  // Compute coarse grid operator (A[i+1] = R * A[i] * P) with R = trans(P).
143  amg_galerkin_prod(list_of_A[i], list_of_P[i], list_of_R[i], list_of_A[i+1]);
144 
145  // send matrices to target context:
146  list_of_A[i].switch_memory_context(tag.get_target_context());
147  list_of_P[i].switch_memory_context(tag.get_target_context());
148  list_of_R[i].switch_memory_context(tag.get_target_context());
149 
150  // If Limit of coarse points is reached then stop. Coarsest level is level i+1.
151  if (tag.get_coarse_levels() == 0 && c_points <= tag.get_coarsening_cutoff())
152  return i+1;
153  }
154 
155  return iterations;
156  }
157 
158 
168  template<typename MatrixT, typename InternalT1, typename InternalT2>
169  void amg_init(MatrixT const & mat, InternalT1 & list_of_A, InternalT1 & list_of_P, InternalT1 & list_of_R, InternalT2 & list_of_amg_level_context, amg_tag & tag)
170  {
171  typedef typename InternalT1::value_type SparseMatrixType;
172 
173  vcl_size_t num_levels = (tag.get_coarse_levels() > 0) ? tag.get_coarse_levels() : VIENNACL_AMG_MAX_LEVELS;
174 
175  list_of_A.resize(num_levels+1, SparseMatrixType(tag.get_setup_context()));
176  list_of_P.resize(num_levels, SparseMatrixType(tag.get_setup_context()));
177  list_of_R.resize(num_levels, SparseMatrixType(tag.get_setup_context()));
178  list_of_amg_level_context.resize(num_levels);
179 
180  // Insert operator matrix as operator for finest level.
181  //SparseMatrixType A0(mat);
182  //A.insert_element(0, A0);
183  list_of_A[0].switch_memory_context(viennacl::traits::context(mat));
184  list_of_A[0] = mat;
185  list_of_A[0].switch_memory_context(tag.get_setup_context());
186  }
187 
198  template<typename InternalVectorT, typename SparseMatrixT>
199  void amg_setup_apply(InternalVectorT & result,
200  InternalVectorT & result_backup,
201  InternalVectorT & rhs,
202  InternalVectorT & residual,
203  SparseMatrixT const & A,
204  vcl_size_t coarse_levels,
205  amg_tag const & tag)
206  {
207  typedef typename InternalVectorT::value_type VectorType;
208 
209  result.resize(coarse_levels + 1);
210  result_backup.resize(coarse_levels + 1);
211  rhs.resize(coarse_levels + 1);
212  residual.resize(coarse_levels);
213 
214  for (vcl_size_t level=0; level <= coarse_levels; ++level)
215  {
216  result[level] = VectorType(A[level].size1(), tag.get_target_context());
217  result_backup[level] = VectorType(A[level].size1(), tag.get_target_context());
218  rhs[level] = VectorType(A[level].size1(), tag.get_target_context());
219  }
220  for (vcl_size_t level=0; level < coarse_levels; ++level)
221  {
222  residual[level] = VectorType(A[level].size1(), tag.get_target_context());
223  }
224  }
225 
226 
235  template<typename NumericT, typename SparseMatrixT>
237  SparseMatrixT const & A,
238  amg_tag const & tag)
239  {
241  op.resize(A.size1(), A.size2(), false);
243 
246  }
247 
248 }
249 
252 template<typename MatrixT>
254 
255 
260 template<typename NumericT, unsigned int AlignmentV>
262 {
266 
267 public:
268 
270 
277  amg_tag const & tag)
278  {
279  tag_ = tag;
280 
281  // Initialize data structures.
282  detail::amg_init(mat, A_list_, P_list_, R_list_, amg_context_list_, tag_);
283  }
284 
287  void setup()
288  {
289  // Start setup phase.
290  vcl_size_t num_coarse_levels = detail::amg_setup(A_list_, P_list_, R_list_, amg_context_list_, tag_);
291 
292  // Setup precondition phase (Data structures).
293  detail::amg_setup_apply(result_list_, result_backup_list_, rhs_list_, residual_list_, A_list_, num_coarse_levels, tag_);
294 
295  // LU factorization for direct solve.
296  detail::amg_lu(coarsest_op_, A_list_[num_coarse_levels], tag_);
297  }
298 
299 
304  template<typename VectorT>
305  void apply(VectorT & vec) const
306  {
307  vcl_size_t level;
308 
309  // Precondition operation (Yang, p.3).
310  rhs_list_[0] = vec;
311 
312  // Part 1: Restrict down to coarsest level
313  for (level=0; level < residual_list_.size(); level++)
314  {
315  result_list_[level].clear();
316 
317  // Apply Smoother presmooth_ times.
318  viennacl::linalg::detail::amg::smooth_jacobi(static_cast<unsigned int>(tag_.get_presmooth_steps()),
319  A_list_[level],
320  result_list_[level],
321  result_backup_list_[level],
322  rhs_list_[level],
323  static_cast<NumericT>(tag_.get_jacobi_weight()));
324 
325  // Compute residual.
326  //residual[level] = rhs_[level] - viennacl::linalg::prod(A_[level], result_[level]);
327  residual_list_[level] = viennacl::linalg::prod(A_list_[level], result_list_[level]);
328  residual_list_[level] = rhs_list_[level] - residual_list_[level];
329 
330  // Restrict to coarse level. Result is RHS of coarse level equation.
331  //residual_coarse[level] = viennacl::linalg::prod(R[level],residual[level]);
332  rhs_list_[level+1] = viennacl::linalg::prod(R_list_[level], residual_list_[level]);
333  }
334 
335  // Part 2: On highest level use direct solve to solve equation (on the CPU)
336  result_list_[level] = rhs_list_[level];
337  viennacl::linalg::lu_substitute(coarsest_op_, result_list_[level]);
338 
339  // Part 3: Prolongation to finest level
340  for (int level2 = static_cast<int>(residual_list_.size()-1); level2 >= 0; level2--)
341  {
342  level = static_cast<vcl_size_t>(level2);
343 
344  // Interpolate error to fine level and correct solution.
345  result_backup_list_[level] = viennacl::linalg::prod(P_list_[level], result_list_[level+1]);
346  result_list_[level] += result_backup_list_[level];
347 
348  // Apply Smoother postsmooth_ times.
349  viennacl::linalg::detail::amg::smooth_jacobi(static_cast<unsigned int>(tag_.get_postsmooth_steps()),
350  A_list_[level],
351  result_list_[level],
352  result_backup_list_[level],
353  rhs_list_[level],
354  static_cast<NumericT>(tag_.get_jacobi_weight()));
355  }
356  vec = result_list_[0];
357  }
358 
360  vcl_size_t levels() const { return residual_list_.size(); }
361 
362 
368  {
369  assert(level < levels() && bool("Level index out of bounds!"));
370  return residual_list_[level].size();
371  }
372 
374  amg_tag const & tag() const { return tag_; }
375 
376 private:
377  std::vector<SparseMatrixType> A_list_;
378  std::vector<SparseMatrixType> P_list_;
379  std::vector<SparseMatrixType> R_list_;
380  std::vector<AMGContextType> amg_context_list_;
381 
382  viennacl::matrix<NumericT> coarsest_op_;
383 
384  mutable std::vector<VectorType> result_list_;
385  mutable std::vector<VectorType> result_backup_list_;
386  mutable std::vector<VectorType> rhs_list_;
387  mutable std::vector<VectorType> residual_list_;
388 
389  amg_tag tag_;
390 };
391 
392 }
393 }
394 
395 
396 
397 #endif
398 
void apply(VectorT &vec) const
Precondition Operation.
Definition: amg.hpp:305
viennacl::context const & get_setup_context() const
Returns the ViennaCL context for the preconditioenr setup.
Definition: amg_base.hpp:142
void amg_interpol(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &P, AMGContextT &amg_context, amg_tag &tag)
vcl_size_t size(vcl_size_t level) const
Returns the problem/operator size at the respective multigrid level.
Definition: amg.hpp:367
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Various little tools used here and there in ViennaCL.
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:163
void lu_substitute(matrix< NumericT, F1, AlignmentV1 > const &A, matrix< NumericT, F2, AlignmentV2 > &B)
LU substitution for the system LU = rhs.
Definition: lu.hpp:201
void amg_lu(viennacl::matrix< NumericT > &op, SparseMatrixT const &A, amg_tag const &tag)
Pre-compute LU factorization for direct solve (ublas library).
Definition: amg.hpp:236
vcl_size_t levels() const
Returns the total number of multigrid levels in the hierarchy including the finest level...
Definition: amg.hpp:360
This file provides the forward declarations for the main types used within ViennaCL.
void amg_coarse(compressed_matrix< NumericT > const &A, AMGContextT &amg_context, amg_tag &tag)
vcl_size_t get_coarsening_cutoff() const
Returns the coarse grid size for which the recursive multigrid scheme is stopped and a direct solver ...
Definition: amg_base.hpp:133
void amg_setup_apply(InternalVectorT &result, InternalVectorT &result_backup, InternalVectorT &rhs, InternalVectorT &residual, SparseMatrixT const &A, vcl_size_t coarse_levels, amg_tag const &tag)
Setup data structures for precondition phase for later use on the GPU.
Definition: amg.hpp:199
float NumericT
Definition: bisect.cpp:40
void amg_init(MatrixT const &mat, InternalT1 &list_of_A, InternalT1 &list_of_P, InternalT1 &list_of_R, InternalT2 &list_of_amg_level_context, amg_tag &tag)
Initialize AMG preconditioner.
Definition: amg.hpp:169
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:34
#define VIENNACL_AMG_MAX_LEVELS
Definition: amg.hpp:48
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
viennacl::context const & get_target_context() const
Returns the ViennaCL context for the solver cycle stage (i.e. preconditioner applications).
Definition: amg_base.hpp:150
void switch_memory_context(viennacl::context new_ctx)
Definition: matrix_def.hpp:249
viennacl::enable_if< viennacl::is_any_sparse_matrix< SparseMatrixType >::value >::type assign_to_dense(SparseMatrixType const &A, viennacl::matrix_base< NumericT > &B)
vcl_size_t coarse_points() const
Returns the number of coarse points for which no further coarsening could be applied.
Definition: amg.hpp:61
void smooth_jacobi(unsigned int iterations, compressed_matrix< NumericT > const &A, vector< NumericT > &x, vector< NumericT > &x_backup, vector< NumericT > const &rhs_smooth, NumericT weight)
Definition: blas3.hpp:36
Implementation of the compressed_matrix class.
Implementations of operations using sparse matrices.
AMG preconditioner class, can be supplied to solve()-routines.
Definition: amg.hpp:253
std::size_t vcl_size_t
Definition: forwards.h:75
void amg_transpose(compressed_matrix< NumericT > &A, compressed_matrix< NumericT > &B)
amg_tag const & tag() const
Returns the associated preconditioner tag containing the configuration for the multigrid precondition...
Definition: amg.hpp:374
Implementations of LU factorization for row-major and column-major dense matrices.
Implementations of dense direct solvers are found here.
amg_precond(compressed_matrix< NumericT, AlignmentV > const &mat, amg_tag const &tag)
The constructor. Builds data structures.
Definition: amg.hpp:276
vcl_size_t get_coarse_levels() const
Returns the number of coarse levels. If zero, then coarse levels are constructed until the cutoff siz...
Definition: amg_base.hpp:128
A simple, yet (mostly) sufficiently accurate timer for benchmarking and profiling.
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
Definition: amg_base.hpp:64
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
void resize(size_type rows, size_type columns, bool preserve=true)
Resizes the matrix. Existing entries can optionally be preserved.
Definition: matrix.hpp:813
void setup()
Start setup phase for this class and copy data structures.
Definition: amg.hpp:287
vcl_size_t amg_setup(std::vector< compressed_matrix< NumericT > > &list_of_A, std::vector< compressed_matrix< NumericT > > &list_of_P, std::vector< compressed_matrix< NumericT > > &list_of_R, AMGContextListT &list_of_amg_level_context, amg_tag &tag)
Setup AMG preconditioner.
Definition: amg.hpp:105
Helper classes and functions for the AMG preconditioner. Experimental.
amg_coarse_problem_too_large_exception(std::string const &msg, vcl_size_t num_points)
Definition: amg.hpp:58
Implementations of operations for algebraic multigrid.
void amg_galerkin_prod(compressed_matrix< NumericT > &A_fine, compressed_matrix< NumericT > &P, compressed_matrix< NumericT > &R, compressed_matrix< NumericT > &A_coarse)
Sparse Galerkin product: Calculates A_coarse = trans(P)*A_fine*P = R*A_fine*P.
Definition: amg.hpp:78
void lu_factorize(matrix< NumericT, viennacl::row_major > &A)
LU factorization of a row-major dense matrix.
Definition: lu.hpp:42