dune-istl
2.1.1
|
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*- 00002 // vi: set et ts=4 sw=2 sts=2: 00003 00004 #ifndef DUNE_BCRSMATRIX_HH 00005 #define DUNE_BCRSMATRIX_HH 00006 00007 #include<cmath> 00008 #include<complex> 00009 #include<set> 00010 #include<iostream> 00011 #include<algorithm> 00012 #include<numeric> 00013 #include<vector> 00014 00015 #include "istlexception.hh" 00016 #include "bvector.hh" 00017 #include <dune/common/shared_ptr.hh> 00018 #include <dune/common/stdstreams.hh> 00019 #include <dune/common/iteratorfacades.hh> 00020 #include <dune/common/typetraits.hh> 00021 #include <dune/common/static_assert.hh> 00022 00027 namespace Dune { 00028 00068 template<typename M> 00069 struct MatrixDimension; 00070 00178 template<class B, class A=std::allocator<B> > 00179 class BCRSMatrix 00180 { 00181 friend struct MatrixDimension<BCRSMatrix>; 00182 00183 private: 00184 enum BuildStage{ 00186 notbuilt=0, 00191 rowSizesBuilt=1, 00193 built=2 00194 }; 00195 00196 public: 00197 00198 //===== type definitions and constants 00199 00201 typedef typename B::field_type field_type; 00202 00204 typedef B block_type; 00205 00207 typedef A allocator_type; 00208 00210 typedef CompressedBlockVectorWindow<B,A> row_type; 00211 00213 typedef typename A::size_type size_type; 00214 00216 enum { 00218 blocklevel = B::blocklevel+1 00219 }; 00220 00222 enum BuildMode { 00233 row_wise, 00242 random, 00246 unknown 00247 }; 00248 00249 00250 //===== random access interface to rows of the matrix 00251 00253 row_type& operator[] (size_type i) 00254 { 00255 #ifdef DUNE_ISTL_WITH_CHECKING 00256 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet"); 00257 if (i>=n) DUNE_THROW(ISTLError,"index out of range"); 00258 if (r[i].getptr()==0) DUNE_THROW(ISTLError,"row not initialized yet"); 00259 #endif 00260 return r[i]; 00261 } 00262 00264 const row_type& operator[] (size_type i) const 00265 { 00266 #ifdef DUNE_ISTL_WITH_CHECKING 00267 if (built!=ready) DUNE_THROW(ISTLError,"row not initialized yet"); 00268 if (i>=n) DUNE_THROW(ISTLError,"index out of range"); 00269 #endif 00270 return r[i]; 00271 } 00272 00273 00274 //===== iterator interface to rows of the matrix 00275 00277 template<class T> 00278 class RealRowIterator 00279 : public RandomAccessIteratorFacade<RealRowIterator<T>, T> 00280 { 00281 00282 public: 00284 typedef typename remove_const<T>::type ValueType; 00285 00286 friend class RandomAccessIteratorFacade<RealRowIterator<const ValueType>, const ValueType>; 00287 friend class RandomAccessIteratorFacade<RealRowIterator<ValueType>, ValueType>; 00288 friend class RealRowIterator<const ValueType>; 00289 friend class RealRowIterator<ValueType>; 00290 00292 RealRowIterator (row_type* _p, size_type _i) 00293 : p(_p), i(_i) 00294 {} 00295 00297 RealRowIterator () 00298 : p(0), i(0) 00299 {} 00300 00301 RealRowIterator(const RealRowIterator<ValueType>& it) 00302 : p(it.p), i(it.i) 00303 {} 00304 00305 00307 size_type index () const 00308 { 00309 return i; 00310 } 00311 00312 std::ptrdiff_t distanceTo(const RealRowIterator<ValueType>& other) const 00313 { 00314 assert(other.p==p); 00315 return (other.i-i); 00316 } 00317 00318 std::ptrdiff_t distanceTo(const RealRowIterator<const ValueType>& other) const 00319 { 00320 assert(other.p==p); 00321 return (other.i-i); 00322 } 00323 00325 bool equals (const RealRowIterator<ValueType>& other) const 00326 { 00327 assert(other.p==p); 00328 return i==other.i; 00329 } 00330 00332 bool equals (const RealRowIterator<const ValueType>& other) const 00333 { 00334 assert(other.p==p); 00335 return i==other.i; 00336 } 00337 00338 private: 00340 void increment() 00341 { 00342 ++i; 00343 } 00344 00346 void decrement() 00347 { 00348 --i; 00349 } 00350 00351 void advance(std::ptrdiff_t diff) 00352 { 00353 i+=diff; 00354 } 00355 00356 T& elementAt(std::ptrdiff_t diff) const 00357 { 00358 return p[i+diff]; 00359 } 00360 00362 row_type& dereference () const 00363 { 00364 return p[i]; 00365 } 00366 00367 row_type* p; 00368 size_type i; 00369 }; 00370 00372 typedef RealRowIterator<row_type> iterator; 00373 typedef RealRowIterator<row_type> Iterator; 00374 00376 Iterator begin () 00377 { 00378 return Iterator(r,0); 00379 } 00380 00382 Iterator end () 00383 { 00384 return Iterator(r,n); 00385 } 00386 00391 Iterator rbegin() DUNE_DEPRECATED 00392 { 00393 return beforeEnd(); 00394 } 00395 00398 Iterator beforeEnd () 00399 { 00400 return Iterator(r,n-1); 00401 } 00402 00407 Iterator rend () DUNE_DEPRECATED 00408 { 00409 return beforeBegin(); 00410 } 00411 00414 Iterator beforeBegin () 00415 { 00416 return Iterator(r,-1); 00417 } 00418 00420 typedef Iterator RowIterator; 00421 00423 typedef typename row_type::Iterator ColIterator; 00424 00426 typedef RealRowIterator<const row_type> const_iterator; 00427 typedef RealRowIterator<const row_type> ConstIterator; 00428 00429 00431 ConstIterator begin () const 00432 { 00433 return ConstIterator(r,0); 00434 } 00435 00437 ConstIterator end () const 00438 { 00439 return ConstIterator(r,n); 00440 } 00441 00446 ConstIterator rbegin() const DUNE_DEPRECATED 00447 { 00448 return beforeEnd(); 00449 } 00450 00453 ConstIterator beforeEnd() const 00454 { 00455 return ConstIterator(r,n-1); 00456 } 00457 00462 ConstIterator rend () const DUNE_DEPRECATED 00463 { 00464 return beforeBegin(); 00465 } 00466 00469 ConstIterator beforeBegin () const 00470 { 00471 return ConstIterator(r,-1); 00472 } 00473 00475 typedef ConstIterator ConstRowIterator; 00476 00478 typedef typename row_type::ConstIterator ConstColIterator; 00479 00480 //===== constructors & resizers 00481 00483 BCRSMatrix () 00484 : build_mode(unknown), ready(notbuilt), n(0), m(0), nnz(0), 00485 r(0), a(0) 00486 {} 00487 00489 BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm) 00490 : build_mode(bm), ready(notbuilt) 00491 { 00492 allocate(_n, _m, _nnz); 00493 } 00494 00496 BCRSMatrix (size_type _n, size_type _m, BuildMode bm) 00497 : build_mode(bm), ready(notbuilt) 00498 { 00499 allocate(_n, _m); 00500 } 00501 00507 BCRSMatrix (const BCRSMatrix& Mat) 00508 : n(Mat.n), nnz(0) 00509 { 00510 // deep copy in global array 00511 size_type _nnz = Mat.nnz; 00512 00513 // in case of row-wise allocation 00514 if (_nnz<=0) 00515 { 00516 _nnz = 0; 00517 for (size_type i=0; i<n; i++) 00518 _nnz += Mat.r[i].getsize(); 00519 } 00520 00521 j = Mat.j; // enable column index sharing, release array in case of row-wise allocation 00522 allocate(Mat.n, Mat.m, _nnz); 00523 00524 // build window structure 00525 copyWindowStructure(Mat); 00526 } 00527 00529 ~BCRSMatrix () 00530 { 00531 deallocate(); 00532 } 00533 00538 void setBuildMode(BuildMode bm) 00539 { 00540 if(ready==notbuilt) 00541 build_mode = bm; 00542 else 00543 DUNE_THROW(InvalidStateException, "Matrix structure is already built (ready="<<ready<<")."); 00544 } 00545 00561 void setSize(size_type rows, size_type columns, size_type nnz=0) 00562 { 00563 // deallocate already setup memory 00564 deallocate(); 00565 00566 // allocate matrix memory 00567 allocate(rows, columns, nnz); 00568 } 00569 00576 BCRSMatrix& operator= (const BCRSMatrix& Mat) 00577 { 00578 // return immediately when self-assignment 00579 if (&Mat==this) return *this; 00580 00581 // make it simple: ALWAYS throw away memory for a and j 00582 deallocate(false); 00583 00584 // reallocate the rows if required 00585 if (n>0 && n!=Mat.n) { 00586 // free rows 00587 int i=n; 00588 while (i) 00589 r[--i].~row_type(); 00590 rowAllocator_.deallocate(r,n); 00591 } 00592 00593 nnz=Mat.nnz; 00594 if (nnz<=0) 00595 { 00596 for (size_type i=0; i<Mat.n; i++) 00597 nnz += Mat.r[i].getsize(); 00598 } 00599 00600 // allocate a, share j 00601 j = Mat.j; 00602 allocate(Mat.n, Mat.m, nnz, n!=Mat.n); 00603 00604 // build window structure 00605 copyWindowStructure(Mat); 00606 return *this; 00607 } 00608 00610 BCRSMatrix& operator= (const field_type& k) 00611 { 00612 for (size_type i=0; i<n; i++) r[i] = k; 00613 return *this; 00614 } 00615 00616 //===== row-wise creation interface 00617 00619 class CreateIterator 00620 { 00621 public: 00623 CreateIterator (BCRSMatrix& _Mat, size_type _i) 00624 : Mat(_Mat), i(_i), nnz(0), current_row(Mat.a, Mat.j.get(), 0) 00625 { 00626 if (i==0 && Mat.ready) 00627 DUNE_THROW(ISTLError,"creation only allowed for uninitialized matrix"); 00628 if(Mat.build_mode!=row_wise) 00629 { 00630 if(Mat.build_mode==unknown) 00631 Mat.build_mode=row_wise; 00632 else 00633 DUNE_THROW(ISTLError,"creation only allowed if row wise allocation was requested in the constructor"); 00634 } 00635 } 00636 00638 CreateIterator& operator++() 00639 { 00640 // this should only be called if matrix is in creation 00641 if (Mat.ready) 00642 DUNE_THROW(ISTLError,"matrix already built up"); 00643 00644 // row i is defined through the pattern 00645 // get memory for the row and initialize the j array 00646 // this depends on the allocation mode 00647 00648 // compute size of the row 00649 size_type s = pattern.size(); 00650 00651 if(s>0){ 00652 // update number of nonzeroes including this row 00653 nnz += s; 00654 00655 // alloc memory / set window 00656 if (Mat.nnz>0) 00657 { 00658 // memory is allocated in one long array 00659 00660 // check if that memory is sufficient 00661 if (nnz>Mat.nnz) 00662 DUNE_THROW(ISTLError,"allocated nnz too small"); 00663 00664 // set row i 00665 Mat.r[i].set(s,current_row.getptr(),current_row.getindexptr()); 00666 current_row.setptr(current_row.getptr()+s); 00667 current_row.setindexptr(current_row.getindexptr()+s); 00668 }else{ 00669 // memory is allocated individually per row 00670 // allocate and set row i 00671 B* a = Mat.allocator_.allocate(s); 00672 new (a) B[s]; 00673 size_type* j = Mat.sizeAllocator_.allocate(s); 00674 new (j) size_type[s]; 00675 Mat.r[i].set(s,a,j); 00676 } 00677 }else 00678 // setup empty row 00679 Mat.r[i].set(0,0,0); 00680 00681 // initialize the j array for row i from pattern 00682 size_type k=0; 00683 size_type *j = Mat.r[i].getindexptr(); 00684 for (typename PatternType::const_iterator it=pattern.begin(); it!=pattern.end(); ++it) 00685 j[k++] = *it; 00686 00687 // now go to next row 00688 i++; 00689 pattern.clear(); 00690 00691 // check if this was last row 00692 if (i==Mat.n) 00693 { 00694 Mat.ready = built; 00695 if(Mat.nnz>0) 00696 // Set nnz to the exact number of nonzero blocks inserted 00697 // as some methods rely on it 00698 Mat.nnz=nnz; 00699 } 00700 // done 00701 return *this; 00702 } 00703 00705 bool operator!= (const CreateIterator& it) const 00706 { 00707 return (i!=it.i) || (&Mat!=&it.Mat); 00708 } 00709 00711 bool operator== (const CreateIterator& it) const 00712 { 00713 return (i==it.i) && (&Mat==&it.Mat); 00714 } 00715 00717 size_type index () const 00718 { 00719 return i; 00720 } 00721 00723 void insert (size_type j) 00724 { 00725 pattern.insert(j); 00726 } 00727 00729 bool contains (size_type j) 00730 { 00731 if (pattern.find(j)!=pattern.end()) 00732 return true; 00733 else 00734 return false; 00735 } 00741 size_type size() const 00742 { 00743 return pattern.size(); 00744 } 00745 00746 private: 00747 BCRSMatrix& Mat; // the matrix we are defining 00748 size_type i; // current row to be defined 00749 size_type nnz; // count total number of nonzeros 00750 typedef std::set<size_type,std::less<size_type> > PatternType; 00751 PatternType pattern; // used to compile entries in a row 00752 row_type current_row; // row pointing to the current row to setup 00753 }; 00754 00756 friend class CreateIterator; 00757 00759 CreateIterator createbegin () 00760 { 00761 return CreateIterator(*this,0); 00762 } 00763 00765 CreateIterator createend () 00766 { 00767 return CreateIterator(*this,n); 00768 } 00769 00770 00771 //===== random creation interface 00772 00774 void setrowsize (size_type i, size_type s) 00775 { 00776 if (build_mode!=random) 00777 DUNE_THROW(ISTLError,"requires random build mode"); 00778 if (ready) 00779 DUNE_THROW(ISTLError,"matrix row sizes already built up"); 00780 00781 r[i].setsize(s); 00782 } 00783 00785 size_type getrowsize (size_type i) const 00786 { 00787 #ifdef DUNE_ISTL_WITH_CHECKING 00788 if (r==0) DUNE_THROW(ISTLError,"row not initialized yet"); 00789 if (i>=n) DUNE_THROW(ISTLError,"index out of range"); 00790 #endif 00791 return r[i].getsize(); 00792 } 00793 00795 void incrementrowsize (size_type i, size_type s = 1) 00796 { 00797 if (build_mode!=random) 00798 DUNE_THROW(ISTLError,"requires random build mode"); 00799 if (ready) 00800 DUNE_THROW(ISTLError,"matrix row sizes already built up"); 00801 00802 r[i].setsize(r[i].getsize()+s); 00803 } 00804 00806 void endrowsizes () 00807 { 00808 if (build_mode!=random) 00809 DUNE_THROW(ISTLError,"requires random build mode"); 00810 if (ready) 00811 DUNE_THROW(ISTLError,"matrix row sizes already built up"); 00812 00813 // compute total size, check positivity 00814 size_type total=0; 00815 for (size_type i=0; i<n; i++) 00816 { 00817 if (r[i].getsize()<0) 00818 DUNE_THROW(ISTLError,"rowsize must be nonnegative"); 00819 total += r[i].getsize(); 00820 } 00821 00822 if(nnz==0) 00823 // allocate/check memory 00824 allocate(n,m,total,false); 00825 else if(nnz<total) 00826 DUNE_THROW(ISTLError,"Specified number of nonzeros ("<<nnz<<") not " 00827 <<"sufficient for calculated nonzeros ("<<total<<"! "); 00828 00829 // set the window pointers correctly 00830 setWindowPointers(begin()); 00831 00832 // initialize j array with m (an invalid column index) 00833 // this indicates an unused entry 00834 for (size_type k=0; k<nnz; k++) 00835 j.get()[k] = m; 00836 ready = rowSizesBuilt; 00837 } 00838 00840 00850 void addindex (size_type row, size_type col) 00851 { 00852 if (build_mode!=random) 00853 DUNE_THROW(ISTLError,"requires random build mode"); 00854 if (ready==built) 00855 DUNE_THROW(ISTLError,"matrix already built up"); 00856 if (ready==notbuilt) 00857 DUNE_THROW(ISTLError,"matrix row sizes not built up yet"); 00858 00859 if (col >= m) 00860 DUNE_THROW(ISTLError,"column index exceeds matrix size"); 00861 00862 // get row range 00863 size_type* const first = r[row].getindexptr(); 00864 size_type* const last = first + r[row].getsize(); 00865 00866 // find correct insertion position for new column index 00867 size_type* pos = std::lower_bound(first,last,col); 00868 00869 // check if index is already in row 00870 if (pos!=last && *pos == col) return; 00871 00872 // find end of already inserted column indices 00873 size_type* end = std::lower_bound(pos,last,m); 00874 if (end==last) 00875 DUNE_THROW(ISTLError,"row is too small"); 00876 00877 // insert new column index at correct position 00878 std::copy_backward(pos,end,end+1); 00879 *pos = col; 00880 00881 } 00882 00884 void endindices () 00885 { 00886 if (build_mode!=random) 00887 DUNE_THROW(ISTLError,"requires random build mode"); 00888 if (ready==built) 00889 DUNE_THROW(ISTLError,"matrix already built up"); 00890 if (ready==notbuilt) 00891 DUNE_THROW(ISTLError,"row sizes are not built up yet"); 00892 00893 // check if there are undefined indices 00894 RowIterator endi=end(); 00895 for (RowIterator i=begin(); i!=endi; ++i) 00896 { 00897 ColIterator endj = (*i).end(); 00898 for (ColIterator j=(*i).begin(); j!=endj; ++j){ 00899 if (j.index()<0) 00900 { 00901 std::cout << "j[" << j.offset() << "]=" << j.index() << std::endl; 00902 DUNE_THROW(ISTLError,"undefined index detected"); 00903 } 00904 if (j.index()>=m){ 00905 dwarn << "WARNING: size of row "<< i.index()<<" is "<<j.offset()<<". But was specified as being "<< (*i).end().offset() 00906 <<". This means you are wasting valuable space and creating additional cache misses!"<<std::endl; 00907 r[i.index()].setsize(j.offset()); 00908 break; 00909 } 00910 } 00911 } 00912 00913 // if not, set matrix to built 00914 ready = built; 00915 } 00916 00917 //===== vector space arithmetic 00918 00920 BCRSMatrix& operator*= (const field_type& k) 00921 { 00922 if (nnz>0) 00923 { 00924 // process 1D array 00925 for (size_type i=0; i<nnz; i++) 00926 a[i] *= k; 00927 } 00928 else 00929 { 00930 RowIterator endi=end(); 00931 for (RowIterator i=begin(); i!=endi; ++i) 00932 { 00933 ColIterator endj = (*i).end(); 00934 for (ColIterator j=(*i).begin(); j!=endj; ++j) 00935 (*j) *= k; 00936 } 00937 } 00938 00939 return *this; 00940 } 00941 00943 BCRSMatrix& operator/= (const field_type& k) 00944 { 00945 if (nnz>0) 00946 { 00947 // process 1D array 00948 for (size_type i=0; i<nnz; i++) 00949 a[i] /= k; 00950 } 00951 else 00952 { 00953 RowIterator endi=end(); 00954 for (RowIterator i=begin(); i!=endi; ++i) 00955 { 00956 ColIterator endj = (*i).end(); 00957 for (ColIterator j=(*i).begin(); j!=endj; ++j) 00958 (*j) /= k; 00959 } 00960 } 00961 00962 return *this; 00963 } 00964 00965 00971 BCRSMatrix& operator+= (const BCRSMatrix& b) 00972 { 00973 #ifdef DUNE_ISTL_WITH_CHECKING 00974 if(N()!=b.N() || M() != b.M()) 00975 DUNE_THROW(RangeError, "Matrix sizes do not match!"); 00976 #endif 00977 RowIterator endi=end(); 00978 ConstRowIterator j=b.begin(); 00979 for (RowIterator i=begin(); i!=endi; ++i, ++j){ 00980 i->operator+=(*j); 00981 } 00982 00983 return *this; 00984 } 00985 00991 BCRSMatrix& operator-= (const BCRSMatrix& b) 00992 { 00993 #ifdef DUNE_ISTL_WITH_CHECKING 00994 if(N()!=b.N() || M() != b.M()) 00995 DUNE_THROW(RangeError, "Matrix sizes do not match!"); 00996 #endif 00997 RowIterator endi=end(); 00998 ConstRowIterator j=b.begin(); 00999 for (RowIterator i=begin(); i!=endi; ++i, ++j){ 01000 i->operator-=(*j); 01001 } 01002 01003 return *this; 01004 } 01005 01014 BCRSMatrix& axpy(field_type alpha, const BCRSMatrix& b) 01015 { 01016 #ifdef DUNE_ISTL_WITH_CHECKING 01017 if(N()!=b.N() || M() != b.M()) 01018 DUNE_THROW(RangeError, "Matrix sizes do not match!"); 01019 #endif 01020 RowIterator endi=end(); 01021 ConstRowIterator j=b.begin(); 01022 for(RowIterator i=begin(); i!=endi; ++i, ++j) 01023 i->axpy(alpha, *j); 01024 01025 return *this; 01026 } 01027 01028 //===== linear maps 01029 01031 template<class X, class Y> 01032 void mv (const X& x, Y& y) const 01033 { 01034 #ifdef DUNE_ISTL_WITH_CHECKING 01035 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01036 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01037 #endif 01038 ConstRowIterator endi=end(); 01039 for (ConstRowIterator i=begin(); i!=endi; ++i) 01040 { 01041 y[i.index()]=0; 01042 ConstColIterator endj = (*i).end(); 01043 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01044 (*j).umv(x[j.index()],y[i.index()]); 01045 } 01046 } 01047 01049 template<class X, class Y> 01050 void umv (const X& x, Y& y) const 01051 { 01052 #ifdef DUNE_ISTL_WITH_CHECKING 01053 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01054 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01055 #endif 01056 ConstRowIterator endi=end(); 01057 for (ConstRowIterator i=begin(); i!=endi; ++i) 01058 { 01059 ConstColIterator endj = (*i).end(); 01060 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01061 (*j).umv(x[j.index()],y[i.index()]); 01062 } 01063 } 01064 01066 template<class X, class Y> 01067 void mmv (const X& x, Y& y) const 01068 { 01069 #ifdef DUNE_ISTL_WITH_CHECKING 01070 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01071 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01072 #endif 01073 ConstRowIterator endi=end(); 01074 for (ConstRowIterator i=begin(); i!=endi; ++i) 01075 { 01076 ConstColIterator endj = (*i).end(); 01077 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01078 (*j).mmv(x[j.index()],y[i.index()]); 01079 } 01080 } 01081 01083 template<class X, class Y> 01084 void usmv (const field_type& alpha, const X& x, Y& y) const 01085 { 01086 #ifdef DUNE_ISTL_WITH_CHECKING 01087 if (x.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01088 if (y.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01089 #endif 01090 ConstRowIterator endi=end(); 01091 for (ConstRowIterator i=begin(); i!=endi; ++i) 01092 { 01093 ConstColIterator endj = (*i).end(); 01094 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01095 (*j).usmv(alpha,x[j.index()],y[i.index()]); 01096 } 01097 } 01098 01100 template<class X, class Y> 01101 void mtv (const X& x, Y& y) const 01102 { 01103 #ifdef DUNE_ISTL_WITH_CHECKING 01104 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01105 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01106 #endif 01107 for(size_type i=0; i<y.N(); ++i) 01108 y[i]=0; 01109 umtv(x,y); 01110 } 01111 01113 template<class X, class Y> 01114 void umtv (const X& x, Y& y) const 01115 { 01116 #ifdef DUNE_ISTL_WITH_CHECKING 01117 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01118 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01119 #endif 01120 ConstRowIterator endi=end(); 01121 for (ConstRowIterator i=begin(); i!=endi; ++i) 01122 { 01123 ConstColIterator endj = (*i).end(); 01124 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01125 (*j).umtv(x[i.index()],y[j.index()]); 01126 } 01127 } 01128 01130 template<class X, class Y> 01131 void mmtv (const X& x, Y& y) const 01132 { 01133 #ifdef DUNE_ISTL_WITH_CHECKING 01134 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01135 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01136 #endif 01137 ConstRowIterator endi=end(); 01138 for (ConstRowIterator i=begin(); i!=endi; ++i) 01139 { 01140 ConstColIterator endj = (*i).end(); 01141 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01142 (*j).mmtv(x[i.index()],y[j.index()]); 01143 } 01144 } 01145 01147 template<class X, class Y> 01148 void usmtv (const field_type& alpha, const X& x, Y& y) const 01149 { 01150 #ifdef DUNE_ISTL_WITH_CHECKING 01151 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01152 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01153 #endif 01154 ConstRowIterator endi=end(); 01155 for (ConstRowIterator i=begin(); i!=endi; ++i) 01156 { 01157 ConstColIterator endj = (*i).end(); 01158 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01159 (*j).usmtv(alpha,x[i.index()],y[j.index()]); 01160 } 01161 } 01162 01164 template<class X, class Y> 01165 void umhv (const X& x, Y& y) const 01166 { 01167 #ifdef DUNE_ISTL_WITH_CHECKING 01168 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01169 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01170 #endif 01171 ConstRowIterator endi=end(); 01172 for (ConstRowIterator i=begin(); i!=endi; ++i) 01173 { 01174 ConstColIterator endj = (*i).end(); 01175 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01176 (*j).umhv(x[i.index()],y[j.index()]); 01177 } 01178 } 01179 01181 template<class X, class Y> 01182 void mmhv (const X& x, Y& y) const 01183 { 01184 #ifdef DUNE_ISTL_WITH_CHECKING 01185 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01186 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01187 #endif 01188 ConstRowIterator endi=end(); 01189 for (ConstRowIterator i=begin(); i!=endi; ++i) 01190 { 01191 ConstColIterator endj = (*i).end(); 01192 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01193 (*j).mmhv(x[i.index()],y[j.index()]); 01194 } 01195 } 01196 01198 template<class X, class Y> 01199 void usmhv (const field_type& alpha, const X& x, Y& y) const 01200 { 01201 #ifdef DUNE_ISTL_WITH_CHECKING 01202 if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range"); 01203 if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range"); 01204 #endif 01205 ConstRowIterator endi=end(); 01206 for (ConstRowIterator i=begin(); i!=endi; ++i) 01207 { 01208 ConstColIterator endj = (*i).end(); 01209 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01210 (*j).usmhv(alpha,x[i.index()],y[j.index()]); 01211 } 01212 } 01213 01214 01215 //===== norms 01216 01218 double frobenius_norm2 () const 01219 { 01220 double sum=0; 01221 01222 ConstRowIterator endi=end(); 01223 for (ConstRowIterator i=begin(); i!=endi; ++i) 01224 { 01225 ConstColIterator endj = (*i).end(); 01226 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01227 sum += (*j).frobenius_norm2(); 01228 } 01229 01230 return sum; 01231 } 01232 01234 double frobenius_norm () const 01235 { 01236 return sqrt(frobenius_norm2()); 01237 } 01238 01240 double infinity_norm () const 01241 { 01242 double max=0; 01243 ConstRowIterator endi=end(); 01244 for (ConstRowIterator i=begin(); i!=endi; ++i) 01245 { 01246 double sum=0; 01247 ConstColIterator endj = (*i).end(); 01248 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01249 sum += (*j).infinity_norm(); 01250 max = std::max(max,sum); 01251 } 01252 return max; 01253 } 01254 01256 double infinity_norm_real () const 01257 { 01258 double max=0; 01259 ConstRowIterator endi=end(); 01260 for (ConstRowIterator i=begin(); i!=endi; ++i) 01261 { 01262 double sum=0; 01263 ConstColIterator endj = (*i).end(); 01264 for (ConstColIterator j=(*i).begin(); j!=endj; ++j) 01265 sum += (*j).infinity_norm_real(); 01266 max = std::max(max,sum); 01267 } 01268 return max; 01269 } 01270 01271 01272 //===== sizes 01273 01275 size_type N () const 01276 { 01277 return n; 01278 } 01279 01281 size_type M () const 01282 { 01283 return m; 01284 } 01285 01287 size_type nonzeroes () const 01288 { 01289 return nnz; 01290 } 01291 01292 01293 //===== query 01294 01296 bool exists (size_type i, size_type j) const 01297 { 01298 #ifdef DUNE_ISTL_WITH_CHECKING 01299 if (i<0 || i>=n) DUNE_THROW(ISTLError,"index out of range"); 01300 if (j<0 || i>=m) DUNE_THROW(ISTLError,"index out of range"); 01301 #endif 01302 if (r[i].size() && r[i].find(j)!=r[i].end()) 01303 return true; 01304 else 01305 return false; 01306 } 01307 01308 01309 private: 01310 // state information 01311 BuildMode build_mode; // row wise or whole matrix 01312 BuildStage ready; // indicate the stage the matrix building is in 01313 01314 // The allocator used for memory management 01315 typename A::template rebind<B>::other allocator_; 01316 01317 typename A::template rebind<row_type>::other rowAllocator_; 01318 01319 typename A::template rebind<size_type>::other sizeAllocator_; 01320 01321 // size of the matrix 01322 size_type n; // number of rows 01323 size_type m; // number of columns 01324 size_type nnz; // number of nonzeros allocated in the a and j array below 01325 // zero means that memory is allocated separately for each row. 01326 01327 // the rows are dynamically allocated 01328 row_type* r; // [n] the individual rows having pointers into a,j arrays 01329 01330 // dynamically allocated memory 01331 B* a; // [nnz] non-zero entries of the matrix in row-wise ordering 01332 // If a single array of column indices is used, it can be shared 01333 // between different matrices with the same sparsity pattern 01334 Dune::shared_ptr<size_type> j; // [nnz] column indices of entries 01335 01336 01337 void setWindowPointers(ConstRowIterator row) 01338 { 01339 row_type current_row(a,j.get(),0); // Pointers to current row data 01340 for (size_type i=0; i<n; i++, ++row){ 01341 // set row i 01342 size_type s = row->getsize(); 01343 01344 if (s>0){ 01345 // setup pointers and size 01346 r[i].set(s,current_row.getptr(), current_row.getindexptr()); 01347 // update pointer for next row 01348 current_row.setptr(current_row.getptr()+s); 01349 current_row.setindexptr(current_row.getindexptr()+s); 01350 } else{ 01351 // empty row 01352 r[i].set(0,0,0); 01353 } 01354 } 01355 } 01356 01358 void copyWindowStructure(const BCRSMatrix& Mat) 01359 { 01360 setWindowPointers(Mat.begin()); 01361 01362 // copy data 01363 for (size_type i=0; i<n; i++) r[i] = Mat.r[i]; 01364 01365 // finish off 01366 build_mode = row_wise; // dummy 01367 ready = built; 01368 } 01369 01375 void deallocate(bool deallocateRows=true) 01376 { 01377 01378 if (nnz>0) 01379 { 01380 // a,j have been allocated as one long vector 01381 j.reset(); 01382 int i=nnz; 01383 while (i) 01384 a[--i].~B(); 01385 allocator_.deallocate(a,n); 01386 } 01387 else 01388 { 01389 // check if memory for rows have been allocated individually 01390 for (size_type i=0; i<n; i++) 01391 if (r[i].getsize()>0) 01392 { 01393 int j=r[i].getsize(); 01394 while (j) { 01395 r[i].getindexptr()[--j].~size_type(); 01396 r[i].getptr()[j].~B(); 01397 } 01398 sizeAllocator_.deallocate(r[i].getindexptr(),1); 01399 allocator_.deallocate(r[i].getptr(),1); 01400 } 01401 } 01402 01403 // deallocate the rows 01404 if (n>0 && deallocateRows) { 01405 int i=n; 01406 while (i) 01407 r[--i].~row_type(); 01408 rowAllocator_.deallocate(r,n); 01409 } 01410 01411 // Mark matrix as not built at all. 01412 ready=notbuilt; 01413 01414 } 01415 01417 class Deallocator 01418 { 01419 typename A::template rebind<size_type>::other& sizeAllocator_; 01420 01421 public: 01422 Deallocator(typename A::template rebind<size_type>::other& sizeAllocator) 01423 : sizeAllocator_(sizeAllocator) 01424 {} 01425 01426 void operator()(size_type* p) { sizeAllocator_.deallocate(p,1); } 01427 }; 01428 01429 01447 void allocate(size_type rows, size_type columns, size_type nnz_=0, bool allocateRows=true) 01448 { 01449 // Store size 01450 n = rows; 01451 m = columns; 01452 nnz = nnz_; 01453 01454 // allocate rows 01455 if(allocateRows){ 01456 if (n>0){ 01457 r = rowAllocator_.allocate(rows); 01458 new (r) row_type[rows]; 01459 }else{ 01460 r = 0; 01461 } 01462 } 01463 01464 01465 // allocate a and j array 01466 if (nnz>0){ 01467 a = allocator_.allocate(nnz); 01468 // allocate column indices only if not yet present (enable sharing) 01469 if (!j.get()) 01470 j.reset(sizeAllocator_.allocate(nnz),Deallocator(sizeAllocator_)); 01471 }else{ 01472 a = 0; 01473 j.reset(); 01474 } 01475 // Mark the matrix as not built. 01476 ready = notbuilt; 01477 } 01478 01479 }; 01480 01481 01484 } // end namespace 01485 01486 #endif