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/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 #ifdef DUNE_EXPRESSIONTEMPLATES
00145 template<class B, class A>
00146 class BCRSMatrix : public ExprTmpl::Matrix< BCRSMatrix<B,A> >
00147 #else
00148 template<class B, class A=ISTLAllocator>
00149 class BCRSMatrix
00150 #endif
00151 {
00152 private:
00153 enum BuildStage{
00155 notbuilt=0,
00160 rowSizesBuilt=1,
00162 built=2
00163 };
00164
00165 public:
00166
00167
00168
00170 typedef typename B::field_type field_type;
00171
00173 typedef B block_type;
00174
00176 typedef A allocator_type;
00177
00179 typedef CompressedBlockVectorWindow<B,A> row_type;
00180
00182 typedef typename A::size_type size_type;
00183
00185 enum {
00187 blocklevel = B::blocklevel+1
00188 };
00189
00191 enum BuildMode {
00202 row_wise,
00211 random,
00215 unknown
00216 };
00217
00218
00219
00220
00222 row_type& operator[] (size_type i)
00223 {
00224 #ifdef DUNE_ISTL_WITH_CHECKING
00225 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00226 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00227 if (r[i].getptr()==0) DUNE_THROW(ISTLError,"row not initialized yet");
00228 #endif
00229 return r[i];
00230 }
00231
00233 const row_type& operator[] (size_type i) const
00234 {
00235 #ifdef DUNE_ISTL_WITH_CHECKING
00236 if (built!=ready) DUNE_THROW(ISTLError,"row not initialized yet");
00237 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00238 #endif
00239 return r[i];
00240 }
00241
00242
00243
00244
00246 template<class T>
00247 class RealRowIterator
00248 : public RandomAccessIteratorFacade<RealRowIterator<T>, T>
00249 {
00250
00251 public:
00253 typedef typename remove_const<T>::type ValueType;
00254
00255 friend class RandomAccessIteratorFacade<RealRowIterator<const ValueType>, const ValueType>;
00256 friend class RandomAccessIteratorFacade<RealRowIterator<ValueType>, ValueType>;
00257 friend class RealRowIterator<const ValueType>;
00258 friend class RealRowIterator<ValueType>;
00259
00261 RealRowIterator (row_type* _p, size_type _i)
00262 : p(_p), i(_i)
00263 {}
00264
00266 RealRowIterator ()
00267 : p(0), i(0)
00268 {}
00269
00270 RealRowIterator(const RealRowIterator<ValueType>& it)
00271 : p(it.p), i(it.i)
00272 {}
00273
00274
00276 size_type index () const
00277 {
00278 return i;
00279 }
00280
00281 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const
00282 {
00283 assert(other.p==p);
00284 return (other.i-i);
00285 }
00286
00287 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const
00288 {
00289 assert(other.p==p);
00290 return (other.i-i);
00291 }
00292
00294 bool equals (const RealRowIterator<ValueType>& other) const
00295 {
00296 assert(other.p==p);
00297 return i==other.i;
00298 }
00299
00301 bool equals (const RealRowIterator<const ValueType>& other) const
00302 {
00303 assert(other.p==p);
00304 return i==other.i;
00305 }
00306
00307 private:
00309 void increment()
00310 {
00311 ++i;
00312 }
00313
00315 void decrement()
00316 {
00317 --i;
00318 }
00319
00320 void advance(std::ptrdiff_t diff)
00321 {
00322 i+=diff;
00323 }
00324
00325 T& elementAt(std::ptrdiff_t diff) const
00326 {
00327 return p[i+diff];
00328 }
00329
00331 row_type& dereference () const
00332 {
00333 return p[i];
00334 }
00335
00336 row_type* p;
00337 size_type i;
00338 };
00339
00341 typedef RealRowIterator<row_type> iterator;
00342 typedef RealRowIterator<row_type> Iterator;
00343
00344
00346 Iterator begin ()
00347 {
00348 return Iterator(r,0);
00349 }
00350
00352 Iterator end ()
00353 {
00354 return Iterator(r,n);
00355 }
00356
00358 Iterator rbegin ()
00359 {
00360 return Iterator(r,n-1);
00361 }
00362
00364 Iterator rend ()
00365 {
00366 return Iterator(r,-1);
00367 }
00368
00370 typedef Iterator RowIterator;
00371
00373 typedef typename row_type::Iterator ColIterator;
00374
00376 typedef RealRowIterator<const row_type> const_iterator;
00377 typedef RealRowIterator<const row_type> ConstIterator;
00378
00379
00381 ConstIterator begin () const
00382 {
00383 return ConstIterator(r,0);
00384 }
00385
00387 ConstIterator end () const
00388 {
00389 return ConstIterator(r,n);
00390 }
00391
00393 ConstIterator rbegin () const
00394 {
00395 return ConstIterator(r,n-1);
00396 }
00397
00399 ConstIterator rend () const
00400 {
00401 return ConstIterator(r,-1);
00402 }
00403
00405 typedef ConstIterator ConstRowIterator;
00406
00408 typedef typename row_type::ConstIterator ConstColIterator;
00409
00410
00411
00413 BCRSMatrix ()
00414 : build_mode(unknown), ready(notbuilt), n(0), m(0), nnz(0),
00415 r(0), a(0), j(0)
00416 {}
00417
00419 BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm)
00420 : build_mode(bm), ready(notbuilt)
00421 {
00422 allocate(_n, _m, _nnz);
00423 }
00424
00426 BCRSMatrix (size_type _n, size_type _m, BuildMode bm)
00427 : build_mode(bm), ready(notbuilt)
00428 {
00429 allocate(_n, _m);
00430 }
00431
00433 BCRSMatrix (const BCRSMatrix& Mat)
00434 : n(Mat.n), nnz(0)
00435 {
00436
00437 int _nnz = Mat.nnz;
00438
00439
00440 if (_nnz<=0)
00441 {
00442 _nnz = 0;
00443 for (size_type i=0; i<n; i++)
00444 _nnz += Mat.r[i].getsize();
00445 }
00446
00447 allocate(Mat.n, Mat.m, _nnz);
00448
00449
00450 copyWindowStructure(Mat);
00451 }
00452
00454 ~BCRSMatrix ()
00455 {
00456 deallocate();
00457 }
00458
00463 void setBuildMode(BuildMode bm)
00464 {
00465 if(ready==notbuilt)
00466 build_mode = bm;
00467 else
00468 DUNE_THROW(InvalidStateException, "Matrix structure is already built (ready="<<ready<<").");
00469 }
00470
00486 void setSize(size_type rows, size_type columns, size_type nnz=0)
00487 {
00488
00489 deallocate();
00490
00491
00492 allocate(rows, columns, nnz);
00493 }
00494
00496 BCRSMatrix& operator= (const BCRSMatrix& Mat)
00497 {
00498
00499 if (&Mat==this) return *this;
00500
00501
00502 deallocate(false);
00503
00504
00505 if (n>0 && n!=Mat.n)
00506
00507 A::template free<row_type>(r);
00508
00509 nnz=Mat.nnz;
00510 if (nnz<=0)
00511 {
00512 for (size_type i=0; i<Mat.n; i++)
00513 nnz += Mat.r[i].getsize();
00514 }
00515
00516
00517 allocate(Mat.n, Mat.m, nnz, n!=Mat.n);
00518
00519
00520 copyWindowStructure(Mat);
00521 return *this;
00522 }
00523
00525 BCRSMatrix& operator= (const field_type& k)
00526 {
00527 for (size_type i=0; i<n; i++) r[i] = k;
00528 return *this;
00529 }
00530
00531
00532
00534 class CreateIterator
00535 {
00536 public:
00538 CreateIterator (BCRSMatrix& _Mat, size_type _i)
00539 : Mat(_Mat), i(_i), nnz(0), current_row(Mat.a, Mat.j, 0)
00540 {
00541 if (i==0 && Mat.ready)
00542 DUNE_THROW(ISTLError,"creation only allowed for uninitialized matrix");
00543 if(Mat.build_mode!=row_wise)
00544 {
00545 if(Mat.build_mode==unknown)
00546 Mat.build_mode=row_wise;
00547 else
00548 DUNE_THROW(ISTLError,"creation only allowed if row wise allocation was requested in the constructor");
00549 }
00550 }
00551
00553 CreateIterator& operator++()
00554 {
00555
00556 if (Mat.ready)
00557 DUNE_THROW(ISTLError,"matrix already built up");
00558
00559
00560
00561
00562
00563
00564 size_type s = pattern.size();
00565
00566 if(s>0){
00567
00568 nnz += s;
00569
00570
00571 if (Mat.nnz>0)
00572 {
00573
00574
00575
00576 if (nnz>Mat.nnz)
00577 DUNE_THROW(ISTLError,"allocated nnz too small");
00578
00579
00580 Mat.r[i].set(s,current_row.getptr(),current_row.getindexptr());
00581 current_row.setptr(current_row.getptr()+s);
00582 current_row.setindexptr(current_row.getindexptr()+s);
00583 }else{
00584
00585
00586 B* a = A::template malloc<B>(s);
00587 size_type* j = A::template malloc<size_type>(s);
00588 Mat.r[i].set(s,a,j);
00589 }
00590 }else
00591
00592 Mat.r[i].set(0,0,0);
00593
00594
00595 size_type k=0;
00596 size_type *j = Mat.r[i].getindexptr();
00597 for (typename std::set<size_type>::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
00598 j[k++] = *it;
00599
00600
00601 i++;
00602 pattern.clear();
00603
00604
00605 if (i==Mat.n)
00606 {
00607 Mat.ready = built;
00608 if(Mat.nnz>0)
00609
00610
00611 Mat.nnz=nnz;
00612 }
00613
00614 return *this;
00615 }
00616
00618 bool operator!= (const CreateIterator& it) const
00619 {
00620 return (i!=it.i) || (&Mat!=&it.Mat);
00621 }
00622
00624 bool operator== (const CreateIterator& it) const
00625 {
00626 return (i==it.i) && (&Mat==&it.Mat);
00627 }
00628
00630 size_type index () const
00631 {
00632 return i;
00633 }
00634
00636 void insert (size_type j)
00637 {
00638 pattern.insert(j);
00639 }
00640
00642 bool contains (size_type j)
00643 {
00644 if (pattern.find(j)!=pattern.end())
00645 return true;
00646 else
00647 return false;
00648 }
00654 typename std::set<size_type>::size_type size() const
00655 {
00656 return pattern.size();
00657 }
00658
00659 private:
00660 BCRSMatrix& Mat;
00661 size_type i;
00662 size_type nnz;
00663 std::set<size_type> pattern;
00664 row_type current_row;
00665 };
00666
00668 friend class CreateIterator;
00669
00671 CreateIterator createbegin ()
00672 {
00673 return CreateIterator(*this,0);
00674 }
00675
00677 CreateIterator createend ()
00678 {
00679 return CreateIterator(*this,n);
00680 }
00681
00682
00683
00684
00686 void setrowsize (size_type i, size_type s)
00687 {
00688 if (build_mode!=random)
00689 DUNE_THROW(ISTLError,"requires random build mode");
00690 if (ready)
00691 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00692
00693 r[i].setsize(s);
00694 }
00695
00697 size_type getrowsize (size_type i) const
00698 {
00699 #ifdef DUNE_ISTL_WITH_CHECKING
00700 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00701 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00702 #endif
00703 return r[i].getsize();
00704 }
00705
00707 void incrementrowsize (size_type i, size_type s = 1)
00708 {
00709 if (build_mode!=random)
00710 DUNE_THROW(ISTLError,"requires random build mode");
00711 if (ready)
00712 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00713
00714 r[i].setsize(r[i].getsize()+s);
00715 }
00716
00718 void endrowsizes ()
00719 {
00720 if (build_mode!=random)
00721 DUNE_THROW(ISTLError,"requires random build mode");
00722 if (ready)
00723 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00724
00725
00726 size_type total=0;
00727 for (size_type i=0; i<n; i++)
00728 {
00729 if (r[i].getsize()<0)
00730 DUNE_THROW(ISTLError,"rowsize must be nonnegative");
00731 total += r[i].getsize();
00732 }
00733
00734
00735 allocate(n,m,total,false);
00736
00737
00738 setWindowPointers(begin());
00739
00740
00741
00742 for (size_type k=0; k<nnz; k++)
00743 j[k] = m;
00744 ready = rowSizesBuilt;
00745 }
00746
00748
00758 void addindex (size_type row, size_type col)
00759 {
00760 if (build_mode!=random)
00761 DUNE_THROW(ISTLError,"requires random build mode");
00762 if (ready==built)
00763 DUNE_THROW(ISTLError,"matrix already built up");
00764 if (ready==notbuilt)
00765 DUNE_THROW(ISTLError,"matrix row sizes not built up yet");
00766
00767 if (col >= m)
00768 DUNE_THROW(ISTLError,"column index exceeds matrix size");
00769
00770
00771 size_type* const first = r[row].getindexptr();
00772 size_type* const last = first + r[row].getsize();
00773
00774
00775 size_type* pos = std::lower_bound(first,last,col);
00776
00777
00778 if (pos!=last && *pos == col) return;
00779
00780
00781 size_type* end = std::lower_bound(pos,last,m);
00782 if (end==last)
00783 DUNE_THROW(ISTLError,"row is too small");
00784
00785
00786 std::copy_backward(pos,end,end+1);
00787 *pos = col;
00788
00789 }
00790
00792 void endindices ()
00793 {
00794 if (build_mode!=random)
00795 DUNE_THROW(ISTLError,"requires random build mode");
00796 if (ready==built)
00797 DUNE_THROW(ISTLError,"matrix already built up");
00798 if (ready==notbuilt)
00799 DUNE_THROW(ISTLError,"row sizes are not built up yet");
00800
00801
00802 RowIterator endi=end();
00803 for (RowIterator i=begin(); i!=endi; ++i)
00804 {
00805 ColIterator endj = (*i).end();
00806 for (ColIterator j=(*i).begin(); j!=endj; ++j){
00807 if (j.index()<0)
00808 {
00809 std::cout << "j[" << j.offset() << "]=" << j.index() << std::endl;
00810 DUNE_THROW(ISTLError,"undefined index detected");
00811 }
00812 if (j.index()>=m){
00813 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
00814 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
00815 r[i.index()].setsize(j.offset());
00816 break;
00817 }
00818 }
00819 }
00820
00821
00822 ready = built;
00823 }
00824
00825
00826
00828 BCRSMatrix& operator*= (const field_type& k)
00829 {
00830 if (nnz>0)
00831 {
00832
00833 for (size_type i=0; i<nnz; i++)
00834 a[i] *= k;
00835 }
00836 else
00837 {
00838 RowIterator endi=end();
00839 for (RowIterator i=begin(); i!=endi; ++i)
00840 {
00841 ColIterator endj = (*i).end();
00842 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00843 (*j) *= k;
00844 }
00845 }
00846
00847 return *this;
00848 }
00849
00851 BCRSMatrix& operator/= (const field_type& k)
00852 {
00853 if (nnz>0)
00854 {
00855
00856 for (size_type i=0; i<nnz; i++)
00857 a[i] /= k;
00858 }
00859 else
00860 {
00861 RowIterator endi=end();
00862 for (RowIterator i=begin(); i!=endi; ++i)
00863 {
00864 ColIterator endj = (*i).end();
00865 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00866 (*j) /= k;
00867 }
00868 }
00869
00870 return *this;
00871 }
00872
00873
00879 BCRSMatrix& operator+= (const BCRSMatrix& b)
00880 {
00881 #ifdef DUNE_ISTL_WITH_CHECKING
00882 if(N()!=b.N() || M() != b.M())
00883 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00884 #endif
00885 RowIterator endi=end();
00886 ConstRowIterator j=b.begin();
00887 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00888 i->operator+=(*j);
00889 }
00890
00891 return *this;
00892 }
00893
00899 BCRSMatrix& operator-= (const BCRSMatrix& b)
00900 {
00901 #ifdef DUNE_ISTL_WITH_CHECKING
00902 if(N()!=b.N() || M() != b.M())
00903 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00904 #endif
00905 RowIterator endi=end();
00906 ConstRowIterator j=b.begin();
00907 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00908 i->operator-=(*j);
00909 }
00910
00911 return *this;
00912 }
00913
00914
00916 template<class X, class Y>
00917 void mv (const X& x, Y& y) const
00918 {
00919 #ifdef DUNE_ISTL_WITH_CHECKING
00920 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00921 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00922 #endif
00923 ConstRowIterator endi=end();
00924 for (ConstRowIterator i=begin(); i!=endi; ++i)
00925 {
00926 y[i.index()]=0;
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 umv (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).umv(x[j.index()],y[i.index()]);
00947 }
00948 }
00949
00951 template<class X, class Y>
00952 void mmv (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).mmv(x[j.index()],y[i.index()]);
00964 }
00965 }
00966
00968 template<class X, class Y>
00969 void usmv (const field_type& alpha, const X& x, Y& y) const
00970 {
00971 #ifdef DUNE_ISTL_WITH_CHECKING
00972 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00973 if (y.N()!=N()) 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).usmv(alpha,x[j.index()],y[i.index()]);
00981 }
00982 }
00983
00985 template<class X, class Y>
00986 void umtv (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).umtv(x[i.index()],y[j.index()]);
00998 }
00999 }
01000
01002 template<class X, class Y>
01003 void mmtv (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).mmtv(x[i.index()],y[j.index()]);
01015 }
01016 }
01017
01019 template<class X, class Y>
01020 void usmtv (const field_type& alpha, 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).usmtv(alpha,x[i.index()],y[j.index()]);
01032 }
01033 }
01034
01036 template<class X, class Y>
01037 void umhv (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).umhv(x[i.index()],y[j.index()]);
01049 }
01050 }
01051
01053 template<class X, class Y>
01054 void mmhv (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).mmhv(x[i.index()],y[j.index()]);
01066 }
01067 }
01068
01070 template<class X, class Y>
01071 void usmhv (const field_type& alpha, const X& x, Y& y) const
01072 {
01073 #ifdef DUNE_ISTL_WITH_CHECKING
01074 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01075 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01076 #endif
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 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
01083 }
01084 }
01085
01086
01087
01088
01090 double frobenius_norm2 () const
01091 {
01092 double sum=0;
01093
01094 ConstRowIterator endi=end();
01095 for (ConstRowIterator i=begin(); i!=endi; ++i)
01096 {
01097 ConstColIterator endj = (*i).end();
01098 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01099 sum += (*j).frobenius_norm2();
01100 }
01101
01102 return sum;
01103 }
01104
01106 double frobenius_norm () const
01107 {
01108 return sqrt(frobenius_norm2());
01109 }
01110
01112 double infinity_norm () const
01113 {
01114 double max=0;
01115 ConstRowIterator endi=end();
01116 for (ConstRowIterator i=begin(); i!=endi; ++i)
01117 {
01118 double sum=0;
01119 ConstColIterator endj = (*i).end();
01120 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01121 sum += (*j).infinity_norm();
01122 max = std::max(max,sum);
01123 }
01124 return max;
01125 }
01126
01128 double infinity_norm_real () const
01129 {
01130 double max=0;
01131 ConstRowIterator endi=end();
01132 for (ConstRowIterator i=begin(); i!=endi; ++i)
01133 {
01134 double sum=0;
01135 ConstColIterator endj = (*i).end();
01136 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01137 sum += (*j).infinity_norm_real();
01138 max = std::max(max,sum);
01139 }
01140 return max;
01141 }
01142
01143
01144
01145
01147 size_type N () const
01148 {
01149 return n;
01150 }
01151
01153 size_type M () const
01154 {
01155 return m;
01156 }
01157
01159 size_type nonzeroes () const
01160 {
01161 return nnz;
01162 }
01163
01167 size_type rowdim (size_type i) const
01168 {
01169 const B* row = r[i].getptr();
01170 if(row)
01171 return row->rowdim();
01172 else
01173 return 0;
01174 }
01175
01179 size_type coldim (size_type c) const
01180 {
01181
01182 if (nnz>0)
01183 {
01184 for (size_type k=0; k<nnz; k++) {
01185 if (j[k]==c) {
01186 return a[k].coldim();
01187 }
01188 }
01189 }
01190 else
01191 {
01192 for (size_type i=0; i<n; i++)
01193 {
01194 size_type* j = r[i].getindexptr();
01195 B* a = r[i].getptr();
01196 for (size_type k=0; k<r[i].getsize(); k++)
01197 if (j[k]==c) {
01198 return a[k].coldim();
01199 }
01200 }
01201 }
01202
01203
01204 return 0;
01205 }
01206
01208 size_type rowdim () const
01209 {
01210 size_type nn=0;
01211 for (size_type i=0; i<n; i++)
01212 nn += rowdim(i);
01213 return nn;
01214 }
01215
01217 size_type coldim () const
01218 {
01219
01220
01221
01222 std::vector<size_type> coldims(M(),-1);
01223 for (ConstRowIterator row=begin(); row!=end(); ++row)
01224 for (ConstColIterator col=row->begin(); col!=row->end(); ++col)
01225
01226 if (coldims[col.index()]==-1)
01227 coldims[col.index()] = col->coldim();
01228
01229 size_type sum = 0;
01230 for (typename std::vector<size_type>::iterator it=coldims.begin(); it!=coldims.end(); ++it)
01231
01232 if ((*it)>=0)
01233 sum += *it;
01234
01235 return sum;
01236 }
01237
01238
01239
01241 bool exists (size_type i, size_type j) const
01242 {
01243 #ifdef DUNE_ISTL_WITH_CHECKING
01244 if (i<0 || i>=n) DUNE_THROW(ISTLError,"index out of range");
01245 if (j<0 || i>=m) DUNE_THROW(ISTLError,"index out of range");
01246 #endif
01247 if (r[i].size() && r[i].find(j)!=r[i].end())
01248 return true;
01249 else
01250 return false;
01251 }
01252
01253
01254 private:
01255
01256 BuildMode build_mode;
01257 BuildStage ready;
01258
01259
01260 size_type n;
01261 size_type m;
01262 size_type nnz;
01263
01264
01265
01266 row_type* r;
01267
01268
01269 B* a;
01270 size_type* j;
01271
01272 void setWindowPointers(ConstRowIterator row)
01273 {
01274 row_type current_row(a,j,0);
01275 for (size_type i=0; i<n; i++, ++row){
01276
01277 size_type s = row->getsize();
01278
01279 if (s>0){
01280
01281 r[i].set(s,current_row.getptr(), current_row.getindexptr());
01282
01283 current_row.setptr(current_row.getptr()+s);
01284 current_row.setindexptr(current_row.getindexptr()+s);
01285 } else{
01286
01287 r[i].set(0,0,0);
01288 }
01289 }
01290 }
01291
01293 void copyWindowStructure(const BCRSMatrix& Mat)
01294 {
01295 setWindowPointers(Mat.begin());
01296
01297
01298 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
01299
01300
01301 build_mode = row_wise;
01302 ready = built;
01303 }
01304
01310 void deallocate(bool deallocateRows=true)
01311 {
01312
01313 if (nnz>0)
01314 {
01315
01316 A::template free<size_type>(j);
01317 A::template free<B>(a);
01318 }
01319 else
01320 {
01321
01322 for (size_type i=0; i<n; i++)
01323 if (r[i].getsize()>0)
01324 {
01325 A::template free<size_type>(r[i].getindexptr());
01326 A::template free<B>(r[i].getptr());
01327 }
01328 }
01329
01330
01331 if (n>0 && deallocateRows) A::template free<row_type>(r);
01332
01333
01334 ready=notbuilt;
01335
01336 }
01337
01355 void allocate(size_type rows, size_type columns, size_type nnz_=0, bool allocateRows=true)
01356 {
01357
01358 n = rows;
01359 m = columns;
01360 nnz = nnz_;
01361
01362
01363 if(allocateRows){
01364 if (n>0){
01365 r = A::template malloc<row_type>(rows);
01366 }else{
01367 r = 0;
01368 }
01369 }
01370
01371
01372
01373 if (nnz>0){
01374 a = A::template malloc<B>(nnz);
01375 j = A::template malloc<size_type>(nnz);
01376 }else{
01377 a = 0;
01378 j = 0;
01379 }
01380
01381 ready = notbuilt;
01382 }
01383
01384 };
01385
01386
01387 #ifdef DUNE_EXPRESSIONTEMPLATES
01388 template <class B, class A>
01389 struct FieldType< BCRSMatrix<B,A> >
01390 {
01391 typedef typename FieldType<B>::type type;
01392 };
01393
01394 template <class B, class A>
01395 struct BlockType< BCRSMatrix<B,A> >
01396 {
01397 typedef B type;
01398 };
01399 template <class B, class A>
01400 struct RowType< BCRSMatrix<B,A> >
01401 {
01402 typedef CompressedBlockVectorWindow<B,A> type;
01403 };
01404 #endif
01405
01408 }
01409
01410 #endif