ViennaCL - The Vienna Computing Library  1.6.1
Free open-source GPU-accelerated linear algebra and solver library.
amg_interpol.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
2 #define VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
3 
4 /* =========================================================================
5  Copyright (c) 2010-2014, 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 PDF manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
25 #include <boost/numeric/ublas/vector.hpp>
26 #include <cmath>
28 
29 #include <map>
30 #ifdef VIENNACL_WITH_OPENMP
31 #include <omp.h>
32 #endif
33 
35 
36 namespace viennacl
37 {
38 namespace linalg
39 {
40 namespace detail
41 {
42 namespace amg
43 {
44 
52 template<typename InternalT1, typename InternalT2>
53 void amg_interpol(unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector, amg_tag & tag)
54 {
55  switch (tag.get_interpol())
56  {
57  case VIENNACL_AMG_INTERPOL_DIRECT: amg_interpol_direct (level, A, P, pointvector, tag); break;
58  case VIENNACL_AMG_INTERPOL_CLASSIC: amg_interpol_classic(level, A, P, pointvector, tag); break;
59  case VIENNACL_AMG_INTERPOL_AG: amg_interpol_ag (level, A, P, pointvector, tag); break;
60  case VIENNACL_AMG_INTERPOL_SA: amg_interpol_sa (level, A, P, pointvector, tag); break;
61  }
62 }
63 
71 template<typename InternalT1, typename InternalT2>
72 void amg_interpol_direct(unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector, amg_tag & tag)
73 {
74  typedef typename InternalT1::value_type SparseMatrixType;
75  //typedef typename InternalType2::value_type PointVectorType;
76  typedef typename SparseMatrixType::value_type ScalarType;
77  typedef typename SparseMatrixType::iterator1 InternalRowIterator;
78  typedef typename SparseMatrixType::iterator2 InternalColIterator;
79 
80  ScalarType temp_res;
81  ScalarType row_sum, c_sum, diag;
82  //int diag_sign;
83  long x, y;
84  amg_point *pointx, *pointy;
85  unsigned int c_points = pointvector[level].get_cpoints();
86 
87  // Setup Prolongation/Interpolation matrix
88  P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points);
89  P[level].clear();
90 
91  // Assign indices to C points
92  pointvector[level].build_index();
93 
94  // Direct Interpolation (Yang, p.14)
95 #ifdef VIENNACL_WITH_OPENMP
96  #pragma omp parallel for private (pointx, pointy, row_sum, c_sum, temp_res, y, x, diag)
97 #endif
98  for (x=0; x < static_cast<long>(pointvector[level].size()); ++x)
99  {
100  pointx = pointvector[level][static_cast<unsigned int>(x)];
101  /*if (A[level](x,x) > 0)
102  diag_sign = 1;
103  else
104  diag_sign = -1;*/
105 
106  // When the current line corresponds to a C point then the diagonal coefficient is 1 and the rest 0
107  if (pointx->is_cpoint())
108  P[level](static_cast<unsigned int>(x),pointx->get_coarse_index()) = 1;
109 
110  // When the current line corresponds to a F point then the diagonal is 0 and the rest has to be computed (Yang, p.14)
111  if (pointx->is_fpoint())
112  {
113  // Jump to row x
114  InternalRowIterator row_iter = A[level].begin1();
115  row_iter += vcl_size_t(x);
116 
117  // Row sum of coefficients (without diagonal) and sum of influencing C point coefficients has to be computed
118  row_sum = c_sum = diag = 0;
119  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
120  {
121  y = static_cast<long>(col_iter.index2());
122  if (x == y)// || *col_iter * diag_sign > 0)
123  {
124  diag += *col_iter;
125  continue;
126  }
127 
128  // Sum all other coefficients in line x
129  row_sum += *col_iter;
130 
131  pointy = pointvector[level][static_cast<unsigned int>(y)];
132  // Sum all coefficients that correspond to a strongly influencing C point
133  if (pointy->is_cpoint())
134  if (pointx->is_influencing(pointy))
135  c_sum += *col_iter;
136  }
137  temp_res = -row_sum/(c_sum*diag);
138 
139  // Iterate over all strongly influencing points of point x
140  for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
141  {
142  pointy = *iter;
143  // The value is only non-zero for columns that correspond to a C point
144  if (pointy->is_cpoint())
145  {
146  if (temp_res > 0 || temp_res < 0)
147  P[level](static_cast<unsigned int>(x), pointy->get_coarse_index()) = temp_res * A[level](static_cast<unsigned int>(x),pointy->get_index());
148  }
149  }
150 
151  //Truncate interpolation if chosen
152  if (tag.get_interpolweight() > 0)
153  amg_truncate_row(P[level], static_cast<unsigned int>(x), tag);
154  }
155  }
156 
157  // P test
158  //test_interpolation(A[level], P[level], Pointvector[level]);
159 
160  #ifdef VIENNACL_AMG_DEBUG
161  std::cout << "Prolongation Matrix:" << std::endl;
162  printmatrix (P[level]);
163  #endif
164 }
165 
173 template<typename InternalT1, typename InternalT2>
174 void amg_interpol_classic(unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector, amg_tag & tag)
175 {
176  typedef typename InternalT1::value_type SparseMatrixType;
177  //typedef typename InternalType2::value_type PointVectorType;
178  typedef typename SparseMatrixType::value_type ScalarType;
179  typedef typename SparseMatrixType::iterator1 InternalRowIterator;
180  typedef typename SparseMatrixType::iterator2 InternalColIterator;
181 
182  ScalarType temp_res;
183  ScalarType weak_sum, strong_sum;
184  int diag_sign;
186  amg_point *pointx, *pointy, *pointk, *pointm;
187  long x, y, k, m;
188 
189  unsigned int c_points = pointvector[level].get_cpoints();
190 
191  // Setup Prolongation/Interpolation matrix
192  P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points);
193  P[level].clear();
194 
195  // Assign indices to C points
196  pointvector[level].build_index();
197 
198  // Classical Interpolation (Yang, p.13-14)
199 #ifdef VIENNACL_WITH_OPENMP
200  #pragma omp parallel for private (pointx,pointy,pointk,pointm,weak_sum,strong_sum,c_sum_row,temp_res,x,y,k,m,diag_sign)
201 #endif
202  for (x=0; x < static_cast<long>(pointvector[level].size()); ++x)
203  {
204  pointx = pointvector[level][static_cast<unsigned int>(x)];
205  if (A[level](static_cast<unsigned int>(x),static_cast<unsigned int>(x)) > 0)
206  diag_sign = 1;
207  else
208  diag_sign = -1;
209 
210  // When the current line corresponds to a C point then the diagonal coefficient is 1 and the rest 0
211  if (pointx->is_cpoint())
212  P[level](static_cast<unsigned int>(x),pointx->get_coarse_index()) = 1;
213 
214  // When the current line corresponds to a F point then the diagonal is 0 and the rest has to be computed (Yang, p.14)
215  if (pointx->is_fpoint())
216  {
217  // Jump to row x
218  InternalRowIterator row_iter = A[level].begin1();
219  row_iter += vcl_size_t(x);
220 
221  weak_sum = 0;
222  c_sum_row = amg_sparsevector<ScalarType>(static_cast<unsigned int>(A[level].size1()));
223  c_sum_row.clear();
224  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
225  {
226  k = static_cast<unsigned int>(col_iter.index2());
227  pointk = pointvector[level][static_cast<unsigned int>(k)];
228 
229  // Sum of weakly influencing neighbors + diagonal coefficient
230  if (x == k || !pointx->is_influencing(pointk))// || *col_iter * diag_sign > 0)
231  {
232  weak_sum += *col_iter;
233  continue;
234  }
235 
236  // Sums of coefficients in row k (strongly influening F neighbors) of C point neighbors of x are calculated
237  if (pointk->is_fpoint() && pointx->is_influencing(pointk))
238  {
239  for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
240  {
241  pointm = *iter;
242  m = pointm->get_index();
243 
244  if (pointm->is_cpoint())
245  // Only use coefficients that have opposite sign of diagonal.
246  if (A[level](static_cast<unsigned int>(k),static_cast<unsigned int>(m)) * ScalarType(diag_sign) < 0)
247  c_sum_row[static_cast<unsigned int>(k)] += A[level](static_cast<unsigned int>(k), static_cast<unsigned int>(m));
248  }
249  continue;
250  }
251  }
252 
253  // Iterate over all strongly influencing points of point x
254  for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
255  {
256  pointy = *iter;
257  y = pointy->get_index();
258 
259  // The value is only non-zero for columns that correspond to a C point
260  if (pointy->is_cpoint())
261  {
262  strong_sum = 0;
263  // Calculate term for strongly influencing F neighbors
264  for (typename amg_sparsevector<ScalarType>::iterator iter2 = c_sum_row.begin(); iter2 != c_sum_row.end(); ++iter2)
265  {
266  k = iter2.index();
267  // Only use coefficients that have opposite sign of diagonal.
268  if (A[level](static_cast<unsigned int>(k),static_cast<unsigned int>(y)) * ScalarType(diag_sign) < 0)
269  strong_sum += (A[level](static_cast<unsigned int>(x),static_cast<unsigned int>(k)) * A[level](static_cast<unsigned int>(k),static_cast<unsigned int>(y))) / (*iter2);
270  }
271 
272  // Calculate coefficient
273  temp_res = - (A[level](static_cast<unsigned int>(x),static_cast<unsigned int>(y)) + strong_sum) / (weak_sum);
274  if (temp_res < 0 || temp_res > 0)
275  P[level](static_cast<unsigned int>(x),pointy->get_coarse_index()) = temp_res;
276  }
277  }
278 
279  //Truncate iteration if chosen
280  if (tag.get_interpolweight() > 0)
281  amg_truncate_row(P[level], static_cast<unsigned int>(x), tag);
282  }
283  }
284 
285  #ifdef VIENNACL_AMG_DEBUG
286  std::cout << "Prolongation Matrix:" << std::endl;
287  printmatrix (P[level]);
288  #endif
289 }
290 
297 template<typename SparseMatrixT>
298 void amg_truncate_row(SparseMatrixT & P, unsigned int row, amg_tag & tag)
299 {
300  typedef typename SparseMatrixT::value_type ScalarType;
301  typedef typename SparseMatrixT::iterator1 InternalRowIterator;
302  typedef typename SparseMatrixT::iterator2 InternalColIterator;
303 
304  ScalarType row_max, row_min, row_sum_pos, row_sum_neg, row_sum_pos_scale, row_sum_neg_scale;
305 
306  InternalRowIterator row_iter = P.begin1();
307  row_iter += row;
308 
309  row_max = 0;
310  row_min = 0;
311  row_sum_pos = 0;
312  row_sum_neg = 0;
313 
314  // Truncate interpolation by making values to zero that are a lot smaller than the biggest value in a row
315  // Determine max entry and sum of row (seperately for negative and positive entries)
316  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
317  {
318  if (*col_iter > row_max)
319  row_max = *col_iter;
320  if (*col_iter < row_min)
321  row_min = *col_iter;
322  if (*col_iter > 0)
323  row_sum_pos += *col_iter;
324  if (*col_iter < 0)
325  row_sum_neg += *col_iter;
326  }
327 
328  row_sum_pos_scale = row_sum_pos;
329  row_sum_neg_scale = row_sum_neg;
330 
331  // Make certain values to zero (seperately for negative and positive entries)
332  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
333  {
334  if (*col_iter > 0 && *col_iter < tag.get_interpolweight() * row_max)
335  {
336  row_sum_pos_scale -= *col_iter;
337  *col_iter = 0;
338  }
339  if (*col_iter < 0 && *col_iter > tag.get_interpolweight() * row_min)
340  {
341  row_sum_pos_scale -= *col_iter;
342  *col_iter = 0;
343  }
344  }
345 
346  // Scale remaining values such that row sum is unchanged
347  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
348  {
349  if (*col_iter > 0)
350  *col_iter = *col_iter *(row_sum_pos/row_sum_pos_scale);
351  if (*col_iter < 0)
352  *col_iter = *col_iter *(row_sum_neg/row_sum_neg_scale);
353  }
354 }
355 
362 template<typename InternalT1, typename InternalT2>
363 void amg_interpol_ag(unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector, amg_tag)
364 {
365  typedef typename InternalT1::value_type SparseMatrixType;
366  //typedef typename InternalType2::value_type PointVectorType;
367  //typedef typename SparseMatrixType::value_type ScalarType;
368  //typedef typename SparseMatrixType::iterator1 InternalRowIterator;
369  //typedef typename SparseMatrixType::iterator2 InternalColIterator;
370 
371  long x;
372  amg_point *pointx, *pointy;
373  unsigned int c_points = pointvector[level].get_cpoints();
374 
375  P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points);
376  P[level].clear();
377 
378  // Assign indices to C points
379  pointvector[level].build_index();
380 
381  // Set prolongation such that F point is interpolated (weight=1) by the aggregate it belongs to (Vanek et al p.6)
382 #ifdef VIENNACL_WITH_OPENMP
383  #pragma omp parallel for private (x,pointx)
384 #endif
385  for (x=0; x<static_cast<long>(pointvector[level].size()); ++x)
386  {
387  pointx = pointvector[level][static_cast<unsigned int>(x)];
388  pointy = pointvector[level][pointx->get_aggregate()];
389  // Point x belongs to aggregate y.
390  P[level](static_cast<unsigned int>(x), pointy->get_coarse_index()) = 1;
391  }
392 
393  #ifdef VIENNACL_AMG_DEBUG
394  std::cout << "Aggregation based Prolongation:" << std::endl;
395  printmatrix(P[level]);
396  #endif
397 }
398 
406 template<typename InternalT1, typename InternalT2>
407 void amg_interpol_sa(unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector, amg_tag & tag)
408 {
409  typedef typename InternalT1::value_type SparseMatrixType;
410  //typedef typename InternalType2::value_type PointVectorType;
411  typedef typename SparseMatrixType::value_type ScalarType;
412  typedef typename SparseMatrixType::iterator1 InternalRowIterator;
413  typedef typename SparseMatrixType::iterator2 InternalColIterator;
414 
415  long x,y;
416  ScalarType diag = 0;
417  unsigned int c_points = pointvector[level].get_cpoints();
418 
419  InternalT1 P_tentative = InternalT1(P.size());
420  SparseMatrixType Jacobi = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), static_cast<unsigned int>(A[level].size2()));
421  Jacobi.clear();
422  P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points);
423  P[level].clear();
424 
425  // Build Jacobi Matrix via filtered A matrix (Vanek et al. p.6)
426 #ifdef VIENNACL_WITH_OPENMP
427  #pragma omp parallel for private (x,y,diag)
428 #endif
429  for (x=0; x<static_cast<long>(A[level].size1()); ++x)
430  {
431  diag = 0;
432  InternalRowIterator row_iter = A[level].begin1();
433  row_iter += vcl_size_t(x);
434  for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
435  {
436  y = static_cast<long>(col_iter.index2());
437  // Determine the structure of the Jacobi matrix by using a filtered matrix of A:
438  // The diagonal consists of the diagonal coefficient minus all coefficients of points not in the neighborhood of x.
439  // All other coefficients are the same as in A.
440  // Already use Jacobi matrix to save filtered A matrix to speed up computation.
441  if (x == y)
442  diag += *col_iter;
443  else if (!pointvector[level][static_cast<unsigned int>(x)]->is_influencing(pointvector[level][static_cast<unsigned int>(y)]))
444  diag += -*col_iter;
445  else
446  Jacobi (static_cast<unsigned int>(x), static_cast<unsigned int>(y)) = *col_iter;
447  }
448  InternalRowIterator row_iter2 = Jacobi.begin1();
449  row_iter2 += vcl_size_t(x);
450  // Traverse through filtered A matrix and compute the Jacobi filtering
451  for (InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2)
452  {
453  *col_iter2 = - static_cast<ScalarType>(tag.get_interpolweight())/diag * *col_iter2;
454  }
455  // Diagonal can be computed seperately.
456  Jacobi (static_cast<unsigned int>(x), static_cast<unsigned int>(x)) = 1 - static_cast<ScalarType>(tag.get_interpolweight());
457  }
458 
459  #ifdef VIENNACL_AMG_DEBUG
460  std::cout << "Jacobi Matrix:" << std::endl;
461  printmatrix(Jacobi);
462  #endif
463 
464  // Use AG interpolation as tentative prolongation
465  amg_interpol_ag(level, A, P_tentative, pointvector, tag);
466 
467  #ifdef VIENNACL_AMG_DEBUG
468  std::cout << "Tentative Prolongation:" << std::endl;
469  printmatrix(P_tentative[level]);
470  #endif
471 
472  // Multiply Jacobi matrix with tentative prolongation to get actual prolongation
473  amg_mat_prod(Jacobi,P_tentative[level],P[level]);
474 
475  #ifdef VIENNACL_AMG_DEBUG
476  std::cout << "Prolongation Matrix:" << std::endl;
477  printmatrix(P[level]);
478  #endif
479 }
480 
481 } //namespace amg
482 } //namespace detail
483 } //namespace linalg
484 } //namespace viennacl
485 
486 #endif
#define VIENNACL_AMG_INTERPOL_SA
Definition: amg_base.hpp:49
Debug functionality for AMG. To be removed.
#define VIENNACL_AMG_INTERPOL_AG
Definition: amg_base.hpp:48
#define VIENNACL_AMG_INTERPOL_DIRECT
Definition: amg_base.hpp:46
A class for the AMG points. Saves point index and influence measure Holds information whether point i...
Definition: amg_base.hpp:815
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:216
A class for the sparse vector type.
Definition: amg_base.hpp:254
void amg_interpol_ag(unsigned int level, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag)
AG (aggregation based) interpolation. Multi-Threaded! (VIENNACL_INTERPOL_SA)
unsigned int get_coarse_index() const
Definition: amg_base.hpp:876
void amg_interpol_sa(unsigned int level, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
SA (smoothed aggregate) interpolation. Multi-Threaded! (VIENNACL_INTERPOL_SA)
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
Definition: size.hpp:245
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:29
bool is_influencing(amg_point *point) const
Definition: amg_base.hpp:864
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
Definition: amg_base.hpp:61
Definition: blas3.hpp:36
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_interpol() const
Definition: amg_base.hpp:93
void amg_truncate_row(SparseMatrixT &P, unsigned int row, amg_tag &tag)
Interpolation truncation (for VIENNACL_AMG_INTERPOL_DIRECT and VIENNACL_AMG_INTERPOL_CLASSIC) ...
std::size_t vcl_size_t
Definition: forwards.h:74
vector_expression< const matrix_base< NumericT >, const int, op_matrix_diag > diag(const matrix_base< NumericT > &A, int k=0)
Definition: matrix.hpp:838
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
Definition: matrix.hpp:853
float ScalarType
Definition: fft_1d.cpp:42
#define VIENNACL_AMG_INTERPOL_CLASSIC
Definition: amg_base.hpp:47
Defines an iterator for the sparse vector type.
Definition: amg_base.hpp:203
void amg_mat_prod(SparseMatrixT &A, SparseMatrixT &B, SparseMatrixT &RES)
Sparse matrix product. Calculates RES = A*B.
Definition: amg_base.hpp:1319
Helper classes and functions for the AMG preconditioner. Experimental.
void printmatrix(MatrixT &, int)
Definition: amg_debug.hpp:77
void amg_interpol_direct(unsigned int level, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
Direct interpolation. Multi-threaded! (VIENNACL_AMG_INTERPOL_DIRECT)
void amg_interpol_classic(unsigned int level, InternalT1 &A, InternalT1 &P, InternalT2 &pointvector, amg_tag &tag)
Classical interpolation. Don't use with onepass classical coarsening or RS0 (Yang, p.14)! Multi-threaded! (VIENNACL_AMG_INTERPOL_CLASSIC)