1 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
2 #define VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
25 #include <boost/numeric/ublas/vector.hpp>
30 #ifdef VIENNACL_WITH_OPENMP
52 template<
typename InternalT1,
typename InternalT2>
53 void amg_interpol(
unsigned int level, InternalT1 & A, InternalT1 & P, InternalT2 & pointvector,
amg_tag & tag)
71 template<
typename InternalT1,
typename InternalT2>
74 typedef typename InternalT1::value_type SparseMatrixType;
76 typedef typename SparseMatrixType::value_type
ScalarType;
77 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
78 typedef typename SparseMatrixType::iterator2 InternalColIterator;
81 ScalarType row_sum, c_sum,
diag;
85 unsigned int c_points = pointvector[level].get_cpoints();
88 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].
size1()), c_points);
92 pointvector[level].build_index();
95 #ifdef VIENNACL_WITH_OPENMP
96 #pragma omp parallel for private (pointx, pointy, row_sum, c_sum, temp_res, y, x, diag)
98 for (x=0; x < static_cast<long>(pointvector[level].size()); ++x)
100 pointx = pointvector[level][
static_cast<unsigned int>(x)];
114 InternalRowIterator row_iter = A[level].begin1();
118 row_sum = c_sum = diag = 0;
119 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
121 y =
static_cast<long>(col_iter.index2());
129 row_sum += *col_iter;
131 pointy = pointvector[level][
static_cast<unsigned int>(y)];
137 temp_res = -row_sum/(c_sum*
diag);
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());
160 #ifdef VIENNACL_AMG_DEBUG
161 std::cout <<
"Prolongation Matrix:" << std::endl;
173 template<
typename InternalT1,
typename InternalT2>
176 typedef typename InternalT1::value_type SparseMatrixType;
178 typedef typename SparseMatrixType::value_type
ScalarType;
179 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
180 typedef typename SparseMatrixType::iterator2 InternalColIterator;
183 ScalarType weak_sum, strong_sum;
186 amg_point *pointx, *pointy, *pointk, *pointm;
189 unsigned int c_points = pointvector[level].get_cpoints();
192 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].
size1()), c_points);
196 pointvector[level].build_index();
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)
202 for (x=0; x < static_cast<long>(pointvector[level].size()); ++x)
204 pointx = pointvector[level][
static_cast<unsigned int>(x)];
205 if (A[level](static_cast<unsigned int>(x),
static_cast<unsigned int>(x)) > 0)
218 InternalRowIterator row_iter = A[level].begin1();
224 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
226 k =
static_cast<unsigned int>(col_iter.index2());
227 pointk = pointvector[level][
static_cast<unsigned int>(k)];
232 weak_sum += *col_iter;
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));
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);
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;
285 #ifdef VIENNACL_AMG_DEBUG
286 std::cout <<
"Prolongation Matrix:" << std::endl;
297 template<
typename SparseMatrixT>
300 typedef typename SparseMatrixT::value_type
ScalarType;
301 typedef typename SparseMatrixT::iterator1 InternalRowIterator;
302 typedef typename SparseMatrixT::iterator2 InternalColIterator;
304 ScalarType row_max, row_min, row_sum_pos, row_sum_neg, row_sum_pos_scale, row_sum_neg_scale;
306 InternalRowIterator row_iter = P.begin1();
316 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
318 if (*col_iter > row_max)
320 if (*col_iter < row_min)
323 row_sum_pos += *col_iter;
325 row_sum_neg += *col_iter;
328 row_sum_pos_scale = row_sum_pos;
329 row_sum_neg_scale = row_sum_neg;
332 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
336 row_sum_pos_scale -= *col_iter;
341 row_sum_pos_scale -= *col_iter;
347 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
350 *col_iter = *col_iter *(row_sum_pos/row_sum_pos_scale);
352 *col_iter = *col_iter *(row_sum_neg/row_sum_neg_scale);
362 template<
typename InternalT1,
typename InternalT2>
365 typedef typename InternalT1::value_type SparseMatrixType;
373 unsigned int c_points = pointvector[level].get_cpoints();
375 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].
size1()), c_points);
379 pointvector[level].build_index();
382 #ifdef VIENNACL_WITH_OPENMP
383 #pragma omp parallel for private (x,pointx)
385 for (x=0; x<static_cast<long>(pointvector[level].size()); ++x)
387 pointx = pointvector[level][
static_cast<unsigned int>(x)];
390 P[level](
static_cast<unsigned int>(x), pointy->get_coarse_index()) = 1;
393 #ifdef VIENNACL_AMG_DEBUG
394 std::cout <<
"Aggregation based Prolongation:" << std::endl;
406 template<
typename InternalT1,
typename InternalT2>
409 typedef typename InternalT1::value_type SparseMatrixType;
411 typedef typename SparseMatrixType::value_type
ScalarType;
412 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
413 typedef typename SparseMatrixType::iterator2 InternalColIterator;
417 unsigned int c_points = pointvector[level].get_cpoints();
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()));
422 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].
size1()), c_points);
426 #ifdef VIENNACL_WITH_OPENMP
427 #pragma omp parallel for private (x,y,diag)
429 for (x=0; x<static_cast<long>(A[level].size1()); ++x)
432 InternalRowIterator row_iter = A[level].begin1();
434 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
436 y =
static_cast<long>(col_iter.index2());
443 else if (!pointvector[level][static_cast<unsigned int>(x)]->is_influencing(pointvector[level][static_cast<unsigned int>(y)]))
446 Jacobi (static_cast<unsigned int>(x), static_cast<unsigned int>(y)) = *col_iter;
448 InternalRowIterator row_iter2 = Jacobi.begin1();
451 for (InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2)
456 Jacobi (static_cast<unsigned int>(x), static_cast<unsigned int>(x)) = 1 -
static_cast<ScalarType
>(tag.
get_interpolweight());
459 #ifdef VIENNACL_AMG_DEBUG
460 std::cout <<
"Jacobi Matrix:" << std::endl;
467 #ifdef VIENNACL_AMG_DEBUG
468 std::cout <<
"Tentative Prolongation:" << std::endl;
475 #ifdef VIENNACL_AMG_DEBUG
476 std::cout <<
"Prolongation Matrix:" << std::endl;
#define VIENNACL_AMG_INTERPOL_SA
Debug functionality for AMG. To be removed.
#define VIENNACL_AMG_INTERPOL_AG
#define VIENNACL_AMG_INTERPOL_DIRECT
A class for the AMG points. Saves point index and influence measure Holds information whether point i...
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
unsigned int get_aggregate()
A class for the sparse vector type.
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
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...
bool is_influencing(amg_point *point) const
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
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
iterator begin_influencing()
void amg_truncate_row(SparseMatrixT &P, unsigned int row, amg_tag &tag)
Interpolation truncation (for VIENNACL_AMG_INTERPOL_DIRECT and VIENNACL_AMG_INTERPOL_CLASSIC) ...
vector_expression< const matrix_base< NumericT >, const int, op_matrix_diag > diag(const matrix_base< NumericT > &A, int k=0)
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
double get_interpolweight() const
#define VIENNACL_AMG_INTERPOL_CLASSIC
Defines an iterator for the sparse vector type.
unsigned int get_index() const
void amg_mat_prod(SparseMatrixT &A, SparseMatrixT &B, SparseMatrixT &RES)
Sparse matrix product. Calculates RES = A*B.
Helper classes and functions for the AMG preconditioner. Experimental.
iterator end_influencing()
void printmatrix(MatrixT &, int)
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)