KASKADE 7 development version
|
#include <constantspace.hh>
A FE mapper class for FE spaces of globally constant functions. These can be interpreted as scalar variables as well, but the interpretation as constant functions allows for easy integration into the usual assembly framework. Since the FE space has dimension 1, there is exactly one degree of freedom.
Definition at line 33 of file constantspace.hh.
Classes | |
class | 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). More... | |
class | Converter |
A class mapping local shape function values and derivatives to global shape function values and derivatives. More... | |
struct | Element |
Public Types | |
typedef ScalarType | Scalar |
typedef GV | GridView |
typedef GridView::Grid | Grid |
typedef GridView::IndexSet | IndexSet |
using | ctype = typename GridView::ctype |
using | Cell = Kaskade::Cell< GridView > |
typedef LagrangeSimplexShapeFunctionSet< ctype, dim, ScalarType > | ShapeFunctionSetImplementation |
typedef LagrangeShapeFunctionSetContainer< ctype, dim, ScalarType >::value_type | ShapeFunctionSet |
typedef RangeView< std::vector< size_t >::const_iterator > | GlobalIndexRange |
typedef std::pair< size_t, int > | IndexPair |
typedef std::vector< IndexPair >::const_iterator | SortedIndexIterator |
typedef RangeView< SortedIndexIterator > | SortedIndexRange |
Public Member Functions | |
ConstantMapper (GridView const &gridView) | |
virtual | ~ConstantMapper () |
ShapeFunctionSet const & | shapefunctions (Cell const &cell) const |
ShapeFunctionSet const & | shapefunctions (size_t) const |
int | maxOrder () const |
Returns the maximal polynomial order of shape functions encountered in any cell. More... | |
GlobalIndexRange | initGlobalIndexRange () const |
Returns an empty range just for initialization purposes, since RangeView is not default constructible. More... | |
GlobalIndexRange | globalIndices (typename Grid::template Codim< 0 >::Entity const &) const |
Returns a const sequence containing the global indices of the shape functions associated to this cell. More... | |
GlobalIndexRange | globalIndices (size_t) const |
Returns a const sequence containing the global indices of the shape functions associated to this cell. More... | |
SortedIndexRange | sortedIndices (Cell const &) const |
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order. More... | |
SortedIndexRange | sortedIndices (size_t) const |
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order. More... | |
size_t | size () const |
Returns the number of global degrees of freedom managed. More... | |
Combiner | combiner (Cell const &cell, size_t index) const |
void | update () |
Recomputes the internal information after grid changes. Since there is no internal state of this trivial mapper, this method does exactly nothing. More... | |
Static Public Member Functions | |
static SortedIndexRange | initSortedIndexRange () |
Returns an empty range just for initialization, since RangeView is not default constructible. More... | |
Static Public Attributes | |
static int const | dim = Grid::dimension |
static int const | maxLocalSize = 1 |
static int const | continuity = std::numeric_limits<int>::max()-1 |
Continuity index. Since constant functions are \( C^\infty \), even analytic, this is infinity (or as close as int comes to infinity). More... | |
static bool const | globalSupport = true |
Whether the ansatz functions have global support (i.e. lead to dense matrices). More... | |
using Kaskade::ConstantMapper< ScalarType, GV >::Cell = Kaskade::Cell<GridView> |
Definition at line 48 of file constantspace.hh.
using Kaskade::ConstantMapper< ScalarType, GV >::ctype = typename GridView::ctype |
Definition at line 43 of file constantspace.hh.
typedef RangeView<std::vector<size_t>::const_iterator> Kaskade::ConstantMapper< ScalarType, GV >::GlobalIndexRange |
Definition at line 51 of file constantspace.hh.
typedef GridView::Grid Kaskade::ConstantMapper< ScalarType, GV >::Grid |
Definition at line 41 of file constantspace.hh.
typedef GV Kaskade::ConstantMapper< ScalarType, GV >::GridView |
Definition at line 40 of file constantspace.hh.
typedef std::pair<size_t,int> Kaskade::ConstantMapper< ScalarType, GV >::IndexPair |
Definition at line 53 of file constantspace.hh.
typedef GridView::IndexSet Kaskade::ConstantMapper< ScalarType, GV >::IndexSet |
Definition at line 42 of file constantspace.hh.
typedef ScalarType Kaskade::ConstantMapper< ScalarType, GV >::Scalar |
Definition at line 38 of file constantspace.hh.
typedef LagrangeShapeFunctionSetContainer<ctype,dim,ScalarType>::value_type Kaskade::ConstantMapper< ScalarType, GV >::ShapeFunctionSet |
Definition at line 50 of file constantspace.hh.
typedef LagrangeSimplexShapeFunctionSet<ctype,dim,ScalarType> Kaskade::ConstantMapper< ScalarType, GV >::ShapeFunctionSetImplementation |
Definition at line 49 of file constantspace.hh.
typedef std::vector<IndexPair>::const_iterator Kaskade::ConstantMapper< ScalarType, GV >::SortedIndexIterator |
Definition at line 54 of file constantspace.hh.
typedef RangeView<SortedIndexIterator> Kaskade::ConstantMapper< ScalarType, GV >::SortedIndexRange |
Definition at line 55 of file constantspace.hh.
|
inline |
Definition at line 74 of file constantspace.hh.
|
inlinevirtual |
Definition at line 82 of file constantspace.hh.
|
inline |
Returns a combiner for the given cell.
Definition at line 246 of file constantspace.hh.
|
inline |
Returns a const sequence containing the global indices of the shape functions associated to this cell.
Global indices start at 0 and are consecutive.
Definition at line 130 of file constantspace.hh.
|
inline |
Returns a const sequence containing the global indices of the shape functions associated to this cell.
Global indices start at 0 and are consecutive.
Definition at line 119 of file constantspace.hh.
|
inline |
Returns an empty range just for initialization purposes, since RangeView is not default constructible.
Definition at line 108 of file constantspace.hh.
|
inlinestatic |
Returns an empty range just for initialization, since RangeView is not default constructible.
Definition at line 138 of file constantspace.hh.
|
inline |
Returns the maximal polynomial order of shape functions encountered in any cell.
Definition at line 99 of file constantspace.hh.
|
inline |
Definition at line 84 of file constantspace.hh.
|
inline |
Definition at line 89 of file constantspace.hh.
|
inline |
Returns the number of global degrees of freedom managed.
Definition at line 163 of file constantspace.hh.
|
inline |
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order.
Definition at line 147 of file constantspace.hh.
|
inline |
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order.
Definition at line 155 of file constantspace.hh.
|
inline |
Recomputes the internal information after grid changes. Since there is no internal state of this trivial mapper, this method does exactly nothing.
Definition at line 255 of file constantspace.hh.
|
static |
Continuity index. Since constant functions are \( C^\infty \), even analytic, this is infinity (or as close as int comes to infinity).
Definition at line 61 of file constantspace.hh.
|
static |
Definition at line 45 of file constantspace.hh.
Referenced by Kaskade::ConstantMapper< ScalarType, GV >::shapefunctions().
|
static |
Whether the ansatz functions have global support (i.e. lead to dense matrices).
Definition at line 66 of file constantspace.hh.
|
static |
Definition at line 46 of file constantspace.hh.