00001 #ifndef DUNE_BCRSMATRIX_HH
00002 #define DUNE_BCRSMATRIX_HH
00003
00004 #include<cmath>
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/shared_ptr.hh>
00016 #include <dune/common/stdstreams.hh>
00017 #include <dune/common/iteratorfacades.hh>
00018 #include <dune/common/typetraits.hh>
00019 #include <dune/common/static_assert.hh>
00020 #include <dune/common/poolallocator.hh>
00021
00026 namespace Dune {
00027
00067 template<typename M>
00068 class MatrixDimension;
00069
00177 template<class B, class A=ISTLAllocator>
00178 class BCRSMatrix
00179 {
00180 friend class MatrixDimension<BCRSMatrix>;
00181
00182 private:
00183 enum BuildStage{
00185 notbuilt=0,
00190 rowSizesBuilt=1,
00192 built=2
00193 };
00194
00195 public:
00196
00197
00198
00200 typedef typename B::field_type field_type;
00201
00203 typedef B block_type;
00204
00206 typedef A allocator_type;
00207
00209 typedef CompressedBlockVectorWindow<B,A> row_type;
00210
00212 typedef typename A::size_type size_type;
00213
00215 enum {
00217 blocklevel = B::blocklevel+1
00218 };
00219
00221 enum BuildMode {
00232 row_wise,
00241 random,
00245 unknown
00246 };
00247
00248
00249
00250
00252 row_type& operator[] (size_type i)
00253 {
00254 #ifdef DUNE_ISTL_WITH_CHECKING
00255 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00256 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00257 if (r[i].getptr()==0) DUNE_THROW(ISTLError,"row not initialized yet");
00258 #endif
00259 return r[i];
00260 }
00261
00263 const row_type& operator[] (size_type i) const
00264 {
00265 #ifdef DUNE_ISTL_WITH_CHECKING
00266 if (built!=ready) DUNE_THROW(ISTLError,"row not initialized yet");
00267 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00268 #endif
00269 return r[i];
00270 }
00271
00272
00273
00274
00276 template<class T>
00277 class RealRowIterator
00278 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
00279 {
00280
00281 public:
00283 typedef typename remove_const<T>::type ValueType;
00284
00285 friend class RandomAccessIteratorFacade<RealRowIterator<const ValueType>, const ValueType>;
00286 friend class RandomAccessIteratorFacade<RealRowIterator<ValueType>, ValueType>;
00287 friend class RealRowIterator<const ValueType>;
00288 friend class RealRowIterator<ValueType>;
00289
00291 RealRowIterator (row_type* _p, size_type _i)
00292 : p(_p), i(_i)
00293 {}
00294
00296 RealRowIterator ()
00297 : p(0), i(0)
00298 {}
00299
00300 RealRowIterator(const RealRowIterator<ValueType>& it)
00301 : p(it.p), i(it.i)
00302 {}
00303
00304
00306 size_type index () const
00307 {
00308 return i;
00309 }
00310
00311 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
00312 {
00313 assert(other.p==p);
00314 return (other.i-i);
00315 }
00316
00317 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
00318 {
00319 assert(other.p==p);
00320 return (other.i-i);
00321 }
00322
00324 bool equals (const RealRowIterator<ValueType>& other) const
00325 {
00326 assert(other.p==p);
00327 return i==other.i;
00328 }
00329
00331 bool equals (const RealRowIterator<const ValueType>& other) const
00332 {
00333 assert(other.p==p);
00334 return i==other.i;
00335 }
00336
00337 private:
00339 void increment()
00340 {
00341 ++i;
00342 }
00343
00345 void decrement()
00346 {
00347 --i;
00348 }
00349
00350 void advance(std::ptrdiff_t diff)
00351 {
00352 i+=diff;
00353 }
00354
00355 T& elementAt(std::ptrdiff_t diff) const
00356 {
00357 return p[i+diff];
00358 }
00359
00361 row_type& dereference () const
00362 {
00363 return p[i];
00364 }
00365
00366 row_type* p;
00367 size_type i;
00368 };
00369
00371 typedef RealRowIterator<row_type> iterator;
00372 typedef RealRowIterator<row_type> Iterator;
00373
00374
00376 Iterator begin ()
00377 {
00378 return Iterator(r,0);
00379 }
00380
00382 Iterator end ()
00383 {
00384 return Iterator(r,n);
00385 }
00386
00388 Iterator rbegin ()
00389 {
00390 return Iterator(r,n-1);
00391 }
00392
00394 Iterator rend ()
00395 {
00396 return Iterator(r,-1);
00397 }
00398
00400 typedef Iterator RowIterator;
00401
00403 typedef typename row_type::Iterator ColIterator;
00404
00406 typedef RealRowIterator<const row_type> const_iterator;
00407 typedef RealRowIterator<const row_type> ConstIterator;
00408
00409
00411 ConstIterator begin () const
00412 {
00413 return ConstIterator(r,0);
00414 }
00415
00417 ConstIterator end () const
00418 {
00419 return ConstIterator(r,n);
00420 }
00421
00423 ConstIterator rbegin () const
00424 {
00425 return ConstIterator(r,n-1);
00426 }
00427
00429 ConstIterator rend () const
00430 {
00431 return ConstIterator(r,-1);
00432 }
00433
00435 typedef ConstIterator ConstRowIterator;
00436
00438 typedef typename row_type::ConstIterator ConstColIterator;
00439
00440
00441
00443 BCRSMatrix ()
00444 : build_mode(unknown), ready(notbuilt), n(0), m(0), nnz(0),
00445 r(0), a(0)
00446 {}
00447
00449 BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm)
00450 : build_mode(bm), ready(notbuilt)
00451 {
00452 allocate(_n, _m, _nnz);
00453 }
00454
00456 BCRSMatrix (size_type _n, size_type _m, BuildMode bm)
00457 : build_mode(bm), ready(notbuilt)
00458 {
00459 allocate(_n, _m);
00460 }
00461
00467 BCRSMatrix (const BCRSMatrix& Mat)
00468 : n(Mat.n), nnz(0)
00469 {
00470
00471 size_type _nnz = Mat.nnz;
00472
00473
00474 if (_nnz<=0)
00475 {
00476 _nnz = 0;
00477 for (size_type i=0; i<n; i++)
00478 _nnz += Mat.r[i].getsize();
00479 }
00480
00481 j = Mat.j;
00482 allocate(Mat.n, Mat.m, _nnz);
00483
00484
00485 copyWindowStructure(Mat);
00486 }
00487
00489 ~BCRSMatrix ()
00490 {
00491 deallocate();
00492 }
00493
00498 void setBuildMode(BuildMode bm)
00499 {
00500 if(ready==notbuilt)
00501 build_mode = bm;
00502 else
00503 DUNE_THROW(InvalidStateException, "Matrix structure is already built (ready="<<ready<<").");
00504 }
00505
00521 void setSize(size_type rows, size_type columns, size_type nnz=0)
00522 {
00523
00524 deallocate();
00525
00526
00527 allocate(rows, columns, nnz);
00528 }
00529
00536 BCRSMatrix& operator= (const BCRSMatrix& Mat)
00537 {
00538
00539 if (&Mat==this) return *this;
00540
00541
00542 deallocate(false);
00543
00544
00545 if (n>0 && n!=Mat.n)
00546
00547 A::template free<row_type>(r);
00548
00549 nnz=Mat.nnz;
00550 if (nnz<=0)
00551 {
00552 for (size_type i=0; i<Mat.n; i++)
00553 nnz += Mat.r[i].getsize();
00554 }
00555
00556
00557 j = Mat.j;
00558 allocate(Mat.n, Mat.m, nnz, n!=Mat.n);
00559
00560
00561 copyWindowStructure(Mat);
00562 return *this;
00563 }
00564
00566 BCRSMatrix& operator= (const field_type& k)
00567 {
00568 for (size_type i=0; i<n; i++) r[i] = k;
00569 return *this;
00570 }
00571
00572
00573
00575 class CreateIterator
00576 {
00577 public:
00579 CreateIterator (BCRSMatrix& _Mat, size_type _i)
00580 : Mat(_Mat), i(_i), nnz(0), current_row(Mat.a, Mat.j.get(), 0)
00581 {
00582 if (i==0 && Mat.ready)
00583 DUNE_THROW(ISTLError,"creation only allowed for uninitialized matrix");
00584 if(Mat.build_mode!=row_wise)
00585 {
00586 if(Mat.build_mode==unknown)
00587 Mat.build_mode=row_wise;
00588 else
00589 DUNE_THROW(ISTLError,"creation only allowed if row wise allocation was requested in the constructor");
00590 }
00591 }
00592
00594 CreateIterator& operator++()
00595 {
00596
00597 if (Mat.ready)
00598 DUNE_THROW(ISTLError,"matrix already built up");
00599
00600
00601
00602
00603
00604
00605 size_type s = pattern.size();
00606
00607 if(s>0){
00608
00609 nnz += s;
00610
00611
00612 if (Mat.nnz>0)
00613 {
00614
00615
00616
00617 if (nnz>Mat.nnz)
00618 DUNE_THROW(ISTLError,"allocated nnz too small");
00619
00620
00621 Mat.r[i].set(s,current_row.getptr(),current_row.getindexptr());
00622 current_row.setptr(current_row.getptr()+s);
00623 current_row.setindexptr(current_row.getindexptr()+s);
00624 }else{
00625
00626
00627 B* a = A::template malloc<B>(s);
00628 size_type* j = A::template malloc<size_type>(s);
00629 Mat.r[i].set(s,a,j);
00630 }
00631 }else
00632
00633 Mat.r[i].set(0,0,0);
00634
00635
00636 size_type k=0;
00637 size_type *j = Mat.r[i].getindexptr();
00638 for (typename PatternType::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
00639 j[k++] = *it;
00640
00641
00642 i++;
00643 pattern.clear();
00644
00645
00646 if (i==Mat.n)
00647 {
00648 Mat.ready = built;
00649 if(Mat.nnz>0)
00650
00651
00652 Mat.nnz=nnz;
00653 }
00654
00655 return *this;
00656 }
00657
00659 bool operator!= (const CreateIterator& it) const
00660 {
00661 return (i!=it.i) || (&Mat!=&it.Mat);
00662 }
00663
00665 bool operator== (const CreateIterator& it) const
00666 {
00667 return (i==it.i) && (&Mat==&it.Mat);
00668 }
00669
00671 size_type index () const
00672 {
00673 return i;
00674 }
00675
00677 void insert (size_type j)
00678 {
00679 pattern.insert(j);
00680 }
00681
00683 bool contains (size_type j)
00684 {
00685 if (pattern.find(j)!=pattern.end())
00686 return true;
00687 else
00688 return false;
00689 }
00695 size_type size() const
00696 {
00697 return pattern.size();
00698 }
00699
00700 private:
00701 BCRSMatrix& Mat;
00702 size_type i;
00703 size_type nnz;
00704 typedef std::set<size_type,std::less<size_type> > PatternType;
00705 PatternType pattern;
00706 row_type current_row;
00707 };
00708
00710 friend class CreateIterator;
00711
00713 CreateIterator createbegin ()
00714 {
00715 return CreateIterator(*this,0);
00716 }
00717
00719 CreateIterator createend ()
00720 {
00721 return CreateIterator(*this,n);
00722 }
00723
00724
00725
00726
00728 void setrowsize (size_type i, size_type s)
00729 {
00730 if (build_mode!=random)
00731 DUNE_THROW(ISTLError,"requires random build mode");
00732 if (ready)
00733 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00734
00735 r[i].setsize(s);
00736 }
00737
00739 size_type getrowsize (size_type i) const
00740 {
00741 #ifdef DUNE_ISTL_WITH_CHECKING
00742 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00743 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00744 #endif
00745 return r[i].getsize();
00746 }
00747
00749 void incrementrowsize (size_type i, size_type s = 1)
00750 {
00751 if (build_mode!=random)
00752 DUNE_THROW(ISTLError,"requires random build mode");
00753 if (ready)
00754 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00755
00756 r[i].setsize(r[i].getsize()+s);
00757 }
00758
00760 void endrowsizes ()
00761 {
00762 if (build_mode!=random)
00763 DUNE_THROW(ISTLError,"requires random build mode");
00764 if (ready)
00765 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00766
00767
00768 size_type total=0;
00769 for (size_type i=0; i<n; i++)
00770 {
00771 if (r[i].getsize()<0)
00772 DUNE_THROW(ISTLError,"rowsize must be nonnegative");
00773 total += r[i].getsize();
00774 }
00775
00776 if(nnz==0)
00777
00778 allocate(n,m,total,false);
00779 else if(nnz<total)
00780 DUNE_THROW(ISTLError,"Specified number of nonzeros ("<<nnz<<") not "
00781 <<"sufficient for calculated nonzeros ("<<total<<"! ");
00782
00783
00784 setWindowPointers(begin());
00785
00786
00787
00788 for (size_type k=0; k<nnz; k++)
00789 j.get()[k] = m;
00790 ready = rowSizesBuilt;
00791 }
00792
00794
00804 void addindex (size_type row, size_type col)
00805 {
00806 if (build_mode!=random)
00807 DUNE_THROW(ISTLError,"requires random build mode");
00808 if (ready==built)
00809 DUNE_THROW(ISTLError,"matrix already built up");
00810 if (ready==notbuilt)
00811 DUNE_THROW(ISTLError,"matrix row sizes not built up yet");
00812
00813 if (col >= m)
00814 DUNE_THROW(ISTLError,"column index exceeds matrix size");
00815
00816
00817 size_type* const first = r[row].getindexptr();
00818 size_type* const last = first + r[row].getsize();
00819
00820
00821 size_type* pos = std::lower_bound(first,last,col);
00822
00823
00824 if (pos!=last && *pos == col) return;
00825
00826
00827 size_type* end = std::lower_bound(pos,last,m);
00828 if (end==last)
00829 DUNE_THROW(ISTLError,"row is too small");
00830
00831
00832 std::copy_backward(pos,end,end+1);
00833 *pos = col;
00834
00835 }
00836
00838 void endindices ()
00839 {
00840 if (build_mode!=random)
00841 DUNE_THROW(ISTLError,"requires random build mode");
00842 if (ready==built)
00843 DUNE_THROW(ISTLError,"matrix already built up");
00844 if (ready==notbuilt)
00845 DUNE_THROW(ISTLError,"row sizes are not built up yet");
00846
00847
00848 RowIterator endi=end();
00849 for (RowIterator i=begin(); i!=endi; ++i)
00850 {
00851 ColIterator endj = (*i).end();
00852 for (ColIterator j=(*i).begin(); j!=endj; ++j){
00853 if (j.index()<0)
00854 {
00855 std::cout << "j[" << j.offset() << "]=" << j.index() << std::endl;
00856 DUNE_THROW(ISTLError,"undefined index detected");
00857 }
00858 if (j.index()>=m){
00859 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
00860 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
00861 r[i.index()].setsize(j.offset());
00862 break;
00863 }
00864 }
00865 }
00866
00867
00868 ready = built;
00869 }
00870
00871
00872
00874 BCRSMatrix& operator*= (const field_type& k)
00875 {
00876 if (nnz>0)
00877 {
00878
00879 for (size_type i=0; i<nnz; i++)
00880 a[i] *= k;
00881 }
00882 else
00883 {
00884 RowIterator endi=end();
00885 for (RowIterator i=begin(); i!=endi; ++i)
00886 {
00887 ColIterator endj = (*i).end();
00888 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00889 (*j) *= k;
00890 }
00891 }
00892
00893 return *this;
00894 }
00895
00897 BCRSMatrix& operator/= (const field_type& k)
00898 {
00899 if (nnz>0)
00900 {
00901
00902 for (size_type i=0; i<nnz; i++)
00903 a[i] /= k;
00904 }
00905 else
00906 {
00907 RowIterator endi=end();
00908 for (RowIterator i=begin(); i!=endi; ++i)
00909 {
00910 ColIterator endj = (*i).end();
00911 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00912 (*j) /= k;
00913 }
00914 }
00915
00916 return *this;
00917 }
00918
00919
00925 BCRSMatrix& operator+= (const BCRSMatrix& b)
00926 {
00927 #ifdef DUNE_ISTL_WITH_CHECKING
00928 if(N()!=b.N() || M() != b.M())
00929 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00930 #endif
00931 RowIterator endi=end();
00932 ConstRowIterator j=b.begin();
00933 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00934 i->operator+=(*j);
00935 }
00936
00937 return *this;
00938 }
00939
00945 BCRSMatrix& operator-= (const BCRSMatrix& b)
00946 {
00947 #ifdef DUNE_ISTL_WITH_CHECKING
00948 if(N()!=b.N() || M() != b.M())
00949 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00950 #endif
00951 RowIterator endi=end();
00952 ConstRowIterator j=b.begin();
00953 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00954 i->operator-=(*j);
00955 }
00956
00957 return *this;
00958 }
00959
00960
00962 template<class X, class Y>
00963 void mv (const X& x, Y& y) const
00964 {
00965 #ifdef DUNE_ISTL_WITH_CHECKING
00966 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00967 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00968 #endif
00969 ConstRowIterator endi=end();
00970 for (ConstRowIterator i=begin(); i!=endi; ++i)
00971 {
00972 y[i.index()]=0;
00973 ConstColIterator endj = (*i).end();
00974 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00975 (*j).umv(x[j.index()],y[i.index()]);
00976 }
00977 }
00978
00980 template<class X, class Y>
00981 void umv (const X& x, Y& y) const
00982 {
00983 #ifdef DUNE_ISTL_WITH_CHECKING
00984 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00985 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00986 #endif
00987 ConstRowIterator endi=end();
00988 for (ConstRowIterator i=begin(); i!=endi; ++i)
00989 {
00990 ConstColIterator endj = (*i).end();
00991 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00992 (*j).umv(x[j.index()],y[i.index()]);
00993 }
00994 }
00995
00997 template<class X, class Y>
00998 void mmv (const X& x, Y& y) const
00999 {
01000 #ifdef DUNE_ISTL_WITH_CHECKING
01001 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01002 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01003 #endif
01004 ConstRowIterator endi=end();
01005 for (ConstRowIterator i=begin(); i!=endi; ++i)
01006 {
01007 ConstColIterator endj = (*i).end();
01008 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01009 (*j).mmv(x[j.index()],y[i.index()]);
01010 }
01011 }
01012
01014 template<class X, class Y>
01015 void usmv (const field_type& alpha, const X& x, Y& y) const
01016 {
01017 #ifdef DUNE_ISTL_WITH_CHECKING
01018 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01019 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01020 #endif
01021 ConstRowIterator endi=end();
01022 for (ConstRowIterator i=begin(); i!=endi; ++i)
01023 {
01024 ConstColIterator endj = (*i).end();
01025 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01026 (*j).usmv(alpha,x[j.index()],y[i.index()]);
01027 }
01028 }
01029
01031 template<class X, class Y>
01032 void mtv (const X& x, Y& y) const
01033 {
01034 #ifdef DUNE_ISTL_WITH_CHECKING
01035 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01036 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01037 #endif
01038 for(size_type i=0; i<y.N(); ++i)
01039 y[i]=0;
01040 umtv(x,y);
01041 }
01042
01044 template<class X, class Y>
01045 void umtv (const X& x, Y& y) const
01046 {
01047 #ifdef DUNE_ISTL_WITH_CHECKING
01048 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01049 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01050 #endif
01051 ConstRowIterator endi=end();
01052 for (ConstRowIterator i=begin(); i!=endi; ++i)
01053 {
01054 ConstColIterator endj = (*i).end();
01055 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01056 (*j).umtv(x[i.index()],y[j.index()]);
01057 }
01058 }
01059
01061 template<class X, class Y>
01062 void mmtv (const X& x, Y& y) const
01063 {
01064 #ifdef DUNE_ISTL_WITH_CHECKING
01065 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01066 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01067 #endif
01068 ConstRowIterator endi=end();
01069 for (ConstRowIterator i=begin(); i!=endi; ++i)
01070 {
01071 ConstColIterator endj = (*i).end();
01072 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01073 (*j).mmtv(x[i.index()],y[j.index()]);
01074 }
01075 }
01076
01078 template<class X, class Y>
01079 void usmtv (const field_type& alpha, const X& x, Y& y) const
01080 {
01081 #ifdef DUNE_ISTL_WITH_CHECKING
01082 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01083 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01084 #endif
01085 ConstRowIterator endi=end();
01086 for (ConstRowIterator i=begin(); i!=endi; ++i)
01087 {
01088 ConstColIterator endj = (*i).end();
01089 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01090 (*j).usmtv(alpha,x[i.index()],y[j.index()]);
01091 }
01092 }
01093
01095 template<class X, class Y>
01096 void umhv (const X& x, Y& y) const
01097 {
01098 #ifdef DUNE_ISTL_WITH_CHECKING
01099 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01100 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01101 #endif
01102 ConstRowIterator endi=end();
01103 for (ConstRowIterator i=begin(); i!=endi; ++i)
01104 {
01105 ConstColIterator endj = (*i).end();
01106 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01107 (*j).umhv(x[i.index()],y[j.index()]);
01108 }
01109 }
01110
01112 template<class X, class Y>
01113 void mmhv (const X& x, Y& y) const
01114 {
01115 #ifdef DUNE_ISTL_WITH_CHECKING
01116 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01117 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01118 #endif
01119 ConstRowIterator endi=end();
01120 for (ConstRowIterator i=begin(); i!=endi; ++i)
01121 {
01122 ConstColIterator endj = (*i).end();
01123 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01124 (*j).mmhv(x[i.index()],y[j.index()]);
01125 }
01126 }
01127
01129 template<class X, class Y>
01130 void usmhv (const field_type& alpha, const X& x, Y& y) const
01131 {
01132 #ifdef DUNE_ISTL_WITH_CHECKING
01133 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01134 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01135 #endif
01136 ConstRowIterator endi=end();
01137 for (ConstRowIterator i=begin(); i!=endi; ++i)
01138 {
01139 ConstColIterator endj = (*i).end();
01140 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01141 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
01142 }
01143 }
01144
01145
01146
01147
01149 double frobenius_norm2 () const
01150 {
01151 double sum=0;
01152
01153 ConstRowIterator endi=end();
01154 for (ConstRowIterator i=begin(); i!=endi; ++i)
01155 {
01156 ConstColIterator endj = (*i).end();
01157 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01158 sum += (*j).frobenius_norm2();
01159 }
01160
01161 return sum;
01162 }
01163
01165 double frobenius_norm () const
01166 {
01167 return sqrt(frobenius_norm2());
01168 }
01169
01171 double infinity_norm () const
01172 {
01173 double max=0;
01174 ConstRowIterator endi=end();
01175 for (ConstRowIterator i=begin(); i!=endi; ++i)
01176 {
01177 double sum=0;
01178 ConstColIterator endj = (*i).end();
01179 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01180 sum += (*j).infinity_norm();
01181 max = std::max(max,sum);
01182 }
01183 return max;
01184 }
01185
01187 double infinity_norm_real () const
01188 {
01189 double max=0;
01190 ConstRowIterator endi=end();
01191 for (ConstRowIterator i=begin(); i!=endi; ++i)
01192 {
01193 double sum=0;
01194 ConstColIterator endj = (*i).end();
01195 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01196 sum += (*j).infinity_norm_real();
01197 max = std::max(max,sum);
01198 }
01199 return max;
01200 }
01201
01202
01203
01204
01206 size_type N () const
01207 {
01208 return n;
01209 }
01210
01212 size_type M () const
01213 {
01214 return m;
01215 }
01216
01218 size_type nonzeroes () const
01219 {
01220 return nnz;
01221 }
01222
01223
01224
01225
01227 bool exists (size_type i, size_type j) const
01228 {
01229 #ifdef DUNE_ISTL_WITH_CHECKING
01230 if (i<0 || i>=n) DUNE_THROW(ISTLError,"index out of range");
01231 if (j<0 || i>=m) DUNE_THROW(ISTLError,"index out of range");
01232 #endif
01233 if (r[i].size() && r[i].find(j)!=r[i].end())
01234 return true;
01235 else
01236 return false;
01237 }
01238
01239
01240 private:
01241
01242 BuildMode build_mode;
01243 BuildStage ready;
01244
01245
01246 size_type n;
01247 size_type m;
01248 size_type nnz;
01249
01250
01251
01252 row_type* r;
01253
01254
01255 B* a;
01256
01257
01258 Dune::shared_ptr<size_type> j;
01259
01260
01261 void setWindowPointers(ConstRowIterator row)
01262 {
01263 row_type current_row(a,j.get(),0);
01264 for (size_type i=0; i<n; i++, ++row){
01265
01266 size_type s = row->getsize();
01267
01268 if (s>0){
01269
01270 r[i].set(s,current_row.getptr(), current_row.getindexptr());
01271
01272 current_row.setptr(current_row.getptr()+s);
01273 current_row.setindexptr(current_row.getindexptr()+s);
01274 } else{
01275
01276 r[i].set(0,0,0);
01277 }
01278 }
01279 }
01280
01282 void copyWindowStructure(const BCRSMatrix& Mat)
01283 {
01284 setWindowPointers(Mat.begin());
01285
01286
01287 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
01288
01289
01290 build_mode = row_wise;
01291 ready = built;
01292 }
01293
01299 void deallocate(bool deallocateRows=true)
01300 {
01301
01302 if (nnz>0)
01303 {
01304
01305 j.reset();
01306 A::template free<B>(a);
01307 }
01308 else
01309 {
01310
01311 for (size_type i=0; i<n; i++)
01312 if (r[i].getsize()>0)
01313 {
01314 A::template free<size_type>(r[i].getindexptr());
01315 A::template free<B>(r[i].getptr());
01316 }
01317 }
01318
01319
01320 if (n>0 && deallocateRows) A::template free<row_type>(r);
01321
01322
01323 ready=notbuilt;
01324
01325 }
01326
01327 struct Deallocator
01328 {
01329 void operator()(size_type* p) const { A::template free<size_type>(p); }
01330 };
01331
01332
01350 void allocate(size_type rows, size_type columns, size_type nnz_=0, bool allocateRows=true)
01351 {
01352
01353 n = rows;
01354 m = columns;
01355 nnz = nnz_;
01356
01357
01358 if(allocateRows){
01359 if (n>0){
01360 r = A::template malloc<row_type>(rows);
01361 }else{
01362 r = 0;
01363 }
01364 }
01365
01366
01367
01368 if (nnz>0){
01369 a = A::template malloc<B>(nnz);
01370
01371 if (!j.get())
01372 j.reset(A::template malloc<size_type>(nnz),Deallocator());
01373 }else{
01374 a = 0;
01375 j.reset();
01376 }
01377
01378 ready = notbuilt;
01379 }
01380
01381 };
01382
01383
01386 }
01387
01388 #endif