ViennaCL - The Vienna Computing Library  1.5.1
viennacl/linalg/detail/amg/amg_coarse.hpp
Go to the documentation of this file.
00001 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_COARSE_HPP
00002 #define VIENNACL_LINALG_DETAIL_AMG_AMG_COARSE_HPP
00003 
00004 /* =========================================================================
00005    Copyright (c) 2010-2014, Institute for Microelectronics,
00006                             Institute for Analysis and Scientific Computing,
00007                             TU Wien.
00008    Portions of this software are copyright by UChicago Argonne, LLC.
00009 
00010                             -----------------
00011                   ViennaCL - The Vienna Computing Library
00012                             -----------------
00013 
00014    Project Head:    Karl Rupp                   rupp@iue.tuwien.ac.at
00015 
00016    (A list of authors and contributors can be found in the PDF manual)
00017 
00018    License:         MIT (X11), see file LICENSE in the base directory
00019 ============================================================================= */
00020 
00025 #include <cmath>
00026 #include "viennacl/linalg/detail/amg/amg_base.hpp"
00027 
00028 #include <map>
00029 #ifdef VIENNACL_WITH_OPENMP
00030 #include <omp.h>
00031 #endif
00032 
00033 #include "viennacl/linalg/detail/amg/amg_debug.hpp"
00034 
00035 namespace viennacl
00036 {
00037   namespace linalg
00038   {
00039     namespace detail
00040     {
00041       namespace amg
00042       {
00043 
00051     template <typename InternalType1, typename InternalType2, typename InternalType3>
00052     void amg_coarse(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, InternalType3 & Slicing, amg_tag & tag)
00053     {
00054       switch (tag.get_coarse())
00055       {
00056         case VIENNACL_AMG_COARSE_RS: amg_coarse_classic (level, A, Pointvector, tag); break;
00057         case VIENNACL_AMG_COARSE_ONEPASS: amg_coarse_classic_onepass (level, A, Pointvector, tag); break;
00058         case VIENNACL_AMG_COARSE_RS0: amg_coarse_rs0 (level, A, Pointvector, Slicing, tag); break;
00059         case VIENNACL_AMG_COARSE_RS3: amg_coarse_rs3 (level, A, Pointvector, Slicing, tag); break;
00060         case VIENNACL_AMG_COARSE_AG:   amg_coarse_ag (level, A, Pointvector, tag); break;
00061       }
00062     }
00063 
00070     template <typename InternalType1, typename InternalType2>
00071     void amg_influence(unsigned int level, InternalType1 const & A, InternalType2 & Pointvector, amg_tag & tag)
00072     {
00073       typedef typename InternalType1::value_type SparseMatrixType;
00074       typedef typename InternalType2::value_type PointVectorType;
00075       typedef typename SparseMatrixType::value_type ScalarType;
00076       typedef typename SparseMatrixType::value_type ScalarType;
00077       typedef typename SparseMatrixType::const_iterator1 ConstRowIterator;
00078       typedef typename SparseMatrixType::const_iterator2 ConstColIterator;
00079 
00080       ScalarType max;
00081       int diag_sign;
00082       //unsigned int i;
00083 
00084 #ifdef VIENNACL_WITH_OPENMP
00085       #pragma omp parallel for private (max,diag_sign)
00086 #endif
00087       for (long i=0; i<static_cast<long>(A[level].size1()); ++i)
00088       {
00089         diag_sign = 1;
00090         if (A[level](i,i) < 0)
00091           diag_sign = -1;
00092 
00093         ConstRowIterator row_iter = A[level].begin1();
00094         row_iter += i;
00095         // Find greatest non-diagonal negative value (positive if diagonal is negative) in row
00096         max = 0;
00097         for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00098         {
00099             if (i == (unsigned int) col_iter.index2()) continue;
00100             if (diag_sign == 1)
00101               if (max > *col_iter)  max = *col_iter;
00102             if (diag_sign == -1)
00103               if (max < *col_iter)  max = *col_iter;
00104         }
00105 
00106         // If maximum is 0 then the row is independent of the others
00107         if (max == 0)
00108           continue;
00109 
00110         // Find all points that strongly influence current point (Yang, p.5)
00111         for (ConstColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00112         {
00113           unsigned int j = static_cast<unsigned int>(col_iter.index2());
00114           if (i == j) continue;
00115           if (diag_sign * (-*col_iter) >= tag.get_threshold() * (diag_sign * (-max)))
00116           {
00117             // Strong influence from j to i found, save information
00118             Pointvector[level][i]->add_influencing_point(Pointvector[level][j]);
00119           }
00120         }
00121       }
00122 
00123       #ifdef VIENNACL_AMG_DEBUG
00124       std::cout << "Influence Matrix: " << std::endl;
00125       boost::numeric::ublas::matrix<bool> mat;
00126       Pointvector[level].get_influence_matrix(mat);
00127       printmatrix (mat);
00128       #endif
00129 
00130       // Save influenced points
00131       for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
00132       {
00133         for (typename amg_point::iterator iter2 = (*iter)->begin_influencing(); iter2 != (*iter)->end_influencing(); ++iter2)
00134         {
00135           (*iter2)->add_influenced_point(*iter);
00136         }
00137       }
00138 
00139       #ifdef VIENNACL_AMG_DEBUG
00140       std::cout << "Influence Measures: " << std::endl;
00141       boost::numeric::ublas::vector<unsigned int> temp;
00142       Pointvector[level].get_influence(temp);
00143       printvector (temp);
00144       std::cout << "Point Sorting: " << std::endl;
00145       Pointvector[level].get_sorting(temp);
00146       printvector (temp);
00147       #endif
00148     }
00149 
00156     template <typename InternalType1, typename InternalType2>
00157     void amg_coarse_classic_onepass(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, amg_tag & tag)
00158     {
00159       amg_point* c_point, *point1, *point2;
00160 
00161       // Check and save all strong influences
00162       amg_influence (level, A, Pointvector, tag);
00163 
00164       // Traverse through points and calculate initial influence measure
00165       long i;
00166 #ifdef VIENNACL_WITH_OPENMP
00167       #pragma omp parallel for private (i)
00168 #endif
00169       for (i=0; i<static_cast<long>(Pointvector[level].size()); ++i)
00170   Pointvector[level][i]->calc_influence();
00171 
00172        // Do initial sorting
00173       Pointvector[level].sort();
00174 
00175       // Get undecided point with highest influence measure
00176       while ((c_point = Pointvector[level].get_nextpoint()) != NULL)
00177       {
00178         // Make this point C point
00179         Pointvector[level].make_cpoint(c_point);
00180 
00181         // All strongly influenced points become F points
00182         for (typename amg_point::iterator iter = c_point->begin_influenced(); iter != c_point->end_influenced(); ++iter)
00183         {
00184           point1 = *iter;
00185           // Found strong influence from C point (c_point influences point1), check whether point is still undecided, otherwise skip
00186           if (!point1->is_undecided()) continue;
00187           // Make this point F point if it is still undecided point
00188           Pointvector[level].make_fpoint(point1);
00189 
00190           // Add +1 to influence measure for all undecided points that strongly influence new F point
00191           for (typename amg_point::iterator iter2 = point1->begin_influencing(); iter2 != point1->end_influencing(); ++iter2)
00192           {
00193             point2 = *iter2;
00194             // Found strong influence to F point (point2 influences point1)
00195             if (point2->is_undecided())
00196               Pointvector[level].add_influence(point2,1);
00197           }
00198         }
00199       }
00200 
00201       // If a point is neither C nor F point but is nevertheless influenced by other points make it F point
00202       // (this situation can happen when this point does not influence other points and the points that influence this point became F points already)
00203       /*#pragma omp parallel for private (i,point1)
00204       for (long i=0; i<static_cast<long>(Pointvector[level].size()); ++i)
00205       {
00206         point1 = Pointvector[level][i];
00207         if (point1->is_undecided())
00208         {
00209           // Undecided point found. Check whether it is influenced by other point and if so: Make it F point.
00210           if (point1->number_influencing() > 0)
00211           {
00212             #pragma omp critical
00213             Pointvector[level].make_fpoint(point1);
00214           }
00215         }
00216       }*/
00217 
00218       #if defined (VIENNACL_AMG_DEBUG)//  or defined (VIENNACL_AMG_DEBUGBENCH)
00219       unsigned int c_points = Pointvector[level].get_cpoints();
00220       unsigned int f_points = Pointvector[level].get_fpoints();
00221       std::cout << "1st pass: Level " << level << ": ";
00222       std::cout << "No of C points = " << c_points << ", ";
00223       std::cout << "No of F points = " << f_points << std::endl;
00224       #endif
00225 
00226       #ifdef VIENNACL_AMG_DEBUG
00227       std::cout << "Coarse Points:" << std::endl;
00228       boost::numeric::ublas::vector<bool> C;
00229       Pointvector[level].get_C(C);
00230       printvector (C);
00231       std::cout << "Fine Points:" << std::endl;
00232       boost::numeric::ublas::vector<bool> F;
00233       Pointvector[level].get_F(F);
00234       printvector (F);
00235       #endif
00236     }
00237 
00244     template <typename InternalType1, typename InternalType2>
00245     void amg_coarse_classic(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, amg_tag & tag)
00246     {
00247       typedef typename InternalType2::value_type PointVectorType;
00248 
00249       bool add_C;
00250       amg_point *c_point, *point1, *point2;
00251 
00252       // Use one-pass-coarsening as first pass.
00253       amg_coarse_classic_onepass(level, A, Pointvector, tag);
00254 
00255       // 2nd pass: Add more C points if F-F connection does not have a common C point.
00256       for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
00257       {
00258         point1 = *iter;
00259         // If point is F point, check for strong connections.
00260         if (point1->is_fpoint())
00261         {
00262           // Check for strong connections from influencing and influenced points.
00263           amg_point::iterator iter2 = point1->begin_influencing();
00264           amg_point::iterator iter3 = point1->begin_influenced();
00265 
00266           // Iterate over both lists at once. This makes sure that points are no checked twice when influence relation is symmetric (which is often the case).
00267           // Note: Only works because influencing and influenced lists are sorted by point-index.
00268           while(iter2 != point1->end_influencing() || iter3 != point1->end_influenced())
00269           {
00270             if (iter2 == point1->end_influencing())
00271             {
00272               point2 = *iter3;
00273               ++iter3;
00274             }
00275             else if (iter3 == point1->end_influenced())
00276             {
00277               point2 = *iter2;
00278               ++iter2;
00279             }
00280             else
00281             {
00282               if ((*iter2)->get_index() == (*iter3)->get_index())
00283               {
00284                 point2 = *iter2;
00285                 ++iter2;
00286                 ++iter3;
00287               }
00288               else if ((*iter2)->get_index() < (*iter3)->get_index())
00289               {
00290                 point2 = *iter2;
00291                 ++iter2;
00292               }
00293               else
00294               {
00295                 point2 = *iter3;
00296                 ++iter3;
00297               }
00298             }
00299             // Only check points with higher index as points with lower index have been checked already.
00300             if (point2->get_index() < point1->get_index())
00301               continue;
00302 
00303             // If there is a strong connection then it has to either be a C point or a F point with common C point.
00304             // C point? Then skip as everything is ok.
00305             if (point2->is_cpoint())
00306               continue;
00307             // F point? Then check whether F points point1 and point2 have a common C point.
00308             if (point2->is_fpoint())
00309             {
00310               add_C = true;
00311               // C point is common for two F points if they are both strongly influenced by that C point.
00312               // Compare strong influences for point1 and point2.
00313               for (amg_point::iterator iter3 = point1->begin_influencing(); iter3 != point1 -> end_influencing(); ++iter3)
00314               {
00315                 c_point = *iter3;
00316                 // Stop search when strong common influence is found via c_point.
00317                 if (c_point->is_cpoint())
00318                 {
00319                   if (point2->is_influencing(c_point))
00320                   {
00321                     add_C = false;
00322                     break;
00323                   }
00324                 }
00325               }
00326               // No common C point found? Then make second F point to C point.
00327               if (add_C == true)
00328                 Pointvector[level].switch_ftoc(point2);
00329             }
00330           }
00331         }
00332       }
00333 
00334       #ifdef VIENNACL_AMG_DEBUG
00335       std::cout << "After 2nd pass:" << std::endl;
00336       std::cout << "Coarse Points:" << std::endl;
00337       boost::numeric::ublas::vector<bool> C;
00338       Pointvector[level].get_C(C);
00339       printvector (C);
00340       std::cout << "Fine Points:" << std::endl;
00341       boost::numeric::ublas::vector<bool> F;
00342       Pointvector[level].get_F(F);
00343       printvector (F);
00344       #endif
00345 
00346       #ifdef VIENNACL_AMG_DEBUG
00347 #ifdef VIENNACL_WITH_OPENMP
00348       #pragma omp critical
00349 #endif
00350       {
00351       std::cout << "No C and no F point: ";
00352       for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
00353   if ((*iter)->is_undecided())
00354     std::cout << (*iter)->get_index() << " ";
00355       std::cout << std::endl;
00356       }
00357       #endif
00358     }
00359 
00367     template <typename InternalType1, typename InternalType2, typename InternalType3>
00368     void amg_coarse_rs0(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, InternalType3 & Slicing, amg_tag & tag)
00369     {
00370       unsigned int total_points;
00371 
00372       // Slice matrix into parts such that points are distributed among threads
00373       Slicing.slice(level, A, Pointvector);
00374 
00375       // Run classical coarsening in parallel
00376       total_points = 0;
00377 #ifdef VIENNACL_WITH_OPENMP
00378       #pragma omp parallel for
00379 #endif
00380       for (long i=0; i<static_cast<long>(Slicing.threads_); ++i)
00381       {
00382         amg_coarse_classic(level,Slicing.A_slice[i],Slicing.Pointvector_slice[i],tag);
00383 
00384         // Save C points (using Slicing.Offset on the next level as temporary memory)
00385         // Note: Number of C points for point i is saved in i+1!! (makes it easier later to compute offset)
00386         Slicing.Offset[i+1][level+1] = Slicing.Pointvector_slice[i][level].get_cpoints();
00387       #ifdef VIENNACL_WITH_OPENMP
00388         #pragma omp critical
00389       #endif
00390         total_points += Slicing.Pointvector_slice[i][level].get_cpoints();
00391       }
00392 
00393       // If no coarser level can be found on any level then resume and coarsening will stop in amg_coarse()
00394       if (total_points != 0)
00395       {
00396       #ifdef VIENNACL_WITH_OPENMP
00397         #pragma omp parallel for
00398       #endif
00399         for (long i=0; i<static_cast<long>(Slicing.threads_); ++i)
00400         {
00401           // If no higher coarse level can be found on slice i (saved in Slicing.Offset[i+1][level+1]) then pull C point(s) to the next level
00402           if (Slicing.Offset[i+1][level+1] == 0)
00403           {
00404             // All points become C points
00405             for (unsigned int j=0; j<Slicing.A_slice[i][level].size1(); ++j)
00406               Slicing.Pointvector_slice[i][level].make_cpoint(Slicing.Pointvector_slice[i][level][j]);
00407             Slicing.Offset[i+1][level+1] = static_cast<unsigned int>(Slicing.A_slice[i][level].size1());
00408           }
00409         }
00410 
00411         // Build slicing offset from number of C points (offset = total sum of C points on threads with lower number)
00412         for (unsigned int i=2; i<=Slicing.threads_; ++i)
00413           Slicing.Offset[i][level+1] += Slicing.Offset[i-1][level+1];
00414 
00415         // Join C and F points
00416         Slicing.join(level, Pointvector);
00417       }
00418 
00419       // Calculate global influence measures for interpolation and/or RS3.
00420       amg_influence(level, A, Pointvector, tag);
00421 
00422       #if defined(VIENNACL_AMG_DEBUG)// or defined (VIENNACL_AMG_DEBUGBENCH)
00423       for (unsigned int i=0; i<Slicing._threads; ++i)
00424       {
00425         unsigned int c_points = Slicing.Pointvector_slice[i][level].get_cpoints();
00426         unsigned int f_points = Slicing.Pointvector_slice[i][level].get_fpoints();
00427         std::cout << "Thread " << i << ": ";
00428         std::cout << "No of C points = " << c_points << ", ";
00429         std::cout << "No of F points = " << f_points << std::endl;
00430       }
00431       #endif
00432     }
00433 
00441     template <typename InternalType1, typename InternalType2, typename InternalType3>
00442     void amg_coarse_rs3(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, InternalType3 & Slicing, amg_tag & tag)
00443     {
00444       amg_point *c_point, *point1, *point2;
00445       bool add_C;
00446       unsigned int i, j;
00447 
00448       // Run RS0 first (parallel).
00449       amg_coarse_rs0(level, A, Pointvector, Slicing, tag);
00450 
00451       // Save slicing offset
00452       boost::numeric::ublas::vector<unsigned int> Offset = boost::numeric::ublas::vector<unsigned int> (Slicing.Offset.size());
00453       for (i=0; i<Slicing.Offset.size(); ++i)
00454         Offset[i] = Slicing.Offset[i][level];
00455 
00456       // Correct the coarsening with a third pass: Don't allow strong F-F connections without common C point
00457       for (i=0; i<Slicing.threads_; ++i)
00458       {
00459       //for (j=Slicing.Offset[i][level]; j<Slicing.Offset[i+1][level]; ++j)
00460       for (j=Offset[i]; j<Offset[i+1]; ++j)
00461       {
00462         point1 = Pointvector[level][j];
00463         // If point is F point, check for strong connections.
00464         if (point1->is_fpoint())
00465         {
00466           // Check for strong connections from influencing and influenced points.
00467           amg_point::iterator iter2 = point1->begin_influencing();
00468           amg_point::iterator iter3 = point1->begin_influenced();
00469 
00470           // Iterate over both lists at once. This makes sure that points are no checked twice when influence relation is symmetric (which is often the case).
00471           // Note: Only works because influencing and influenced lists are sorted by point-index.
00472           while(iter2 != point1->end_influencing() || iter3 != point1->end_influenced())
00473           {
00474             if (iter2 == point1->end_influencing())
00475             {
00476               point2 = *iter3;
00477               ++iter3;
00478             }
00479             else if (iter3 == point1->end_influenced())
00480             {
00481               point2 = *iter2;
00482               ++iter2;
00483             }
00484             else
00485             {
00486               if ((*iter2)->get_index() == (*iter3)->get_index())
00487               {
00488                 point2 = *iter2;
00489                 ++iter2;
00490                 ++iter3;
00491               }
00492               else if ((*iter2)->get_index() < (*iter3)->get_index())
00493               {
00494                 point2 = *iter2;
00495                 ++iter2;
00496               }
00497               else
00498               {
00499                 point2 = *iter3;
00500                 ++iter3;
00501               }
00502             }
00503 
00504             // Only check points with higher index as points with lower index have been checked already.
00505             if (point2->get_index() < point1->get_index())
00506               continue;
00507 
00508             // Only check points that are outside the slicing boundaries (interior F-F connections have already been checked in second pass)
00509             //if (point2->get_index() >= Slicing.Offset[i][level] || point2->get_index() < Slicing.Offset[i+1][level])
00510             if (point2->get_index() >= Offset[i] && point2->get_index() < Offset[i+1])
00511               continue;
00512 
00513             // If there is a strong connection then it has to either be a C point or a F point with common C point.
00514             // C point? Then skip as everything is ok.
00515             if (point2->is_cpoint())
00516               continue;
00517             // F point? Then check whether F points point1 and point2 have a common C point.
00518             if (point2->is_fpoint())
00519             {
00520               add_C = true;
00521               // C point is common for two F points if they are both strongly influenced by that C point.
00522               // Compare strong influences for point1 and point2.
00523               for (amg_point::iterator iter3 = point1->begin_influencing(); iter3 != point1 -> end_influencing(); ++iter3)
00524               {
00525                 c_point = *iter3;
00526                 // Stop search when strong common influence is found via c_point.
00527                 if (c_point->is_cpoint())
00528                 {
00529                   if (point2->is_influencing(c_point))
00530                   {
00531                     add_C = false;
00532                     break;
00533                   }
00534                 }
00535               }
00536               // No common C point found? Then make second F point to C point.
00537               if (add_C == true)
00538               {
00539                 Pointvector[level].switch_ftoc(point2);
00540                 // Add +1 to offsets as one C point has been added.
00541                 for (unsigned int j=i+1; j<=Slicing.threads_; ++j)
00542                   Slicing.Offset[j][level+1]++;
00543               }
00544                   }
00545                 }
00546               }
00547             }
00548           }
00549 
00550           #ifdef VIENNACL_AMG_DEBUG
00551           std::cout << "After 3rd pass:" << std::endl;
00552           std::cout << "Coarse Points:" << std::endl;
00553           boost::numeric::ublas::vector<bool> C;
00554           Pointvector[level].get_C(C);
00555           printvector (C);
00556           std::cout << "Fine Points:" << std::endl;
00557           boost::numeric::ublas::vector<bool> F;
00558           Pointvector[level].get_F(F);
00559           printvector (F);
00560           #endif
00561         }
00562 
00570         template <typename InternalType1, typename InternalType2>
00571         void amg_coarse_ag(unsigned int level, InternalType1 & A, InternalType2 & Pointvector, amg_tag & tag)
00572         {
00573           typedef typename InternalType1::value_type SparseMatrixType;
00574           typedef typename InternalType2::value_type PointVectorType;
00575           typedef typename SparseMatrixType::value_type ScalarType;
00576 
00577           typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00578           typedef typename SparseMatrixType::iterator2 InternalColIterator;
00579 
00580           long x,y;
00581           ScalarType diag;
00582           amg_point *pointx, *pointy;
00583 
00584           // Cannot determine aggregates if size == 1 as then a new aggregate would always consist of this point (infinite loop)
00585           if (A[level].size1() == 1) return;
00586 
00587           // SA algorithm (Vanek et al. p.6)
00588           // Build neighborhoods
00589     #ifdef VIENNACL_WITH_OPENMP
00590           #pragma omp parallel for private (x,y,diag)
00591     #endif
00592           for (x=0; x<static_cast<long>(A[level].size1()); ++x)
00593           {
00594             InternalRowIterator row_iter = A[level].begin1();
00595             row_iter += x;
00596             diag = A[level](x,x);
00597             for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00598             {
00599               y = static_cast<long>(col_iter.index2());
00600               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](y,y)))))
00601               {
00602                 // Neighborhood x includes point y
00603                 Pointvector[level][x]->add_influencing_point(Pointvector[level][y]);
00604               }
00605             }
00606           }
00607 
00608           #ifdef VIENNACL_AMG_DEBUG
00609           std::cout << "Neighborhoods:" << std::endl;
00610           boost::numeric::ublas::matrix<bool> mat;
00611           Pointvector[level].get_influence_matrix(mat);
00612           printmatrix (mat);
00613           #endif
00614 
00615           // Build aggregates from neighborhoods
00616           for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
00617           {
00618             pointx = (*iter);
00619 
00620             if (pointx->is_undecided())
00621             {
00622               // Make center of aggregate to C point and include it to aggregate x.
00623               Pointvector[level].make_cpoint(pointx);
00624               pointx->set_aggregate (pointx->get_index());
00625               for (amg_point::iterator iter2 = pointx->begin_influencing(); iter2 != pointx->end_influencing(); ++iter2)
00626               {
00627               pointy = (*iter2);
00628 
00629                 if (pointy->is_undecided())
00630                 {
00631                   // Make neighbor y to F point and include it to aggregate x.
00632                   Pointvector[level].make_fpoint(pointy);
00633                   pointy->set_aggregate (pointx->get_index());
00634                 }
00635               }
00636             }
00637           }
00638 
00639           #ifdef VIENNACL_AMG_DEBUG
00640           std::cout << "After aggregation:" << std::endl;
00641           std::cout << "Coarse Points:" << std::endl;
00642           boost::numeric::ublas::vector<bool> C;
00643           Pointvector[level].get_C(C);
00644           printvector (C);
00645           std::cout << "Fine Points:" << std::endl;
00646           boost::numeric::ublas::vector<bool> F;
00647           Pointvector[level].get_F(F);
00648           printvector (F);
00649           std::cout << "Aggregates:" << std::endl;
00650           printvector (Aggregates[level]);
00651           #endif
00652         }
00653       } //namespace amg
00654     }
00655   }
00656 }
00657 
00658 #endif