5#ifndef DUNE_ISTL_MATRIX_HH
6#define DUNE_ISTL_MATRIX_HH
15#include <dune/common/ftraits.hh>
16#include <dune/common/typetraits.hh>
17#include <dune/common/scalarvectorview.hh>
18#include <dune/common/scalarmatrixview.hh>
39 template<
class B,
class A=std::allocator<B> >
50 using field_type =
typename Imp::BlockTraits<B>::field_type;
98 this->n = rows*columns;
102 this->p = allocator_.allocate(this->n);
103 new (this->p)B[this->n];
120 columns_ = a.columns_;
124 this->p = allocator_.allocate(this->n);
125 new (this->p)B[this->n];
148 allocator_.deallocate(this->p,this->n);
160 allocator_.deallocate(this->p,this->n);
164 this->n = rows*columns;
167 this->p = allocator_.allocate(this->n);
168 new (this->p)B[this->n];
186 columns_ = a.columns_;
189 if (this->n!=a.n || rows_!=a.rows_)
196 allocator_.deallocate(this->p,this->n);
204 this->p = allocator_.allocate(this->n);
205 new (this->p)B[this->n];
231 (
static_cast<Imp::block_vector_unmanaged<B,A>&
>(*this)) = k;
243#ifdef DUNE_ISTL_WITH_CHECKING
244 if (i>=rows_) DUNE_THROW(
ISTLError,
"index out of range");
246 return window_type(this->p + i*columns_, columns_);
252#ifdef DUNE_ISTL_WITH_CHECKING
253 if (i<0 || i>=rows_) DUNE_THROW(
ISTLError,
"index out of range");
255 return window_type(this->p + i*columns_, columns_);
278 window_(data + _i*columns, columns)
286 window_.set(other.window_.getsize(),other.window_.getptr());
295 window_.set(other.window_.getsize(),other.window_.getptr());
303 window_.setptr(window_.getptr()+window_.getsize());
311 window_.setptr(window_.getptr()-window_.getsize());
318 return window_.getptr() == it.window_.getptr();
324 return window_.getptr() != it.window_.getptr();
330 return window_.getptr() == it.window_.getptr();
336 return window_.getptr() != it.window_.getptr();
367 return Iterator(this->p, columns_, 0);
373 return Iterator(this->p, columns_, rows_);
380 return Iterator(this->p, columns_, rows_-1);
387 return Iterator(this->p, columns_, -1);
393 return Iterator(this->p, columns_, std::min(i,rows_));
416 window_(const_cast<B*>(data + _i * columns), columns)
421 : i(it.i), window_(it.window_.getptr(),it.window_.getsize())
428 window_.set(other.window_.getsize(),other.window_.getptr());
436 window_.set(other.window_.getsize(),other.window_.getptr());
444 window_.setptr(window_.getptr()+window_.getsize());
452 window_.setptr(window_.getptr()-window_.getsize());
459 return window_.getptr() == it.window_.getptr();
465 return window_.getptr() != it.window_.getptr();
471 return window_.getptr() == it.window_.getptr();
477 return window_.getptr() != it.window_.getptr();
559 template<
class T,
class A=std::allocator<T> >
592 [[deprecated(
"Use free function blockLevel(). Will be removed after 2.8.")]]
609 data_.resize(rows,cols);
616 return data_.begin();
629 return data_.beforeEnd();
636 return data_.beforeBegin();
642 return data_.begin();
655 return data_.beforeEnd();
662 return data_.beforeBegin();
674#ifdef DUNE_ISTL_WITH_CHECKING
676 DUNE_THROW(
ISTLError,
"Can't access negative rows!");
678 DUNE_THROW(
ISTLError,
"Row index out of range!");
685#ifdef DUNE_ISTL_WITH_CHECKING
687 DUNE_THROW(
ISTLError,
"Can't access negative rows!");
689 DUNE_THROW(
ISTLError,
"Row index out of range!");
722#ifdef DUNE_ISTL_WITH_CHECKING
723 if(
N()!=b.
N() ||
M() != b.
M())
724 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
736#ifdef DUNE_ISTL_WITH_CHECKING
737 if(
N()!=b.
N() ||
M() != b.
M())
738 DUNE_THROW(RangeError,
"Matrix sizes do not match!");
749 out[j][i] = (*
this)[i][j];
762 out[i][j] += m1[i][k]*m2[k][j];
769 template <
class X,
class Y>
771#ifdef DUNE_ISTL_WITH_CHECKING
772 if (m.
M()!=vec.size())
773 DUNE_THROW(
ISTLError,
"Vector size doesn't match the number of matrix columns!");
778 for (
size_type i=0; i<out.size(); i++ ) {
780 out[i] += m[i][j]*vec[j];
787 template <
class X,
class Y>
788 void mv(
const X& x, Y& y)
const
790#ifdef DUNE_ISTL_WITH_CHECKING
791 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
792 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
798 auto&& xj = Impl::asVector(x[j]);
799 auto&& yi = Impl::asVector(y[i]);
800 Impl::asMatrix((*
this)[i][j]).umv(xj, yi);
806 template<
class X,
class Y>
807 void mtv (
const X& x, Y& y)
const
809#ifdef DUNE_ISTL_WITH_CHECKING
810 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"index out of range");
811 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"index out of range");
819 template <
class X,
class Y>
820 void umv(
const X& x, Y& y)
const
822#ifdef DUNE_ISTL_WITH_CHECKING
823 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
824 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
829 auto&& xj = Impl::asVector(x[j]);
830 auto&& yi = Impl::asVector(y[i]);
831 Impl::asMatrix((*
this)[i][j]).umv(xj, yi);
836 template<
class X,
class Y>
837 void mmv (
const X& x, Y& y)
const
839#ifdef DUNE_ISTL_WITH_CHECKING
840 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
841 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
846 auto&& xj = Impl::asVector(x[j]);
847 auto&& yi = Impl::asVector(y[i]);
848 Impl::asMatrix((*
this)[i][j]).mmv(xj, yi);
853 template <
class X,
class Y>
856#ifdef DUNE_ISTL_WITH_CHECKING
857 if (x.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
858 if (y.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
863 auto&& xj = Impl::asVector(x[j]);
864 auto&& yi = Impl::asVector(y[i]);
865 Impl::asMatrix((*
this)[i][j]).usmv(alpha, xj, yi);
870 template<
class X,
class Y>
871 void umtv (
const X& x, Y& y)
const
873#ifdef DUNE_ISTL_WITH_CHECKING
874 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
875 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
880 auto&& xi = Impl::asVector(x[i]);
881 auto&& yj = Impl::asVector(y[j]);
882 Impl::asMatrix((*
this)[i][j]).umtv(xi, yj);
887 template<
class X,
class Y>
888 void mmtv (
const X& x, Y& y)
const
890#ifdef DUNE_ISTL_WITH_CHECKING
891 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
892 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
897 auto&& xi = Impl::asVector(x[i]);
898 auto&& yj = Impl::asVector(y[j]);
899 Impl::asMatrix((*
this)[i][j]).mmtv(xi, yj);
904 template<
class X,
class Y>
907#ifdef DUNE_ISTL_WITH_CHECKING
908 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
909 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
914 auto&& xi = Impl::asVector(x[i]);
915 auto&& yj = Impl::asVector(y[j]);
916 Impl::asMatrix((*
this)[i][j]).usmtv(alpha, xi, yj);
921 template<
class X,
class Y>
922 void umhv (
const X& x, Y& y)
const
924#ifdef DUNE_ISTL_WITH_CHECKING
925 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
926 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
931 auto&& xi = Impl::asVector(x[i]);
932 auto&& yj = Impl::asVector(y[j]);
933 Impl::asMatrix((*
this)[i][j]).umhv(xi,yj);
938 template<
class X,
class Y>
939 void mmhv (
const X& x, Y& y)
const
941#ifdef DUNE_ISTL_WITH_CHECKING
942 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
943 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
948 auto&& xi = Impl::asVector(x[i]);
949 auto&& yj = Impl::asVector(y[j]);
950 Impl::asMatrix((*
this)[i][j]).mmhv(xi,yj);
955 template<
class X,
class Y>
958#ifdef DUNE_ISTL_WITH_CHECKING
959 if (x.N()!=
N()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
960 if (y.N()!=
M()) DUNE_THROW(
ISTLError,
"vector/matrix size mismatch!");
965 auto&& xi = Impl::asVector(x[i]);
966 auto&& yj = Impl::asVector(y[j]);
967 Impl::asMatrix((*
this)[i][j]).usmhv(alpha,xi,yj);
982 typename FieldTraits<field_type>::real_type sum=0;
985 sum += Impl::asMatrix(
data_[i][j]).frobenius_norm2();
991 typename std::enable_if<!HasNaN<ft>::value,
int>::type = 0>
993 using real_type =
typename FieldTraits<ft>::real_type;
997 for (
auto const &x : *
this) {
999 for (
auto const &y : x)
1000 sum += Impl::asMatrix(y).infinity_norm();
1001 norm = max(sum, norm);
1010 typename std::enable_if<!HasNaN<ft>::value,
int>::type = 0>
1012 using real_type =
typename FieldTraits<ft>::real_type;
1016 for (
auto const &x : *
this) {
1018 for (
auto const &y : x)
1019 sum += Impl::asMatrix(y).infinity_norm_real();
1020 norm = max(sum, norm);
1027 typename std::enable_if<HasNaN<ft>::value,
int>::type = 0>
1029 using real_type =
typename FieldTraits<ft>::real_type;
1033 real_type isNaN = 1;
1034 for (
auto const &x : *
this) {
1036 for (
auto const &y : x)
1037 sum += Impl::asMatrix(y).infinity_norm();
1038 norm = max(sum, norm);
1042 return norm * (isNaN / isNaN);
1047 typename std::enable_if<HasNaN<ft>::value,
int>::type = 0>
1049 using real_type =
typename FieldTraits<ft>::real_type;
1053 real_type isNaN = 1;
1054 for (
auto const &x : *
this) {
1056 for (
auto const &y : x)
1057 sum += Impl::asMatrix(y).infinity_norm_real();
1058 norm = max(sum, norm);
1062 return norm * (isNaN / isNaN);
1070#ifdef DUNE_ISTL_WITH_CHECKING
1071 if (i<0 || i>=
N()) DUNE_THROW(
ISTLError,
"row index out of range");
1072 if (j<0 || i>=
M()) DUNE_THROW(
ISTLError,
"column index out of range");
1091 template<
class T,
class A>
1095 using real_type =
typename FieldTraits<field_type>::real_type;
Helper functions for determining the vector/matrix block level.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Definition allocator.hh:11
derive error class from the base class in common
Definition istlexception.hh:19
A Vector of blocks with different blocksizes.
Definition matrix.hh:44
Imp::BlockVectorWindow< B, A > window_type
Definition matrix.hh:70
BlockVector< B, A > block_type
Same as value_type, here for historical reasons.
Definition matrix.hh:67
DenseMatrixBase & operator=(const DenseMatrixBase &a)
assignment
Definition matrix.hh:182
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition matrix.hh:50
Iterator beforeBegin() const
Definition matrix.hh:385
void resize(size_type rows, size_type columns)
same effect as constructor with same argument
Definition matrix.hh:153
const window_type const_reference
Definition matrix.hh:74
DenseMatrixBase(size_type rows, size_type columns)
Definition matrix.hh:95
Iterator end()
end Iterator
Definition matrix.hh:371
reference operator[](size_type i)
random access to blocks
Definition matrix.hh:241
BlockVector< B, A > value_type
Type of the elements of the outer vector, i.e., dynamic vectors of B.
Definition matrix.hh:63
Iterator find(size_type i)
random access returning iterator (end if not contained)
Definition matrix.hh:391
size_type N() const
number of blocks in the vector (are of variable size here)
Definition matrix.hh:539
ConstIterator beforeEnd() const
Definition matrix.hh:525
ConstIterator end() const
end ConstIterator
Definition matrix.hh:518
window_type reference
Definition matrix.hh:72
ConstIterator rend() const
end ConstIterator
Definition matrix.hh:531
Iterator beforeEnd()
Definition matrix.hh:378
ConstIterator begin() const
begin ConstIterator
Definition matrix.hh:512
A allocator_type
export the allocator type
Definition matrix.hh:53
DenseMatrixBase()
Definition matrix.hh:82
ConstIterator find(size_type i) const
random access returning iterator (end if not contained)
Definition matrix.hh:397
A::size_type size_type
The size type for the index access.
Definition matrix.hh:56
~DenseMatrixBase()
free dynamic memory
Definition matrix.hh:142
Iterator begin()
begin Iterator
Definition matrix.hh:365
DenseMatrixBase(const DenseMatrixBase &a)
copy constructor, has copy semantics
Definition matrix.hh:116
Iterator class for sequential access.
Definition matrix.hh:263
Iterator & operator--()
prefix decrement
Definition matrix.hh:308
size_type index() const
Definition matrix.hh:352
Iterator(Iterator &other)=default
Iterator(Iterator &&other)=default
bool operator!=(const Iterator &it) const
inequality
Definition matrix.hh:322
Iterator & operator=(Iterator &&other)
Move assignment.
Definition matrix.hh:282
Iterator & operator++()
prefix increment
Definition matrix.hh:300
Iterator()
constructor, no arguments
Definition matrix.hh:266
window_type & operator*() const
dereferencing
Definition matrix.hh:340
bool operator==(const Iterator &it) const
equality
Definition matrix.hh:316
Iterator & operator=(Iterator &other)
Copy assignment.
Definition matrix.hh:291
Iterator(B *data, size_type columns, size_type _i)
constructor
Definition matrix.hh:276
window_type * operator->() const
arrow
Definition matrix.hh:346
ConstIterator class for sequential access.
Definition matrix.hh:404
const window_type * operator->() const
arrow
Definition matrix.hh:487
const window_type & operator*() const
dereferencing
Definition matrix.hh:481
ConstIterator & operator++()
prefix increment
Definition matrix.hh:441
ConstIterator(const B *data, size_type columns, size_type _i)
constructor from pointer
Definition matrix.hh:414
ConstIterator & operator--()
prefix decrement
Definition matrix.hh:449
ConstIterator(const Iterator &it)
constructor from non_const iterator
Definition matrix.hh:420
bool operator!=(const ConstIterator &it) const
inequality
Definition matrix.hh:463
ConstIterator()
constructor
Definition matrix.hh:407
bool operator==(const ConstIterator &it) const
equality
Definition matrix.hh:457
size_type index() const
Definition matrix.hh:493
ConstIterator & operator=(Iterator &&other)
Definition matrix.hh:424
ConstIterator & operator=(Iterator &other)
Definition matrix.hh:432
A generic dynamic dense matrix.
Definition matrix.hh:561
size_type cols_
Number of columns of the matrix.
Definition matrix.hh:1088
FieldTraits< ft >::real_type infinity_norm() const
infinity norm (row sum norm, how to generalize for blocks?)
Definition matrix.hh:992
A allocator_type
Export the allocator.
Definition matrix.hh:571
FieldTraits< ft >::real_type infinity_norm_real() const
simplified infinity norm (uses Manhattan norm for complex values)
Definition matrix.hh:1011
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition matrix.hh:956
void usmv(const field_type &alpha, const X &x, Y &y) const
Definition matrix.hh:854
A::size_type size_type
Type for indices and sizes.
Definition matrix.hh:577
MatrixImp::DenseMatrixBase< T, A > data_
Abuse DenseMatrixBase as an engine for a 2d array ISTL-style.
Definition matrix.hh:1081
Matrix transpose() const
Return the transpose of the matrix.
Definition matrix.hh:745
void mtv(const X &x, Y &y) const
y = A^T x
Definition matrix.hh:807
void umv(const X &x, Y &y) const
y += A x
Definition matrix.hh:820
void mv(const X &x, Y &y) const
y = A x
Definition matrix.hh:788
MatrixImp::DenseMatrixBase< T, A >::ConstIterator ConstRowIterator
Const iterator over the matrix rows.
Definition matrix.hh:586
void setSize(size_type rows, size_type cols)
Change the matrix size.
Definition matrix.hh:608
RowIterator beforeBegin()
Definition matrix.hh:634
RowIterator beforeEnd()
Definition matrix.hh:627
Matrix()
Create empty matrix.
Definition matrix.hh:596
Matrix & operator-=(const Matrix &b)
Subtract the entries of another matrix from this one.
Definition matrix.hh:735
FieldTraits< field_type >::real_type frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition matrix.hh:980
row_type::iterator ColIterator
Iterator for the entries of each row.
Definition matrix.hh:583
ConstRowIterator beforeEnd() const
Definition matrix.hh:653
Matrix & operator=(const field_type &t)
Assignment from scalar.
Definition matrix.hh:666
RowIterator end()
Get iterator to one beyond last row.
Definition matrix.hh:620
const row_type operator[](size_type row) const
The const index operator.
Definition matrix.hh:684
ConstRowIterator end() const
Get const iterator to one beyond last row.
Definition matrix.hh:646
friend Y operator*(const Matrix< T > &m, const X &vec)
Generic matrix-vector multiplication.
Definition matrix.hh:770
Matrix< T > & operator*=(const field_type &scalar)
Multiplication with a scalar.
Definition matrix.hh:705
row_type operator[](size_type row)
The index operator.
Definition matrix.hh:673
void mmv(const X &x, Y &y) const
y -= A x
Definition matrix.hh:837
Matrix & operator+=(const Matrix &b)
Add the entries of another matrix to this one.
Definition matrix.hh:721
ConstRowIterator begin() const
Get const iterator to first row.
Definition matrix.hh:640
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition matrix.hh:939
RowIterator begin()
Get iterator to first row.
Definition matrix.hh:614
typename Imp::BlockTraits< T >::field_type field_type
Export the type representing the underlying field.
Definition matrix.hh:565
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition matrix.hh:589
static constexpr auto blocklevel
The number of nesting levels the matrix contains.
Definition matrix.hh:593
size_type M() const
Return the number of columns.
Definition matrix.hh:700
T block_type
Export the type representing the components.
Definition matrix.hh:568
bool exists(size_type i, size_type j) const
return true if (i,j) is in pattern
Definition matrix.hh:1068
Matrix< T > & operator/=(const field_type &scalar)
Division by a scalar.
Definition matrix.hh:711
friend Matrix< T > operator*(const Matrix< T > &m1, const Matrix< T > &m2)
Generic matrix multiplication.
Definition matrix.hh:755
void umtv(const X &x, Y &y) const
y += A^T x
Definition matrix.hh:871
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition matrix.hh:888
MatrixImp::DenseMatrixBase< T, A >::window_type row_type
The type implementing a matrix row.
Definition matrix.hh:574
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition matrix.hh:974
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition matrix.hh:905
void umhv(const X &x, Y &y) const
y += A^H x
Definition matrix.hh:922
size_type N() const
Return the number of rows.
Definition matrix.hh:695
ConstRowIterator beforeBegin() const
Definition matrix.hh:660
Matrix(size_type rows, size_type cols)
Create uninitialized matrix of size rows x cols.
Definition matrix.hh:601
MatrixImp::DenseMatrixBase< T, A >::Iterator RowIterator
Iterator over the matrix rows.
Definition matrix.hh:580
typename Matrix< T, A >::field_type field_type
Definition matrix.hh:1094
typename FieldTraits< field_type >::real_type real_type
Definition matrix.hh:1095