• Main Page
  • Namespaces
  • Data Structures
  • Files
  • File List
  • Globals

/data/development/ViennaCL/dev/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-2011, Institute for Microelectronics,
00006                             Institute for Analysis and Scientific Computing,
00007                             TU Wien.
00008 
00009                             -----------------
00010                   ViennaCL - The Vienna Computing Library
00011                             -----------------
00012 
00013    Project Head:    Karl Rupp                   rupp@iue.tuwien.ac.at
00014                
00015    (A list of authors and contributors can be found in the PDF manual)
00016 
00017    License:         MIT (X11), see file LICENSE in the base directory
00018 ============================================================================= */
00019 
00026 #include <boost/numeric/ublas/vector.hpp>
00027 #include <cmath>
00028 #include <set>
00029 #include <list>
00030 #include <algorithm>
00031 
00032 #include <map>
00033 #ifdef _OPENMP
00034 #include <omp.h>
00035 #endif
00036 
00037 #include "amg_debug.hpp"
00038 
00039 #define VIENNACL_AMG_COARSE_RS 1
00040 #define VIENNACL_AMG_COARSE_ONEPASS 2
00041 #define VIENNACL_AMG_COARSE_RS0 3
00042 #define VIENNACL_AMG_COARSE_RS3 4
00043 #define VIENNACL_AMG_COARSE_AG 5
00044 #define VIENNACL_AMG_INTERPOL_DIRECT 1
00045 #define VIENNACL_AMG_INTERPOL_CLASSIC 2
00046 #define VIENNACL_AMG_INTERPOL_AG 3
00047 #define VIENNACL_AMG_INTERPOL_SA 4
00048 
00049 namespace viennacl
00050 {
00051   namespace linalg
00052   {
00053     namespace detail
00054     {
00055       namespace amg
00056       {
00059         class amg_tag
00060         {
00061           public:
00074             amg_tag(unsigned int coarse = 1,
00075                     unsigned int interpol = 1,
00076                     double threshold = 0.25,
00077                     double interpolweight = 0.2,
00078                     double jacobiweight = 1,
00079                     unsigned int presmooth = 1,
00080                     unsigned int postsmooth = 1,
00081                     unsigned int coarselevels = 0)
00082             : _coarse(coarse), _interpol(interpol),
00083               _threshold(threshold), _interpolweight(interpolweight), _jacobiweight(jacobiweight), 
00084               _presmooth(presmooth), _postsmooth(postsmooth), _coarselevels(coarselevels) {}; 
00085 
00086             // Getter-/Setter-Functions
00087             void set_coarse(unsigned int coarse) { if (coarse > 0) _coarse = coarse; }
00088             unsigned int get_coarse() const { return _coarse; }
00089             
00090             void set_interpol(unsigned int interpol) { if (interpol > 0) _interpol = interpol; }
00091             unsigned int get_interpol() const { return _interpol; }
00092             
00093             void set_threshold(double threshold) { if (threshold > 0 && threshold <= 1) _threshold = threshold; }
00094             double get_threshold() const{ return _threshold; }
00095             
00096             void set_as(double jacobiweight) { if (jacobiweight > 0 && jacobiweight <= 2) _jacobiweight = jacobiweight; }
00097             double get_interpolweight() const { return _interpolweight; }
00098             
00099             void set_interpolweight(double interpolweight) { if (interpolweight > 0 && interpolweight <= 2) _interpolweight = interpolweight; }
00100             double get_jacobiweight() const { return _jacobiweight; }
00101             
00102             void set_presmooth(int presmooth) { if (presmooth >= 0) _presmooth = presmooth; }
00103             unsigned int get_presmooth() const { return _presmooth; }
00104             
00105             void set_postsmooth(int postsmooth) { if (postsmooth >= 0) _postsmooth = postsmooth; }
00106             unsigned int get_postsmooth() const { return _postsmooth; }
00107             
00108             void set_coarselevels(int coarselevels)  { if (coarselevels >= 0) _coarselevels = coarselevels; }
00109             unsigned int get_coarselevels() const { return _coarselevels; }
00110 
00111           private:
00112             unsigned int _coarse, _interpol;
00113             double _threshold, _interpolweight, _jacobiweight;
00114             unsigned int _presmooth, _postsmooth, _coarselevels;
00115         };
00116         
00121         template <typename InternalType, typename IteratorType, typename ScalarType>
00122         class amg_nonzero_scalar
00123         {
00124           private:
00125             InternalType *_m;
00126             IteratorType _iter;
00127             unsigned int _i,_j;
00128             ScalarType _s;
00129 
00130           public:
00131             amg_nonzero_scalar();
00132             
00140             amg_nonzero_scalar(InternalType *m,
00141                               IteratorType & iter,
00142                               unsigned int i,
00143                               unsigned int j,
00144                               ScalarType s = 0): _m(m), _iter(iter), _i(i), _j(j), _s(s) {}
00145             
00149             ScalarType operator = (const ScalarType value)
00150             {
00151               _s = value;
00152               // Only write if scalar is nonzero
00153               if (_s == 0) return _s;
00154               // Write to _m using iterator _iter or indices (_i,_j)
00155               _m->addscalar (_iter,_i,_j,_s);
00156               return _s;
00157             }
00158             
00162             ScalarType operator += (const ScalarType value)
00163             {
00164               // If zero is added, then no change necessary
00165               if (value == 0)
00166                 return _s;
00167               
00168               _s += value;
00169               // Remove entry if resulting scalar is zero
00170               if (_s == 0)
00171               {
00172                 _m->removescalar(_iter,_i);
00173                 return _s;
00174               }
00175               //Write to _m using iterator _iter or indices (_i,_j)
00176               _m->addscalar (_iter,_i,_j,_s);
00177               return _s;
00178             }
00179             ScalarType operator ++ (int)
00180             {
00181               _s++;
00182               if (_s == 0)
00183                 _m->removescalar(_iter,_i);
00184               _m->addscalar (_iter,_i,_j,_s);
00185               return _s;
00186             }
00187             ScalarType operator ++ ()
00188             {
00189               _s++;
00190               if (_s == 0)
00191                 _m->removescalar(_iter,_i);
00192               _m->addscalar (_iter,_i,_j,_s);
00193               return _s;
00194             }
00195             operator ScalarType (void) { return _s;  }
00196         };
00197     
00200         template <typename InternalType>
00201         class amg_sparsevector_iterator
00202         {
00203           private:
00204             typedef amg_sparsevector_iterator<InternalType> self_type;  
00205             typedef typename InternalType::mapped_type ScalarType;
00206             
00207             InternalType & internal_vec;
00208             typename InternalType::iterator iter;
00209               
00210           public:
00211             
00216             amg_sparsevector_iterator(InternalType & vec, bool begin=true): internal_vec(vec)
00217             {
00218               if (begin)
00219                 iter = internal_vec.begin();
00220               else
00221                 iter = internal_vec.end();
00222             }
00223             
00224             bool operator == (self_type other)
00225             {
00226               if (iter == other.iter)
00227                 return true;
00228               else
00229                 return false;
00230             }
00231             bool operator != (self_type other)
00232             {
00233               if (iter != other.iter)
00234                 return true;
00235               else
00236                 return false;
00237             }
00238             
00239             self_type & operator ++ () const { iter++; return *this; }
00240             self_type & operator ++ () { iter++; return *this; }
00241             self_type & operator -- () const { iter--; return *this; }
00242             self_type & operator -- () { iter--; return *this; }  
00243             ScalarType & operator * () const { return (*iter).second; }
00244             ScalarType & operator * () { return (*iter).second; }
00245             unsigned int index() const { return (*iter).first; }
00246             unsigned int index() { return (*iter).first; }
00247         };
00248     
00251         template <typename ScalarType>
00252         class amg_sparsevector
00253         {
00254           public:
00255             typedef ScalarType value_type;
00256       
00257           private:
00258             // A map is used internally which saves all non-zero elements with pairs of (index,value)
00259             typedef std::map<unsigned int,ScalarType> InternalType;
00260             typedef amg_sparsevector<ScalarType> self_type;
00261             typedef amg_nonzero_scalar<self_type,typename InternalType::iterator,ScalarType> NonzeroScalarType;
00262               
00263             // Size is only a dummy variable. Not needed for internal map structure but for compatible vector interface.
00264             unsigned int _size;
00265             InternalType internal_vector;
00266       
00267           public:
00268             typedef amg_sparsevector_iterator<InternalType> iterator;
00269             typedef typename InternalType::const_iterator const_iterator;
00270       
00271           public:
00275             amg_sparsevector(unsigned int size = 0): _size(size)
00276             {
00277               internal_vector = InternalType();
00278             }
00279             
00280             void resize(unsigned int size) { _size = size; }
00281             unsigned int size() const { return _size;}
00282             
00283             // Returns number of non-zero entries in vector equal to the size of the underlying map.
00284             unsigned int internal_size() const { return internal_vector.size(); }
00285             // Delete underlying map.
00286             void clear() { internal_vector.clear();  }
00287             // Remove entry at position i.
00288             void remove(unsigned int i) { internal_vector.erase(i); }
00289             
00290             // Add s to the entry at position i
00291             void add (unsigned int i, ScalarType s)
00292             {
00293               typename InternalType::iterator iter = internal_vector.find(i);
00294               // If there is no entry at position i, add new entry at that position
00295               if (iter == internal_vector.end())
00296                 addscalar(iter,i,i,s);
00297               else
00298               {
00299                 s += (*iter).second;
00300                 // If new value is zero, then erase the entry, otherwise write new value
00301                 if (s != 0)
00302                   (*iter).second = s;
00303                 else
00304                   internal_vector.erase(iter);
00305               }
00306             }
00307             
00308             // Write to the map. Is called from non-zero scalar type.
00309             template <typename IteratorType>
00310             void addscalar(IteratorType & iter, unsigned int i, unsigned int j, ScalarType s)
00311             {
00312               // Don't write if value is zero
00313               if (s == 0)
00314                 return;
00315               
00316               // If entry is already present, overwrite value, otherwise make new entry
00317               if (iter != internal_vector.end())  
00318                 (*iter).second = s;
00319               else
00320                 internal_vector[i] = s;
00321             }
00322             
00323             // Remove value from the map. Is called from non-zero scalar type.
00324             template <typename IteratorType>
00325             void removescalar(IteratorType & iter, unsigned int i) { internal_vector.erase(iter); }   
00326             
00327             // Bracket operator. Returns non-zero scalar type with actual values of the respective entry which calls addscalar/removescalar after value is altered.
00328             NonzeroScalarType operator [] (unsigned int i)
00329             {
00330               typename InternalType::iterator it = internal_vector.find(i);
00331               // If value is already present then build non-zero scalar with actual value, otherwise 0.
00332               if (it != internal_vector.end())
00333                 return NonzeroScalarType (this,it,i,i,(*it).second);
00334               else
00335                 return NonzeroScalarType (this,it,i,i,0);  
00336             }
00337             
00338             // Use internal data structure directly for read-only access. No need to use non-zero scalar as no write access possible.
00339             ScalarType operator [] (unsigned int i) const
00340             {
00341               const_iterator it = internal_vector.find(i);
00342               
00343               if (it != internal_vector.end())
00344                 return (*it).second;
00345               else
00346                 return 0;
00347             }
00348             
00349             // Iterator functions.
00350             iterator begin() { return iterator(internal_vector); }
00351             const_iterator begin() const { return internal_vector.begin(); }
00352             iterator end() { return iterator(internal_vector,false); }
00353             const_iterator end() const { return internal_vector.end(); }
00354             
00355             // checks whether value at index i is nonzero. More efficient than doing [] == 0.
00356             bool isnonzero(unsigned int i) const { return internal_vector.find(i) != internal_vector.end();  }
00357             
00358             // Copies data into a ublas vector type.
00359             operator boost::numeric::ublas::vector<ScalarType> (void)
00360             {
00361               boost::numeric::ublas::vector<ScalarType> vec (_size);    
00362               for (iterator iter = begin(); iter != end(); ++iter)
00363                 vec [iter.index()] = *iter;        
00364               return vec;
00365             } 
00366          };
00367     
00373         template <typename ScalarType>
00374         class amg_sparsematrix
00375         {
00376           public:
00377             typedef ScalarType value_type;
00378           private:
00379             typedef std::map<unsigned int,ScalarType> RowType;
00380             typedef std::vector<RowType> InternalType;
00381             typedef amg_sparsematrix<ScalarType> self_type;
00382             
00383             // Adapter is used for certain functionality, especially iterators.
00384             typedef typename viennacl::tools::sparse_matrix_adapter<ScalarType> AdapterType;
00385             typedef typename viennacl::tools::const_sparse_matrix_adapter<ScalarType> ConstAdapterType;
00386             
00387             // Non-zero scalar is used to write to the matrix.
00388             typedef amg_nonzero_scalar<self_type,typename RowType::iterator,ScalarType> NonzeroScalarType;
00389 
00390             // Holds matrix coefficients.
00391             InternalType internal_mat;
00392             // Holds matrix coefficient of transposed matrix if built. 
00393             // Note: Only internal_mat is written using operators and methods while internal_mat_trans is built from internal_mat using do_trans().
00394             InternalType internal_mat_trans;
00395             // Saves sizes.
00396             size_t s1, s2;
00397             
00398             // True if the transposed of the matrix is used (for calculations, iteration, etc.).
00399             bool transposed_mode;
00400             // True if the transposed is already built (saved in internal_mat_trans) and also up to date (no changes to internal_mat).
00401             bool transposed;
00402             
00403           public:          
00404             typedef typename AdapterType::iterator1 iterator1;
00405             typedef typename AdapterType::iterator2 iterator2;
00406             typedef typename ConstAdapterType::const_iterator1 const_iterator1;
00407             typedef typename ConstAdapterType::const_iterator2 const_iterator2;
00408             
00410             amg_sparsematrix ()
00411             {
00412               transposed_mode = false;
00413               transposed = false;
00414             }
00415             
00420             amg_sparsematrix (unsigned int i, unsigned int j)
00421             {
00422               AdapterType a (internal_mat);
00423               a.resize (i,j,false);
00424               AdapterType a_trans (internal_mat_trans);
00425               a_trans.resize (j,i,false);
00426               s1 = i;
00427               s2 = j;
00428               a.clear();
00429               a_trans.clear();
00430               transposed_mode = false;
00431               transposed = false;
00432             }
00433             
00438             amg_sparsematrix (std::vector<std::map<unsigned int, ScalarType> > const & mat)
00439             {  
00440               AdapterType a (internal_mat);
00441               AdapterType a_trans (internal_mat_trans);
00442               a.resize(mat.size(), mat.size());
00443               a_trans.resize(mat.size(), mat.size());
00444               
00445               internal_mat = mat;  
00446               s1 = s2 = mat.size();
00447               
00448               transposed_mode = false;
00449               transposed = false;
00450             }
00451             
00456             template <typename MatrixType>
00457             amg_sparsematrix (MatrixType const & mat)
00458             {  
00459               AdapterType a (internal_mat);
00460               AdapterType a_trans (internal_mat_trans);
00461               a.resize(mat.size1(), mat.size2());
00462               a_trans.resize (mat.size2(), mat.size1());
00463               s1 = mat.size1();
00464               s2 = mat.size2();
00465               a.clear();
00466               a_trans.clear();
00467               
00468               for (typename MatrixType::const_iterator1 row_iter = mat.begin1(); row_iter != mat.end1(); ++row_iter)
00469               {
00470                 for (typename MatrixType::const_iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00471                 {
00472                   if (*col_iter != 0)
00473                   {
00474                     unsigned int x = col_iter.index1();
00475                     unsigned int y = col_iter.index2();
00476                     a (x,y) = *col_iter;
00477                     a_trans (y,x) = *col_iter;
00478                   }
00479                 }
00480               }
00481               transposed_mode = false;
00482               transposed = true;
00483             }
00484                   
00485             // Build transposed of the current matrix.
00486             void do_trans()
00487             {
00488               // Do it only once if called in a parallel section
00489             #ifdef _OPENMP
00490               #pragma omp critical
00491             #endif
00492               { 
00493                 // Only build transposed if it is not built or not up to date
00494                 if (!transposed)
00495                 {
00496                   // Mode has to be set to standard mode temporarily
00497                   bool save_mode = transposed_mode;
00498                   transposed_mode = false;
00499                   
00500                   for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00501                 for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00502                   internal_mat_trans[col_iter.index2()][col_iter.index1()] = *col_iter;
00503                 
00504                   transposed_mode = save_mode;
00505                   transposed = true;
00506                 }
00507               }
00508             } //do_trans()
00509             
00510             // Set transposed mode (true=transposed, false=regular)
00511             void set_trans(bool mode)
00512             {
00513               transposed_mode = mode;
00514               if (mode)
00515                 do_trans();
00516             }
00517             
00518             bool get_trans() const { return transposed_mode; }     
00519                   
00520             // Checks whether coefficient (i,j) is non-zero. More efficient than using (i,j) == 0.
00521             bool isnonzero (unsigned int i, unsigned int j) const
00522             {
00523               if (!transposed_mode)
00524               {
00525                 if (internal_mat[i].find(j) != internal_mat[i].end())
00526                   return true;
00527                 else
00528                   return false;
00529               }
00530               else
00531               {
00532                 if (internal_mat[j].find(i) != internal_mat[j].end())
00533                   return true;
00534                 else
00535                   return false;
00536               }
00537             } //isnonzero()
00538                 
00539             // Add s to value at (i,j)
00540             void add (unsigned int i, unsigned int j, ScalarType s)
00541             {
00542               // If zero is added then do nothing.
00543               if (s == 0)
00544                 return;
00545               
00546               typename RowType::iterator col_iter = internal_mat[i].find(j);
00547               // If there is no entry at position (i,j), then make new entry.
00548               if (col_iter == internal_mat[i].end())
00549                 addscalar(col_iter,i,j,s);
00550               else
00551               {
00552                 s += (*col_iter).second;
00553                 // Update value and erase entry if value is zero.
00554                 if (s != 0)
00555                   (*col_iter).second = s;
00556                 else
00557                   internal_mat[i].erase(col_iter);
00558               }
00559               transposed = false;
00560             } //add()
00561             
00562             // Write to the internal data structure. Is called from non-zero scalar type.
00563             template <typename IteratorType>
00564             void addscalar(IteratorType & iter, unsigned int i, unsigned int j, ScalarType s)
00565             {
00566               // Don't write if value is zero
00567               if (s == 0)
00568                 return;
00569               
00570               if (iter != internal_mat[i].end())  
00571                 (*iter).second = s;
00572               else
00573                 internal_mat[i][j] = s;
00574               
00575               transposed = false;
00576             }
00577             
00578             // Remove entry from internal data structure. Is called from non-zero scalar type.
00579             template <typename IteratorType>
00580             void removescalar(IteratorType & iter, unsigned int i)
00581             {
00582               internal_mat[i].erase(iter);
00583               transposed = false;
00584             }   
00585             
00586             // Return non-zero scalar at position (i,j). Value is written to the non-zero scalar and updated via addscalar()/removescalar().
00587             NonzeroScalarType operator()(unsigned int i, unsigned int j)
00588             {
00589               typename RowType::iterator iter;
00590               
00591               if (!transposed_mode)
00592               {
00593                 iter = internal_mat[i].find(j);
00594                 if (iter != internal_mat[i].end())
00595                   return NonzeroScalarType (this,iter,i,j,(*iter).second);
00596                 else
00597                   return NonzeroScalarType (this,iter,i,j,0);
00598               }
00599               else
00600               {
00601                 iter = internal_mat[j].find(i);
00602                 if (iter != internal_mat[j].end())
00603                   return NonzeroScalarType (this,iter,j,i,(*iter).second);
00604                 else
00605                   return NonzeroScalarType (this,iter,j,i,0);
00606               }
00607             }
00608             
00609             // For read-only access return the actual value directly. Non-zero datatype not needed as no write access possible.
00610             ScalarType operator()(unsigned int i, unsigned int j) const
00611             {
00612               typename RowType::const_iterator iter;
00613               
00614               if (!transposed_mode)
00615               {
00616                 iter = internal_mat[i].find(j);
00617                 if (iter != internal_mat[i].end())
00618                   return (*iter).second;
00619                 else
00620                   return 0;
00621               }
00622               else
00623               {
00624                 iter = internal_mat[j].find(i);
00625                 if (iter != internal_mat[j].end())
00626                   return (*iter).second;
00627                 else
00628                   return 0;
00629               }
00630             }
00631               
00632             void resize(unsigned int i, unsigned int j, bool preserve = true)
00633             {
00634               AdapterType a (internal_mat);
00635               a.resize(i,j,preserve);
00636               AdapterType a_trans (internal_mat_trans);
00637               a_trans.resize(j,i,preserve);
00638               s1 = i;
00639               s2 = j;
00640             }
00641             
00642             void clear() 
00643             {
00644               AdapterType a (internal_mat);
00645               a.clear();
00646               AdapterType a_trans (internal_mat_trans);
00647               a_trans.clear();
00648               transposed = true;
00649             }
00650 
00651             size_t size1()
00652             {
00653               if (!transposed_mode)
00654                 return s1;
00655               else
00656                 return s2;
00657             }
00658             
00659             size_t size1() const
00660             {
00661               if (!transposed_mode)
00662                 return s1;
00663               else
00664                 return s2;
00665             }
00666             
00667             
00668             size_t size2()
00669             {
00670               if (!transposed_mode)
00671                 return s2;
00672               else
00673                 return s1;
00674             }
00675             size_t size2() const
00676             {
00677               if (!transposed_mode)
00678                 return s2;
00679               else
00680                 return s1;
00681             }
00682             
00683             iterator1 begin1(bool trans = false)
00684             {
00685               if (!trans && !transposed_mode)
00686               {
00687                 AdapterType a (internal_mat);
00688                 return a.begin1();
00689               }
00690               else
00691               {
00692                 do_trans();
00693                 AdapterType a_trans (internal_mat_trans);
00694                 return a_trans.begin1();
00695               }
00696             }
00697             
00698             iterator1 end1(bool trans = false)
00699             {
00700               if (!trans && !transposed_mode)
00701               {
00702                 AdapterType a (internal_mat);
00703                 return a.end1();
00704               }
00705               else
00706               {
00707                 //do_trans();
00708                 AdapterType a_trans (internal_mat_trans);
00709                 return a_trans.end1();
00710               }
00711             }
00712             
00713             iterator2 begin2(bool trans = false)
00714             {
00715               if (!trans && !transposed_mode)
00716               {
00717                 AdapterType a (internal_mat);
00718                 return a.begin2();
00719               }
00720               else
00721               {
00722                 do_trans();
00723                 AdapterType a_trans (internal_mat_trans);
00724                 return a_trans.begin2();
00725               }
00726             }
00727             
00728             iterator2 end2(bool trans = false)
00729             {
00730               if (!trans && !transposed_mode)
00731               {
00732                 AdapterType a (internal_mat);
00733                 return a.end2();
00734               }
00735               else
00736               {
00737                 //do_trans();
00738                 AdapterType a_trans (internal_mat_trans);
00739                 return a_trans.end2();
00740               }
00741             }
00742             
00743             const_iterator1 begin1() const
00744             {
00745               // Const_iterator of transposed can only be used if transposed matrix is already built and up to date.
00746               assert((!transposed_mode || (transposed_mode && transposed)) && "Error: Cannot build const_iterator when transposed has not been built yet!");
00747                     ConstAdapterType a_const (internal_mat);
00748               return a_const.begin1();
00749             }
00750             
00751             const_iterator1 end1(bool trans = false) const
00752             {
00753               assert((!transposed_mode || (transposed_mode && transposed)) && "Error: Cannot build const_iterator when transposed has not been built yet!");
00754               ConstAdapterType a_const (internal_mat);
00755               return a_const.end1();
00756             }
00757             
00758             const_iterator2 begin2(bool trans = false) const
00759             {
00760               assert((!transposed_mode || (transposed_mode && transposed)) && "Error: Cannot build const_iterator when transposed has not been built yet!");
00761               ConstAdapterType a_const (internal_mat);
00762               return a_const.begin2();
00763             }
00764             
00765             const_iterator2 end2(bool trans = false) const
00766             {
00767               assert((!transposed_mode || (transposed_mode && transposed)) && "Error: Cannot build const_iterator when transposed has not been built yet!");
00768               ConstAdapterType a_const (internal_mat);
00769               return a_const.end2();
00770             }
00771             
00772             // Returns pointer to the internal data structure. Improves performance of copy operation to GPU.
00773             std::vector<std::map<unsigned int, ScalarType> > * get_internal_pointer()
00774             {    
00775               if (!transposed_mode)
00776                 return &internal_mat;
00777               
00778               if (!transposed)
00779                 do_trans();
00780               return &internal_mat_trans;
00781             }
00782             operator boost::numeric::ublas::compressed_matrix<ScalarType> (void)
00783             {
00784               boost::numeric::ublas::compressed_matrix<ScalarType> mat;
00785               mat.resize(size1(),size2(),false);
00786               mat.clear();
00787               
00788               for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00789                   for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00790                     mat (col_iter.index1(), col_iter.index2()) = *col_iter;
00791                   
00792               return mat;
00793             } 
00794             operator boost::numeric::ublas::matrix<ScalarType> (void)
00795             {
00796               boost::numeric::ublas::matrix<ScalarType> mat;
00797               mat.resize(size1(),size2(),false);
00798               mat.clear();
00799               
00800               for (iterator1 row_iter = begin1(); row_iter != end1(); ++row_iter)
00801                   for (iterator2 col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00802                     mat (col_iter.index1(), col_iter.index2()) = *col_iter;
00803                   
00804               return mat;
00805             }
00806         };
00807           
00813         class amg_point
00814         {
00815           private:
00816             typedef amg_sparsevector<amg_point*> ListType;
00817             
00818             unsigned int _index;
00819             unsigned int _influence;
00820             // Determines whether point is undecided.
00821             bool _undecided;
00822             // Determines wheter point is C point (true) or F point (false).
00823             bool _cpoint;
00824             unsigned int _coarse_index;
00825             // 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
00826             unsigned int _offset;
00827             // Aggregate the point belongs to.
00828             unsigned int _aggregate;
00829             
00830             // Holds all points influencing this point.
00831             ListType influencing_points;
00832             // Holds all points that are influenced by this point.
00833             ListType influenced_points;
00834       
00835           public:
00836             typedef ListType::iterator iterator;
00837             typedef ListType::const_iterator const_iterator;
00838             
00841             amg_point (unsigned int index, unsigned int size): _index(index), _influence(0), _undecided(true), _cpoint(false), _coarse_index(0), _offset(0), _aggregate(0)
00842             {
00843               influencing_points = ListType(size);
00844               influenced_points = ListType(size);
00845             }
00846             
00847             void set_offset(unsigned int offset) { _offset = offset; }
00848             unsigned int get_offset() { return _offset; }
00849             void set_index(unsigned int index) { _index = index+_offset; }
00850             unsigned int get_index() const { return _index-_offset;  }
00851             unsigned int get_influence() const { return _influence;  }
00852             void set_aggregate(unsigned int aggregate) { _aggregate = aggregate; }
00853             unsigned int get_aggregate () { return _aggregate; }
00854             
00855             bool is_cpoint() const { return _cpoint && !_undecided;  }
00856             bool is_fpoint() const { return !_cpoint && !_undecided; }
00857             bool is_undecided() const { return _undecided; }
00858             
00859             // Returns number of influencing points
00860             unsigned int number_influencing() const  { return influencing_points.internal_size(); }
00861             // Returns true if *point is influencing this point
00862             bool is_influencing(amg_point* point) const { return influencing_points.isnonzero(point->get_index()+_offset); }
00863             // Add *point to influencing points
00864             void add_influencing_point(amg_point* point) { influencing_points[point->get_index()+_offset] = point;  }
00865             // Add *point to influenced points
00866             void add_influenced_point(amg_point* point) { influenced_points[point->get_index()+_offset] = point; }
00867             
00868             // Clear influencing points
00869             void clear_influencing() { influencing_points.clear(); }
00870             // Clear influenced points
00871             void clear_influenced() {influenced_points.clear(); }
00872             
00873             
00874             unsigned int get_coarse_index() const { return _coarse_index; }
00875             void set_coarse_index(unsigned int index) { _coarse_index = index; }
00876             
00877             // Calculates the initial influence measure equal to the number of influenced points.
00878             void calc_influence() { _influence = influenced_points.internal_size();  }
00879             
00880             // Add to influence measure.
00881             unsigned int add_influence(int add)
00882             {
00883               _influence += add;
00884               return _influence;
00885             }
00886             // Make this point C point. Only call via amg_pointvector.
00887             void make_cpoint() 
00888             { 
00889               _undecided = false;
00890               _cpoint = true; 
00891               _influence = 0;
00892             }
00893             // Make this point F point. Only call via amg_pointvector.
00894             void make_fpoint()
00895             {
00896               _undecided = false;
00897               _cpoint = false;
00898               _influence = 0;
00899             }
00900             // Switch point from F to C point. Only call via amg_pointvector.
00901             void switch_ftoc() { _cpoint = true; }  
00902             
00903             // Iterator handling for influencing and influenced points.
00904             iterator begin_influencing() { return influencing_points.begin(); }
00905             iterator end_influencing() { return influencing_points.end(); }
00906             const_iterator begin_influencing() const { return influencing_points.begin(); }
00907             const_iterator end_influencing() const { return influencing_points.end(); }
00908             iterator begin_influenced() { return influenced_points.begin();  }
00909             iterator end_influenced() { return influenced_points.end(); }
00910             const_iterator begin_influenced() const { return influenced_points.begin(); }
00911             const_iterator end_influenced() const { return influenced_points.end(); }
00912         };
00913         
00916         struct classcomp
00917         {
00918           // Function returns true if l comes before r in the ordering.
00919           bool operator() (amg_point* l, amg_point* r) const
00920           {
00921             // Map is sorted by influence number starting with the highest
00922             // If influence number is the same then lowest point index comes first
00923             return (l->get_influence() < r->get_influence() || (l->get_influence() == r->get_influence() && l->get_index() > r->get_index()));
00924           }
00925         };
00926       
00932         class amg_pointvector
00933         {
00934           private:
00935             // Type for the sorted list
00936             typedef std::set<amg_point*,classcomp> ListType;
00937             // Type for the vector of pointers
00938             typedef std::vector<amg_point*> VectorType;
00939             VectorType pointvector;
00940             ListType pointlist;
00941             unsigned int _size;
00942             unsigned int c_points, f_points;
00943       
00944           public:
00945             typedef VectorType::iterator iterator;
00946             typedef VectorType::const_iterator const_iterator;
00947             
00951             amg_pointvector(unsigned int size = 0): _size(size)
00952             {
00953               pointvector = VectorType(size);
00954               c_points = f_points = 0;
00955             }
00956             
00957             // Construct all the points dynamically and save pointers into vector.
00958             void init_points()
00959             {  
00960               for (unsigned int i=0; i<size(); ++i)
00961                 pointvector[i] = new amg_point(i,size());
00962             }
00963             // Delete all the points.
00964             void delete_points()
00965             {  
00966               for (unsigned int i=0; i<size(); ++i)
00967                 delete pointvector[i];
00968             }
00969             // 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!
00970             void add_point(amg_point *point)
00971             {
00972               pointvector[point->get_index()] = point;
00973               if (point->is_cpoint()) c_points++;
00974               else if (point->is_fpoint()) f_points++;
00975             }
00976 
00977             // Update C and F count for point *point.
00978             // Necessary if C and F points were constructed outside this data structure (e.g. by parallel coarsening RS0 or RS3).
00979             void update_cf(amg_point *point) 
00980             {
00981               if (point->is_cpoint()) c_points++;
00982               else if (point->is_fpoint()) f_points++;
00983             }
00984             // Clear the C and F point count.
00985             void clear_cf() { c_points = f_points = 0; }
00986             
00987             // Clear both point lists.
00988             void clear_influencelists()
00989             {
00990               for (iterator iter = pointvector.begin(); iter != pointvector.end(); ++iter)
00991               {
00992                 (*iter)->clear_influencing();
00993                 (*iter)->clear_influenced();
00994               }
00995             }
00996             
00997             amg_point* operator [] (unsigned int i) const { return pointvector[i]; }
00998             iterator begin() { return pointvector.begin(); }
00999             iterator end() { return pointvector.end(); }
01000             const_iterator begin() const { return pointvector.begin(); }
01001             const_iterator end() const { return pointvector.end(); }
01002             
01003             void resize(unsigned int size)
01004             {
01005               _size = size;
01006               pointvector = VectorType(size);
01007             }
01008             unsigned int size() const { return _size; }
01009             
01010             // Returns number of C points
01011             unsigned int get_cpoints() const { return c_points; }
01012             // Returns number of F points
01013             unsigned int get_fpoints() const { return f_points; }
01014             
01015             // Does the initial sorting of points into the list. Sorting is automatically done by the std::set data type.
01016             void sort()
01017             {
01018               for (iterator iter = begin(); iter != end(); ++iter)
01019                 pointlist.insert(*iter);
01020             }
01021             // Returns the point with the highest influence measure
01022             amg_point* get_nextpoint()
01023             {
01024               // No points remaining? Return NULL such that coarsening will stop.
01025               if (pointlist.size() == 0)
01026                 return NULL;
01027               // If point with highest influence measure (end of the list) has measure of zero, then no further C points can be constructed. Return NULL.
01028               if ((*(--pointlist.end()))->get_influence() == 0)
01029                 return NULL;
01030               // Otherwise, return the point with highest influence measure located at the end of the list.
01031               else
01032                 return *(--pointlist.end());
01033             }
01034             // Add "add" to influence measure for point *point in the sorted list.
01035             void add_influence(amg_point* point, unsigned int add)
01036             {
01037               ListType::iterator iter = pointlist.find(point);
01038               // If point is not in the list then stop.
01039               if (iter == pointlist.end()) return;
01040               
01041               // Save iterator and decrement
01042               ListType::iterator iter2 = iter;
01043               iter2--;
01044               
01045               // Point has to be erased first as changing the value does not re-order the std::set
01046               pointlist.erase(iter);
01047               point->add_influence(add);
01048               
01049               // 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.
01050               pointlist.insert(iter2,point);
01051             }
01052             // Make *point to C point and remove from sorted list
01053             void make_cpoint(amg_point* point)
01054             {
01055               pointlist.erase(point);
01056               point->make_cpoint();
01057               c_points++;
01058             }
01059             // Make *point to F point and remove from sorted list
01060             void make_fpoint(amg_point* point)
01061             {
01062               pointlist.erase(point);
01063               point->make_fpoint();
01064               f_points++;
01065             }
01066             // Swich *point from F to C point
01067             void switch_ftoc(amg_point* point)
01068             { 
01069               point->switch_ftoc();
01070               c_points++;
01071               f_points--;
01072             }
01073 
01074             // Build vector of indices for C point on the coarse level.
01075             void build_index()
01076             {
01077               unsigned int count = 0;
01078               // Use simple counter for index creation.
01079               for (iterator iter = pointvector.begin(); iter != pointvector.end(); ++iter)
01080               {
01081                 // Set index on coarse level using counter variable
01082                 if ((*iter)->is_cpoint())
01083                 {
01084                   (*iter)->set_coarse_index(count);
01085                   count++;
01086                 }
01087               }
01088             }
01089             
01090             // Return information for debugging purposes
01091             template <typename MatrixType>
01092             void get_influence_matrix(MatrixType & mat) const
01093             {
01094               mat = MatrixType(size(),size());
01095               mat.clear();
01096               
01097               for (const_iterator row_iter = begin(); row_iter != end(); ++row_iter)
01098                 for (amg_point::iterator col_iter = (*row_iter)->begin_influencing(); col_iter != (*row_iter)->end_influencing(); ++col_iter)
01099                   mat((*row_iter)->get_index(),(*col_iter)->get_index()) = true;  
01100             }
01101             template <typename VectorType>
01102             void get_influence(VectorType & vec) const
01103             {
01104               vec = VectorType(_size);
01105               vec.clear();
01106               
01107               for (const_iterator iter = begin(); iter != end(); ++iter)
01108                 vec[(*iter)->get_index()] = (*iter)->get_influence();
01109             }
01110             template <typename VectorType>
01111             void get_sorting(VectorType & vec) const
01112             {
01113               vec = VectorType(pointlist.size());
01114               vec.clear();
01115               unsigned int i=0;
01116               
01117               for (ListType::const_iterator iter = pointlist.begin(); iter != pointlist.end(); ++iter)
01118               {
01119                 vec[i] = (*iter)->get_index();
01120                 i++;
01121               }
01122             }
01123             template <typename VectorType>
01124             void get_C(VectorType & vec) const
01125             {
01126               vec = VectorType(_size);
01127               vec.clear();
01128               
01129               for (const_iterator iter = begin(); iter != end(); ++iter)
01130               {
01131                 if ((*iter)->is_cpoint())
01132                   vec[(*iter)->get_index()] = true;
01133               }
01134             }
01135             template <typename VectorType>
01136             void get_F(VectorType & vec) const
01137             {
01138               vec = VectorType(_size);
01139               vec.clear();
01140               
01141               for (const_iterator iter = begin(); iter != end(); ++iter)
01142               {
01143                 if ((*iter)->is_fpoint())
01144                   vec[(*iter)->get_index()] = true;
01145               }
01146             }
01147             template <typename MatrixType>
01148             void get_Aggregates(MatrixType & mat) const
01149             {
01150               mat = MatrixType(_size,_size);
01151               mat.clear();
01152               
01153               for (const_iterator iter = begin(); iter != end(); ++iter)
01154               {
01155                 if (!(*iter)->is_undecided())
01156                   mat((*iter)->get_aggregate(),(*iter)->get_index()) = true;
01157               }
01158             }
01159         };
01160         
01164         template <typename InternalType1, typename InternalType2>
01165         class amg_slicing
01166         {
01167             typedef typename InternalType1::value_type SparseMatrixType;
01168             typedef typename InternalType2::value_type PointVectorType;    
01169             
01170           public:
01171             // Data structures on a per-processor basis.
01172             boost::numeric::ublas::vector<InternalType1> A_slice;
01173             boost::numeric::ublas::vector<InternalType2> Pointvector_slice;
01174             // Holds the offsets showing the indices for which a new slice begins.
01175             boost::numeric::ublas::vector<boost::numeric::ublas::vector<unsigned int> > Offset;
01176             
01177             unsigned int _threads;
01178             unsigned int _levels;
01179             
01180             void init(unsigned int levels, unsigned int threads = 0)
01181             {
01182               // Either use the number of threads chosen by the user or the maximum number of threads available on the processor.
01183               if (threads == 0)
01184             #ifdef _OPENMP
01185                 _threads = omp_get_num_procs();
01186             #else
01187               _threads = 1;
01188             #endif   
01189               else 
01190                 _threads = threads;
01191               
01192               _levels = levels;
01193               
01194               A_slice.resize(_threads);
01195               Pointvector_slice.resize(_threads);
01196               // Offset has _threads+1 entries to also hold the total size
01197               Offset.resize(_threads+1);
01198               
01199               for (unsigned int i=0; i<_threads; ++i)
01200               {
01201                 A_slice[i].resize(_levels);
01202                 Pointvector_slice[i].resize(_levels);
01203                 // Offset needs one more level for the build-up of the next offset
01204                 Offset[i].resize(_levels+1);
01205               }
01206               Offset[_threads].resize(_levels+1);
01207             } //init()
01208             
01209             // Slice matrix A into as many parts as threads are used.
01210             void slice (unsigned int level, InternalType1 const & A, InternalType2 const & Pointvector)
01211             {
01212               // On the finest level, build a new slicing first.
01213               if (level == 0)
01214                 slice_new (level, A);
01215               
01216               // On coarser levels use the same slicing as on the finest level (Points stay together on the same thread on all levels).
01217               // This is necessary as due to interpolation and galerkin product there only exist connections between points on the same thread on coarser levels.
01218               // Note: Offset is determined in amg_coarse_rs0() after fine level was built.
01219               slice_build (level, A, Pointvector);
01220             }
01221 
01222             // Join point data structure into Pointvector
01223             void join (unsigned int level, InternalType2 & Pointvector) const
01224             {
01225               typedef typename InternalType2::value_type PointVectorType;
01226               
01227               // Reset index offset of all points and update overall C and F point count
01228               Pointvector[level].clear_cf();
01229               for (typename PointVectorType::iterator iter = Pointvector[level].begin(); iter != Pointvector[level].end(); ++iter)
01230               {
01231                 (*iter)->set_offset(0);
01232                 Pointvector[level].update_cf(*iter);
01233               }
01234             }
01235               
01236           private:     
01241             void slice_new (unsigned int level, InternalType1 const & A)
01242             {  
01243               typedef typename SparseMatrixType::const_iterator1 ConstRowIterator;
01244               typedef typename SparseMatrixType::const_iterator2 ConstColIterator;
01245               
01246               // Determine index offset of all the slices (index of A[level] when the respective slice starts).
01247             #ifdef _OPENMP
01248               #pragma omp parallel for 
01249             #endif
01250               for (unsigned int i=0; i<=_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] = A[level].size1();
01255                 else Offset[i][level] = 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 _OPENMP
01273               #pragma omp parallel for private (x,y,point)
01274             #endif
01275               for (unsigned int i=0; i<_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 = 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 = 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 = col_iter.index2();
01304                   }
01305                   
01306                   ++row_iter;
01307                   x = 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           unsigned int x,y,z;
01326           ScalarType prod;
01327           RES = SparseMatrixType(A.size1(), B.size2());
01328           RES.clear();
01329           
01330     #ifdef _OPENMP
01331           #pragma omp parallel for private (x,y,z,prod) shared (A,B,RES)
01332     #endif
01333           for (x=0; x<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 = 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 = 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           unsigned int x,y1,y2,z;
01366           amg_sparsevector<ScalarType> row;
01367           RES = SparseMatrixType(P.size2(), P.size2());
01368           RES.clear();
01369           
01370     #ifdef _OPENMP
01371           #pragma omp parallel for private (x,y1,y2,z,row) shared (A,P,RES)
01372     #endif      
01373           for (x=0; x<P.size2(); ++x)
01374           {
01375             row = amg_sparsevector<ScalarType>(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 = 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 = 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 = col_iter3.index2();
01400                 RES.add (x, z, *col_iter3 * *iter);
01401               }
01402             }
01403           }
01404           
01405           #ifdef 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 (abs((ScalarType)RAP(x,y) - (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

Generated on Fri Dec 30 2011 23:20:43 for ViennaCL - The Vienna Computing Library by  doxygen 1.7.1