DUNE PDELab (git)

Dune::PDELab::LocalMatrix< T, W > Class Template Reference

A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpaces. More...

#include <dune/pdelab/gridoperator/common/localmatrix.hh>

Public Types

typedef std::vector< T > BaseContainer
 The type of the underlying storage container. More...
 
typedef BaseContainer::value_type value_type
 The value type of this container.
 
typedef BaseContainer::size_type size_type
 The size type of this container.
 
typedef BaseContainer::reference reference
 The reference type of this container.
 
typedef BaseContainer::const_reference const_reference
 The const reference type of this container.
 
typedef W weight_type
 The weight type of this container. More...
 
typedef WeightedMatrixAccumulationView< LocalMatrixWeightedAccumulationView
 An accumulate-only view of this container that automatically applies a weight to all contributions.
 

Public Member Functions

 LocalMatrix ()
 Default constructor.
 
 LocalMatrix (size_type r, size_type c)
 Construct a LocalMatrix with r rows and c columns.
 
 LocalMatrix (size_type r, size_type c, const T &t)
 Construct a LocalMatrix with r rows and c columns and initialize its entries with t.
 
void resize (size_type r, size_type c)
 Resize the matrix.
 
LocalMatrixoperator= (const T &t)
 Assign t to all entries of the matrix.
 
void assign (size_type r, size_type c, const T &t)
 Resize the matrix and assign t to all entries.
 
template<typename LFSU , typename LFSV >
T & operator() (const LFSV &lfsv, size_type i, const LFSU &lfsu, size_type j)
 Access the value associated with the i-th DOF of lfsv and the j-th DOF of lfsu. More...
 
template<typename LFSU , typename LFSV >
const T & operator() (const LFSV &lfsv, size_type i, const LFSU &lfsu, size_type j) const
 Access the value associated with the i-th DOF of lfsv and the j-th DOF of lfsu (const version). More...
 
LocalMatrixoperator*= (const T &x)
 Multiplies all entries of the matrix with x.
 
size_type nrows () const
 Returns the number of rows.
 
size_type ncols () const
 Returns the number of columns.
 
template<class X , class R >
void umv (const X &x, R &y) const
 y += A x
 
template<class X , class R >
void usmv (const value_type &alpha, const X &x, R &y) const
 y += alpha A x
 
WeightedAccumulationView weightedAccumulationView (weight_type weight)
 Returns a weighted accumulate-only view of this matrix with the given weight.
 
BaseContainerbase ()
 Returns the underlying storage container. More...
 
const BaseContainerbase () const
 Returns the underlying storage container (const version). More...
 
value_typegetEntry (size_type i, size_type j)
 Direct (unmapped) access to the (i,j)-th entry of the matrix. More...
 
const value_typegetEntry (size_type i, size_type j) const
 Direct (unmapped) access to the (i,j)-th entry of the matrix (const version). More...
 

Detailed Description

template<typename T, typename W = T>
class Dune::PDELab::LocalMatrix< T, W >

A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpaces.

This container represents a dense matrix based on a std::vector-like storage container and supports accessing its entries indexed by pairs of (LocalFunctionSpace,DOF of LocalFunctionSpace). If the LocalFunctionSpaces contain tags indicating whether they are trial or test spaces, the access methods will also assert that the first space is a test space and the second space is a trial space.

Template Parameters
TThe type of values to store in the matrix.
WThe type of weight applied in a WeightedAccumulationView.

Member Typedef Documentation

◆ BaseContainer

template<typename T , typename W = T>
typedef std::vector<T> Dune::PDELab::LocalMatrix< T, W >::BaseContainer

The type of the underlying storage container.

Warning
This is not a matrix-like container anymore, but a std::vector-like one!

◆ weight_type

template<typename T , typename W = T>
typedef W Dune::PDELab::LocalMatrix< T, W >::weight_type

The weight type of this container.

A value of this type will be used to assign a weight to contributions in a WeightedAccumulationView.

Member Function Documentation

◆ base() [1/2]

template<typename T , typename W = T>
BaseContainer& Dune::PDELab::LocalMatrix< T, W >::base ( )
inline

Returns the underlying storage container.

Warning
This is not a matrix-like container anymore, but a std::vector-like one!

◆ base() [2/2]

template<typename T , typename W = T>
const BaseContainer& Dune::PDELab::LocalMatrix< T, W >::base ( ) const
inline

Returns the underlying storage container (const version).

Warning
This is not a matrix-like container anymore, but a std::vector-like one!

◆ getEntry() [1/2]

template<typename T , typename W = T>
value_type& Dune::PDELab::LocalMatrix< T, W >::getEntry ( size_type  i,
size_type  j 
)
inline

Direct (unmapped) access to the (i,j)-th entry of the matrix.

Note
This method is just a convenience method for helper function doing things like printing the matrix, as the BaseContainer is a linear std::vector-like object and does not know about the shape of the matrix anymore.
Warning
This method should normally not be called by normal users. Only use it if you really know what you are doing!

Referenced by Dune::PDELab::LocalMatrix< T, W >::operator()().

◆ getEntry() [2/2]

template<typename T , typename W = T>
const value_type& Dune::PDELab::LocalMatrix< T, W >::getEntry ( size_type  i,
size_type  j 
) const
inline

Direct (unmapped) access to the (i,j)-th entry of the matrix (const version).

Note
This method is just a convenience method for helper function doing things like printing the matrix, as the BaseContainer is a linear std::vector-like object and does not know about the shape of the matrix anymore.
Warning
This method should normally not be called by normal users. Only use it if you really know what you are doing!

◆ operator()() [1/2]

template<typename T , typename W = T>
template<typename LFSU , typename LFSV >
T& Dune::PDELab::LocalMatrix< T, W >::operator() ( const LFSV &  lfsv,
size_type  i,
const LFSU &  lfsu,
size_type  j 
)
inline

Access the value associated with the i-th DOF of lfsv and the j-th DOF of lfsu.

Parameters
lfsvThe LocalFunctionSpace whose DOF will determine the row index of the entry.
iThe index of the DOF of lfsv that will determine the row index of the entry.
lfsuThe LocalFunctionSpace whose DOF will determine the column index of the entry.
jThe index of the DOF of lfsu that will determine the column index of the entry.

References Dune::PDELab::LocalMatrix< T, W >::getEntry().

◆ operator()() [2/2]

template<typename T , typename W = T>
template<typename LFSU , typename LFSV >
const T& Dune::PDELab::LocalMatrix< T, W >::operator() ( const LFSV &  lfsv,
size_type  i,
const LFSU &  lfsu,
size_type  j 
) const
inline

Access the value associated with the i-th DOF of lfsv and the j-th DOF of lfsu (const version).

Parameters
lfsvThe LocalFunctionSpace whose DOF will determine the row index of the entry.
iThe index of the DOF of lfsv that will determine the row index of the entry.
lfsuThe LocalFunctionSpace whose DOF will determine the column index of the entry.
jThe index of the DOF of lfsu that will determine the column index of the entry.

References Dune::PDELab::LocalMatrix< T, W >::getEntry().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 8, 22:30, 2024)