1 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_COARSE_HPP
2 #define VIENNACL_LINALG_DETAIL_AMG_AMG_COARSE_HPP
29 #ifdef VIENNACL_WITH_OPENMP
51 template<
typename InternalT1,
typename InternalT2,
typename InternalT3>
52 void amg_coarse(
unsigned int level, InternalT1 & A, InternalT2 & pointvector, InternalT3 & slicing,
amg_tag & tag)
70 template<
typename InternalT1,
typename InternalT2>
73 typedef typename InternalT1::value_type SparseMatrixType;
74 typedef typename InternalT2::value_type PointVectorType;
75 typedef typename SparseMatrixType::value_type
ScalarType;
76 typedef typename SparseMatrixType::const_iterator1 ConstRowIterator;
77 typedef typename SparseMatrixType::const_iterator2 ConstColIterator;
82 #ifdef VIENNACL_WITH_OPENMP
83 #pragma omp parallel for private (max,diag_sign)
85 for (
long i=0; i<static_cast<long>(A[level].size1()); ++i)
88 if (A[level](static_cast<unsigned int>(i),
static_cast<unsigned int>(i)) < 0)
91 ConstRowIterator row_iter = A[level].begin1();
92 row_iter +=
static_cast<unsigned int>(i);
95 for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
97 if (i == (
unsigned int) col_iter.index2())
continue;
99 if (max > *col_iter) max = *col_iter;
101 if (max < *col_iter) max = *col_iter;
105 if (std::fabs(max) <= 0)
109 for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
111 unsigned int j =
static_cast<unsigned int>(col_iter.index2());
119 pointvector[level][
static_cast<unsigned int>(i)]->add_influencing_point(pointvector[level][static_cast<unsigned int>(j)]);
124 #ifdef VIENNACL_AMG_DEBUG
125 std::cout <<
"Influence Matrix: " << std::endl;
126 boost::numeric::ublas::matrix<bool> mat;
127 Pointvector[level].get_influence_matrix(mat);
132 for (
typename PointVectorType::iterator iter = pointvector[level].begin(); iter != pointvector[level].end(); ++iter)
133 for (
typename amg_point::iterator iter2 = (*iter)->begin_influencing(); iter2 != (*iter)->end_influencing(); ++iter2)
134 (*iter2)->add_influenced_point(*iter);
136 #ifdef VIENNACL_AMG_DEBUG
137 std::cout <<
"Influence Measures: " << std::endl;
138 boost::numeric::ublas::vector<unsigned int> temp;
139 pointvector[level].get_influence(temp);
141 std::cout <<
"Point Sorting: " << std::endl;
142 Pointvector[level].get_sorting(temp);
155 template<
typename InternalT1,
typename InternalT2>
165 #ifdef VIENNACL_WITH_OPENMP
166 #pragma omp parallel for private (i)
168 for (i=0; i<static_cast<long>(pointvector[level].size()); ++i)
169 pointvector[level][static_cast<unsigned int>(i)]->calc_influence();
172 pointvector[level].sort();
175 while ((c_point = pointvector[level].get_nextpoint()) != NULL)
178 pointvector[level].make_cpoint(c_point);
195 pointvector[level].add_influence(point2,1);
217 #if defined (VIENNACL_AMG_DEBUG)// or defined (VIENNACL_AMG_DEBUGBENCH)
218 unsigned int c_points = pointvector[level].get_cpoints();
219 unsigned int f_points = pointvector[level].get_fpoints();
220 std::cout <<
"1st pass: Level " << level <<
": ";
221 std::cout <<
"No of C points = " << c_points <<
", ";
222 std::cout <<
"No of F points = " << f_points << std::endl;
225 #ifdef VIENNACL_AMG_DEBUG
226 std::cout <<
"Coarse Points:" << std::endl;
227 boost::numeric::ublas::vector<bool> C;
228 pointvector[level].get_C(C);
230 std::cout <<
"Fine Points:" << std::endl;
231 boost::numeric::ublas::vector<bool> F;
232 pointvector[level].get_F(F);
243 template<
typename InternalT1,
typename InternalT2>
246 typedef typename InternalT2::value_type PointVectorType;
255 for (
typename PointVectorType::iterator iter = pointvector[level].begin(); iter != pointvector[level].end(); ++iter)
281 if ((*iter2)->get_index() == (*iter3)->get_index())
287 else if ((*iter2)->get_index() < (*iter3)->get_index())
333 #ifdef VIENNACL_AMG_DEBUG
334 std::cout <<
"After 2nd pass:" << std::endl;
335 std::cout <<
"Coarse Points:" << std::endl;
336 boost::numeric::ublas::vector<bool> C;
337 pointvector[level].get_C(C);
339 std::cout <<
"Fine Points:" << std::endl;
340 boost::numeric::ublas::vector<bool> F;
341 pointvector[level].get_F(F);
345 #ifdef VIENNACL_AMG_DEBUG
346 #ifdef VIENNACL_WITH_OPENMP
350 std::cout <<
"No C and no F point: ";
351 for (
typename PointVectorType::iterator iter = pointvector[level].begin(); iter != pointvector[level].end(); ++iter)
352 if ((*iter)->is_undecided())
353 std::cout << (*iter)->get_index() <<
" ";
354 std::cout << std::endl;
366 template<
typename InternalT1,
typename InternalT2,
typename InternalT3>
367 void amg_coarse_rs0(
unsigned int level, InternalT1 & A, InternalT2 & pointvector, InternalT3 & slicing,
amg_tag & tag)
369 unsigned int total_points;
372 slicing.slice(level, A, pointvector);
376 #ifdef VIENNACL_WITH_OPENMP
377 #pragma omp parallel for
379 for (
long i2=0; i2<static_cast<long>(slicing.threads_); ++i2)
381 std::size_t i =
static_cast<std::size_t
>(i2);
386 slicing.offset_[i+1][level+1] = slicing.pointvector_slice_[i][level].get_cpoints();
387 #ifdef VIENNACL_WITH_OPENMP
390 total_points += slicing.pointvector_slice_[i][level].get_cpoints();
394 if (total_points != 0)
396 #ifdef VIENNACL_WITH_OPENMP
397 #pragma omp parallel for
399 for (
long i2=0; i2<static_cast<long>(slicing.threads_); ++i2)
401 std::size_t i =
static_cast<std::size_t
>(i2);
404 if (slicing.offset_[i+1][level+1] == 0)
407 for (
unsigned int j=0; j<slicing.A_slice_[i][level].size1(); ++j)
408 slicing.pointvector_slice_[i][level].make_cpoint(slicing.pointvector_slice_[i][level][j]);
409 slicing.offset_[i+1][level+1] =
static_cast<unsigned int>(slicing.A_slice_[i][level].size1());
414 for (
unsigned int i2=2; i2<=slicing.threads_; ++i2)
416 std::size_t i =
static_cast<std::size_t
>(i2);
417 slicing.offset_[i][level+1] += slicing.offset_[i-1][level+1];
421 slicing.join(level, pointvector);
427 #if defined(VIENNACL_AMG_DEBUG)// or defined (VIENNACL_AMG_DEBUGBENCH)
428 for (
unsigned int i=0; i<slicing._threads; ++i)
430 unsigned int c_points = slicing.pointvector_slice_[i][level].get_cpoints();
431 unsigned int f_points = slicing.pointvector_slice_[i][level].get_fpoints();
432 std::cout <<
"Thread " << i <<
": ";
433 std::cout <<
"No of C points = " << c_points <<
", ";
434 std::cout <<
"No of F points = " << f_points << std::endl;
446 template<
typename InternalT1,
typename InternalT2,
typename InternalT3>
447 void amg_coarse_rs3(
unsigned int level, InternalT1 & A, InternalT2 & pointvector, InternalT3 & slicing,
amg_tag & tag)
457 boost::numeric::ublas::vector<unsigned int> offset = boost::numeric::ublas::vector<unsigned int> (slicing.offset_.size());
458 for (i=0; i<slicing.offset_.size(); ++i)
459 offset[i] = slicing.offset_[i][level];
462 for (i=0; i<slicing.threads_; ++i)
465 for (j=offset[i]; j<offset[i+1]; ++j)
467 point1 = pointvector[level][j];
491 if ((*iter2)->get_index() == (*iter3)->get_index())
497 else if ((*iter2)->get_index() < (*iter3)->get_index())
546 for (
unsigned int k=i+1; k<=slicing.threads_; ++k)
547 slicing.offset_[k][level+1]++;
555 #ifdef VIENNACL_AMG_DEBUG
556 std::cout <<
"After 3rd pass:" << std::endl;
557 std::cout <<
"Coarse Points:" << std::endl;
558 boost::numeric::ublas::vector<bool> C;
559 pointvector[level].get_C(C);
561 std::cout <<
"Fine Points:" << std::endl;
562 boost::numeric::ublas::vector<bool> F;
563 pointvector[level].get_F(F);
575 template<
typename InternalT1,
typename InternalT2>
578 typedef typename InternalT1::value_type SparseMatrixType;
579 typedef typename InternalT2::value_type PointVectorType;
580 typedef typename SparseMatrixType::value_type
ScalarType;
582 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
583 typedef typename SparseMatrixType::iterator2 InternalColIterator;
590 if (A[level].
size1() == 1)
595 #ifdef VIENNACL_WITH_OPENMP
596 #pragma omp parallel for private (x,y,diag)
598 for (x=0; x<static_cast<long>(A[level].size1()); ++x)
600 InternalRowIterator row_iter = A[level].begin1();
601 row_iter += std::size_t(x);
602 diag = A[level](
static_cast<unsigned int>(x),static_cast<unsigned int>(x));
603 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
605 y =
static_cast<long>(col_iter.index2());
606 if (y == x || (std::fabs(*col_iter) >= tag.
get_threshold()*pow(0.5, static_cast<double>(level-1)) * std::sqrt(std::fabs(diag*A[level](static_cast<unsigned int>(y),static_cast<unsigned int>(y))))))
609 pointvector[level][
static_cast<unsigned int>(x)]->add_influencing_point(pointvector[level][static_cast<unsigned int>(y)]);
614 #ifdef VIENNACL_AMG_DEBUG
615 std::cout <<
"Neighborhoods:" << std::endl;
616 boost::numeric::ublas::matrix<bool> mat;
617 pointvector[level].get_influence_matrix(mat);
622 for (
typename PointVectorType::iterator iter = pointvector[level].begin(); iter != pointvector[level].end(); ++iter)
645 #ifdef VIENNACL_AMG_DEBUG
646 std::cout <<
"After aggregation:" << std::endl;
647 std::cout <<
"Coarse Points:" << std::endl;
648 boost::numeric::ublas::vector<bool> C;
649 pointvector[level].get_C(C);
651 std::cout <<
"Fine Points:" << std::endl;
652 boost::numeric::ublas::vector<bool> F;
653 pointvector[level].get_F(F);
655 std::cout <<
"Aggregates:" << std::endl;
void amg_coarse_rs3(unsigned int level, InternalT1 &A, InternalT2 &pointvector, InternalT3 &slicing, amg_tag &tag)
RS3 coarsening. Single-Threaded! (VIENNACL_AMG_COARSE_RS3)
Debug functionality for AMG. To be removed.
void set_aggregate(unsigned int aggregate)
double get_threshold() const
void amg_coarse_classic(unsigned int level, InternalT1 &A, InternalT2 &pointvector, amg_tag &tag)
Classical (RS) two-pass coarsening. Single-Threaded! (VIENNACL_AMG_COARSE_CLASSIC) ...
A class for the AMG points. Saves point index and influence measure Holds information whether point i...
#define VIENNACL_AMG_COARSE_RS
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
#define VIENNACL_AMG_COARSE_RS3
void amg_coarse_classic_onepass(unsigned int level, InternalT1 &A, InternalT2 &pointvector, amg_tag &tag)
Classical (RS) one-pass coarsening. Single-Threaded! (VIENNACL_AMG_COARSE_CLASSIC_ONEPASS) ...
T max(const T &lhs, const T &rhs)
Maximum.
#define VIENNACL_AMG_COARSE_RS0
#define VIENNACL_AMG_COARSE_AG
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
iterator begin_influenced()
bool is_influencing(amg_point *point) const
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
iterator begin_influencing()
vector_expression< const matrix_base< NumericT >, const int, op_matrix_diag > diag(const matrix_base< NumericT > &A, int k=0)
void printvector(VectorT const &)
bool is_undecided() const
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 amg_coarse_ag(unsigned int level, InternalT1 &A, InternalT2 &pointvector, amg_tag &tag)
AG (aggregation based) coarsening. Single-Threaded! (VIENNACL_AMG_COARSE_SA)
#define VIENNACL_AMG_COARSE_ONEPASS
void amg_coarse_rs0(unsigned int level, InternalT1 &A, InternalT2 &pointvector, InternalT3 &slicing, amg_tag &tag)
Parallel classical RS0 coarsening. Multi-Threaded! (VIENNACL_AMG_COARSE_RS0 || VIENNACL_AMG_COARSE_RS...
Defines an iterator for the sparse vector type.
iterator end_influenced()
unsigned int get_index() const
Helper classes and functions for the AMG preconditioner. Experimental.
iterator end_influencing()
void printmatrix(MatrixT &, int)
void amg_influence(unsigned int level, InternalT1 const &A, InternalT2 &pointvector, amg_tag &tag)
Determines strong influences in system matrix, classical approach (RS). Multithreaded! ...