4 #ifndef DUNE_BCRSMATRIX_HH
5 #define DUNE_BCRSMATRIX_HH
17 #include <dune/common/shared_ptr.hh>
18 #include <dune/common/stdstreams.hh>
19 #include <dune/common/iteratorfacades.hh>
20 #include <dune/common/typetraits.hh>
21 #include <dune/common/static_assert.hh>
69 struct MatrixDimension;
178 template<
class B,
class A=std::allocator<B> >
255 #ifdef DUNE_ISTL_WITH_CHECKING
256 if (r==0) DUNE_THROW(
ISTLError,
"row not initialized yet");
257 if (i>=n) DUNE_THROW(
ISTLError,
"index out of range");
258 if (r[i].getptr()==0) DUNE_THROW(
ISTLError,
"row not initialized yet");
266 #ifdef DUNE_ISTL_WITH_CHECKING
267 if (built!=ready) DUNE_THROW(
ISTLError,
"row not initialized yet");
268 if (i>=n) DUNE_THROW(
ISTLError,
"index out of range");
279 :
public RandomAccessIteratorFacade<RealRowIterator<T>, T>
351 void advance(std::ptrdiff_t diff)
356 T& elementAt(std::ptrdiff_t diff)
const
448 : build_mode(
unknown), ready(notbuilt), n(0), m(0), nnz(0),
454 : build_mode(bm), ready(notbuilt)
456 allocate(_n, _m, _nnz);
461 : build_mode(bm), ready(notbuilt)
486 allocate(Mat.n, Mat.m, _nnz);
489 copyWindowStructure(Mat);
507 DUNE_THROW(InvalidStateException,
"Matrix structure is already built (ready="<<ready<<
").");
531 allocate(rows, columns, nnz);
543 if (&Mat==
this)
return *
this;
549 if (n>0 && n!=Mat.n) {
551 for(
row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
552 rowAllocator_.destroy(riter);
553 rowAllocator_.deallocate(r,n);
565 allocate(Mat.n, Mat.m, nnz, n!=Mat.n);
568 copyWindowStructure(Mat);
587 : Mat(_Mat), i(_i), nnz(0), current_row(Mat.a, Mat.j.
get(), 0)
589 if (i==0 && Mat.ready)
590 DUNE_THROW(
ISTLError,
"creation only allowed for uninitialized matrix");
596 DUNE_THROW(
ISTLError,
"creation only allowed if row wise allocation was requested in the constructor");
605 DUNE_THROW(
ISTLError,
"matrix already built up");
625 DUNE_THROW(
ISTLError,
"allocated nnz too small");
634 B* a = Mat.allocator_.allocate(s);
636 size_type* j = Mat.sizeAllocator_.allocate(s);
647 for (
typename PatternType::const_iterator it=pattern.begin(); it!=pattern.end(); ++it)
670 return (i!=it.i) || (&Mat!=&it.Mat);
676 return (i==it.i) && (&Mat==&it.Mat);
694 if (pattern.find(j)!=pattern.end())
706 return pattern.size();
713 typedef std::set<size_type,std::less<size_type> > PatternType;
740 DUNE_THROW(
ISTLError,
"requires random build mode");
742 DUNE_THROW(
ISTLError,
"matrix row sizes already built up");
750 #ifdef DUNE_ISTL_WITH_CHECKING
751 if (r==0) DUNE_THROW(
ISTLError,
"row not initialized yet");
752 if (i>=n) DUNE_THROW(
ISTLError,
"index out of range");
761 DUNE_THROW(
ISTLError,
"requires random build mode");
763 DUNE_THROW(
ISTLError,
"matrix row sizes already built up");
765 r[i].
setsize(r[i].getsize()+s);
772 DUNE_THROW(
ISTLError,
"requires random build mode");
774 DUNE_THROW(
ISTLError,
"matrix row sizes already built up");
780 if (r[i].getsize()<0)
781 DUNE_THROW(
ISTLError,
"rowsize must be nonnegative");
787 allocate(n,m,total,
false);
789 DUNE_THROW(
ISTLError,
"Specified number of nonzeros ("<<nnz<<
") not "
790 <<
"sufficient for calculated nonzeros ("<<total<<
"! ");
793 setWindowPointers(
begin());
799 ready = rowSizesBuilt;
816 DUNE_THROW(
ISTLError,
"requires random build mode");
818 DUNE_THROW(
ISTLError,
"matrix already built up");
820 DUNE_THROW(
ISTLError,
"matrix row sizes not built up yet");
823 DUNE_THROW(
ISTLError,
"column index exceeds matrix size");
830 size_type* pos = std::lower_bound(first,last,col);
833 if (pos!=last && *pos == col)
return;
838 DUNE_THROW(
ISTLError,
"row is too small");
841 std::copy_backward(pos,end,end+1);
850 DUNE_THROW(
ISTLError,
"requires random build mode");
852 DUNE_THROW(
ISTLError,
"matrix already built up");
854 DUNE_THROW(
ISTLError,
"row sizes are not built up yet");
864 std::cout <<
"j[" << j.offset() <<
"]=" << j.index() << std::endl;
865 DUNE_THROW(
ISTLError,
"undefined index detected");
868 dwarn <<
"WARNING: size of row "<< i.index()<<
" is "<<j.offset()<<
". But was specified as being "<< (*i).end().offset()
869 <<
". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
870 r[i.index()].
setsize(j.offset());
936 #ifdef DUNE_ISTL_WITH_CHECKING
937 if(
N()!=b.
N() ||
M() != b.
M())
938 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
956 #ifdef DUNE_ISTL_WITH_CHECKING
957 if(
N()!=b.
N() ||
M() != b.
M())
958 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
979 #ifdef DUNE_ISTL_WITH_CHECKING
980 if(
N()!=b.
N() ||
M() != b.
M())
981 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
994 template<
class X,
class Y>
995 void mv (
const X& x, Y& y)
const
997 #ifdef DUNE_ISTL_WITH_CHECKING
999 "Size mismatch: M: " <<
N() <<
"x" <<
M() <<
" x: " << x.N());
1001 "Size mismatch: M: " <<
N() <<
"x" <<
M() <<
" y: " << y.N());
1009 (*j).umv(x[j.index()],y[i.index()]);
1014 template<
class X,
class Y>
1015 void umv (
const X& x, Y& y)
const
1017 #ifdef DUNE_ISTL_WITH_CHECKING
1018 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1019 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1026 (*j).umv(x[j.index()],y[i.index()]);
1031 template<
class X,
class Y>
1032 void mmv (
const X& x, Y& y)
const
1034 #ifdef DUNE_ISTL_WITH_CHECKING
1035 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1036 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1043 (*j).mmv(x[j.index()],y[i.index()]);
1048 template<
class X,
class Y>
1051 #ifdef DUNE_ISTL_WITH_CHECKING
1052 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1053 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1060 (*j).usmv(alpha,x[j.index()],y[i.index()]);
1065 template<
class X,
class Y>
1066 void mtv (
const X& x, Y& y)
const
1068 #ifdef DUNE_ISTL_WITH_CHECKING
1069 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1070 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1078 template<
class X,
class Y>
1081 #ifdef DUNE_ISTL_WITH_CHECKING
1082 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1083 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1090 (*j).umtv(x[i.index()],y[j.index()]);
1095 template<
class X,
class Y>
1098 #ifdef DUNE_ISTL_WITH_CHECKING
1099 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1100 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1107 (*j).mmtv(x[i.index()],y[j.index()]);
1112 template<
class X,
class Y>
1115 #ifdef DUNE_ISTL_WITH_CHECKING
1116 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1117 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1124 (*j).usmtv(alpha,x[i.index()],y[j.index()]);
1129 template<
class X,
class Y>
1132 #ifdef DUNE_ISTL_WITH_CHECKING
1133 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1134 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1141 (*j).umhv(x[i.index()],y[j.index()]);
1146 template<
class X,
class Y>
1149 #ifdef DUNE_ISTL_WITH_CHECKING
1150 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1151 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1158 (*j).mmhv(x[i.index()],y[j.index()]);
1163 template<
class X,
class Y>
1166 #ifdef DUNE_ISTL_WITH_CHECKING
1167 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
1168 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
1175 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
1192 sum += (*j).frobenius_norm2();
1214 sum += (*j).infinity_norm();
1215 max = std::max(max,sum);
1230 sum += (*j).infinity_norm_real();
1231 max = std::max(max,sum);
1263 #ifdef DUNE_ISTL_WITH_CHECKING
1264 if (i<0 || i>=n) DUNE_THROW(
ISTLError,
"row index out of range");
1265 if (j<0 || j>=m) DUNE_THROW(
ISTLError,
"column index out of range");
1267 if (r[i].size() && r[i].find(j)!=r[i].
end())
1280 typename A::template rebind<B>::other allocator_;
1282 typename A::template rebind<row_type>::other rowAllocator_;
1284 typename A::template rebind<size_type>::other sizeAllocator_;
1299 Dune::shared_ptr<size_type> j;
1311 r[i].
set(s,current_row.getptr(), current_row.getindexptr());
1313 current_row.setptr(current_row.getptr()+s);
1314 current_row.setindexptr(current_row.getindexptr()+s);
1323 void copyWindowStructure(
const BCRSMatrix& Mat)
1325 setWindowPointers(Mat.begin());
1328 for (
size_type i=0; i<n; i++) r[i] = Mat.r[i];
1340 void deallocate(
bool deallocateRows=
true)
1347 for(B *aiter=a+(nnz-1), *aend=a-1; aiter!=aend; --aiter)
1348 allocator_.destroy(aiter);
1349 allocator_.deallocate(a,n);
1355 if (r[i].getsize()>0)
1357 for (B *
col=r[i].getptr()+(r[i].getsize()-1),
1358 *colend = r[i].getptr()-1;
col!=colend; --
col) {
1359 allocator_.destroy(
col);
1361 sizeAllocator_.deallocate(r[i].getindexptr(),1);
1362 allocator_.deallocate(r[i].getptr(),1);
1367 if (n>0 && deallocateRows) {
1368 for(
row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
1369 rowAllocator_.destroy(riter);
1370 rowAllocator_.deallocate(r,n);
1381 typename A::template rebind<size_type>::other& sizeAllocator_;
1384 Deallocator(
typename A::template rebind<size_type>::other& sizeAllocator)
1385 : sizeAllocator_(sizeAllocator)
1388 void operator()(
size_type* p) { sizeAllocator_.deallocate(p,1); }
1419 r = rowAllocator_.allocate(rows);
1429 a = allocator_.allocate(nnz);
1432 j.reset(sizeAllocator_.allocate(nnz),Deallocator(sizeAllocator_));
1436 for(
row_type* ri=r; ri!=r+rows;++ri)
1437 rowAllocator_.construct(ri,
row_type());