KASKADE 7 development version
Public Member Functions | List of all members
Kaskade::NedelecMapper< ScalarType, GV >::Combiner Class Reference

A class implementing a matrix \( K \) mapping a subset of global degrees of freedom (those given by globalIndices()) to local degrees of freedom (shape functions). More...

#include <nedelecspace.hh>

Detailed Description

template<class ScalarType, class GV>
class Kaskade::NedelecMapper< ScalarType, GV >::Combiner

A class implementing a matrix \( K \) mapping a subset of global degrees of freedom (those given by globalIndices()) to local degrees of freedom (shape functions).

For edge elements, this is a diagonal matrix with entries 1 or -1 depending on the global vs. local orientation of the edge.

Definition at line 657 of file nedelecspace.hh.

Public Member Functions

 Combiner (Cell const &cell, IndexSet const &is, std::vector< std::pair< int, int > > const &vertexIndex)
 
template<class Matrix >
void rightTransform (Matrix &A) const
 In-place computation of \( A \leftarrow A K \). More...
 
template<int n, int m>
void rightTransform (std::vector< VariationalArg< Scalar, n, m > > &v) const
 In-place computation of row vectors \( v \leftarrow v K \). More...
 
template<class Matrix >
void leftPseudoInverse (Matrix &A) const
 In-place computation of \( A \leftarrow K^+ A \). More...
 
 operator Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > > () const
 

Constructor & Destructor Documentation

◆ Combiner()

template<class ScalarType , class GV >
Kaskade::NedelecMapper< ScalarType, GV >::Combiner::Combiner ( Cell const &  cell,
IndexSet const &  is,
std::vector< std::pair< int, int > > const &  vertexIndex 
)
inline
Parameters
cellthe cell for which the combiner is to construct
isthe index set used for global orientation of edges by global vertex indices
vertexIndexfor each shape function, i.e. for each edge, a pair of local starting vertex index and local end vertex index

Definition at line 662 of file nedelecspace.hh.

Member Function Documentation

◆ leftPseudoInverse()

template<class ScalarType , class GV >
template<class Matrix >
void Kaskade::NedelecMapper< ScalarType, GV >::Combiner::leftPseudoInverse ( Matrix &  A) const
inline

In-place computation of \( A \leftarrow K^+ A \).

Template Parameters
MatrixA matrix class satisfying the Dune::DenseMatrix interface.

Definition at line 708 of file nedelecspace.hh.

◆ operator Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > >()

template<class ScalarType , class GV >
Kaskade::NedelecMapper< ScalarType, GV >::Combiner::operator Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > > ( ) const
inline

Implicit conversion to a sparse matrix. This is just the diagonal.

Definition at line 719 of file nedelecspace.hh.

◆ rightTransform() [1/2]

template<class ScalarType , class GV >
template<class Matrix >
void Kaskade::NedelecMapper< ScalarType, GV >::Combiner::rightTransform ( Matrix &  A) const
inline

In-place computation of \( A \leftarrow A K \).

Template Parameters
MatrixA matrix class satisfying the Dune::DenseMatrix interface.

Definition at line 685 of file nedelecspace.hh.

◆ rightTransform() [2/2]

template<class ScalarType , class GV >
template<int n, int m>
void Kaskade::NedelecMapper< ScalarType, GV >::Combiner::rightTransform ( std::vector< VariationalArg< Scalar, n, m > > &  v) const
inline

In-place computation of row vectors \( v \leftarrow v K \).

Definition at line 695 of file nedelecspace.hh.


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