6#ifndef DUNE_ISTL_BCRSMATRIX_HH
7#define DUNE_ISTL_BCRSMATRIX_HH
22#include <dune/common/stdstreams.hh>
23#include <dune/common/iteratorfacades.hh>
24#include <dune/common/typetraits.hh>
25#include <dune/common/ftraits.hh>
26#include <dune/common/scalarvectorview.hh>
27#include <dune/common/scalarmatrixview.hh>
86 template<
typename size_type>
144 return _m.entry(_i,j);
173 if (m.buildMode() != Matrix::implicit)
174 DUNE_THROW(
BCRSMatrixError,
"You can only create an ImplicitBuilder for a matrix in implicit build mode");
175 if (m.buildStage() != Matrix::building)
176 DUNE_THROW(
BCRSMatrixError,
"You can only create an ImplicitBuilder for a matrix with set size that has not been compressed() yet");
197 if (m.buildStage() != Matrix::notAllocated)
198 DUNE_THROW(
BCRSMatrixError,
"You can only set up a matrix for this ImplicitBuilder if it has no memory allocated yet");
199 m.setBuildMode(Matrix::implicit);
200 m.setImplicitBuildModeParameters(avg_cols_per_row,overflow_fraction);
201 m.setSize(rows,cols);
464 template<
class B,
class A=std::allocator<B> >
467 friend struct MatrixDimension<BCRSMatrix>;
500 typedef Imp::CompressedBlockVectorWindow<B,size_type>
row_type;
547#ifdef DUNE_ISTL_WITH_CHECKING
549 DUNE_THROW(
BCRSMatrixError,
"You cannot use operator[] in implicit build mode before calling compress()");
559#ifdef DUNE_ISTL_WITH_CHECKING
561 DUNE_THROW(BCRSMatrixError,
"You cannot use operator[] in implicit build mode before calling compress()");
562 if (
built!=ready) DUNE_THROW(BCRSMatrixError,
"row not initialized yet");
563 if (i>=n) DUNE_THROW(BCRSMatrixError,
"index out of range");
574 :
public RandomAccessIteratorFacade<RealRowIterator<T>, T>
646 void advance(std::ptrdiff_t diff)
651 T& elementAt(std::ptrdiff_t diff)
const
667 typedef RealRowIterator<row_type> iterator;
668 typedef RealRowIterator<row_type> Iterator;
673 return Iterator(r,0);
679 return Iterator(r,n);
684 Iterator beforeEnd ()
686 return Iterator(r,n-1);
691 Iterator beforeBegin ()
693 return Iterator(r,-1);
697 typedef Iterator RowIterator;
700 typedef typename row_type::Iterator ColIterator;
703 typedef RealRowIterator<const row_type> const_iterator;
704 typedef RealRowIterator<const row_type> ConstIterator;
708 ConstIterator begin ()
const
710 return ConstIterator(r,0);
714 ConstIterator end ()
const
716 return ConstIterator(r,n);
721 ConstIterator beforeEnd()
const
723 return ConstIterator(r,n-1);
728 ConstIterator beforeBegin ()
const
730 return ConstIterator(r,-1);
734 typedef ConstIterator ConstRowIterator;
737 typedef typename row_type::ConstIterator ConstColIterator;
747 allocationSize_(0), r(0), a(0),
748 avg(0), compressionBufferSize_(-1.0)
753 : build_mode(bm), ready(
notAllocated), n(0), m(0), nnz_(0),
754 allocationSize_(0), r(0), a(0),
755 avg(0), compressionBufferSize_(-1.0)
757 allocate(_n, _m, _nnz,
true,
false);
762 : build_mode(bm), ready(
notAllocated), n(0), m(0), nnz_(0),
763 allocationSize_(0), r(0), a(0),
764 avg(0), compressionBufferSize_(-1.0)
766 allocate(_n, _m,0,
true,
false);
782 : build_mode(bm), ready(
notAllocated), n(0), m(0), nnz_(0),
783 allocationSize_(0), r(0), a(0),
784 avg(_avg), compressionBufferSize_(compressionBufferSize)
787 DUNE_THROW(BCRSMatrixError,
"Only call this constructor when using the implicit build mode");
792 if (compressionBufferSize_ < 0.0)
793 DUNE_THROW(BCRSMatrixError,
"You cannot set a negative overflow fraction");
794 implicit_allocate(_n,_m);
802 BCRSMatrix (
const BCRSMatrix& Mat)
803 : build_mode(Mat.build_mode), ready(
notAllocated), n(0), m(0), nnz_(0),
804 allocationSize_(0), r(0), a(0),
805 avg(Mat.avg), compressionBufferSize_(Mat.compressionBufferSize_)
808 DUNE_THROW(InvalidStateException,
"BCRSMatrix can only be copy-constructed when source matrix is completely empty (size not set) or fully built)");
814 allocate(Mat.n, Mat.m, _nnz,
true,
true);
817 copyWindowStructure(Mat);
840 DUNE_THROW(InvalidStateException,
"Matrix structure cannot be changed at this stage anymore (ready == "<<ready<<
").");
866 DUNE_THROW(Dune::BCRSMatrixError,
"number of non-zeroes may not be set in implicit mode, use setImplicitBuildModeParameters() instead");
869 implicit_allocate(rows,columns);
874 allocate(rows, columns, nnz,
true,
false);
886 void setImplicitBuildModeParameters(
size_type _avg,
double compressionBufferSize)
892 if (compressionBufferSize < 0.0)
893 DUNE_THROW(BCRSMatrixError,
"You cannot set a negative compressionBufferSize value");
897 DUNE_THROW(InvalidStateException,
"You cannot modify build mode parameters at this stage anymore");
899 compressionBufferSize_ = compressionBufferSize;
908 BCRSMatrix& operator= (
const BCRSMatrix& Mat)
911 if (&Mat==
this)
return *
this;
914 DUNE_THROW(InvalidStateException,
"BCRSMatrix can only be copied when both target and source are empty or fully built)");
918 deallocate(n!=Mat.n);
921 if (n>0 && n!=Mat.n) {
923 for(
row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
924 std::allocator_traits<
decltype(rowAllocator_)>::destroy(rowAllocator_, riter);
925 rowAllocator_.deallocate(r,n);
928 nnz_ = Mat.nonzeroes();
932 allocate(Mat.n, Mat.m, nnz_, n!=Mat.n,
true);
935 copyWindowStructure(Mat);
944 DUNE_THROW(InvalidStateException,
"Scalar assignment only works on fully built BCRSMatrix)");
957 CreateIterator (BCRSMatrix& _Mat, size_type _i)
958 : Mat(_Mat), i(_i), nnz(0), current_row(nullptr, Mat.j_.
get(), 0)
960 if (Mat.build_mode == unknown && Mat.ready == building)
962 Mat.build_mode = row_wise;
964 if (i==0 && Mat.ready != building)
965 DUNE_THROW(BCRSMatrixError,
"creation only allowed for uninitialized matrix");
966 if(Mat.build_mode!=row_wise)
967 DUNE_THROW(BCRSMatrixError,
"creation only allowed if row wise allocation was requested in the constructor");
968 if(i==0 && _Mat.N()==0)
974 CreateIterator& operator++()
977 if (Mat.ready != building)
978 DUNE_THROW(BCRSMatrixError,
"matrix already built up");
998 DUNE_THROW(BCRSMatrixError,
"allocated nnz too small");
1001 Mat.r[i].set(s,
nullptr,current_row.getindexptr());
1002 current_row.setindexptr(current_row.getindexptr()+s);
1006 B* b = Mat.allocator_.allocate(s);
1010 size_type* j = Mat.sizeAllocator_.allocate(s);
1011 Mat.r[i].set(s,b,j);
1015 Mat.r[i].set(0,
nullptr,
nullptr);
1018 std::copy(
pattern.cbegin(),
pattern.cend(), Mat.r[i].getindexptr());
1035 Mat.setDataPointers();
1043 bool operator!= (
const CreateIterator& it)
const
1045 return (i!=it.i) || (&Mat!=&it.Mat);
1049 bool operator== (
const CreateIterator& it)
const
1051 return (i==it.i) && (&Mat==&it.Mat);
1061 void insert (size_type j)
1067 bool contains (size_type j)
1085 typedef std::set<size_type,std::less<size_type> > PatternType;
1091 friend class CreateIterator;
1094 CreateIterator createbegin ()
1096 return CreateIterator(*
this,0);
1100 CreateIterator createend ()
1102 return CreateIterator(*
this,n);
1117 DUNE_THROW(BCRSMatrixError,
"requires random build mode");
1119 DUNE_THROW(BCRSMatrixError,
"matrix row sizes already built up");
1127#ifdef DUNE_ISTL_WITH_CHECKING
1128 if (r==0) DUNE_THROW(BCRSMatrixError,
"row not initialized yet");
1129 if (i>=n) DUNE_THROW(BCRSMatrixError,
"index out of range");
1131 return r[i].getsize();
1138 DUNE_THROW(BCRSMatrixError,
"requires random build mode");
1140 DUNE_THROW(BCRSMatrixError,
"matrix row sizes already built up");
1142 r[i].setsize(r[i].getsize()+s);
1149 DUNE_THROW(BCRSMatrixError,
"requires random build mode");
1151 DUNE_THROW(BCRSMatrixError,
"matrix row sizes already built up");
1157 total += r[i].getsize();
1162 allocate(n,m,total,
false,
false);
1163 else if(nnz_ < total)
1164 DUNE_THROW(BCRSMatrixError,
"Specified number of nonzeros ("<<nnz_<<
") not "
1165 <<
"sufficient for calculated nonzeros ("<<total<<
"! ");
1168 setColumnPointers(begin());
1191 DUNE_THROW(BCRSMatrixError,
"requires random build mode");
1193 DUNE_THROW(BCRSMatrixError,
"matrix already built up");
1195 DUNE_THROW(BCRSMatrixError,
"matrix row sizes not built up yet");
1197 DUNE_THROW(BCRSMatrixError,
"matrix size not set and no memory allocated yet");
1200 DUNE_THROW(BCRSMatrixError,
"column index exceeds matrix size");
1203 size_type*
const first = r[row].getindexptr();
1204 size_type*
const last = first + r[row].getsize();
1210 if (pos!=last && *pos ==
col)
return;
1213 size_type* end = std::lower_bound(pos,last,m);
1215 DUNE_THROW(BCRSMatrixError,
"row is too small");
1218 std::copy_backward(pos,end,end+1);
1231 template<
typename It>
1232 void setIndicesNoSort(
size_type row, It begin, It end)
1235 size_type* col_begin = r[row].getindexptr();
1238 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1239 DUNE_THROW(BCRSMatrixError,
"Given size of row " << row
1241 <<
") does not match number of passed entries (" << (col_end - col_begin) <<
")");
1254 template<
typename It>
1255 void setIndices(
size_type row, It begin, It end)
1258 size_type* col_begin = r[row].getindexptr();
1261 if ((col_end = std::copy(begin,end,r[row].getindexptr())) != col_begin + row_size)
1262 DUNE_THROW(BCRSMatrixError,
"Given size of row " << row
1264 <<
") does not match number of passed entries (" << (col_end - col_begin) <<
")");
1265 std::sort(col_begin,col_end);
1272 DUNE_THROW(BCRSMatrixError,
"requires random build mode");
1274 DUNE_THROW(BCRSMatrixError,
"matrix already built up");
1276 DUNE_THROW(BCRSMatrixError,
"row sizes are not built up yet");
1278 DUNE_THROW(BCRSMatrixError,
"matrix size not set and no memory allocated yet");
1281 RowIterator endi=end();
1282 for (RowIterator i=begin(); i!=endi; ++i)
1284 ColIterator endj = (*i).end();
1285 for (ColIterator j=(*i).begin(); j!=endj; ++j) {
1286 if (j.index() >= m) {
1287 dwarn <<
"WARNING: size of row "<< i.index()<<
" is "<<j.offset()<<
". But was specified as being "<< (*i).end().offset()
1288 <<
". This means you are wasting valuable space and creating additional cache misses!"<<std::endl;
1289 nnz_ -= ((*i).end().offset() - j.offset());
1290 r[i.index()].setsize(j.offset());
1319#ifdef DUNE_ISTL_WITH_CHECKING
1321 DUNE_THROW(BCRSMatrixError,
"requires implicit build mode");
1323 DUNE_THROW(BCRSMatrixError,
"matrix already built up, use operator[] for entry access now");
1325 DUNE_THROW(BCRSMatrixError,
"matrix size not set and no memory allocated yet");
1327 DUNE_THROW(InvalidStateException,
"You may only use entry() during the 'building' stage");
1330 DUNE_THROW(BCRSMatrixError,
"row index exceeds matrix size");
1332 DUNE_THROW(BCRSMatrixError,
"column index exceeds matrix size");
1335 size_type* begin = r[row].getindexptr();
1336 size_type* end = begin + r[row].getsize();
1344 std::ptrdiff_t offset = pos - r[row].getindexptr();
1345 B* aptr = r[row].getptr() + offset;
1351 if (r[row].getsize() == avg)
1352 return overflow[std::make_pair(row,
col)];
1359 std::ptrdiff_t offset = end - r[row].getindexptr();
1360 B* apos = r[row].getptr() + offset;
1363 r[row].setsize(r[row].getsize()+1);
1381 CompressionStatistics compress()
1384 DUNE_THROW(BCRSMatrixError,
"requires implicit build mode");
1386 DUNE_THROW(BCRSMatrixError,
"matrix already built up, no more need for compression");
1388 DUNE_THROW(BCRSMatrixError,
"matrix size not set and no memory allocated yet");
1390 DUNE_THROW(InvalidStateException,
"You may only call compress() at the end of the 'building' stage");
1393 CompressionStatistics stats;
1394 stats.overflow_total = overflow.size();
1402 typename OverflowType::iterator oit = overflow.begin();
1405 std::vector<size_type*> perm;
1417 typename std::vector<size_type*>::iterator it = perm.begin();
1418 for (
size_type* iit = begin; iit < begin + size; ++iit, ++it)
1422 std::sort(perm.begin(),perm.end(),PointerCompare<size_type>());
1425 r[i].setindexptr(jiit);
1428 for (it = perm.begin(); it != perm.end(); ++it)
1431 while ((oit!=overflow.end()) && (oit->first < std::make_pair(i,**it)))
1435 DUNE_THROW(Dune::ImplicitModeCompressionBufferExhausted,
1436 "Allocated memory for BCRSMatrix exhausted during compress()!"
1437 "Please increase either the average number of entries per row or the compressionBufferSize value."
1440 *jiit = oit->first.second;
1442 *aiit = oit->second;
1445 r[i].setsize(r[i].getsize()+1);
1450 DUNE_THROW(Dune::ImplicitModeCompressionBufferExhausted,
1451 "Allocated memory for BCRSMatrix exhausted during compress()!"
1452 "Please increase either the average number of entries per row or the compressionBufferSize value."
1458 B* apos = *it - j_.get() + a;
1464 while ((oit!=overflow.end()) && (oit->first.first == i))
1468 DUNE_THROW(Dune::ImplicitModeCompressionBufferExhausted,
1469 "Allocated memory for BCRSMatrix exhausted during compress()!"
1470 "Please increase either the average number of entries per row or the compressionBufferSize value."
1474 *jiit = oit->first.second;
1476 *aiit = oit->second;
1479 r[i].setsize(r[i].getsize()+1);
1483 if (r[i].getsize()>stats.maximum)
1484 stats.maximum = r[i].getsize();
1494 stats.mem_ratio = 1;
1498 std::ptrdiff_t diff = (r[n-1].getindexptr() + r[n-1].getsize() - j_.get());
1500 stats.avg = (double) (nnz_) / (double) n;
1501 stats.mem_ratio = (double) (nnz_) / (double) allocationSize_;
1515#ifdef DUNE_ISTL_WITH_CHECKING
1517 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1528 RowIterator endi=end();
1529 for (RowIterator i=begin(); i!=endi; ++i)
1531 ColIterator endj = (*i).end();
1532 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1543#ifdef DUNE_ISTL_WITH_CHECKING
1545 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1556 RowIterator endi=end();
1557 for (RowIterator i=begin(); i!=endi; ++i)
1559 ColIterator endj = (*i).end();
1560 for (ColIterator j=(*i).begin(); j!=endj; ++j)
1574 BCRSMatrix& operator+= (
const BCRSMatrix& b)
1576#ifdef DUNE_ISTL_WITH_CHECKING
1578 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1579 if(N()!=b.N() || M() != b.M())
1580 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
1582 RowIterator endi=end();
1583 ConstRowIterator j=b.begin();
1584 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1596 BCRSMatrix& operator-= (
const BCRSMatrix& b)
1598#ifdef DUNE_ISTL_WITH_CHECKING
1600 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1601 if(N()!=b.N() || M() != b.M())
1602 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
1604 RowIterator endi=end();
1605 ConstRowIterator j=b.begin();
1606 for (RowIterator i=begin(); i!=endi; ++i, ++j) {
1621 BCRSMatrix& axpy(
field_type alpha,
const BCRSMatrix& b)
1623#ifdef DUNE_ISTL_WITH_CHECKING
1625 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1626 if(N()!=b.N() || M() != b.M())
1627 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
1629 RowIterator endi=end();
1630 ConstRowIterator j=b.begin();
1631 for(RowIterator i=begin(); i!=endi; ++i, ++j)
1640 template<
class X,
class Y>
1641 void mv (
const X& x, Y& y)
const
1643#ifdef DUNE_ISTL_WITH_CHECKING
1645 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1646 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
1647 "Size mismatch: M: " << N() <<
"x" << M() <<
" x: " << x.N());
1648 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
1649 "Size mismatch: M: " << N() <<
"x" << M() <<
" y: " << y.N());
1651 ConstRowIterator endi=end();
1652 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1655 ConstColIterator endj = (*i).end();
1656 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1658 auto&& xj = Impl::asVector(x[j.index()]);
1659 auto&& yi = Impl::asVector(y[i.index()]);
1660 Impl::asMatrix(*j).umv(xj, yi);
1666 template<
class X,
class Y>
1667 void umv (
const X& x, Y& y)
const
1669#ifdef DUNE_ISTL_WITH_CHECKING
1671 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1672 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1673 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1675 ConstRowIterator endi=end();
1676 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1678 ConstColIterator endj = (*i).end();
1679 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1681 auto&& xj = Impl::asVector(x[j.index()]);
1682 auto&& yi = Impl::asVector(y[i.index()]);
1683 Impl::asMatrix(*j).umv(xj,yi);
1689 template<
class X,
class Y>
1690 void mmv (
const X& x, Y& y)
const
1692#ifdef DUNE_ISTL_WITH_CHECKING
1694 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1695 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1696 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1698 ConstRowIterator endi=end();
1699 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1701 ConstColIterator endj = (*i).end();
1702 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1704 auto&& xj = Impl::asVector(x[j.index()]);
1705 auto&& yi = Impl::asVector(y[i.index()]);
1706 Impl::asMatrix(*j).mmv(xj,yi);
1712 template<
class X,
class Y,
class F>
1713 void usmv (F&& alpha,
const X& x, Y& y)
const
1715#ifdef DUNE_ISTL_WITH_CHECKING
1717 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1718 if (x.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1719 if (y.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1721 ConstRowIterator endi=end();
1722 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1724 ConstColIterator endj = (*i).end();
1725 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1727 auto&& xj = Impl::asVector(x[j.index()]);
1728 auto&& yi = Impl::asVector(y[i.index()]);
1729 Impl::asMatrix(*j).usmv(alpha,xj,yi);
1735 template<
class X,
class Y>
1736 void mtv (
const X& x, Y& y)
const
1738#ifdef DUNE_ISTL_WITH_CHECKING
1740 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1741 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1742 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1750 template<
class X,
class Y>
1751 void umtv (
const X& x, Y& y)
const
1753#ifdef DUNE_ISTL_WITH_CHECKING
1755 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1756 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1757 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1759 ConstRowIterator endi=end();
1760 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1762 ConstColIterator endj = (*i).end();
1763 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1765 auto&& xi = Impl::asVector(x[i.index()]);
1766 auto&& yj = Impl::asVector(y[j.index()]);
1767 Impl::asMatrix(*j).umtv(xi,yj);
1773 template<
class X,
class Y>
1774 void mmtv (
const X& x, Y& y)
const
1776#ifdef DUNE_ISTL_WITH_CHECKING
1777 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1778 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1780 ConstRowIterator endi=end();
1781 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1783 ConstColIterator endj = (*i).end();
1784 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1786 auto&& xi = Impl::asVector(x[i.index()]);
1787 auto&& yj = Impl::asVector(y[j.index()]);
1788 Impl::asMatrix(*j).mmtv(xi,yj);
1794 template<
class X,
class Y>
1795 void usmtv (
const field_type& alpha,
const X& x, Y& y)
const
1797#ifdef DUNE_ISTL_WITH_CHECKING
1799 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1800 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1801 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1803 ConstRowIterator endi=end();
1804 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1806 ConstColIterator endj = (*i).end();
1807 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1809 auto&& xi = Impl::asVector(x[i.index()]);
1810 auto&& yj = Impl::asVector(y[j.index()]);
1811 Impl::asMatrix(*j).usmtv(alpha,xi,yj);
1817 template<
class X,
class Y>
1818 void umhv (
const X& x, Y& y)
const
1820#ifdef DUNE_ISTL_WITH_CHECKING
1822 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1823 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1824 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1826 ConstRowIterator endi=end();
1827 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1829 ConstColIterator endj = (*i).end();
1830 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1832 auto&& xi = Impl::asVector(x[i.index()]);
1833 auto&& yj = Impl::asVector(y[j.index()]);
1834 Impl::asMatrix(*j).umhv(xi,yj);
1840 template<
class X,
class Y>
1841 void mmhv (
const X& x, Y& y)
const
1843#ifdef DUNE_ISTL_WITH_CHECKING
1845 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1846 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1847 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1849 ConstRowIterator endi=end();
1850 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1852 ConstColIterator endj = (*i).end();
1853 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1855 auto&& xi = Impl::asVector(x[i.index()]);
1856 auto&& yj = Impl::asVector(y[j.index()]);
1857 Impl::asMatrix(*j).mmhv(xi,yj);
1863 template<
class X,
class Y>
1864 void usmhv (
const field_type& alpha,
const X& x, Y& y)
const
1866#ifdef DUNE_ISTL_WITH_CHECKING
1868 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1869 if (x.N()!=N()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1870 if (y.N()!=M()) DUNE_THROW(BCRSMatrixError,
"index out of range");
1872 ConstRowIterator endi=end();
1873 for (ConstRowIterator i=this->begin(); i!=endi; ++i)
1875 ConstColIterator endj = (*i).end();
1876 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
1878 auto&& xi = Impl::asVector(x[i.index()]);
1879 auto&& yj = Impl::asVector(y[j.index()]);
1880 Impl::asMatrix(*j).usmhv(alpha,xi,yj);
1889 typename FieldTraits<field_type>::real_type frobenius_norm2 ()
const
1891#ifdef DUNE_ISTL_WITH_CHECKING
1893 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1896 typename FieldTraits<field_type>::real_type sum=0;
1898 for (
auto&& row : (*
this))
1899 for (
auto&& entry : row)
1900 sum += Impl::asMatrix(entry).frobenius_norm2();
1906 typename FieldTraits<field_type>::real_type frobenius_norm ()
const
1908 return sqrt(frobenius_norm2());
1913 typename std::enable_if<!HasNaN<ft>::value,
int>
::type = 0>
1914 typename FieldTraits<ft>::real_type infinity_norm()
const {
1916 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1918 using real_type =
typename FieldTraits<ft>::real_type;
1922 for (
auto const &x : *
this) {
1924 for (
auto const &y : x)
1925 sum += Impl::asMatrix(y).infinity_norm();
1926 norm = max(sum, norm);
1933 typename std::enable_if<!HasNaN<ft>::value,
int>
::type = 0>
1934 typename FieldTraits<ft>::real_type infinity_norm_real()
const {
1936 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1938 using real_type =
typename FieldTraits<ft>::real_type;
1942 for (
auto const &x : *
this) {
1944 for (
auto const &y : x)
1945 sum += Impl::asMatrix(y).infinity_norm_real();
1946 norm = max(sum, norm);
1953 typename std::enable_if<HasNaN<ft>::value,
int>
::type = 0>
1954 typename FieldTraits<ft>::real_type infinity_norm()
const {
1956 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1958 using real_type =
typename FieldTraits<ft>::real_type;
1962 real_type isNaN = 1;
1963 for (
auto const &x : *
this) {
1965 for (
auto const &y : x)
1966 sum += Impl::asMatrix(y).infinity_norm();
1967 norm = max(sum, norm);
1971 return norm * (isNaN / isNaN);
1976 typename std::enable_if<HasNaN<ft>::value,
int>
::type = 0>
1977 typename FieldTraits<ft>::real_type infinity_norm_real()
const {
1979 DUNE_THROW(BCRSMatrixError,
"You can only call arithmetic operations on fully built BCRSMatrix instances");
1981 using real_type =
typename FieldTraits<ft>::real_type;
1985 real_type isNaN = 1;
1987 for (
auto const &x : *
this) {
1989 for (
auto const &y : x)
1990 sum += Impl::asMatrix(y).infinity_norm_real();
1991 norm = max(sum, norm);
1995 return norm * (isNaN / isNaN);
2017 nnz_ = std::accumulate( begin(), end(),
size_type( 0 ), [] (
size_type s,
const row_type &row ) {
return s+row.getsize(); } );
2022 BuildStage buildStage()
const
2038#ifdef DUNE_ISTL_WITH_CHECKING
2039 if (i<0 || i>=n) DUNE_THROW(BCRSMatrixError,
"row index out of range");
2040 if (j<0 || j>=m) DUNE_THROW(BCRSMatrixError,
"column index out of range");
2042 return (r[i].size() && r[i].find(j) != r[i].end());
2052 typename std::allocator_traits<A>::template rebind_alloc<B> allocator_;
2054 typename std::allocator_traits<A>::template rebind_alloc<row_type> rowAllocator_;
2056 typename std::allocator_traits<A>::template rebind_alloc<size_type> sizeAllocator_;
2072 std::shared_ptr<size_type> j_;
2076 double compressionBufferSize_;
2078 typedef std::map<std::pair<size_type,size_type>, B> OverflowType;
2079 OverflowType overflow;
2081 void setWindowPointers(ConstRowIterator row)
2083 row_type current_row(a,j_.get(),0);
2090 r[i].set(s,current_row.getptr(), current_row.getindexptr());
2092 current_row.setptr(current_row.getptr()+s);
2093 current_row.setindexptr(current_row.getindexptr()+s);
2096 r[i].set(0,
nullptr,
nullptr);
2106 void setColumnPointers(ConstRowIterator row)
2116 r[i].setindexptr(jptr);
2119 r[i].set(0,
nullptr,
nullptr);
2132 void setDataPointers()
2137 if (r[i].getsize() > 0) {
2142 r[i].set(0,
nullptr,
nullptr);
2146 aptr += r[i].getsize();
2151 void copyWindowStructure(
const BCRSMatrix& Mat)
2153 setWindowPointers(Mat.begin());
2156 for (
size_type i=0; i<n; i++) r[i] = Mat.r[i];
2168 void deallocate(
bool deallocateRows=
true)
2174 if (allocationSize_>0)
2180 for(B *aiter=a+(allocationSize_-1), *aend=a-1; aiter!=aend; --aiter)
2181 std::allocator_traits<
decltype(allocator_)>::destroy(allocator_, aiter);
2182 allocator_.deallocate(a,allocationSize_);
2190 if (r[i].getsize()>0)
2192 for (B *
col=r[i].getptr()+(r[i].getsize()-1),
2193 *colend = r[i].getptr()-1;
col!=colend; --
col) {
2194 std::allocator_traits<
decltype(allocator_)>::destroy(allocator_,
col);
2196 sizeAllocator_.deallocate(r[i].getindexptr(),1);
2197 allocator_.deallocate(r[i].getptr(),1);
2200 r[i].set(0,
nullptr,
nullptr);
2205 if (n>0 && deallocateRows && r) {
2206 for(
row_type *riter=r+(n-1), *rend=r-1; riter!=rend; --riter)
2207 std::allocator_traits<
decltype(rowAllocator_)>::destroy(rowAllocator_, riter);
2208 rowAllocator_.deallocate(r,n);
2240 nnz_ = allocationSize;
2241 allocationSize_ = allocationSize;
2247 DUNE_THROW(InvalidStateException,
"Rows have already been allocated, cannot allocate a second time");
2248 r = rowAllocator_.allocate(rows);
2250 for(
row_type* ri=r; ri!=r+rows; ++ri)
2251 std::allocator_traits<
decltype(rowAllocator_)>::construct(rowAllocator_, ri,
row_type());
2261 if (allocationSize_>0) {
2264 j_.reset(sizeAllocator_.allocate(allocationSize_),
2265 [alloc = sizeAllocator_, size = allocationSize_](
auto ptr)
mutable {
2266 alloc.deallocate(ptr, size);
2279 DUNE_THROW(InvalidStateException,
"Cannot allocate data array (already allocated)");
2280 if (allocationSize_>0) {
2281 a = allocator_.allocate(allocationSize_);
2284 new (a) B[allocationSize_];
2298 DUNE_THROW(InvalidStateException,
"implicit_allocate() may only be called in implicit build mode");
2300 DUNE_THROW(InvalidStateException,
"memory has already been allocated");
2303 if (compressionBufferSize_ < 0)
2304 DUNE_THROW(InvalidStateException,
"You have to set the implicit build mode parameters before starting to build the matrix");
2307 allocationSize_ = _n*avg + osize;
2309 allocate(_n, _m, allocationSize_,
true,
true);
2313 B* aptr = a + osize;
2316 r[i].set(0,aptr,jptr);
2326 template<
class B,
class A>
2329 using field_type =
typename BCRSMatrix<B, A>::field_type;
2330 using real_type =
typename FieldTraits<field_type>::real_type;
Some handy generic functions for ISTL matrices.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Helper functions for determining the vector/matrix block level.
*random access to the rows row_type & operator[](size_type i)
Definition bcrsmatrix.hh:545
*The type for the index access and the size typedef A::size_type size_type
Definition bcrsmatrix.hh:497
Matrix::block_type block_type
The block_type of the underlying matrix.
Definition bcrsmatrix.hh:125
built
The matrix structure is fully built.
Definition bcrsmatrix.hh:482
*equality bool equals(const RealRowIterator< ValueType > &other) const
Definition bcrsmatrix.hh:620
ImplicitMatrixBuilder(Matrix &m)
Creates an ImplicitMatrixBuilder for matrix m.
Definition bcrsmatrix.hh:170
std::ptrdiff_t distanceTo(const RealRowIterator< const ValueType > &other) const
Definition bcrsmatrix.hh:613
M_ Matrix
The underlying matrix.
Definition bcrsmatrix.hh:122
*the allocator type typedef A allocator_type
Definition bcrsmatrix.hh:494
*constructor RealRowIterator(row_type *_p, size_type _i)
Definition bcrsmatrix.hh:587
RealRowIterator(const RealRowIterator< ValueType > &it)
Definition bcrsmatrix.hh:596
ImplicitMatrixBuilder(Matrix &m, size_type rows, size_type cols, size_type avg_cols_per_row, double overflow_fraction)
Sets up matrix m for implicit construction using the given parameters and creates an ImplicitBmatrixu...
Definition bcrsmatrix.hh:194
block_type & operator[](size_type j) const
Returns entry in column j.
Definition bcrsmatrix.hh:142
size_type overflow_total
total number of elements written to the overflow area during construction.
Definition bcrsmatrix.hh:94
BuildMode
Definition bcrsmatrix.hh:506
std::ptrdiff_t distanceTo(const RealRowIterator< ValueType > &other) const
Definition bcrsmatrix.hh:607
*return index size_type index() const
Definition bcrsmatrix.hh:602
*equality bool equals(const RealRowIterator< const ValueType > &other) const
Definition bcrsmatrix.hh:627
rowSizesBuilt
The row sizes of the matrix are known.
Definition bcrsmatrix.hh:480
size_type M() const
The number of columns in the matrix.
Definition bcrsmatrix.hh:217
Matrix::size_type size_type
The size_type of the underlying matrix.
Definition bcrsmatrix.hh:128
row_object operator[](size_type i) const
Returns a proxy for entries in row i.
Definition bcrsmatrix.hh:205
building
Matrix is currently being built, some memory has been allocated, build mode and size are fixed.
Definition bcrsmatrix.hh:475
typename Imp::BlockTraits< B >::field_type field_type
Definition bcrsmatrix.hh:488
size_type maximum
maximum number of non-zeroes per row.
Definition bcrsmatrix.hh:92
double avg
average number of non-zeroes per row.
Definition bcrsmatrix.hh:90
notAllocated
Matrix is not built at all, no memory has been allocated, build mode and size can still be set.
Definition bcrsmatrix.hh:473
double mem_ratio
fraction of wasted memory resulting from non-used overflow area.
Definition bcrsmatrix.hh:99
notbuilt
A sparse block matrix with compressed row storage.
Definition bcrsmatrix.hh:471
*implement row_type with compressed vector typedef Imp::CompressedBlockVectorWindow< B, size_type > row_type
Definition bcrsmatrix.hh:500
size_type N() const
The number of rows in the matrix.
Definition bcrsmatrix.hh:211
*the type representing the components typedef B block_type
Definition bcrsmatrix.hh:491
*brief The unqualified value type typedef std::remove_const< T >::type ValueType
Definition bcrsmatrix.hh:579
@ implicit
Build entries randomly with an educated guess for the number of entries per row.
Definition bcrsmatrix.hh:535
@ random
Build entries randomly.
Definition bcrsmatrix.hh:526
@ unknown
Build mode not set!
Definition bcrsmatrix.hh:539
@ row_wise
Build in a row-wise manner.
Definition bcrsmatrix.hh:517
Definition allocator.hh:11
PropertyMapTypeSelector< Amg::VertexVisitedTag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > >::Type get(const Amg::VertexVisitedTag &tag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > &graph)
Definition dependency.hh:293
MatrixSparsityPatternGatherScatter< M, I >::ColIter MatrixSparsityPatternGatherScatter< M, I >::col
Definition matrixredistribute.hh:591
const Dtype_t BaseGetSuperLUType< T >::type
Definition supermatrix.hh:135
@ pattern
Definition matrixmarket.hh:301
Definition matrixutils.hh:211
Statistics about compression achieved in implicit mode.
Definition bcrsmatrix.hh:88
Proxy row object for entry access.
Definition bcrsmatrix.hh:137
Error specific to BCRSMatrix.
Definition istlexception.hh:24
A generic dynamic dense matrix.
Definition matrix.hh:561
A::size_type size_type
Type for indices and sizes.
Definition matrix.hh:577
T block_type
Export the type representing the components.
Definition matrix.hh:568
Definition matrixutils.hh:24