KASKADE 7 development version
Public Types | Public Member Functions | Related Functions | List of all members
Kaskade::DynamicMatrix< K > Class Template Reference

A LAPACK-compatible dense matrix class with shape specified at runtime. More...

#include <dynamicMatrix.hh>

Detailed Description

template<class K>
class Kaskade::DynamicMatrix< K >

A LAPACK-compatible dense matrix class with shape specified at runtime.

Template Parameters
Kthe entry type, either a scalar type or Dune::FieldMatrix or Dune::FieldVector

The memory layout is BLAS-like column major. For K = Dune::FieldMatrix<double,n,1>, classical LAPACK and BLAS routines can be applied.

In contrast to Dune::DynamicMatrix this works with contiguous memory. In contrast to Dune::Matrix this does not reallocate on resize (unless the number of entries grows).

Note that the inherited "solve" method does not work for K=Dune::FieldMatrix<??,1,1>, probably due to a Dune bug around common/densematrix.hh:808 (wrong argument type - pivots are matrix entries, not vector entries...).

Definition at line 377 of file dynamicMatrix.hh.

Inheritance diagram for Kaskade::DynamicMatrix< K >:

Public Types

using size_type = typename Base::size_type
 
using value_type = typename Base::value_type
 
using field_type = typename Dune::FieldTraits< K >::field_type
 
typedef Base::row_type row_type
 

Public Member Functions

DynamicMatrixoperator= (DynamicMatrix const &a)=default
 Copy assignment. More...
 
DynamicMatrixoperator= (DynamicMatrix &&a)=default
 Move assignment. More...
 
template<class OtherK >
DynamicMatrixoperator= (DynamicMatrix< OtherK > const &A)
 Assignment from matrix with other scalar type. More...
 
void resize (size_type r, size_type c)
 Resizes the matrix to r x c, leaving the entries in an undefined state. More...
 
void setSize (size_type r, size_type c)
 Resizes the matrix to r x c, leaving the entries in an undefined state. More...
 
void fill (K val)
 
template<class RRange , class CRange >
DynamicMatrix operator() (RRange const &ridx, CRange const &cidx) const
 Submatrix extraction. More...
 
Linear algebra operations
template<class X , class Y >
void mv (const X &x, Y &y) const
 matrix-vector product \( y \leftarrow A x \) More...
 
template<class X , class Y >
void mtv (const X &x, Y &y) const
 transpose matrix vector product \( y \leftarrow A^T x \) More...
 
template<class X , class Y >
void umtv (const X &x, Y &y) const
 update transpose matrix vector product \( y \leftarrow y + A^T x \). More...
 
template<class X , class Y >
void usmv (typename Dune::FieldTraits< Y >::field_type alpha, X const &x, Y &y) const
 update and scale matrix vector product \( y \leftarrow y + \alpha A x \). More...
 
template<class X , class Y >
void usmtv (typename Dune::FieldTraits< Y >::field_type alpha, X const &x, Y &y) const
 update and scale transposed matrix vector product \( y \leftarrow y + \alpha A^T x \). More...
 
DynamicMatrixoperator*= (field_type s)
 Scales the matrix by a scalar: \( A \gets sA \). More...
 
Low level access
size_type lda () const
 The leading dimension (for calling BLAS/LAPACK). More...
 
field_typedata ()
 Raw access to data (for calling BLAS/LAPACK). More...
 
field_type const * data () const
 Raw access to data (for calling BLAS/LAPACK). More...
 
size_type mat_rows () const
 
size_type mat_cols () const
 
row_reference mat_access (size_type i)
 
const_row_reference mat_access (size_type i) const
 

Related Functions

(Note that these are not member functions.)

template<class Entry >
std::ostream & operator<< (std::ostream &out, DynamicMatrix< Entry > const &A)
 pretty output of dense matrices More...
 
template<class Entry >
DynamicMatrix< typename EntryTraits< Entry >::transpose_type > transpose (DynamicMatrix< Entry > const &A)
 Computes the transpose of a matrix. More...
 
template<class Entry >
std::tuple< DynamicMatrix< Dune::FieldMatrix< typename EntryTraits< Entry >::field_type, 1, 1 > >, std::vector< typename EntryTraits< Entry >::real_type >, DynamicMatrix< Dune::FieldMatrix< typename EntryTraits< Entry >::field_type, 1, 1 > > > svd (DynamicMatrix< Entry > const &A)
 Computes the singular value decomposition \( A = U \Sigma V^T \). More...
 
template<class Entry >
EntryTraits< Entry >::real_type two_norm (DynamicMatrix< Entry > const &A)
 Computes the spectral norm of the given matrix. More...
 
template<class Entry >
EntryTraits< Entry >::real_type condition (DynamicMatrix< Entry > const &A)
 Computes the condition number of the given square matrix. More...
 
template<class MEntry , class Vector >
Vector gesv (DynamicMatrix< MEntry > const &A, Vector const &b)
 Computes the solution of \( Ax = b \) for general invertible matrices. More...
 
template<class MEntry , class Vector >
Vector posv (DynamicMatrix< MEntry > const &A, Vector const &b)
 Computes the solution of \( Ax = b \) if \( A \) is symmetric positive definite. More...
 
template<class Scalar >
void invert (DynamicMatrix< Scalar > &A)
 Computes the inverse \( A^{-1} \). More...
 
template<class Entry >
DynamicMatrix< Entry > & invertSpd (DynamicMatrix< Entry > &A, bool symmetrize=true)
 In place computation of the inverse \( A^{-1} \). More...
 
template<class Entry >
DynamicMatrix< typename EntryTraits< Entry >::transpose_type > pseudoinverse (DynamicMatrix< Entry > A)
 Computes the right pseudoinverse \( A^+ \) of a wide matrix such that \( A A^+ = I \). More...
 
template<class Entry >
std::vector< typename EntryTraits< Entry >::real_type > eigenvalues (DynamicMatrix< Entry > A)
 Computes the eigenvalues of a symmetric matrix \( A \). More...
 
template<class Entry >
std::vector< typename EntryTraits< Entry >::real_type > eig (DynamicMatrix< Entry > &A, bool computeEigenvectors=true)
 Computes the eigenvalues and eigenvectors of a symmetric matrix \( A \). More...
 
template<class Scalar , int d>
bool makePositiveSemiDefinite (Dune::FieldMatrix< Scalar, d, d > &A)
 Computes a positive semi-definite approximation to \( A \) by cutting off negative eigenvalues. More...
 
template<class Entry >
bool makePositiveSemiDefinite (DynamicMatrix< Entry > &A, typename ScalarTraits< typename EntryTraits< Entry >::field_type >::Real threshold=0)
 Computes a positive semi-definite approximation to \( A \) by cutting off negative eigenvalues. More...
 
template<class Entry >
std::tuple< std::vector< typename EntryTraits< Entry >::field_type >, DynamicMatrix< typename EntryTraits< Entry >::field_type >, std::vector< int > > qr (DynamicMatrix< Entry > const &A)
 Computes a rank-revealing QR decomposition of \( A \) using column pivoting. More...
 
template<class MEntry , class VEntry >
auto operator* (DynamicMatrix< MEntry > const &A, Dune::DynamicVector< VEntry > const &x)
 matrix-vector product \( Ax \) More...
 
template<class MEntry , class VEntry >
auto operator* (DynamicMatrix< MEntry > const &A, Dune::BlockVector< VEntry > const &x)
 matrix-vector product \( Ax \) More...
 
template<class MEntry >
auto operator* (typename Dune::FieldTraits< DynamicMatrix< MEntry > >::field_type x, DynamicMatrix< MEntry > A)
 multiplication with scalar More...
 
template<class MEntry >
auto operator- (DynamicMatrix< MEntry > const &A)
 Multiplication with -1. More...
 
template<class MEntry >
auto operator- (DynamicMatrix< MEntry > A, DynamicMatrix< MEntry > const &B)
 matrix subtraction \( A-B \) More...
 
template<class MEntry >
auto operator+ (DynamicMatrix< MEntry > A, DynamicMatrix< MEntry > const &B)
 matrix addition \( A+B \) More...
 
template<class EntryA , class EntryB >
void matrixMatrixProduct (DynamicMatrix< EntryA > const &A, DynamicMatrix< EntryB > const &B, DynamicMatrix< typename DynamicMatrixDetail::ProductTraits< EntryA, EntryB >::type > &C)
 matrix-matrix product \( A B \) More...
 
template<class EntryA , class EntryB >
DynamicMatrix< typename DynamicMatrixDetail::ProductTraits< EntryA, EntryB >::type > operator* (DynamicMatrix< EntryA > const &A, DynamicMatrix< EntryB > const &B)
 matrix-matrix product \( A B \) More...
 
template<class MEntry >
DynamicMatrix< MEntry > submatrix (DynamicMatrix< MEntry > const &A, int rowstart, int rowend, int colstart, int colend)
 extracts contiguous submatrices, copying the data More...
 
template<class EntryA , class EntryB >
auto kron (DynamicMatrix< EntryA > const &A, DynamicMatrix< EntryB > const &B)
 Computes the Kronecker product of two matrices. More...
 
template<class EntryA , class EntryB >
auto horzcat (DynamicMatrix< EntryA > const &A, DynamicMatrix< EntryB > const &B)
 Concatenates two matrices horizontally. More...
 
template<class EntryX , class EntryY >
auto vertcat (Dune::DynamicVector< EntryX > const &x, Dune::DynamicVector< EntryY > const &y)
 Concatenes two (column) vectors vertically. More...
 
template<class Entry >
DynamicMatrix< Entry > vertcat (DynamicMatrix< Entry > const &A, DynamicMatrix< Entry > const &B)
 Concatenates two matrices vertically. More...
 
template<class Entry >
auto dynamicUnitMatrix (Entry const &d, int n)
 Returns the diagonal matrix with the given entry replicated on the diagonal. More...
 

Constructors

template<class OtherK >
class DynamicMatrix
 
 DynamicMatrix ()
 Creates a 0 x 0 matrix. More...
 
 DynamicMatrix (Self const &a)=default
 Copy constructor. More...
 
template<class OtherK >
 DynamicMatrix (DynamicMatrix< OtherK > const &A)
 Copies from a matrix with different entry type. More...
 
 DynamicMatrix (Self &&a)=default
 Move constructor. More...
 
 DynamicMatrix (size_type r, size_type c, value_type v=value_type())
 Creates a r x c matrix, initialized with given value. More...
 

Member Typedef Documentation

◆ field_type

template<class K >
using Kaskade::DynamicMatrix< K >::field_type = typename Dune::FieldTraits<K>::field_type

Definition at line 390 of file dynamicMatrix.hh.

◆ row_type

template<class K >
typedef Base::row_type Kaskade::DynamicMatrix< K >::row_type

Definition at line 391 of file dynamicMatrix.hh.

◆ size_type

template<class K >
using Kaskade::DynamicMatrix< K >::size_type = typename Base::size_type

Definition at line 388 of file dynamicMatrix.hh.

◆ value_type

template<class K >
using Kaskade::DynamicMatrix< K >::value_type = typename Base::value_type

Definition at line 389 of file dynamicMatrix.hh.

Constructor & Destructor Documentation

◆ DynamicMatrix() [1/5]

template<class K >
Kaskade::DynamicMatrix< K >::DynamicMatrix ( )
inline

Creates a 0 x 0 matrix.

Definition at line 401 of file dynamicMatrix.hh.

◆ DynamicMatrix() [2/5]

template<class K >
Kaskade::DynamicMatrix< K >::DynamicMatrix ( Self const &  a)
default

Copy constructor.

◆ DynamicMatrix() [3/5]

template<class K >
template<class OtherK >
Kaskade::DynamicMatrix< K >::DynamicMatrix ( DynamicMatrix< OtherK > const &  A)
inline

Copies from a matrix with different entry type.

If the entry types have different size, the "flattend" matrices are copied scalar by scalar.

Definition at line 417 of file dynamicMatrix.hh.

◆ DynamicMatrix() [4/5]

template<class K >
Kaskade::DynamicMatrix< K >::DynamicMatrix ( Self &&  a)
default

Move constructor.

◆ DynamicMatrix() [5/5]

template<class K >
Kaskade::DynamicMatrix< K >::DynamicMatrix ( size_type  r,
size_type  c,
value_type  v = value_type() 
)
inline

Creates a r x c matrix, initialized with given value.

Definition at line 445 of file dynamicMatrix.hh.

Member Function Documentation

◆ data() [1/2]

template<class K >
field_type * Kaskade::DynamicMatrix< K >::data ( )
inline

Raw access to data (for calling BLAS/LAPACK).

Warning
Unless the entry type is scalar or Dune::FieldMatrix<Scalar,n,1>, the order of matrix entries is NOT column-major.

Definition at line 691 of file dynamicMatrix.hh.

Referenced by Kaskade::DynamicMatrix< K >::matrixMatrixProduct(), Kaskade::DynamicMatrix< K >::mv(), Kaskade::DynamicMatrix< K >::usmtv(), and Kaskade::DynamicMatrix< K >::usmv().

◆ data() [2/2]

template<class K >
field_type const * Kaskade::DynamicMatrix< K >::data ( ) const
inline

Raw access to data (for calling BLAS/LAPACK).

Warning
Unless the entry type is scalar or Dune::FieldMatrix<Scalar,n,1>, the order of matrix entries is NOT column-major.

Definition at line 697 of file dynamicMatrix.hh.

◆ fill()

template<class K >
void Kaskade::DynamicMatrix< K >::fill ( val)
inline

◆ lda()

template<class K >
size_type Kaskade::DynamicMatrix< K >::lda ( ) const
inline

The leading dimension (for calling BLAS/LAPACK).

Since entries are stored gap-less, this is always the number of matrix rows times the number of scalar rows in the entries.

Definition at line 682 of file dynamicMatrix.hh.

Referenced by Kaskade::DynamicMatrix< K >::matrixMatrixProduct(), Kaskade::DynamicMatrix< K >::mv(), Kaskade::DynamicMatrix< K >::usmtv(), and Kaskade::DynamicMatrix< K >::usmv().

◆ mat_access() [1/2]

template<class K >
row_reference Kaskade::DynamicMatrix< K >::mat_access ( size_type  i)
inline

Definition at line 702 of file dynamicMatrix.hh.

◆ mat_access() [2/2]

template<class K >
const_row_reference Kaskade::DynamicMatrix< K >::mat_access ( size_type  i) const
inline

Definition at line 703 of file dynamicMatrix.hh.

◆ mat_cols()

template<class K >
size_type Kaskade::DynamicMatrix< K >::mat_cols ( ) const
inline

Definition at line 701 of file dynamicMatrix.hh.

◆ mat_rows()

template<class K >
size_type Kaskade::DynamicMatrix< K >::mat_rows ( ) const
inline

Definition at line 700 of file dynamicMatrix.hh.

◆ mtv()

template<class K >
template<class X , class Y >
void Kaskade::DynamicMatrix< K >::mtv ( const X &  x,
Y &  y 
) const
inline

transpose matrix vector product \( y \leftarrow A^T x \)

This computes \( y \leftarrow A^T x \). We cannot rely on the Dune implementation in DenseMatrix, as it does not transpose the individual entries. This leads to errors when the entries are not square.

Definition at line 578 of file dynamicMatrix.hh.

◆ mv()

template<class K >
template<class X , class Y >
void Kaskade::DynamicMatrix< K >::mv ( const X &  x,
Y &  y 
) const
inline

matrix-vector product \( y \leftarrow A x \)

Template Parameters
Xa vector type with subscript operator[]. Its value type must be multiplicable by the matrix entries.
Ya vector type with subscript operator[]. Matrix entries times x entries need to be convertible to y's value type.

This implementation calls BLAS gemv if possible, and falls back to a simple double loop otherwise. The BLAS call is possible if the entries are either scalars or Dune FieldMatrix of size n x 1, with arbitrary n.

Definition at line 549 of file dynamicMatrix.hh.

Referenced by Kaskade::DynamicMatrix< K >::operator*().

◆ operator()()

template<class K >
template<class RRange , class CRange >
DynamicMatrix Kaskade::DynamicMatrix< K >::operator() ( RRange const &  ridx,
CRange const &  cidx 
) const
inline

Submatrix extraction.

This extracts a submatrix, as in Matlab A(ridx,cidx).

Template Parameters
RRangean STL random access range type
CRangean STL random access range type
Todo:
change the implementation to work with forward ranges

Definition at line 522 of file dynamicMatrix.hh.

◆ operator*=()

template<class K >
DynamicMatrix & Kaskade::DynamicMatrix< K >::operator*= ( field_type  s)
inline

Scales the matrix by a scalar: \( A \gets sA \).

This should be provided by the base class Dune::DenseMatrix, but as of Dune 2.5.2, the definition of field_type as value_type in the base class prevents this operator to work for matrix-valued entries. Thus we provide a working implementation.

Definition at line 660 of file dynamicMatrix.hh.

◆ operator=() [1/3]

template<class K >
DynamicMatrix & Kaskade::DynamicMatrix< K >::operator= ( DynamicMatrix< K > &&  a)
default

Move assignment.

◆ operator=() [2/3]

template<class K >
DynamicMatrix & Kaskade::DynamicMatrix< K >::operator= ( DynamicMatrix< K > const &  a)
default

Copy assignment.

The matrix size is adjusted as needed.

◆ operator=() [3/3]

template<class K >
template<class OtherK >
DynamicMatrix & Kaskade::DynamicMatrix< K >::operator= ( DynamicMatrix< OtherK > const &  A)
inline

Assignment from matrix with other scalar type.

Definition at line 475 of file dynamicMatrix.hh.

◆ resize()

template<class K >
void Kaskade::DynamicMatrix< K >::resize ( size_type  r,
size_type  c 
)
inline

Resizes the matrix to r x c, leaving the entries in an undefined state.

On resize, all row objects and iterators are invalidated. This method adheres to the Dune::DynamicMatrix interface.

Definition at line 489 of file dynamicMatrix.hh.

Referenced by Kaskade::EmptyShapeFunctionSet< ctype, dimension, T, comp >::interpolate(), Kaskade::DynamicMatrix< K >::operator=(), Kaskade::DynamicMatrix< K >::setSize(), and Kaskade::DynamicMatrixDetail::unflatten().

◆ setSize()

template<class K >
void Kaskade::DynamicMatrix< K >::setSize ( size_type  r,
size_type  c 
)
inline

Resizes the matrix to r x c, leaving the entries in an undefined state.

On resize, all row objects and iterators are invalidated. This method adheres to the Dune::Matrix interface.

Definition at line 501 of file dynamicMatrix.hh.

Referenced by Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::evaluate(), Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::initHierarchicalProjection(), Kaskade::interpolateGloballyWeak(), and Kaskade::CoarseningDetail::GetLocalTransferProjection< Cell >::operator()().

◆ umtv()

template<class K >
template<class X , class Y >
void Kaskade::DynamicMatrix< K >::umtv ( const X &  x,
Y &  y 
) const
inline

update transpose matrix vector product \( y \leftarrow y + A^T x \).

We cannot rely on the Dune implementation in DenseMatrix, as it does not transpose the individual entries. This leads to errors when the entries are not square.

Definition at line 592 of file dynamicMatrix.hh.

Referenced by Kaskade::DynamicMatrix< K >::mtv().

◆ usmtv()

template<class K >
template<class X , class Y >
void Kaskade::DynamicMatrix< K >::usmtv ( typename Dune::FieldTraits< Y >::field_type  alpha,
X const &  x,
Y &  y 
) const
inline

update and scale transposed matrix vector product \( y \leftarrow y + \alpha A^T x \).

Template Parameters
Xa vector type with subscript operator[]. Its value type must be multiplicable by the matrix entries.
Ya vector type with subscript operator[]. Matrix entries times x entries need to be convertible to y's value type.

Dune 2.6 has a bug (the matrix entries are not transposed for multiplication), such that we provide a corrected version here.

Definition at line 639 of file dynamicMatrix.hh.

Referenced by Kaskade::DynamicMatrix< K >::umtv().

◆ usmv()

template<class K >
template<class X , class Y >
void Kaskade::DynamicMatrix< K >::usmv ( typename Dune::FieldTraits< Y >::field_type  alpha,
X const &  x,
Y &  y 
) const
inline

update and scale matrix vector product \( y \leftarrow y + \alpha A x \).

Template Parameters
Xa vector type with subscript operator[]. Its value type must be multiplicable by the matrix entries.
Ya vector type with subscript operator[]. Matrix entries times x entries need to be convertible to y's value type.

Dune::DenseMatrix::usmv appears to run into a compile time bug, because DenseMatrix::field_type is defined as the same as value_type, which may not be scalar (as of Dune 2.5).

Definition at line 609 of file dynamicMatrix.hh.

Friends And Related Function Documentation

◆ DynamicMatrix

template<class K >
template<class OtherK >
friend class DynamicMatrix
friend

Definition at line 409 of file dynamicMatrix.hh.

◆ operator<<()

template<class Entry >
std::ostream & operator<< ( std::ostream &  out,
DynamicMatrix< Entry > const &  A 
)
related

pretty output of dense matrices

Definition at line 732 of file dynamicMatrix.hh.


The documentation for this class was generated from the following files: