bcrsmatrix.hh

Go to the documentation of this file.
00001 #ifndef DUNE_BCRSMATRIX_HH
00002 #define DUNE_BCRSMATRIX_HH
00003 
00004 #include<math.h>
00005 #include<complex>
00006 #include<set>
00007 #include<iostream>
00008 #include<algorithm>
00009 #include<numeric>
00010 #include<vector>
00011 
00012 #include "istlexception.hh"
00013 #include "allocator.hh"
00014 #include "bvector.hh"
00015 #include <dune/common/stdstreams.hh>
00016 #include <dune/common/iteratorfacades.hh>
00017 #include <dune/common/typetraits.hh>
00018 #include <dune/common/helpertemplates.hh>
00019 
00024 namespace Dune {
00025   
00144   template<class B, class A=ISTLAllocator>
00145   class BCRSMatrix
00146   {
00147   private:
00148     enum BuildStage{
00150       notbuilt=0, 
00155       rowSizesBuilt=1, 
00157       built=2
00158     };
00159 
00160   public:
00161 
00162         //===== type definitions and constants
00163 
00165         typedef typename B::field_type field_type;
00166 
00168         typedef B block_type;
00169 
00171         typedef A allocator_type;
00172 
00174         typedef CompressedBlockVectorWindow<B,A> row_type;
00175         
00177     typedef typename A::size_type size_type;
00178     
00180         enum {
00182           blocklevel = B::blocklevel+1
00183         };
00184 
00186         enum BuildMode {
00197           row_wise, 
00206           random,
00210           unknown
00211         };
00212         
00213 
00214         //===== random access interface to rows of the matrix
00215 
00217         row_type& operator[] (size_type i)
00218         {
00219 #ifdef DUNE_ISTL_WITH_CHECKING
00220           if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00221           if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00222           if (r[i].getptr()==0) DUNE_THROW(ISTLError,"row not initialized yet");
00223 #endif
00224           return r[i];
00225         }
00226 
00228         const row_type& operator[] (size_type i) const
00229         {
00230 #ifdef DUNE_ISTL_WITH_CHECKING
00231           if (built!=ready) DUNE_THROW(ISTLError,"row not initialized yet");
00232           if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00233 #endif
00234           return r[i];
00235         }
00236 
00237 
00238         //===== iterator interface to rows of the matrix
00239 
00241     template<class T>
00242     class RealRowIterator
00243       : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
00244     {
00245 
00246     public:
00248       typedef typename Dune::RemoveConst<T>::Type ValueType;
00249 
00250       friend class RandomAccessIteratorFacade<RealRowIterator<const ValueType>, const ValueType>;
00251       friend class RandomAccessIteratorFacade<RealRowIterator<ValueType>, ValueType>;
00252       friend class RealRowIterator<const ValueType>;
00253       friend class RealRowIterator<ValueType>;
00254       
00256       RealRowIterator (row_type* _p, size_type _i)
00257         : p(_p), i(_i)
00258       {}
00259 
00261       RealRowIterator ()
00262         : p(0), i(0)
00263       {}
00264       
00265       RealRowIterator(const RealRowIterator<ValueType>& it)
00266         : p(it.p), i(it.i)
00267       {}
00268 
00269       
00271       size_type index () const
00272       {
00273         return i;
00274       }
00275       
00276       std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
00277       {
00278         assert(other.p==p);
00279         return (other.i-i);
00280       }
00281       
00282       std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
00283       {
00284         assert(other.p==p);
00285         return (other.i-i);
00286       }
00287       
00289       bool equals (const RealRowIterator<ValueType>& other) const
00290       {
00291         assert(other.p==p);
00292         return i==other.i;
00293       }
00294       
00296       bool equals (const RealRowIterator<const ValueType>& other) const
00297       {
00298         assert(other.p==p);
00299         return i==other.i;
00300       }
00301 
00302     private:
00304       void increment()
00305       {
00306         ++i;
00307       }
00308       
00310       void decrement()
00311       {
00312         --i;
00313       }
00314           
00315       void advance(std::ptrdiff_t diff)
00316       {
00317         i+=diff;
00318       }
00319       
00320       T& elementAt(std::ptrdiff_t diff) const
00321       {
00322         return p[i+diff];
00323       }
00324 
00326       row_type& dereference () const
00327       {
00328         return p[i];
00329       }
00330       
00331       row_type* p;
00332       size_type i;
00333     };
00334 
00336     typedef RealRowIterator<row_type> iterator;
00337     typedef RealRowIterator<row_type> Iterator;
00338     
00339 
00341         Iterator begin ()
00342         {
00343           return Iterator(r,0);
00344         }
00345           
00347         Iterator end ()
00348         {
00349           return Iterator(r,n);
00350         }
00351 
00353         Iterator rbegin ()
00354         {
00355           return Iterator(r,n-1);
00356         }
00357           
00359         Iterator rend ()
00360         {
00361           return Iterator(r,-1);
00362         }
00363 
00365         typedef Iterator RowIterator;
00366 
00368         typedef typename row_type::Iterator ColIterator;
00369 
00371     typedef RealRowIterator<const row_type> const_iterator;
00372     typedef RealRowIterator<const row_type> ConstIterator;
00373 
00374 
00376         ConstIterator begin () const
00377         {
00378           return ConstIterator(r,0);
00379         }
00380           
00382         ConstIterator end () const
00383         {
00384           return ConstIterator(r,n);
00385         }
00386 
00388         ConstIterator rbegin () const
00389         {
00390           return ConstIterator(r,n-1);
00391         }
00392           
00394         ConstIterator rend () const
00395         {
00396           return ConstIterator(r,-1);
00397         }
00398 
00400         typedef ConstIterator ConstRowIterator;
00401 
00403         typedef typename row_type::ConstIterator ConstColIterator;
00404 
00405         //===== constructors & resizers
00406 
00408         BCRSMatrix () 
00409           : build_mode(unknown), ready(notbuilt), n(0), m(0), nnz(0),
00410             r(0), a(0), j(0)
00411         {}
00412 
00414         BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm)
00415           : build_mode(bm), ready(notbuilt)
00416         {
00417           allocate(_n, _m, _nnz);
00418         }
00419 
00421         BCRSMatrix (size_type _n, size_type _m, BuildMode bm)
00422           : build_mode(bm), ready(notbuilt)
00423         {
00424           allocate(_n, _m);
00425         }
00426 
00428         BCRSMatrix (const BCRSMatrix& Mat)
00429           : n(Mat.n), nnz(0)
00430         {
00431           // deep copy in global array
00432           int _nnz = Mat.nnz;
00433 
00434           // in case of row-wise allocation
00435           if (_nnz<=0)
00436                 {
00437                   _nnz = 0;
00438                   for (size_type i=0; i<n; i++)
00439                         _nnz += Mat.r[i].getsize();
00440                 }
00441 
00442           allocate(Mat.n, Mat.m, _nnz);
00443 
00444           // build window structure
00445           copyWindowStructure(Mat);
00446         }
00447 
00449         ~BCRSMatrix ()
00450         {
00451           deallocate();
00452         }
00453 
00458     void setBuildMode(BuildMode bm)
00459     {
00460       if(ready==notbuilt)
00461         build_mode = bm;
00462       else
00463         DUNE_THROW(InvalidStateException, "Matrix structure is already built (ready="<<ready<<").");
00464     }
00465     
00481     void setSize(size_type rows, size_type columns, size_type nnz=0)
00482     {
00483       // deallocate already setup memory
00484       deallocate();
00485       
00486       // allocate matrix memory
00487       allocate(rows, columns, nnz);
00488     }
00489     
00491     BCRSMatrix& operator= (const BCRSMatrix& Mat)
00492     {
00493       // return immediately when self-assignment
00494       if (&Mat==this) return *this;
00495       
00496       // make it simple: ALWAYS throw away memory for a and j
00497       deallocate(false);
00498       
00499       // reallocate the rows if required
00500       if (n>0 && n!=Mat.n)
00501         // free rows
00502         A::template free<row_type>(r); 
00503       
00504       nnz=Mat.nnz;
00505       if (nnz<=0)
00506         {
00507           for (size_type i=0; i<Mat.n; i++)
00508             nnz += Mat.r[i].getsize();
00509         }
00510 
00511       // allocate a,j
00512       allocate(Mat.n, Mat.m, nnz, n!=Mat.n);
00513       
00514       // build window structure
00515       copyWindowStructure(Mat);
00516       return *this;
00517     }
00518 
00520         BCRSMatrix& operator= (const field_type& k)
00521         {
00522             for (size_type i=0; i<n; i++) r[i] = k;
00523             return *this;
00524         }
00525 
00526         //===== row-wise creation interface 
00527 
00529         class CreateIterator 
00530         {
00531         public:
00533           CreateIterator (BCRSMatrix& _Mat, size_type _i) 
00534             : Mat(_Mat), i(_i), nnz(0), current_row(Mat.a, Mat.j, 0)
00535           {
00536                 if (i==0 && Mat.ready)
00537                   DUNE_THROW(ISTLError,"creation only allowed for uninitialized matrix");
00538                 if(Mat.build_mode!=row_wise)
00539                   if(Mat.build_mode==unknown)
00540                     Mat.build_mode=row_wise;
00541                   else
00542                     DUNE_THROW(ISTLError,"creation only allowed if row wise allocation was requested in the constructor");
00543           }
00544 
00546           CreateIterator& operator++()
00547           {
00548                 // this should only be called if matrix is in creation
00549                 if (Mat.ready)
00550                   DUNE_THROW(ISTLError,"matrix already built up");
00551 
00552                 // row i is defined through the pattern
00553                 // get memory for the row and initialize the j array
00554                 // this depends on the allocation mode
00555 
00556                 // compute size of the row
00557                 size_type s = pattern.size();
00558 
00559                 if(s>0){                  
00560                   // update number of nonzeroes including this row
00561                   nnz += s;
00562                   
00563                   // alloc memory / set window
00564                   if (Mat.nnz>0)
00565                     {
00566                       // memory is allocated in one long array
00567 
00568                       // check if that memory is sufficient
00569                       if (nnz>Mat.nnz) 
00570                         DUNE_THROW(ISTLError,"allocated nnz too small");
00571                       
00572                       // set row i
00573                       Mat.r[i].set(s,current_row.getptr(),current_row.getindexptr());
00574                       current_row.setptr(current_row.getptr()+s);
00575                       current_row.setindexptr(current_row.getindexptr()+s);
00576                     }else{
00577                       // memory is allocated individually per row
00578                       // allocate and set row i
00579                       B*   a = A::template malloc<B>(s);
00580                       size_type* j = A::template malloc<size_type>(s);
00581                       Mat.r[i].set(s,a,j);
00582                     }
00583                 }else
00584                   // setup empty row
00585                   Mat.r[i].set(0,0,0);
00586 
00587                 // initialize the j array for row i from pattern
00588                 size_type k=0;
00589                 size_type *j =  Mat.r[i].getindexptr();
00590                 for (typename std::set<size_type>::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
00591                   j[k++] = *it;
00592 
00593                 // now go to next row
00594         i++;
00595                 pattern.clear();
00596 
00597                 // check if this was last row
00598                 if (i==Mat.n)
00599                   {
00600                         Mat.ready = built;
00601                         if(Mat.nnz>0)
00602                           // Set nnz to the exact number of nonzero blocks inserted
00603                           // as some methods rely on it
00604                           Mat.nnz=nnz;
00605                   }
00606                 // done
00607                 return *this;
00608           }
00609           
00611           bool operator!= (const CreateIterator& it) const
00612           {
00613         return (i!=it.i) || (&Mat!=&it.Mat);
00614           }
00615 
00617           bool operator== (const CreateIterator& it) const
00618           {
00619         return (i==it.i) && (&Mat==&it.Mat);
00620           }
00621 
00623           size_type index () const
00624           {
00625         return i;
00626           }
00627 
00629           void insert (size_type j)
00630           {
00631                 pattern.insert(j);
00632           }
00633 
00635           bool contains (size_type j)
00636           {
00637                 if (pattern.find(j)!=pattern.end())
00638                   return true;
00639                 else
00640                   return false;
00641           }
00647           typename std::set<size_type>::size_type size() const
00648           {
00649             return pattern.size();
00650           }
00651           
00652         private:
00653           BCRSMatrix& Mat; // the matrix we are defining
00654           size_type i;           // current row to be defined 
00655           size_type nnz;         // count total number of nonzeros
00656           std::set<size_type> pattern; // used to compile entries in a row
00657           row_type current_row; // row poiting to the current row to setup
00658         };
00659 
00661         friend class CreateIterator;
00662 
00664         CreateIterator createbegin ()
00665         {
00666           return CreateIterator(*this,0);
00667         }
00668 
00670         CreateIterator createend ()
00671         {
00672           return CreateIterator(*this,n);
00673         }
00674 
00675 
00676         //===== random creation interface 
00677 
00679         void setrowsize (size_type i, size_type s)
00680         {
00681           if (build_mode!=random)
00682                 DUNE_THROW(ISTLError,"requires random build mode");       
00683           if (ready)
00684                 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00685 
00686           r[i].setsize(s);
00687         }
00688 
00690         void incrementrowsize (size_type i)
00691         {
00692           if (build_mode!=random)
00693                 DUNE_THROW(ISTLError,"requires random build mode");       
00694           if (ready)
00695                 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00696 
00697           r[i].setsize(r[i].getsize()+1);
00698         }
00699 
00701         void endrowsizes ()
00702         {
00703           if (build_mode!=random)
00704                 DUNE_THROW(ISTLError,"requires random build mode");       
00705           if (ready)
00706                 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00707 
00708           // compute total size, check positivity
00709           size_type total=0;
00710           for (size_type i=0; i<n; i++)
00711                 {
00712                     if (r[i].getsize()<0)
00713                         DUNE_THROW(ISTLError,"rowsize must be nonnegative");      
00714                   total += r[i].getsize();
00715                 }
00716 
00717           // allocate/check memory
00718           allocate(n,m,total,false);
00719           
00720           // set the window pointers correctly
00721           setWindowPointers(begin());
00722           
00723           // initialize j array with m (an invalid column index)
00724           // this indicates an unused entry
00725           for (size_type k=0; k<nnz; k++)
00726                 j[k] = m;
00727           ready = rowSizesBuilt;
00728         }
00729 
00731 
00741         void addindex (size_type row, size_type col)
00742         {
00743           if (build_mode!=random)
00744                 DUNE_THROW(ISTLError,"requires random build mode");       
00745           if (ready==built)
00746                 DUNE_THROW(ISTLError,"matrix already built up");          
00747           if (ready==notbuilt)
00748                 DUNE_THROW(ISTLError,"matrix row sizes not built up yet");
00749 
00750           if (col >= m)
00751             DUNE_THROW(ISTLError,"column index exceeds matrix size");
00752 
00753           // get row range
00754           size_type* const first = r[row].getindexptr();
00755           size_type* const last = first + r[row].getsize();
00756 
00757           // find correct insertion position for new column index
00758           size_type* pos = std::lower_bound(first,last,col);
00759 
00760           // check if index is already in row
00761           if (pos!=last && *pos == col) return;
00762 
00763           // find end of already inserted column indices
00764           size_type* end = std::lower_bound(pos,last,m);
00765           if (end==last)
00766             DUNE_THROW(ISTLError,"row is too small");
00767 
00768           // insert new column index at correct position
00769           std::copy_backward(pos,end,end+1);
00770           *pos = col;
00771 
00772         }
00773 
00775         void endindices ()
00776         {
00777           if (build_mode!=random)
00778                 DUNE_THROW(ISTLError,"requires random build mode");       
00779           if (ready==built)
00780                 DUNE_THROW(ISTLError,"matrix already built up");
00781           if (ready==notbuilt)
00782             DUNE_THROW(ISTLError,"row sizes are not built up yet");
00783 
00784           // check if there are undefined indices
00785           RowIterator endi=end();
00786           for (RowIterator i=begin(); i!=endi; ++i)
00787                 {
00788                   ColIterator endj = (*i).end();
00789                   for (ColIterator j=(*i).begin(); j!=endj; ++j){
00790                         if (j.index()<0)
00791                           {
00792                                 std::cout << "j[" << j.offset() << "]=" << j.index() << std::endl;
00793                                 DUNE_THROW(ISTLError,"undefined index detected");
00794                           }
00795                         if (j.index()>=m){
00796                           dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
00797                           <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
00798                           r[i.index()].setsize(j.offset());
00799                           break;
00800                         }
00801                   }
00802           }
00803 
00804           // if not, set matrix to built
00805           ready = built;
00806         }
00807 
00808         //===== vector space arithmetic
00809 
00811         BCRSMatrix& operator*= (const field_type& k)
00812         {
00813           if (nnz>0)
00814                 {
00815                   // process 1D array
00816                   for (size_type i=0; i<nnz; i++)
00817                         a[i] *= k;
00818                 }
00819           else
00820                 {
00821                   RowIterator endi=end();
00822                   for (RowIterator i=begin(); i!=endi; ++i)
00823                         {
00824                           ColIterator endj = (*i).end();
00825                           for (ColIterator j=(*i).begin(); j!=endj; ++j)
00826                                 (*j) *= k;
00827                         }
00828                 }
00829 
00830           return *this;
00831         }
00832 
00834         BCRSMatrix& operator/= (const field_type& k)
00835         {
00836           if (nnz>0)
00837                 {
00838                   // process 1D array
00839                   for (size_type i=0; i<nnz; i++)
00840                         a[i] /= k;
00841                 }
00842           else
00843                 {
00844                   RowIterator endi=end();
00845                   for (RowIterator i=begin(); i!=endi; ++i)
00846                         {
00847                           ColIterator endj = (*i).end();
00848                           for (ColIterator j=(*i).begin(); j!=endj; ++j)
00849                                 (*j) /= k;
00850                         }
00851                 }
00852 
00853           return *this;
00854         }
00855 
00856 
00862         BCRSMatrix& operator+= (const BCRSMatrix& b)
00863         {
00864 #ifdef DUNE_ISTL_WITH_CHECKING
00865           if(N()!=b.N() || M() != b.M())
00866             DUNE_THROW(RangeError, "Matrix sizes do not match!");
00867 #endif
00868           RowIterator endi=end();
00869           ConstRowIterator j=b.begin();
00870           for (RowIterator i=begin(); i!=endi; ++i, ++j){
00871             i->operator+=(*j);
00872           }
00873         
00874           return *this;
00875         }
00876 
00882         BCRSMatrix& operator-= (const BCRSMatrix& b)
00883         {
00884 #ifdef DUNE_ISTL_WITH_CHECKING
00885           if(N()!=b.N() || M() != b.M())
00886             DUNE_THROW(RangeError, "Matrix sizes do not match!");
00887 #endif
00888           RowIterator endi=end();
00889           ConstRowIterator j=b.begin();
00890           for (RowIterator i=begin(); i!=endi; ++i, ++j){
00891             i->operator-=(*j);
00892           }
00893         
00894           return *this;
00895         }
00896         //===== linear maps
00897    
00899         template<class X, class Y>
00900         void mv (const X& x, Y& y) const
00901         {
00902 #ifdef DUNE_ISTL_WITH_CHECKING
00903           if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00904           if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00905 #endif
00906           ConstRowIterator endi=end();
00907           for (ConstRowIterator i=begin(); i!=endi; ++i)
00908                 {
00909                   y[i.index()]=0;
00910                   ConstColIterator endj = (*i).end();
00911                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00912                         (*j).umv(x[j.index()],y[i.index()]);
00913                 }
00914         }
00915 
00917         template<class X, class Y>
00918         void umv (const X& x, Y& y) const
00919         {
00920 #ifdef DUNE_ISTL_WITH_CHECKING
00921           if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00922           if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00923 #endif
00924           ConstRowIterator endi=end();
00925           for (ConstRowIterator i=begin(); i!=endi; ++i)
00926                 {
00927                   ConstColIterator endj = (*i).end();
00928                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00929                         (*j).umv(x[j.index()],y[i.index()]);
00930                 }
00931         }
00932 
00934         template<class X, class Y>
00935         void mmv (const X& x, Y& y) const
00936         {
00937 #ifdef DUNE_ISTL_WITH_CHECKING
00938           if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00939           if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00940 #endif
00941           ConstRowIterator endi=end();
00942           for (ConstRowIterator i=begin(); i!=endi; ++i)
00943                 {
00944                   ConstColIterator endj = (*i).end();
00945                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00946                         (*j).mmv(x[j.index()],y[i.index()]);
00947                 }
00948         }
00949 
00951         template<class X, class Y>
00952         void usmv (const field_type& alpha, const X& x, Y& y) const
00953         {
00954 #ifdef DUNE_ISTL_WITH_CHECKING
00955           if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00956           if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00957 #endif
00958           ConstRowIterator endi=end();
00959           for (ConstRowIterator i=begin(); i!=endi; ++i)
00960                 {
00961                   ConstColIterator endj = (*i).end();
00962                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00963                         (*j).usmv(alpha,x[j.index()],y[i.index()]);
00964                 }
00965         }
00966 
00968         template<class X, class Y>
00969         void umtv (const X& x, Y& y) const
00970         {
00971 #ifdef DUNE_ISTL_WITH_CHECKING
00972           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00973           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00974 #endif
00975           ConstRowIterator endi=end();
00976           for (ConstRowIterator i=begin(); i!=endi; ++i)
00977                 {
00978                   ConstColIterator endj = (*i).end();
00979                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00980                         (*j).umtv(x[i.index()],y[j.index()]);
00981                 }
00982         }
00983 
00985         template<class X, class Y>
00986         void mmtv (const X& x, Y& y) const
00987         {
00988 #ifdef DUNE_ISTL_WITH_CHECKING
00989           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00990           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00991 #endif
00992           ConstRowIterator endi=end();
00993           for (ConstRowIterator i=begin(); i!=endi; ++i)
00994                 {
00995                   ConstColIterator endj = (*i).end();
00996                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00997                         (*j).mmtv(x[i.index()],y[j.index()]);
00998                 }
00999         }
01000 
01002         template<class X, class Y>
01003         void usmtv (const field_type& alpha, const X& x, Y& y) const
01004         {
01005 #ifdef DUNE_ISTL_WITH_CHECKING
01006           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01007           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01008 #endif
01009           ConstRowIterator endi=end();
01010           for (ConstRowIterator i=begin(); i!=endi; ++i)
01011                 {
01012                   ConstColIterator endj = (*i).end();
01013                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01014                         (*j).usmtv(alpha,x[i.index()],y[j.index()]);
01015                 }
01016         }
01017           
01019         template<class X, class Y>
01020         void umhv (const X& x, Y& y) const
01021         {
01022 #ifdef DUNE_ISTL_WITH_CHECKING
01023           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01024           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01025 #endif
01026           ConstRowIterator endi=end();
01027           for (ConstRowIterator i=begin(); i!=endi; ++i)
01028                 {
01029                   ConstColIterator endj = (*i).end();
01030                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01031                         (*j).umhv(x[i.index()],y[j.index()]);
01032                 }
01033         }
01034 
01036         template<class X, class Y>
01037         void mmhv (const X& x, Y& y) const
01038         {
01039 #ifdef DUNE_ISTL_WITH_CHECKING
01040           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01041           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01042 #endif
01043           ConstRowIterator endi=end();
01044           for (ConstRowIterator i=begin(); i!=endi; ++i)
01045                 {
01046                   ConstColIterator endj = (*i).end();
01047                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01048                         (*j).mmhv(x[i.index()],y[j.index()]);
01049                 }
01050         }
01051 
01053         template<class X, class Y>
01054         void usmhv (const field_type& alpha, const X& x, Y& y) const
01055         {
01056 #ifdef DUNE_ISTL_WITH_CHECKING
01057           if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01058           if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01059 #endif
01060           ConstRowIterator endi=end();
01061           for (ConstRowIterator i=begin(); i!=endi; ++i)
01062                 {
01063                   ConstColIterator endj = (*i).end();
01064                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01065                         (*j).usmhv(alpha,x[i.index()],y[j.index()]);
01066                 }
01067         }
01068           
01069 
01070         //===== norms
01071 
01073     double frobenius_norm2 () const
01074         {
01075           double sum=0;
01076 
01077           ConstRowIterator endi=end();
01078           for (ConstRowIterator i=begin(); i!=endi; ++i)
01079                 {
01080                   ConstColIterator endj = (*i).end();
01081                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01082                         sum += (*j).frobenius_norm2();
01083                 }
01084 
01085           return sum;
01086         }
01087 
01089     double frobenius_norm () const
01090         {
01091           return sqrt(frobenius_norm2());
01092         }
01093 
01095     double infinity_norm () const
01096         {
01097           double max=0;
01098           ConstRowIterator endi=end();
01099           for (ConstRowIterator i=begin(); i!=endi; ++i)
01100                 {
01101                   double sum=0;
01102                   ConstColIterator endj = (*i).end();
01103                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01104                         sum += (*j).infinity_norm();
01105                   max = std::max(max,sum);
01106                 }
01107           return max;
01108         }
01109 
01111         double infinity_norm_real () const
01112         {
01113           double max=0;
01114           ConstRowIterator endi=end();
01115           for (ConstRowIterator i=begin(); i!=endi; ++i)
01116                 {
01117                   double sum=0;
01118                   ConstColIterator endj = (*i).end();
01119                   for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01120                         sum += (*j).infinity_norm_real();
01121                   max = std::max(max,sum);
01122                 }
01123           return max;
01124         }
01125 
01126 
01127         //===== sizes
01128 
01130         size_type N () const
01131         {
01132           return n;
01133         }
01134 
01136         size_type M () const
01137         {
01138           return m;
01139         }
01140 
01144         size_type rowdim (size_type i) const
01145         {
01146           const B* row = r[i].getptr();
01147           if(row)
01148             return row->rowdim();
01149           else
01150             return 0;
01151         }
01152 
01156         size_type coldim (size_type c) const
01157         {
01158           // find an entry in column j
01159           if (nnz>0)
01160                 {
01161                   for (size_type k=0; k<nnz; k++) {
01162                         if (j[k]==c) {
01163                           return a[k].coldim();
01164                         }
01165                   }
01166                 }
01167           else
01168                 {
01169                   for (size_type i=0; i<n; i++)
01170                         {
01171                           size_type* j = r[i].getindexptr();
01172                           B*   a = r[i].getptr();
01173                           for (size_type k=0; k<r[i].getsize(); k++)
01174                                 if (j[k]==c) {
01175                                   return a[k].coldim();
01176                                 }
01177                         }
01178                 }
01179 
01180           // not found
01181           return 0;
01182         }
01183 
01185         size_type rowdim () const
01186         {
01187           size_type nn=0;
01188           for (size_type i=0; i<n; i++)
01189                 nn += rowdim(i);
01190           return nn;
01191         }
01192 
01194         size_type coldim () const
01195         {
01196           // The following code has a complexity of nnz, and
01197           // typically a very small constant.
01198           //
01199           std::vector<size_type> coldims(M(),-1);
01200           for (ConstRowIterator row=begin(); row!=end(); ++row) 
01201             for (ConstColIterator col=row->begin(); col!=row->end(); ++col)
01202               // only compute blocksizes we don't already have
01203               if (coldims[col.index()]==-1) 
01204                 coldims[col.index()] = col->coldim();
01205 
01206           size_type sum = 0;
01207           for (typename std::vector<size_type>::iterator it=coldims.begin(); it!=coldims.end(); ++it)
01208             // skip rows for which no coldim could be determined
01209             if ((*it)>=0)
01210               sum += *it;
01211 
01212           return sum;
01213         }
01214 
01215         //===== query
01216         
01218         bool exists (size_type i, size_type j) const
01219         {
01220 #ifdef DUNE_ISTL_WITH_CHECKING
01221           if (i<0 || i>=n) DUNE_THROW(ISTLError,"index out of range");
01222           if (j<0 || i>=m) DUNE_THROW(ISTLError,"index out of range");
01223 #endif
01224           if (r[i].size() && r[i].find(j)!=r[i].end()) 
01225                 return true;
01226           else
01227                 return false;
01228         }
01229 
01230         
01231   private:
01232         // state information
01233         BuildMode build_mode; // row wise or whole matrix
01234         BuildStage ready;           // indicate the stage the matrix building is in
01235 
01236         // size of the matrix
01237         size_type  n;  // number of rows
01238         size_type  m;  // number of columns
01239         size_type nnz; // number of nonzeros allocated in the a and j array below
01240              // zero means that memory is allocated separately for each row.
01241 
01242         // the rows are dynamically allocated
01243         row_type* r; // [n] the individual rows having pointers into a,j arrays
01244 
01245         // dynamically allocated memory
01246         B*   a;  // [nnz] non-zero entries of the matrix in row-wise ordering
01247         size_type* j;  // [nnz] column indices of entries
01248 
01249     void setWindowPointers(ConstRowIterator row)
01250     {
01251       row_type current_row(a,j,0); // Pointers to current row data
01252       for (size_type i=0; i<n; i++, ++row){               
01253         // set row i
01254         size_type s = row->getsize();
01255         
01256         if (s>0){
01257           // setup pointers and size
01258           r[i].set(s,current_row.getptr(), current_row.getindexptr());
01259           // update pointer for next row
01260           current_row.setptr(current_row.getptr()+s);
01261           current_row.setindexptr(current_row.getindexptr()+s);
01262         } else{
01263           // empty row
01264           r[i].set(0,0,0);
01265         }
01266       }
01267     }
01268     
01270     void copyWindowStructure(const BCRSMatrix& Mat)
01271     {
01272       setWindowPointers(Mat.begin());
01273             
01274       // copy data
01275       for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
01276 
01277       // finish off
01278       build_mode = row_wise; // dummy
01279       ready = built;
01280     }
01281   
01287     void deallocate(bool deallocateRows=true)
01288     {
01289       
01290       if (nnz>0)
01291         {
01292           // a,j have been allocated as one long vector
01293           A::template free<size_type>(j); 
01294           A::template free<B>(a); 
01295         }
01296       else
01297         {
01298           // check if memory for rows have been allocated individually
01299           for (size_type i=0; i<n; i++)
01300             if (r[i].getsize()>0) 
01301               {
01302                 A::template free<size_type>(r[i].getindexptr());
01303                 A::template free<B>(r[i].getptr());
01304               }
01305         }
01306       
01307       // deallocate the rows
01308       if (n>0 && deallocateRows) A::template free<row_type>(r); 
01309       
01310       // Mark matrix as not built at all.
01311       ready=notbuilt;
01312       
01313     }
01314     
01332     void allocate(size_type rows, size_type columns, size_type nnz_=0, bool allocateRows=true)
01333     {
01334       // Store size
01335       n = rows;
01336       m = columns;
01337       nnz = nnz_;
01338 
01339       // allocate rows
01340       if(allocateRows){
01341         if (n>0){ 
01342           r = A::template malloc<row_type>(rows);
01343         }else{
01344           r = 0;
01345         }
01346       }
01347       
01348 
01349       // allocate a and j array
01350       if (nnz>0){ 
01351         a = A::template malloc<B>(nnz);
01352         j = A::template malloc<size_type>(nnz);
01353       }else{
01354         a = 0;
01355         j = 0;
01356       }
01357       // Mark the matrix as not built.
01358       ready = notbuilt;
01359     }
01360     
01361   };
01362 
01363 
01366 } // end namespace
01367 
01368 #endif

Generated on 12 Dec 2007 with Doxygen (ver 1.5.1)