ViennaCL - The Vienna Computing Library  1.5.0
viennacl/linalg/detail/amg/amg_base.hpp
Go to the documentation of this file.
00001 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_BASE_HPP_
00002 #define VIENNACL_LINALG_DETAIL_AMG_AMG_BASE_HPP_
00003 
00004 /* =========================================================================
00005    Copyright (c) 2010-2013, 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 
00027 #include <boost/numeric/ublas/operation.hpp>
00028 #include <boost/numeric/ublas/vector.hpp>
00029 #include <cmath>
00030 #include <set>
00031 #include <list>
00032 #include <algorithm>
00033 
00034 #include <map>
00035 #ifdef VIENNACL_WITH_OPENMP
00036 #include <omp.h>
00037 #endif
00038 
00039 #include "amg_debug.hpp"
00040 
00041 #define VIENNACL_AMG_COARSE_RS 1
00042 #define VIENNACL_AMG_COARSE_ONEPASS 2
00043 #define VIENNACL_AMG_COARSE_RS0 3
00044 #define VIENNACL_AMG_COARSE_RS3 4
00045 #define VIENNACL_AMG_COARSE_AG 5
00046 #define VIENNACL_AMG_INTERPOL_DIRECT 1
00047 #define VIENNACL_AMG_INTERPOL_CLASSIC 2
00048 #define VIENNACL_AMG_INTERPOL_AG 3
00049 #define VIENNACL_AMG_INTERPOL_SA 4
00050 
00051 namespace viennacl
00052 {
00053   namespace linalg
00054   {
00055     namespace detail
00056     {
00057       namespace amg
00058       {
00061         class amg_tag
00062         {
00063           public:
00076             amg_tag(unsigned int coarse = 1,
00077                     unsigned int interpol = 1,
00078                     double threshold = 0.25,
00079                     double interpolweight = 0.2,
00080                     double jacobiweight = 1,
00081                     unsigned int presmooth = 1,
00082                     unsigned int postsmooth = 1,
00083                     unsigned int coarselevels = 0)
00084             : coarse_(coarse), interpol_(interpol),
00085               threshold_(threshold), interpolweight_(interpolweight), jacobiweight_(jacobiweight),
00086               presmooth_(presmooth), postsmooth_(postsmooth), coarselevels_(coarselevels) {}
00087 
00088             // Getter-/Setter-Functions
00089             void set_coarse(unsigned int coarse) { if (coarse > 0) coarse_ = coarse; }
00090             unsigned int get_coarse() const { return coarse_; }
00091 
00092             void set_interpol(unsigned int interpol) { if (interpol > 0) interpol_ = interpol; }
00093             unsigned int get_interpol() const { return interpol_; }
00094 
00095             void set_threshold(double threshold) { if (threshold > 0 && threshold <= 1) threshold_ = threshold; }
00096             double get_threshold() const{ return threshold_; }
00097 
00098             void set_as(double jacobiweight) { if (jacobiweight > 0 && jacobiweight <= 2) jacobiweight_ = jacobiweight; }
00099             double get_interpolweight() const { return interpolweight_; }
00100 
00101             void set_interpolweight(double interpolweight) { if (interpolweight > 0 && interpolweight <= 2) interpolweight_ = interpolweight; }
00102             double get_jacobiweight() const { return jacobiweight_; }
00103 
00104             void set_presmooth(int presmooth) { if (presmooth >= 0) presmooth_ = presmooth; }
00105             unsigned int get_presmooth() const { return presmooth_; }
00106 
00107             void set_postsmooth(int postsmooth) { if (postsmooth >= 0) postsmooth_ = postsmooth; }
00108             unsigned int get_postsmooth() const { return postsmooth_; }
00109 
00110             void set_coarselevels(int coarselevels)  { if (coarselevels >= 0) coarselevels_ = coarselevels; }
00111             unsigned int get_coarselevels() const { return coarselevels_; }
00112 
00113           private:
00114             unsigned int coarse_, interpol_;
00115             double threshold_, interpolweight_, jacobiweight_;
00116             unsigned int presmooth_, postsmooth_, coarselevels_;
00117         };
00118 
00123         template <typename InternalType, typename IteratorType, typename ScalarType>
00124         class amg_nonzero_scalar
00125         {
00126           private:
00127             InternalType *m_;
00128             IteratorType iter_;
00129             unsigned int i_,j_;
00130             ScalarType s_;
00131 
00132           public:
00133             amg_nonzero_scalar();
00134 
00142             amg_nonzero_scalar(InternalType *m,
00143                               IteratorType & iter,
00144                               unsigned int i,
00145                               unsigned int j,
00146                               ScalarType s = 0): m_(m), iter_(iter), i_(i), j_(j), s_(s) {}
00147 
00151             ScalarType operator = (const ScalarType value)
00152             {
00153               s_ = value;
00154               // Only write if scalar is nonzero
00155               if (s_ == 0) return s_;
00156               // Write to m_ using iterator iter_ or indices (i_,j_)
00157               m_->addscalar (iter_,i_,j_,s_);
00158               return s_;
00159             }
00160 
00164             ScalarType operator += (const ScalarType value)
00165             {
00166               // If zero is added, then no change necessary
00167               if (value == 0)
00168                 return s_;
00169 
00170               s_ += value;
00171               // Remove entry if resulting scalar is zero
00172               if (s_ == 0)
00173               {
00174                 m_->removescalar(iter_,i_);
00175                 return s_;
00176               }
00177               //Write to m_ using iterator iter_ or indices (i_,j_)
00178               m_->addscalar (iter_,i_,j_,s_);
00179               return s_;
00180             }
00181             ScalarType operator ++ (int)
00182             {
00183               s_++;
00184               if (s_ == 0)
00185                 m_->removescalar(iter_,i_);
00186               m_->addscalar (iter_,i_,j_,s_);
00187               return s_;
00188             }
00189             ScalarType operator ++ ()
00190             {
00191               s_++;
00192               if (s_ == 0)
00193                 m_->removescalar(iter_,i_);
00194               m_->addscalar (iter_,i_,j_,s_);
00195               return s_;
00196             }
00197             operator ScalarType (void) { return s_;  }
00198         };
00199 
00202         template <typename InternalType>
00203         class amg_sparsevector_iterator
00204         {
00205           private:
00206             typedef amg_sparsevector_iterator<InternalType> self_type;
00207             typedef typename InternalType::mapped_type ScalarType;
00208 
00209             InternalType & internal_vec;
00210             typename InternalType::iterator iter;
00211 
00212           public:
00213 
00218             amg_sparsevector_iterator(InternalType & vec, bool begin=true): internal_vec(vec)
00219             {
00220               if (begin)
00221                 iter = internal_vec.begin();
00222               else
00223                 iter = internal_vec.end();
00224             }
00225 
00226             bool operator == (self_type other)
00227             {
00228               if (iter == other.iter)
00229                 return true;
00230               else
00231                 return false;
00232             }
00233             bool operator != (self_type other)
00234             {
00235               if (iter != other.iter)
00236                 return true;
00237               else
00238                 return false;
00239             }
00240 
00241             self_type & operator ++ () const { iter++; return *this; }
00242             self_type & operator ++ () { iter++; return *this; }
00243             self_type & operator -- () const { iter--; return *this; }
00244             self_type & operator -- () { iter--; return *this; }
00245             ScalarType & operator * () const { return (*iter).second; }
00246             ScalarType & operator * () { return (*iter).second; }
00247             unsigned int index() const { return (*iter).first; }
00248             unsigned int index() { return (*iter).first; }
00249         };
00250 
00253         template <typename ScalarType>
00254         class amg_sparsevector
00255         {
00256           public:
00257             typedef ScalarType value_type;
00258 
00259           private:
00260             // A map is used internally which saves all non-zero elements with pairs of (index,value)
00261             typedef std::map<unsigned int,ScalarType> InternalType;
00262             typedef amg_sparsevector<ScalarType> self_type;
00263             typedef amg_nonzero_scalar<self_type,typename InternalType::iterator,ScalarType> NonzeroScalarType;
00264 
00265             // Size is only a dummy variable. Not needed for internal map structure but for compatible vector interface.
00266             unsigned int size_;
00267             InternalType internal_vector;
00268 
00269           public:
00270             typedef amg_sparsevector_iterator<InternalType> iterator;
00271             typedef typename InternalType::const_iterator const_iterator;
00272 
00273           public:
00277             amg_sparsevector(unsigned int size = 0): size_(size)
00278             {
00279               internal_vector = InternalType();
00280             }
00281 
00282             void resize(unsigned int size) { size_ = size; }
00283             unsigned int size() const { return size_;}
00284 
00285             // Returns number of non-zero entries in vector equal to the size of the underlying map.
00286             unsigned int internal_size() const { return static_cast<unsigned int>(internal_vector.size()); }
00287             // Delete underlying map.
00288             void clear() { internal_vector.clear();  }
00289             // Remove entry at position i.
00290             void remove(unsigned int i) { internal_vector.erase(i); }
00291 
00292             // Add s to the entry at position i
00293             void add (unsigned int i, ScalarType s)
00294             {
00295               typename InternalType::iterator iter = internal_vector.find(i);
00296               // If there is no entry at position i, add new entry at that position
00297               if (iter == internal_vector.end())
00298                 addscalar(iter,i,i,s);
00299               else
00300               {
00301                 s += (*iter).second;
00302                 // If new value is zero, then erase the entry, otherwise write new value
00303                 if (s != 0)
00304                   (*iter).second = s;
00305                 else
00306                   internal_vector.erase(iter);
00307               }
00308             }
00309 
00310             // Write to the map. Is called from non-zero scalar type.
00311             template <typename IteratorType>
00312             void addscalar(IteratorType & iter, unsigned int i, unsigned int /* j */, ScalarType s)
00313             {
00314               // Don't write if value is zero
00315               if (s == 0)
00316                 return;
00317 
00318               // If entry is already present, overwrite value, otherwise make new entry
00319               if (iter != internal_vector.end())
00320                 (*iter).second = s;
00321               else
00322                 internal_vector[i] = s;
00323             }
00324 
00325             // Remove value from the map. Is called from non-zero scalar type.
00326             template <typename IteratorType>
00327             void removescalar(IteratorType & iter, unsigned int /* i */) { internal_vector.erase(iter); }
00328 
00329             // Bracket operator. Returns non-zero scalar type with actual values of the respective entry which calls addscalar/removescalar after value is altered.
00330             NonzeroScalarType operator [] (unsigned int i)
00331             {
00332               typename InternalType::iterator it = internal_vector.find(i);
00333               // If value is already present then build non-zero scalar with actual value, otherwise 0.
00334               if (it != internal_vector.end())
00335                 return NonzeroScalarType (this,it,i,i,(*it).second);
00336               else
00337                 return NonzeroScalarType (this,it,i,i,0);
00338             }
00339 
00340             // Use internal data structure directly for read-only access. No need to use non-zero scalar as no write access possible.
00341             ScalarType operator [] (unsigned int i) const
00342             {
00343               const_iterator it = internal_vector.find(i);
00344 
00345               if (it != internal_vector.end())
00346                 return (*it).second;
00347               else
00348                 return 0;
00349             }
00350 
00351             // Iterator functions.
00352             iterator begin() { return iterator(internal_vector); }
00353             const_iterator begin() const { return internal_vector.begin(); }
00354             iterator end() { return iterator(internal_vector,false); }
00355             const_iterator end() const { return internal_vector.end(); }
00356 
00357             // checks whether value at index i is nonzero. More efficient than doing [] == 0.
00358             bool isnonzero(unsigned int i) const { return internal_vector.find(i) != internal_vector.end();  }
00359 
00360             // Copies data into a ublas vector type.
00361             operator boost::numeric::ublas::vector<ScalarType> (void)
00362             {
00363               boost::numeric::ublas::vector<ScalarType> vec (size_);
00364               for (iterator iter = begin(); iter != end(); ++iter)
00365                 vec [iter.index()] = *iter;
00366               return vec;
00367             }
00368          };
00369 
00375         template <typename ScalarType>
00376         class amg_sparsematrix
00377         {
00378           public:
00379             typedef ScalarType value_type;
00380           private:
00381             typedef std::map<unsigned int,ScalarType> RowType;
00382             typedef std::vector<RowType> InternalType;
00383             typedef amg_sparsematrix<ScalarType> self_type;
00384 
00385             // Adapter is used for certain functionality, especially iterators.
00386             typedef typename viennacl::tools::sparse_matrix_adapter<ScalarType> AdapterType;
00387             typedef typename viennacl::tools::const_sparse_matrix_adapter<ScalarType> ConstAdapterType;
00388 
00389             // Non-zero scalar is used to write to the matrix.
00390             typedef amg_nonzero_scalar<self_type,typename RowType::iterator,ScalarType> NonzeroScalarType;
00391 
00392             // Holds matrix coefficients.
00393             InternalType internal_mat;
00394             // Holds matrix coefficient of transposed matrix if built.
00395             // Note: Only internal_mat is written using operators and methods while internal_mat_trans is built from internal_mat using do_trans().
00396             InternalType internal_mat_trans;
00397             // Saves sizes.
00398             vcl_size_t s1, s2;
00399 
00400             // True if the transposed of the matrix is used (for calculations, iteration, etc.).
00401             bool transposed_mode;
00402             // True if the transposed is already built (saved in internal_mat_trans) and also up to date (no changes to internal_mat).
00403             bool transposed;
00404 
00405           public:
00406             typedef typename AdapterType::iterator1 iterator1;
00407             typedef typename AdapterType::iterator2 iterator2;
00408             typedef typename ConstAdapterType::const_iterator1 const_iterator1;
00409             typedef typename ConstAdapterType::const_iterator2 const_iterator2;
00410 
00412             amg_sparsematrix ()
00413             {
00414               transposed_mode = false;
00415               transposed = false;
00416             }
00417 
00422             amg_sparsematrix (unsigned int i, unsigned int j)
00423             {
00424               AdapterType a (internal_mat, i, j);
00425               a.resize (i,j,false);
00426               AdapterType a_trans (internal_mat_trans, j, i);
00427               a_trans.resize (j,i,false);
00428               s1 = i;
00429               s2 = j;
00430               a.clear();
00431               a_trans.clear();
00432               transposed_mode = false;
00433               transposed = false;
00434             }
00435 
00440             amg_sparsematrix (std::vector<std::map<unsigned int, ScalarType> > const & mat)
00441             {
00442               AdapterType a (internal_mat, mat.size(), mat.size());
00443               AdapterType a_trans (internal_mat_trans, mat.size(), mat.size());
00444               a.resize(mat.size(), mat.size());
00445               a_trans.resize(mat.size(), mat.size());
00446 
00447               internal_mat = mat;
00448               s1 = s2 = mat.size();
00449 
00450               transposed_mode = false;
00451               transposed = false;
00452             }
00453 
00458             template <typename MatrixType>
00459             amg_sparsematrix (MatrixType const & mat)
00460             {
00461               AdapterType a (internal_mat, mat.size1(), mat.size2());
00462               AdapterType a_trans (internal_mat_trans, mat.size2(), mat.size1());
00463               a.resize(mat.size1(), mat.size2());
00464               a_trans.resize(mat.size2(), mat.size1());
00465               s1 = mat.size1();
00466               s2 = mat.size2();
00467               a.clear();
00468               a_trans.clear();
00469 
00470               for (typename MatrixType::const_iterator1 row_iter = mat.begin1(); row_iter != mat.end1(); ++row_iter)
00471               {
00472                 for (typename MatrixType::const_iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00473                 {
00474                   if (*col_iter != 0)
00475                   {
00476                     unsigned int x = static_cast<unsigned int>(col_iter.index1());
00477                     unsigned int y = static_cast<unsigned int>(col_iter.index2());
00478                     a (x,y) = *col_iter;
00479                     a_trans (y,x) = *col_iter;
00480                   }
00481                 }
00482               }
00483               transposed_mode = false;
00484               transposed = true;
00485             }
00486 
00487             // Build transposed of the current matrix.
00488             void do_trans()
00489             {
00490               // Do it only once if called in a parallel section
00491             #ifdef VIENNACL_WITH_OPENMP
00492               #pragma omp critical
00493             #endif
00494               {
00495                 // Only build transposed if it is not built or not up to date
00496                 if (!transposed)
00497                 {
00498                   // Mode has to be set to standard mode temporarily
00499                   bool save_mode = transposed_mode;
00500                   transposed_mode = false;
00501 
00502                   for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00503                     for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00504                       internal_mat_trans[col_iter.index2()][static_cast<unsigned int>(col_iter.index1())] = *col_iter;
00505 
00506                   transposed_mode = save_mode;
00507                   transposed = true;
00508                 }
00509               }
00510             } //do_trans()
00511 
00512             // Set transposed mode (true=transposed, false=regular)
00513             void set_trans(bool mode)
00514             {
00515               transposed_mode = mode;
00516               if (mode)
00517                 do_trans();
00518             }
00519 
00520             bool get_trans() const { return transposed_mode; }
00521 
00522             // Checks whether coefficient (i,j) is non-zero. More efficient than using (i,j) == 0.
00523             bool isnonzero (unsigned int i, unsigned int j) const
00524             {
00525               if (!transposed_mode)
00526               {
00527                 if (internal_mat[i].find(j) != internal_mat[i].end())
00528                   return true;
00529                 else
00530                   return false;
00531               }
00532               else
00533               {
00534                 if (internal_mat[j].find(i) != internal_mat[j].end())
00535                   return true;
00536                 else
00537                   return false;
00538               }
00539             } //isnonzero()
00540 
00541             // Add s to value at (i,j)
00542             void add (unsigned int i, unsigned int j, ScalarType s)
00543             {
00544               // If zero is added then do nothing.
00545               if (s == 0)
00546                 return;
00547 
00548               typename RowType::iterator col_iter = internal_mat[i].find(j);
00549               // If there is no entry at position (i,j), then make new entry.
00550               if (col_iter == internal_mat[i].end())
00551                 addscalar(col_iter,i,j,s);
00552               else
00553               {
00554                 s += (*col_iter).second;
00555                 // Update value and erase entry if value is zero.
00556                 if (s != 0)
00557                   (*col_iter).second = s;
00558                 else
00559                   internal_mat[i].erase(col_iter);
00560               }
00561               transposed = false;
00562             } //add()
00563 
00564             // Write to the internal data structure. Is called from non-zero scalar type.
00565             template <typename IteratorType>
00566             void addscalar(IteratorType & iter, unsigned int i, unsigned int j, ScalarType s)
00567             {
00568               // Don't write if value is zero
00569               if (s == 0)
00570                 return;
00571 
00572               if (iter != internal_mat[i].end())
00573                 (*iter).second = s;
00574               else
00575                 internal_mat[i][j] = s;
00576 
00577               transposed = false;
00578             }
00579 
00580             // Remove entry from internal data structure. Is called from non-zero scalar type.
00581             template <typename IteratorType>
00582             void removescalar(IteratorType & iter, unsigned int i)
00583             {
00584               internal_mat[i].erase(iter);
00585               transposed = false;
00586             }
00587 
00588             // Return non-zero scalar at position (i,j). Value is written to the non-zero scalar and updated via addscalar()/removescalar().
00589             NonzeroScalarType operator()(unsigned int i, unsigned int j)
00590             {
00591               typename RowType::iterator iter;
00592 
00593               if (!transposed_mode)
00594               {
00595                 iter = internal_mat[i].find(j);
00596                 if (iter != internal_mat[i].end())
00597                   return NonzeroScalarType (this,iter,i,j,(*iter).second);
00598                 else
00599                   return NonzeroScalarType (this,iter,i,j,0);
00600               }
00601               else
00602               {
00603                 iter = internal_mat[j].find(i);
00604                 if (iter != internal_mat[j].end())
00605                   return NonzeroScalarType (this,iter,j,i,(*iter).second);
00606                 else
00607                   return NonzeroScalarType (this,iter,j,i,0);
00608               }
00609             }
00610 
00611             // For read-only access return the actual value directly. Non-zero datatype not needed as no write access possible.
00612             ScalarType operator()(unsigned int i, unsigned int j) const
00613             {
00614               typename RowType::const_iterator iter;
00615 
00616               if (!transposed_mode)
00617               {
00618                 iter = internal_mat[i].find(j);
00619                 if (iter != internal_mat[i].end())
00620                   return (*iter).second;
00621                 else
00622                   return 0;
00623               }
00624               else
00625               {
00626                 iter = internal_mat[j].find(i);
00627                 if (iter != internal_mat[j].end())
00628                   return (*iter).second;
00629                 else
00630                   return 0;
00631               }
00632             }
00633 
00634             void resize(unsigned int i, unsigned int j, bool preserve = true)
00635             {
00636               AdapterType a (internal_mat);
00637               a.resize(i,j,preserve);
00638               AdapterType a_trans (internal_mat_trans);
00639               a_trans.resize(j,i,preserve);
00640               s1 = i;
00641               s2 = j;
00642             }
00643 
00644             void clear()
00645             {
00646               AdapterType a (internal_mat, s1, s2);
00647               a.clear();
00648               AdapterType a_trans (internal_mat_trans, s2, s1);
00649               a_trans.clear();
00650               transposed = true;
00651             }
00652 
00653             vcl_size_t size1()
00654             {
00655               if (!transposed_mode)
00656                 return s1;
00657               else
00658                 return s2;
00659             }
00660 
00661             vcl_size_t size1() const
00662             {
00663               if (!transposed_mode)
00664                 return s1;
00665               else
00666                 return s2;
00667             }
00668 
00669 
00670             vcl_size_t size2()
00671             {
00672               if (!transposed_mode)
00673                 return s2;
00674               else
00675                 return s1;
00676             }
00677 
00678             vcl_size_t size2() const
00679             {
00680               if (!transposed_mode)
00681                 return s2;
00682               else
00683                 return s1;
00684             }
00685 
00686             iterator1 begin1(bool trans = false)
00687             {
00688               if (!trans && !transposed_mode)
00689               {
00690                 AdapterType a (internal_mat, s1, s2);
00691                 return a.begin1();
00692               }
00693               else
00694               {
00695                 do_trans();
00696                 AdapterType a_trans (internal_mat_trans, s2, s1);
00697                 return a_trans.begin1();
00698               }
00699             }
00700 
00701             iterator1 end1(bool trans = false)
00702             {
00703               if (!trans && !transposed_mode)
00704               {
00705                 AdapterType a (internal_mat, s1, s2);
00706                 return a.end1();
00707               }
00708               else
00709               {
00710                 //do_trans();
00711                 AdapterType a_trans (internal_mat_trans, s2, s1);
00712                 return a_trans.end1();
00713               }
00714             }
00715 
00716             iterator2 begin2(bool trans = false)
00717             {
00718               if (!trans && !transposed_mode)
00719               {
00720                 AdapterType a (internal_mat, s1, s2);
00721                 return a.begin2();
00722               }
00723               else
00724               {
00725                 do_trans();
00726                 AdapterType a_trans (internal_mat_trans, s2, s1);
00727                 return a_trans.begin2();
00728               }
00729             }
00730 
00731             iterator2 end2(bool trans = false)
00732             {
00733               if (!trans && !transposed_mode)
00734               {
00735                 AdapterType a (internal_mat, s1, s2);
00736                 return a.end2();
00737               }
00738               else
00739               {
00740                 //do_trans();
00741                 AdapterType a_trans (internal_mat_trans, s2, s1);
00742                 return a_trans.end2();
00743               }
00744             }
00745 
00746             const_iterator1 begin1() const
00747             {
00748               // Const_iterator of transposed can only be used if transposed matrix is already built and up to date.
00749               assert((!transposed_mode || (transposed_mode && transposed)) && bool("Error: Cannot build const_iterator when transposed has not been built yet!"));
00750                     ConstAdapterType a_const (internal_mat, s1, s2);
00751               return a_const.begin1();
00752             }
00753 
00754             const_iterator1 end1(bool trans = false) const
00755             {
00756               assert((!transposed_mode || (transposed_mode && transposed)) && bool("Error: Cannot build const_iterator when transposed has not been built yet!"));
00757               ConstAdapterType a_const (internal_mat, trans ? s2 : s1, trans ? s1 : s2);
00758               return a_const.end1();
00759             }
00760 
00761             const_iterator2 begin2(bool trans = false) const
00762             {
00763               assert((!transposed_mode || (transposed_mode && transposed)) && bool("Error: Cannot build const_iterator when transposed has not been built yet!"));
00764               ConstAdapterType a_const (internal_mat, trans ? s2 : s1, trans ? s1 : s2);
00765               return a_const.begin2();
00766             }
00767 
00768             const_iterator2 end2(bool trans = false) const
00769             {
00770               assert((!transposed_mode || (transposed_mode && transposed)) && bool("Error: Cannot build const_iterator when transposed has not been built yet!"));
00771               ConstAdapterType a_const (internal_mat, trans ? s2 : s1, trans ? s1 : s2);
00772               return a_const.end2();
00773             }
00774 
00775             // Returns pointer to the internal data structure. Improves performance of copy operation to GPU.
00776             std::vector<std::map<unsigned int, ScalarType> > * get_internal_pointer()
00777             {
00778               if (!transposed_mode)
00779                 return &internal_mat;
00780 
00781               if (!transposed)
00782                 do_trans();
00783               return &internal_mat_trans;
00784             }
00785             operator boost::numeric::ublas::compressed_matrix<ScalarType> (void)
00786             {
00787               boost::numeric::ublas::compressed_matrix<ScalarType> mat;
00788               mat.resize(size1(),size2(),false);
00789               mat.clear();
00790 
00791               for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00792                   for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00793                     mat (col_iter.index1(), col_iter.index2()) = *col_iter;
00794 
00795               return mat;
00796             }
00797             operator boost::numeric::ublas::matrix<ScalarType> (void)
00798             {
00799               boost::numeric::ublas::matrix<ScalarType> mat;
00800               mat.resize(size1(),size2(),false);
00801               mat.clear();
00802 
00803               for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00804                   for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00805                     mat (col_iter.index1(), col_iter.index2()) = *col_iter;
00806 
00807               return mat;
00808             }
00809         };
00810 
00816         class amg_point
00817         {
00818           private:
00819             typedef amg_sparsevector<amg_point*> ListType;
00820 
00821             unsigned int index_;
00822             unsigned int influence_;
00823             // Determines whether point is undecided.
00824             bool undecided_;
00825             // Determines wheter point is C point (true) or F point (false).
00826             bool cpoint_;
00827             unsigned int coarse_index_;
00828             // Index offset of parallel coarsening. In that case a point acts as if it had an index of index_-offset_ and treats other points as if they had an index of index+offset_
00829             unsigned int offset_;
00830             // Aggregate the point belongs to.
00831             unsigned int aggregate_;
00832 
00833             // Holds all points influencing this point.
00834             ListType influencing_points;
00835             // Holds all points that are influenced by this point.
00836             ListType influenced_points;
00837 
00838           public:
00839             typedef ListType::iterator iterator;
00840             typedef ListType::const_iterator const_iterator;
00841 
00844             amg_point (unsigned int index, unsigned int size): index_(index), influence_(0), undecided_(true), cpoint_(false), coarse_index_(0), offset_(0), aggregate_(0)
00845             {
00846               influencing_points = ListType(size);
00847               influenced_points = ListType(size);
00848             }
00849 
00850             void set_offset(unsigned int offset) { offset_ = offset; }
00851             unsigned int get_offset() { return offset_; }
00852             void set_index(unsigned int index) { index_ = index+offset_; }
00853             unsigned int get_index() const { return index_-offset_;  }
00854             unsigned int get_influence() const { return influence_;  }
00855             void set_aggregate(unsigned int aggregate) { aggregate_ = aggregate; }
00856             unsigned int get_aggregate () { return aggregate_; }
00857 
00858             bool is_cpoint() const { return cpoint_ && !undecided_;  }
00859             bool is_fpoint() const { return !cpoint_ && !undecided_; }
00860             bool is_undecided() const { return undecided_; }
00861 
00862             // Returns number of influencing points
00863             unsigned int number_influencing() const  { return influencing_points.internal_size(); }
00864             // Returns true if *point is influencing this point
00865             bool is_influencing(amg_point* point) const { return influencing_points.isnonzero(point->get_index()+offset_); }
00866             // Add *point to influencing points
00867             void add_influencing_point(amg_point* point) { influencing_points[point->get_index()+offset_] = point;  }
00868             // Add *point to influenced points
00869             void add_influenced_point(amg_point* point) { influenced_points[point->get_index()+offset_] = point; }
00870 
00871             // Clear influencing points
00872             void clear_influencing() { influencing_points.clear(); }
00873             // Clear influenced points
00874             void clear_influenced() {influenced_points.clear(); }
00875 
00876 
00877             unsigned int get_coarse_index() const { return coarse_index_; }
00878             void set_coarse_index(unsigned int index) { coarse_index_ = index; }
00879 
00880             // Calculates the initial influence measure equal to the number of influenced points.
00881             void calc_influence() { influence_ = influenced_points.internal_size();  }
00882 
00883             // Add to influence measure.
00884             unsigned int add_influence(int add)
00885             {
00886               influence_ += add;
00887               return influence_;
00888             }
00889             // Make this point C point. Only call via amg_pointvector.
00890             void make_cpoint()
00891             {
00892               undecided_ = false;
00893               cpoint_ = true;
00894               influence_ = 0;
00895             }
00896             // Make this point F point. Only call via amg_pointvector.
00897             void make_fpoint()
00898             {
00899               undecided_ = false;
00900               cpoint_ = false;
00901               influence_ = 0;
00902             }
00903             // Switch point from F to C point. Only call via amg_pointvector.
00904             void switch_ftoc() { cpoint_ = true; }
00905 
00906             // Iterator handling for influencing and influenced points.
00907             iterator begin_influencing() { return influencing_points.begin(); }
00908             iterator end_influencing() { return influencing_points.end(); }
00909             const_iterator begin_influencing() const { return influencing_points.begin(); }
00910             const_iterator end_influencing() const { return influencing_points.end(); }
00911             iterator begin_influenced() { return influenced_points.begin();  }
00912             iterator end_influenced() { return influenced_points.end(); }
00913             const_iterator begin_influenced() const { return influenced_points.begin(); }
00914             const_iterator end_influenced() const { return influenced_points.end(); }
00915         };
00916 
00919         struct classcomp
00920         {
00921           // Function returns true if l comes before r in the ordering.
00922           bool operator() (amg_point* l, amg_point* r) const
00923           {
00924             // Map is sorted by influence number starting with the highest
00925             // If influence number is the same then lowest point index comes first
00926             return (l->get_influence() < r->get_influence() || (l->get_influence() == r->get_influence() && l->get_index() > r->get_index()));
00927           }
00928         };
00929 
00935         class amg_pointvector
00936         {
00937           private:
00938             // Type for the sorted list
00939             typedef std::set<amg_point*,classcomp> ListType;
00940             // Type for the vector of pointers
00941             typedef std::vector<amg_point*> VectorType;
00942             VectorType pointvector;
00943             ListType pointlist;
00944             unsigned int size_;
00945             unsigned int c_points, f_points;
00946 
00947           public:
00948             typedef VectorType::iterator iterator;
00949             typedef VectorType::const_iterator const_iterator;
00950 
00954             amg_pointvector(unsigned int size = 0): size_(size)
00955             {
00956               pointvector = VectorType(size);
00957               c_points = f_points = 0;
00958             }
00959 
00960             // Construct all the points dynamically and save pointers into vector.
00961             void init_points()
00962             {
00963               for (unsigned int i=0; i<size(); ++i)
00964                 pointvector[i] = new amg_point(i,size());
00965             }
00966             // Delete all the points.
00967             void delete_points()
00968             {
00969               for (unsigned int i=0; i<size(); ++i)
00970                 delete pointvector[i];
00971             }
00972             // Add point to the vector. Note: User has to make sure that point at point->get_index() does not exist yet, otherwise it will be overwritten!
00973             void add_point(amg_point *point)
00974             {
00975               pointvector[point->get_index()] = point;
00976               if (point->is_cpoint()) c_points++;
00977               else if (point->is_fpoint()) f_points++;
00978             }
00979 
00980             // Update C and F count for point *point.
00981             // Necessary if C and F points were constructed outside this data structure (e.g. by parallel coarsening RS0 or RS3).
00982             void update_cf(amg_point *point)
00983             {
00984               if (point->is_cpoint()) c_points++;
00985               else if (point->is_fpoint()) f_points++;
00986             }
00987             // Clear the C and F point count.
00988             void clear_cf() { c_points = f_points = 0; }
00989 
00990             // Clear both point lists.
00991             void clear_influencelists()
00992             {
00993               for (iterator iter = pointvector.begin(); iter != pointvector.end(); ++iter)
00994               {
00995                 (*iter)->clear_influencing();
00996                 (*iter)->clear_influenced();
00997               }
00998             }
00999 
01000             amg_point* operator [] (unsigned int i) const { return pointvector[i]; }
01001             iterator begin() { return pointvector.begin(); }
01002             iterator end() { return pointvector.end(); }
01003             const_iterator begin() const { return pointvector.begin(); }
01004             const_iterator end() const { return pointvector.end(); }
01005 
01006             void resize(unsigned int size)
01007             {
01008               size_ = size;
01009               pointvector = VectorType(size);
01010             }
01011             unsigned int size() const { return size_; }
01012 
01013             // Returns number of C points
01014             unsigned int get_cpoints() const { return c_points; }
01015             // Returns number of F points
01016             unsigned int get_fpoints() const { return f_points; }
01017 
01018             // Does the initial sorting of points into the list. Sorting is automatically done by the std::set data type.
01019             void sort()
01020             {
01021               for (iterator iter = begin(); iter != end(); ++iter)
01022                 pointlist.insert(*iter);
01023             }
01024             // Returns the point with the highest influence measure
01025             amg_point* get_nextpoint()
01026             {
01027               // No points remaining? Return NULL such that coarsening will stop.
01028               if (pointlist.size() == 0)
01029                 return NULL;
01030               // If point with highest influence measure (end of the list) has measure of zero, then no further C points can be constructed. Return NULL.
01031               if ((*(--pointlist.end()))->get_influence() == 0)
01032                 return NULL;
01033               // Otherwise, return the point with highest influence measure located at the end of the list.
01034               else
01035                 return *(--pointlist.end());
01036             }
01037             // Add "add" to influence measure for point *point in the sorted list.
01038             void add_influence(amg_point* point, unsigned int add)
01039             {
01040               ListType::iterator iter = pointlist.find(point);
01041               // If point is not in the list then stop.
01042               if (iter == pointlist.end()) return;
01043 
01044               // Save iterator and decrement
01045               ListType::iterator iter2 = iter;
01046               iter2--;
01047 
01048               // Point has to be erased first as changing the value does not re-order the std::set
01049               pointlist.erase(iter);
01050               point->add_influence(add);
01051 
01052               // Insert point back into the list. Using the iterator improves performance. The new position has to be at the same position or to the right of the old.
01053               pointlist.insert(iter2,point);
01054             }
01055             // Make *point to C point and remove from sorted list
01056             void make_cpoint(amg_point* point)
01057             {
01058               pointlist.erase(point);
01059               point->make_cpoint();
01060               c_points++;
01061             }
01062             // Make *point to F point and remove from sorted list
01063             void make_fpoint(amg_point* point)
01064             {
01065               pointlist.erase(point);
01066               point->make_fpoint();
01067               f_points++;
01068             }
01069             // Swich *point from F to C point
01070             void switch_ftoc(amg_point* point)
01071             {
01072               point->switch_ftoc();
01073               c_points++;
01074               f_points--;
01075             }
01076 
01077             // Build vector of indices for C point on the coarse level.
01078             void build_index()
01079             {
01080               unsigned int count = 0;
01081               // Use simple counter for index creation.
01082               for (iterator iter = pointvector.begin(); iter != pointvector.end(); ++iter)
01083               {
01084                 // Set index on coarse level using counter variable
01085                 if ((*iter)->is_cpoint())
01086                 {
01087                   (*iter)->set_coarse_index(count);
01088                   count++;
01089                 }
01090               }
01091             }
01092 
01093             // Return information for debugging purposes
01094             template <typename MatrixType>
01095             void get_influence_matrix(MatrixType & mat) const
01096             {
01097               mat = MatrixType(size(),size());
01098               mat.clear();
01099 
01100               for (const_iterator row_iter = begin(); row_iter != end(); ++row_iter)
01101                 for (amg_point::iterator col_iter = (*row_iter)->begin_influencing(); col_iter != (*row_iter)->end_influencing(); ++col_iter)
01102                   mat((*row_iter)->get_index(),(*col_iter)->get_index()) = true;
01103             }
01104             template <typename VectorType>
01105             void get_influence(VectorType & vec) const
01106             {
01107               vec = VectorType(size_);
01108               vec.clear();
01109 
01110               for (const_iterator iter = begin(); iter != end(); ++iter)
01111                 vec[(*iter)->get_index()] = (*iter)->get_influence();
01112             }
01113             template <typename VectorType>
01114             void get_sorting(VectorType & vec) const
01115             {
01116               vec = VectorType(pointlist.size());
01117               vec.clear();
01118               unsigned int i=0;
01119 
01120               for (ListType::const_iterator iter = pointlist.begin(); iter != pointlist.end(); ++iter)
01121               {
01122                 vec[i] = (*iter)->get_index();
01123                 i++;
01124               }
01125             }
01126             template <typename VectorType>
01127             void get_C(VectorType & vec) const
01128             {
01129               vec = VectorType(size_);
01130               vec.clear();
01131 
01132               for (const_iterator iter = begin(); iter != end(); ++iter)
01133               {
01134                 if ((*iter)->is_cpoint())
01135                   vec[(*iter)->get_index()] = true;
01136               }
01137             }
01138             template <typename VectorType>
01139             void get_F(VectorType & vec) const
01140             {
01141               vec = VectorType(size_);
01142               vec.clear();
01143 
01144               for (const_iterator iter = begin(); iter != end(); ++iter)
01145               {
01146                 if ((*iter)->is_fpoint())
01147                   vec[(*iter)->get_index()] = true;
01148               }
01149             }
01150             template <typename MatrixType>
01151             void get_Aggregates(MatrixType & mat) const
01152             {
01153               mat = MatrixType(size_,size_);
01154               mat.clear();
01155 
01156               for (const_iterator iter = begin(); iter != end(); ++iter)
01157               {
01158                 if (!(*iter)->is_undecided())
01159                   mat((*iter)->get_aggregate(),(*iter)->get_index()) = true;
01160               }
01161             }
01162         };
01163 
01167         template <typename InternalType1, typename InternalType2>
01168         class amg_slicing
01169         {
01170             typedef typename InternalType1::value_type SparseMatrixType;
01171             typedef typename InternalType2::value_type PointVectorType;
01172 
01173           public:
01174             // Data structures on a per-processor basis.
01175             boost::numeric::ublas::vector<InternalType1> A_slice;
01176             boost::numeric::ublas::vector<InternalType2> Pointvector_slice;
01177             // Holds the offsets showing the indices for which a new slice begins.
01178             boost::numeric::ublas::vector<boost::numeric::ublas::vector<unsigned int> > Offset;
01179 
01180             unsigned int threads_;
01181             unsigned int levels_;
01182 
01183             void init(unsigned int levels, unsigned int threads = 0)
01184             {
01185               // Either use the number of threads chosen by the user or the maximum number of threads available on the processor.
01186               if (threads == 0)
01187             #ifdef VIENNACL_WITH_OPENMP
01188                 threads_ = omp_get_num_procs();
01189             #else
01190               threads_ = 1;
01191             #endif
01192               else
01193                 threads_ = threads;
01194 
01195               levels_ = levels;
01196 
01197               A_slice.resize(threads_);
01198               Pointvector_slice.resize(threads_);
01199               // Offset has threads_+1 entries to also hold the total size
01200               Offset.resize(threads_+1);
01201 
01202               for (unsigned int i=0; i<threads_; ++i)
01203               {
01204                 A_slice[i].resize(levels_);
01205                 Pointvector_slice[i].resize(levels_);
01206                 // Offset needs one more level for the build-up of the next offset
01207                 Offset[i].resize(levels_+1);
01208               }
01209               Offset[threads_].resize(levels_+1);
01210             } //init()
01211 
01212             // Slice matrix A into as many parts as threads are used.
01213             void slice (unsigned int level, InternalType1 const & A, InternalType2 const & Pointvector)
01214             {
01215               // On the finest level, build a new slicing first.
01216               if (level == 0)
01217                 slice_new (level, A);
01218 
01219               // On coarser levels use the same slicing as on the finest level (Points stay together on the same thread on all levels).
01220               // This is necessary as due to interpolation and galerkin product there only exist connections between points on the same thread on coarser levels.
01221               // Note: Offset is determined in amg_coarse_rs0() after fine level was built.
01222               slice_build (level, A, Pointvector);
01223             }
01224 
01225             // Join point data structure into Pointvector
01226             void join (unsigned int level, InternalType2 & Pointvector) const
01227             {
01228               typedef typename InternalType2::value_type PointVectorType;
01229 
01230               // Reset index offset of all points and update overall C and F point count
01231               Pointvector[level].clear_cf();
01232               for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
01233               {
01234                 (*iter)->set_offset(0);
01235                 Pointvector[level].update_cf(*iter);
01236               }
01237             }
01238 
01239           private:
01244             void slice_new (unsigned int level, InternalType1 const & A)
01245             {
01246               // Determine index offset of all the slices (index of A[level] when the respective slice starts).
01247             #ifdef VIENNACL_WITH_OPENMP
01248               #pragma omp parallel for
01249             #endif
01250               for (long i=0; i<=static_cast<long>(threads_); ++i)
01251               {
01252                 // Offset of first piece is zero. Pieces 1,...,threads-1 have equal size while the last one might be greater.
01253                 if (i == 0) Offset[i][level] = 0;
01254                 else if (i == threads_) Offset[i][level] = static_cast<unsigned int>(A[level].size1());
01255                 else Offset[i][level] = static_cast<unsigned int>(i*(A[level].size1()/threads_));
01256               }
01257             }
01258 
01264             void slice_build (unsigned int level, InternalType1 const & A, InternalType2 const & Pointvector)
01265             {
01266               typedef typename SparseMatrixType::const_iterator1 ConstRowIterator;
01267               typedef typename SparseMatrixType::const_iterator2 ConstColIterator;
01268 
01269               unsigned int x, y;
01270               amg_point *point;
01271 
01272             #ifdef VIENNACL_WITH_OPENMP
01273               #pragma omp parallel for private (x,y,point)
01274             #endif
01275               for (long i=0; i<static_cast<long>(threads_); ++i)
01276               {
01277                 // Allocate space for the matrix slice and the pointvector.
01278                 A_slice[i][level] = SparseMatrixType(Offset[i+1][level]-Offset[i][level],Offset[i+1][level]-Offset[i][level]);
01279                 Pointvector_slice[i][level] = PointVectorType(Offset[i+1][level]-Offset[i][level]);
01280 
01281                 // Iterate over the part that belongs to thread i (from Offset[i][level] to Offset[i+1][level]).
01282                 ConstRowIterator row_iter = A[level].begin1();
01283                 row_iter += Offset[i][level];
01284                 x = static_cast<unsigned int>(row_iter.index1());
01285 
01286                 while (x < Offset[i+1][level] && row_iter != A[level].end1())
01287                 {
01288                   // Set offset for point index and save point for the respective thread
01289                   point = Pointvector[level][x];
01290                   point->set_offset(Offset[i][level]);
01291                   Pointvector_slice[i][level].add_point(point);
01292 
01293                   ConstColIterator col_iter = row_iter.begin();
01294                   y = static_cast<unsigned int>(col_iter.index2());
01295 
01296                   // Save all coefficients from the matrix slice
01297                   while (y < Offset[i+1][level] && col_iter != row_iter.end())
01298                   {
01299                     if (y >= Offset[i][level])
01300                 A_slice[i][level](x-Offset[i][level],y-Offset[i][level]) = *col_iter;
01301 
01302                     ++col_iter;
01303                     y = static_cast<unsigned int>(col_iter.index2());
01304                   }
01305 
01306                   ++row_iter;
01307                   x = static_cast<unsigned int>(row_iter.index1());
01308                 }
01309               }
01310             }
01311         };
01312 
01318         template <typename SparseMatrixType>
01319         void amg_mat_prod (SparseMatrixType & A, SparseMatrixType & B, SparseMatrixType & RES)
01320         {
01321           typedef typename SparseMatrixType::value_type ScalarType;
01322           typedef typename SparseMatrixType::iterator1 InternalRowIterator;
01323           typedef typename SparseMatrixType::iterator2 InternalColIterator;
01324 
01325           long x,y,z;
01326           ScalarType prod;
01327           RES = SparseMatrixType(static_cast<unsigned int>(A.size1()), static_cast<unsigned int>(B.size2()));
01328           RES.clear();
01329 
01330     #ifdef VIENNACL_WITH_OPENMP
01331           #pragma omp parallel for private (x,y,z,prod)
01332     #endif
01333           for (x=0; x<static_cast<long>(A.size1()); ++x)
01334           {
01335             InternalRowIterator row_iter = A.begin1();
01336             row_iter += x;
01337             for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
01338             {
01339               y = static_cast<unsigned int>(col_iter.index2());
01340               InternalRowIterator row_iter2 = B.begin1();
01341               row_iter2 += y;
01342 
01343               for(InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2)
01344               {
01345                 z = static_cast<unsigned int>(col_iter2.index2());
01346                 prod = *col_iter * *col_iter2;
01347                 RES.add(x,z,prod);
01348               }
01349             }
01350           }
01351         }
01352 
01358         template <typename SparseMatrixType>
01359         void amg_galerkin_prod (SparseMatrixType & A, SparseMatrixType & P, SparseMatrixType & RES)
01360         {
01361           typedef typename SparseMatrixType::value_type ScalarType;
01362           typedef typename SparseMatrixType::iterator1 InternalRowIterator;
01363           typedef typename SparseMatrixType::iterator2 InternalColIterator;
01364 
01365           long x,y1,y2,z;
01366           amg_sparsevector<ScalarType> row;
01367           RES = SparseMatrixType(static_cast<unsigned int>(P.size2()), static_cast<unsigned int>(P.size2()));
01368           RES.clear();
01369 
01370     #ifdef VIENNACL_WITH_OPENMP
01371           #pragma omp parallel for private (x,y1,y2,z,row)
01372     #endif
01373           for (x=0; x<static_cast<long>(P.size2()); ++x)
01374           {
01375             row = amg_sparsevector<ScalarType>(static_cast<unsigned int>(A.size2()));
01376             InternalRowIterator row_iter = P.begin1(true);
01377             row_iter += x;
01378 
01379             for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
01380             {
01381               y1 = static_cast<long>(col_iter.index2());
01382               InternalRowIterator row_iter2 = A.begin1();
01383               row_iter2 += y1;
01384 
01385               for(InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2)
01386               {
01387                 y2 = static_cast<long>(col_iter2.index2());
01388                 row.add (y2, *col_iter * *col_iter2);
01389               }
01390             }
01391             for (typename amg_sparsevector<ScalarType>::iterator iter = row.begin(); iter != row.end(); ++iter)
01392             {
01393               y2 = iter.index();
01394               InternalRowIterator row_iter3 = P.begin1();
01395               row_iter3 += y2;
01396 
01397               for (InternalColIterator col_iter3 = row_iter3.begin(); col_iter3 != row_iter3.end(); ++col_iter3)
01398               {
01399                 z = static_cast<long>(col_iter3.index2());
01400                 RES.add (x, z, *col_iter3 * *iter);
01401               }
01402             }
01403           }
01404 
01405           #ifdef VIENNACL_AMG_DEBUG
01406           std::cout << "Galerkin Operator: " << std::endl;
01407           printmatrix (RES);
01408           #endif
01409         }
01410 
01416         template <typename SparseMatrixType>
01417         void test_triplematprod(SparseMatrixType & A, SparseMatrixType & P, SparseMatrixType  & A_i1)
01418         {
01419           typedef typename SparseMatrixType::value_type ScalarType;
01420 
01421           boost::numeric::ublas::compressed_matrix<ScalarType> A_temp (A.size1(), A.size2());
01422           A_temp = A;
01423           boost::numeric::ublas::compressed_matrix<ScalarType> P_temp (P.size1(), P.size2());
01424           P_temp = P;
01425           P.set_trans(true);
01426           boost::numeric::ublas::compressed_matrix<ScalarType> R_temp (P.size1(), P.size2());
01427           R_temp = P;
01428           P.set_trans(false);
01429 
01430           boost::numeric::ublas::compressed_matrix<ScalarType> RA (R_temp.size1(),A_temp.size2());
01431           RA = boost::numeric::ublas::prod(R_temp,A_temp);
01432           boost::numeric::ublas::compressed_matrix<ScalarType> RAP (RA.size1(),P_temp.size2());
01433           RAP = boost::numeric::ublas::prod(RA,P_temp);
01434 
01435           for (unsigned int x=0; x<RAP.size1(); ++x)
01436           {
01437             for (unsigned int y=0; y<RAP.size2(); ++y)
01438             {
01439               if (std::fabs(static_cast<ScalarType>(RAP(x,y)) - static_cast<ScalarType>(A_i1(x,y))) > 0.0001)
01440                 std::cout << x << " " << y << " " << RAP(x,y) << " " << A_i1(x,y) << std::endl;
01441             }
01442           }
01443         }
01444 
01450         template <typename SparseMatrixType, typename PointVectorType>
01451         void test_interpolation(SparseMatrixType & A, SparseMatrixType & P, PointVectorType & Pointvector)
01452         {
01453           for (unsigned int i=0; i<P.size1(); ++i)
01454           {
01455             if (Pointvector.is_cpoint(i))
01456             {
01457               bool set = false;
01458               for (unsigned int j=0; j<P.size2(); ++j)
01459               {
01460                 if (P.isnonzero(i,j))
01461                 {
01462                   if (P(i,j) != 1)
01463                     std::cout << "Error 1 in row " << i << std::endl;
01464                   if (P(i,j) == 1 && set)
01465                     std::cout << "Error 2 in row " << i << std::endl;
01466                   if (P(i,j) == 1 && !set)
01467                     set = true;
01468                 }
01469               }
01470             }
01471 
01472             if (Pointvector.is_fpoint(i))
01473               for (unsigned int j=0; j<P.size2(); ++j)
01474               {
01475                 if (P.isnonzero(i,j) && j> Pointvector.get_cpoints()-1)
01476                   std::cout << "Error 3 in row " << i << std::endl;
01477                 if (P.isnonzero(i,j))
01478                 {
01479                   bool set = false;
01480                   for (unsigned int k=0; k<P.size1(); ++k)
01481                   {
01482                     if (P.isnonzero(k,j))
01483                     {
01484                       if (Pointvector.is_cpoint(k) && P(k,j) == 1 && A.isnonzero(i,k))
01485                         set = true;
01486                     }
01487                   }
01488                   if (!set)
01489                     std::cout << "Error 4 in row " << i << std::endl;
01490                 }
01491               }
01492             }
01493         }
01494 
01495 
01496       } //namespace amg
01497     }
01498   }
01499 }
01500 
01501 #endif