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 if(Mat.build_mode==unknown)
00545 Mat.build_mode=row_wise;
00546 else
00547 DUNE_THROW(ISTLError,"creation only allowed if row wise allocation was requested in the constructor");
00548 }
00549
00551 CreateIterator& operator++()
00552 {
00553
00554 if (Mat.ready)
00555 DUNE_THROW(ISTLError,"matrix already built up");
00556
00557
00558
00559
00560
00561
00562 size_type s = pattern.size();
00563
00564 if(s>0){
00565
00566 nnz += s;
00567
00568
00569 if (Mat.nnz>0)
00570 {
00571
00572
00573
00574 if (nnz>Mat.nnz)
00575 DUNE_THROW(ISTLError,"allocated nnz too small");
00576
00577
00578 Mat.r[i].set(s,current_row.getptr(),current_row.getindexptr());
00579 current_row.setptr(current_row.getptr()+s);
00580 current_row.setindexptr(current_row.getindexptr()+s);
00581 }else{
00582
00583
00584 B* a = A::template malloc<B>(s);
00585 size_type* j = A::template malloc<size_type>(s);
00586 Mat.r[i].set(s,a,j);
00587 }
00588 }else
00589
00590 Mat.r[i].set(0,0,0);
00591
00592
00593 size_type k=0;
00594 size_type *j = Mat.r[i].getindexptr();
00595 for (typename std::set<size_type>::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
00596 j[k++] = *it;
00597
00598
00599 i++;
00600 pattern.clear();
00601
00602
00603 if (i==Mat.n)
00604 {
00605 Mat.ready = built;
00606 if(Mat.nnz>0)
00607
00608
00609 Mat.nnz=nnz;
00610 }
00611
00612 return *this;
00613 }
00614
00616 bool operator!= (const CreateIterator& it) const
00617 {
00618 return (i!=it.i) || (&Mat!=&it.Mat);
00619 }
00620
00622 bool operator== (const CreateIterator& it) const
00623 {
00624 return (i==it.i) && (&Mat==&it.Mat);
00625 }
00626
00628 size_type index () const
00629 {
00630 return i;
00631 }
00632
00634 void insert (size_type j)
00635 {
00636 pattern.insert(j);
00637 }
00638
00640 bool contains (size_type j)
00641 {
00642 if (pattern.find(j)!=pattern.end())
00643 return true;
00644 else
00645 return false;
00646 }
00652 typename std::set<size_type>::size_type size() const
00653 {
00654 return pattern.size();
00655 }
00656
00657 private:
00658 BCRSMatrix& Mat;
00659 size_type i;
00660 size_type nnz;
00661 std::set<size_type> pattern;
00662 row_type current_row;
00663 };
00664
00666 friend class CreateIterator;
00667
00669 CreateIterator createbegin ()
00670 {
00671 return CreateIterator(*this,0);
00672 }
00673
00675 CreateIterator createend ()
00676 {
00677 return CreateIterator(*this,n);
00678 }
00679
00680
00681
00682
00684 void setrowsize (size_type i, size_type s)
00685 {
00686 if (build_mode!=random)
00687 DUNE_THROW(ISTLError,"requires random build mode");
00688 if (ready)
00689 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00690
00691 r[i].setsize(s);
00692 }
00693
00695 size_type getrowsize (size_type i) const
00696 {
00697 #ifdef DUNE_ISTL_WITH_CHECKING
00698 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet");
00699 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
00700 #endif
00701 return r[i].getsize();
00702 }
00703
00705 void incrementrowsize (size_type i, size_type s = 1)
00706 {
00707 if (build_mode!=random)
00708 DUNE_THROW(ISTLError,"requires random build mode");
00709 if (ready)
00710 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00711
00712 r[i].setsize(r[i].getsize()+s);
00713 }
00714
00716 void endrowsizes ()
00717 {
00718 if (build_mode!=random)
00719 DUNE_THROW(ISTLError,"requires random build mode");
00720 if (ready)
00721 DUNE_THROW(ISTLError,"matrix row sizes already built up");
00722
00723
00724 size_type total=0;
00725 for (size_type i=0; i<n; i++)
00726 {
00727 if (r[i].getsize()<0)
00728 DUNE_THROW(ISTLError,"rowsize must be nonnegative");
00729 total += r[i].getsize();
00730 }
00731
00732
00733 allocate(n,m,total,false);
00734
00735
00736 setWindowPointers(begin());
00737
00738
00739
00740 for (size_type k=0; k<nnz; k++)
00741 j[k] = m;
00742 ready = rowSizesBuilt;
00743 }
00744
00746
00756 void addindex (size_type row, size_type col)
00757 {
00758 if (build_mode!=random)
00759 DUNE_THROW(ISTLError,"requires random build mode");
00760 if (ready==built)
00761 DUNE_THROW(ISTLError,"matrix already built up");
00762 if (ready==notbuilt)
00763 DUNE_THROW(ISTLError,"matrix row sizes not built up yet");
00764
00765 if (col >= m)
00766 DUNE_THROW(ISTLError,"column index exceeds matrix size");
00767
00768
00769 size_type* const first = r[row].getindexptr();
00770 size_type* const last = first + r[row].getsize();
00771
00772
00773 size_type* pos = std::lower_bound(first,last,col);
00774
00775
00776 if (pos!=last && *pos == col) return;
00777
00778
00779 size_type* end = std::lower_bound(pos,last,m);
00780 if (end==last)
00781 DUNE_THROW(ISTLError,"row is too small");
00782
00783
00784 std::copy_backward(pos,end,end+1);
00785 *pos = col;
00786
00787 }
00788
00790 void endindices ()
00791 {
00792 if (build_mode!=random)
00793 DUNE_THROW(ISTLError,"requires random build mode");
00794 if (ready==built)
00795 DUNE_THROW(ISTLError,"matrix already built up");
00796 if (ready==notbuilt)
00797 DUNE_THROW(ISTLError,"row sizes are not built up yet");
00798
00799
00800 RowIterator endi=end();
00801 for (RowIterator i=begin(); i!=endi; ++i)
00802 {
00803 ColIterator endj = (*i).end();
00804 for (ColIterator j=(*i).begin(); j!=endj; ++j){
00805 if (j.index()<0)
00806 {
00807 std::cout << "j[" << j.offset() << "]=" << j.index() << std::endl;
00808 DUNE_THROW(ISTLError,"undefined index detected");
00809 }
00810 if (j.index()>=m){
00811 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset()
00812 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
00813 r[i.index()].setsize(j.offset());
00814 break;
00815 }
00816 }
00817 }
00818
00819
00820 ready = built;
00821 }
00822
00823
00824
00826 BCRSMatrix& operator*= (const field_type& k)
00827 {
00828 if (nnz>0)
00829 {
00830
00831 for (size_type i=0; i<nnz; i++)
00832 a[i] *= k;
00833 }
00834 else
00835 {
00836 RowIterator endi=end();
00837 for (RowIterator i=begin(); i!=endi; ++i)
00838 {
00839 ColIterator endj = (*i).end();
00840 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00841 (*j) *= k;
00842 }
00843 }
00844
00845 return *this;
00846 }
00847
00849 BCRSMatrix& operator/= (const field_type& k)
00850 {
00851 if (nnz>0)
00852 {
00853
00854 for (size_type i=0; i<nnz; i++)
00855 a[i] /= k;
00856 }
00857 else
00858 {
00859 RowIterator endi=end();
00860 for (RowIterator i=begin(); i!=endi; ++i)
00861 {
00862 ColIterator endj = (*i).end();
00863 for (ColIterator j=(*i).begin(); j!=endj; ++j)
00864 (*j) /= k;
00865 }
00866 }
00867
00868 return *this;
00869 }
00870
00871
00877 BCRSMatrix& operator+= (const BCRSMatrix& b)
00878 {
00879 #ifdef DUNE_ISTL_WITH_CHECKING
00880 if(N()!=b.N() || M() != b.M())
00881 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00882 #endif
00883 RowIterator endi=end();
00884 ConstRowIterator j=b.begin();
00885 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00886 i->operator+=(*j);
00887 }
00888
00889 return *this;
00890 }
00891
00897 BCRSMatrix& operator-= (const BCRSMatrix& b)
00898 {
00899 #ifdef DUNE_ISTL_WITH_CHECKING
00900 if(N()!=b.N() || M() != b.M())
00901 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00902 #endif
00903 RowIterator endi=end();
00904 ConstRowIterator j=b.begin();
00905 for (RowIterator i=begin(); i!=endi; ++i, ++j){
00906 i->operator-=(*j);
00907 }
00908
00909 return *this;
00910 }
00911
00912
00914 template<class X, class Y>
00915 void mv (const X& x, Y& y) const
00916 {
00917 #ifdef DUNE_ISTL_WITH_CHECKING
00918 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00919 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00920 #endif
00921 ConstRowIterator endi=end();
00922 for (ConstRowIterator i=begin(); i!=endi; ++i)
00923 {
00924 y[i.index()]=0;
00925 ConstColIterator endj = (*i).end();
00926 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00927 (*j).umv(x[j.index()],y[i.index()]);
00928 }
00929 }
00930
00932 template<class X, class Y>
00933 void umv (const X& x, Y& y) const
00934 {
00935 #ifdef DUNE_ISTL_WITH_CHECKING
00936 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00937 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00938 #endif
00939 ConstRowIterator endi=end();
00940 for (ConstRowIterator i=begin(); i!=endi; ++i)
00941 {
00942 ConstColIterator endj = (*i).end();
00943 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00944 (*j).umv(x[j.index()],y[i.index()]);
00945 }
00946 }
00947
00949 template<class X, class Y>
00950 void mmv (const X& x, Y& y) const
00951 {
00952 #ifdef DUNE_ISTL_WITH_CHECKING
00953 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00954 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00955 #endif
00956 ConstRowIterator endi=end();
00957 for (ConstRowIterator i=begin(); i!=endi; ++i)
00958 {
00959 ConstColIterator endj = (*i).end();
00960 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00961 (*j).mmv(x[j.index()],y[i.index()]);
00962 }
00963 }
00964
00966 template<class X, class Y>
00967 void usmv (const field_type& alpha, const X& x, Y& y) const
00968 {
00969 #ifdef DUNE_ISTL_WITH_CHECKING
00970 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00971 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00972 #endif
00973 ConstRowIterator endi=end();
00974 for (ConstRowIterator i=begin(); i!=endi; ++i)
00975 {
00976 ConstColIterator endj = (*i).end();
00977 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00978 (*j).usmv(alpha,x[j.index()],y[i.index()]);
00979 }
00980 }
00981
00983 template<class X, class Y>
00984 void umtv (const X& x, Y& y) const
00985 {
00986 #ifdef DUNE_ISTL_WITH_CHECKING
00987 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00988 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00989 #endif
00990 ConstRowIterator endi=end();
00991 for (ConstRowIterator i=begin(); i!=endi; ++i)
00992 {
00993 ConstColIterator endj = (*i).end();
00994 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00995 (*j).umtv(x[i.index()],y[j.index()]);
00996 }
00997 }
00998
01000 template<class X, class Y>
01001 void mmtv (const X& x, Y& y) const
01002 {
01003 #ifdef DUNE_ISTL_WITH_CHECKING
01004 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01005 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01006 #endif
01007 ConstRowIterator endi=end();
01008 for (ConstRowIterator i=begin(); i!=endi; ++i)
01009 {
01010 ConstColIterator endj = (*i).end();
01011 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01012 (*j).mmtv(x[i.index()],y[j.index()]);
01013 }
01014 }
01015
01017 template<class X, class Y>
01018 void usmtv (const field_type& alpha, const X& x, Y& y) const
01019 {
01020 #ifdef DUNE_ISTL_WITH_CHECKING
01021 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01022 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01023 #endif
01024 ConstRowIterator endi=end();
01025 for (ConstRowIterator i=begin(); i!=endi; ++i)
01026 {
01027 ConstColIterator endj = (*i).end();
01028 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01029 (*j).usmtv(alpha,x[i.index()],y[j.index()]);
01030 }
01031 }
01032
01034 template<class X, class Y>
01035 void umhv (const X& x, Y& y) const
01036 {
01037 #ifdef DUNE_ISTL_WITH_CHECKING
01038 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01039 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01040 #endif
01041 ConstRowIterator endi=end();
01042 for (ConstRowIterator i=begin(); i!=endi; ++i)
01043 {
01044 ConstColIterator endj = (*i).end();
01045 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01046 (*j).umhv(x[i.index()],y[j.index()]);
01047 }
01048 }
01049
01051 template<class X, class Y>
01052 void mmhv (const X& x, Y& y) const
01053 {
01054 #ifdef DUNE_ISTL_WITH_CHECKING
01055 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01056 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01057 #endif
01058 ConstRowIterator endi=end();
01059 for (ConstRowIterator i=begin(); i!=endi; ++i)
01060 {
01061 ConstColIterator endj = (*i).end();
01062 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01063 (*j).mmhv(x[i.index()],y[j.index()]);
01064 }
01065 }
01066
01068 template<class X, class Y>
01069 void usmhv (const field_type& alpha, const X& x, Y& y) const
01070 {
01071 #ifdef DUNE_ISTL_WITH_CHECKING
01072 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
01073 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
01074 #endif
01075 ConstRowIterator endi=end();
01076 for (ConstRowIterator i=begin(); i!=endi; ++i)
01077 {
01078 ConstColIterator endj = (*i).end();
01079 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01080 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
01081 }
01082 }
01083
01084
01085
01086
01088 double frobenius_norm2 () const
01089 {
01090 double sum=0;
01091
01092 ConstRowIterator endi=end();
01093 for (ConstRowIterator i=begin(); i!=endi; ++i)
01094 {
01095 ConstColIterator endj = (*i).end();
01096 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01097 sum += (*j).frobenius_norm2();
01098 }
01099
01100 return sum;
01101 }
01102
01104 double frobenius_norm () const
01105 {
01106 return sqrt(frobenius_norm2());
01107 }
01108
01110 double infinity_norm () const
01111 {
01112 double max=0;
01113 ConstRowIterator endi=end();
01114 for (ConstRowIterator i=begin(); i!=endi; ++i)
01115 {
01116 double sum=0;
01117 ConstColIterator endj = (*i).end();
01118 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01119 sum += (*j).infinity_norm();
01120 max = std::max(max,sum);
01121 }
01122 return max;
01123 }
01124
01126 double infinity_norm_real () const
01127 {
01128 double max=0;
01129 ConstRowIterator endi=end();
01130 for (ConstRowIterator i=begin(); i!=endi; ++i)
01131 {
01132 double sum=0;
01133 ConstColIterator endj = (*i).end();
01134 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
01135 sum += (*j).infinity_norm_real();
01136 max = std::max(max,sum);
01137 }
01138 return max;
01139 }
01140
01141
01142
01143
01145 size_type N () const
01146 {
01147 return n;
01148 }
01149
01151 size_type M () const
01152 {
01153 return m;
01154 }
01155
01157 size_type nonzeroes () const
01158 {
01159 return nnz;
01160 }
01161
01165 size_type rowdim (size_type i) const
01166 {
01167 const B* row = r[i].getptr();
01168 if(row)
01169 return row->rowdim();
01170 else
01171 return 0;
01172 }
01173
01177 size_type coldim (size_type c) const
01178 {
01179
01180 if (nnz>0)
01181 {
01182 for (size_type k=0; k<nnz; k++) {
01183 if (j[k]==c) {
01184 return a[k].coldim();
01185 }
01186 }
01187 }
01188 else
01189 {
01190 for (size_type i=0; i<n; i++)
01191 {
01192 size_type* j = r[i].getindexptr();
01193 B* a = r[i].getptr();
01194 for (size_type k=0; k<r[i].getsize(); k++)
01195 if (j[k]==c) {
01196 return a[k].coldim();
01197 }
01198 }
01199 }
01200
01201
01202 return 0;
01203 }
01204
01206 size_type rowdim () const
01207 {
01208 size_type nn=0;
01209 for (size_type i=0; i<n; i++)
01210 nn += rowdim(i);
01211 return nn;
01212 }
01213
01215 size_type coldim () const
01216 {
01217
01218
01219
01220 std::vector<size_type> coldims(M(),-1);
01221 for (ConstRowIterator row=begin(); row!=end(); ++row)
01222 for (ConstColIterator col=row->begin(); col!=row->end(); ++col)
01223
01224 if (coldims[col.index()]==-1)
01225 coldims[col.index()] = col->coldim();
01226
01227 size_type sum = 0;
01228 for (typename std::vector<size_type>::iterator it=coldims.begin(); it!=coldims.end(); ++it)
01229
01230 if ((*it)>=0)
01231 sum += *it;
01232
01233 return sum;
01234 }
01235
01236
01237
01239 bool exists (size_type i, size_type j) const
01240 {
01241 #ifdef DUNE_ISTL_WITH_CHECKING
01242 if (i<0 || i>=n) DUNE_THROW(ISTLError,"index out of range");
01243 if (j<0 || i>=m) DUNE_THROW(ISTLError,"index out of range");
01244 #endif
01245 if (r[i].size() && r[i].find(j)!=r[i].end())
01246 return true;
01247 else
01248 return false;
01249 }
01250
01251
01252 private:
01253
01254 BuildMode build_mode;
01255 BuildStage ready;
01256
01257
01258 size_type n;
01259 size_type m;
01260 size_type nnz;
01261
01262
01263
01264 row_type* r;
01265
01266
01267 B* a;
01268 size_type* j;
01269
01270 void setWindowPointers(ConstRowIterator row)
01271 {
01272 row_type current_row(a,j,0);
01273 for (size_type i=0; i<n; i++, ++row){
01274
01275 size_type s = row->getsize();
01276
01277 if (s>0){
01278
01279 r[i].set(s,current_row.getptr(), current_row.getindexptr());
01280
01281 current_row.setptr(current_row.getptr()+s);
01282 current_row.setindexptr(current_row.getindexptr()+s);
01283 } else{
01284
01285 r[i].set(0,0,0);
01286 }
01287 }
01288 }
01289
01291 void copyWindowStructure(const BCRSMatrix& Mat)
01292 {
01293 setWindowPointers(Mat.begin());
01294
01295
01296 for (size_type i=0; i<n; i++) r[i] = Mat.r[i];
01297
01298
01299 build_mode = row_wise;
01300 ready = built;
01301 }
01302
01308 void deallocate(bool deallocateRows=true)
01309 {
01310
01311 if (nnz>0)
01312 {
01313
01314 A::template free<size_type>(j);
01315 A::template free<B>(a);
01316 }
01317 else
01318 {
01319
01320 for (size_type i=0; i<n; i++)
01321 if (r[i].getsize()>0)
01322 {
01323 A::template free<size_type>(r[i].getindexptr());
01324 A::template free<B>(r[i].getptr());
01325 }
01326 }
01327
01328
01329 if (n>0 && deallocateRows) A::template free<row_type>(r);
01330
01331
01332 ready=notbuilt;
01333
01334 }
01335
01353 void allocate(size_type rows, size_type columns, size_type nnz_=0, bool allocateRows=true)
01354 {
01355
01356 n = rows;
01357 m = columns;
01358 nnz = nnz_;
01359
01360
01361 if(allocateRows){
01362 if (n>0){
01363 r = A::template malloc<row_type>(rows);
01364 }else{
01365 r = 0;
01366 }
01367 }
01368
01369
01370
01371 if (nnz>0){
01372 a = A::template malloc<B>(nnz);
01373 j = A::template malloc<size_type>(nnz);
01374 }else{
01375 a = 0;
01376 j = 0;
01377 }
01378
01379 ready = notbuilt;
01380 }
01381
01382 };
01383
01384
01385 #ifdef DUNE_EXPRESSIONTEMPLATES
01386 template <class B, class A>
01387 struct FieldType< BCRSMatrix<B,A> >
01388 {
01389 typedef typename FieldType<B>::type type;
01390 };
01391
01392 template <class B, class A>
01393 struct BlockType< BCRSMatrix<B,A> >
01394 {
01395 typedef B type;
01396 };
01397 template <class B, class A>
01398 struct RowType< BCRSMatrix<B,A> >
01399 {
01400 typedef CompressedBlockVectorWindow<B,A> type;
01401 };
01402 #endif
01403
01406 }
01407
01408 #endif