KASKADE 7 development version
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Attributes | List of all members
Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType > Class Template Reference

A degrees of freedom manager for piecewise continuous FEFunctionSpace s of local polynomials. More...

#include <partitionedspace.hh>

Detailed Description

template<class Tagger, class GV, class ScalarType = double>
class Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >

A degrees of freedom manager for piecewise continuous FEFunctionSpace s of local polynomials.

auto tagger = [](auto const& cell)
{
auto x = cell.geometry().center();
if (x[0]+x[1] > 1) return 0;
if (x[0]>0.5) return 1;
return 2;
};
using Tagger = decltype(tagger);
BrokenH1 brokenH1 = (gridManager,gridManager.grid().leafGridView(),order,tagger);
A representation of a finite element function space defined over a domain covered by a grid.
Template Parameters
Taggera callable type mapping cells to subdomain tags
ScalarTypea scalar arithmetic type to be used as shape function value type
GVgrid view

Definition at line 683 of file partitionedspace.hh.

Inheritance diagram for Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >:
Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV >, Tagger >

Classes

struct  Element
 

Public Types

typedef ScalarType Scalar
 
typedef GV GridView
 
typedef Base::ShapeFunctionSet ShapeFunctionSet
 
typedef int ConstructorArgument
 
template<int m>
using Element_t = FunctionSpaceElement< FEFunctionSpace< PiecewiseContinuousLagrangeMapper >, m >
 Type of the FunctionSpaceElement, associated to the FEFunctionSpace. More...
 
typedef Implementation::Grid Grid
 
using Cell = Kaskade::Cell< Grid >
 
typedef Implementation::Converter Converter
 
typedef Implementation::Combiner Combiner
 
typedef GridView::IndexSet IndexSet
 
typedef std::pair< size_t, int > IndexPair
 Indexing information for localized ansatz functions. More...
 
typedef RangeView< GlobalIndexIterator > GlobalIndexRange
 
typedef RangeView< SortedIndexIterator > SortedIndexRange
 

Public Member Functions

 PiecewiseContinuousLagrangeMapper (GridView const &gridView, int order, Tagger const &tagger)
 Constructor. More...
 
GridView const & gridView () const
 Returns the grid view used. More...
 
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 (Cell const &cell) const
 Returns an immutable sequence containing the global indices of the ansatz functions with support intersecting this cell. More...
 
GlobalIndexRange globalIndices (size_t n) const
 Returns an immutable sequence containing the global indices of the ansatz functions with support intersecting this cell. More...
 
SortedIndexRange sortedIndices (Cell const &cell) const
 Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order. More...
 
SortedIndexRange sortedIndices (size_t n) 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 scalar degrees of freedom managed. More...
 
ShapeFunctionSet const & shapefunctions (Cell const &cell, bool contained=false) const
 Returns the set of shape functions defined on this cell. More...
 
ShapeFunctionSet const & shapefunctions (size_t n) const
 Returns the set of shape functions defined on this cell. More...
 
ShapeFunctionSetshapefunctions_non_const (Cell const &cell)
 
ShapeFunctionSetshapefunctions_non_const (size_t n)
 
ShapeFunctionSet const & lowerShapeFunctions (Cell const &cell) const
 Returns a shape function set for a lower ansatz order (or an empty shape function set if there is no lower order). More...
 
Combiner combiner (Cell const &cell, size_t index) const
 Returns a combiner for the given leaf cell. More...
 
std::pair< bool, size_t > unrestrictedToRestrictedIndex (size_t unrestrictedIndex)
 
void update ()
 (Re)computes the internal data. 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 continuity = -1
 
static bool const globalSupport
 Whether the ansatz functions have global support (i.e. lead to dense matrices). More...
 

Protected Attributes

ContinuousLagrangeMapperImplementation< double, GV > implementation
 

Member Typedef Documentation

◆ Cell

Definition at line 238 of file partitionedspace.hh.

◆ Combiner

typedef Implementation::Combiner Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::Combiner
inherited

Definition at line 241 of file partitionedspace.hh.

◆ ConstructorArgument

template<class Tagger , class GV , class ScalarType = double>
typedef int Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::ConstructorArgument

Definition at line 693 of file partitionedspace.hh.

◆ Converter

Definition at line 240 of file partitionedspace.hh.

◆ Element_t

template<class Tagger , class GV , class ScalarType = double>
template<int m>
using Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::Element_t = FunctionSpaceElement<FEFunctionSpace<PiecewiseContinuousLagrangeMapper>,m>

Type of the FunctionSpaceElement, associated to the FEFunctionSpace.

Template Parameters
mnumber of components

Definition at line 702 of file partitionedspace.hh.

◆ GlobalIndexRange

typedef RangeView<GlobalIndexIterator> Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::GlobalIndexRange
inherited

Definition at line 276 of file partitionedspace.hh.

◆ Grid

Definition at line 237 of file partitionedspace.hh.

◆ GridView

template<class Tagger , class GV , class ScalarType = double>
typedef GV Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::GridView

Definition at line 691 of file partitionedspace.hh.

◆ IndexPair

typedef std::pair<size_t,int> Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::IndexPair
inherited

Indexing information for localized ansatz functions.

Localized ansatz functions on a cell have global indices \( I \in {\bf N}^m \), compare the LocalToGlobalMapper concept. This index vector need not be sorted, in fact, as the extracted dofs are to be multiplied with the combiner matrix, their order matters. Consequently, there is a one-to-one mapping of cell-local indices (second component) and affected global indices (first component).

Definition at line 255 of file partitionedspace.hh.

◆ IndexSet

typedef GridView::IndexSet Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::IndexSet
inherited

Definition at line 244 of file partitionedspace.hh.

◆ Scalar

template<class Tagger , class GV , class ScalarType = double>
typedef ScalarType Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::Scalar

Definition at line 690 of file partitionedspace.hh.

◆ ShapeFunctionSet

template<class Tagger , class GV , class ScalarType = double>
typedef Base::ShapeFunctionSet Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::ShapeFunctionSet

Definition at line 692 of file partitionedspace.hh.

◆ SortedIndexRange

typedef RangeView<SortedIndexIterator> Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::SortedIndexRange
inherited

Definition at line 277 of file partitionedspace.hh.

Constructor & Destructor Documentation

◆ PiecewiseContinuousLagrangeMapper()

template<class Tagger , class GV , class ScalarType = double>
Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::PiecewiseContinuousLagrangeMapper ( GridView const &  gridView,
int  order,
Tagger const &  tagger 
)
inline

Constructor.

Parameters
gridViewthe grid view on which to define the space, usually a leaf grid view
orderpolynomial ansatz order of shape functions (> 0)
taggera callable mapping cells to subdomain tags. The tagger is copied and stored in the space - make copying cheap.

Definition at line 717 of file partitionedspace.hh.

Member Function Documentation

◆ combiner()

Combiner Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::combiner ( Cell const &  cell,
size_t  index 
) const
inlineinherited

Returns a combiner for the given leaf cell.

Parameters
cellthe grid cell for which the combiner is requested
indexthe index of the cell

Definition at line 444 of file partitionedspace.hh.

◆ globalIndices() [1/2]

GlobalIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::globalIndices ( Cell const &  cell) const
inlineinherited

Returns an immutable sequence containing the global indices of the ansatz functions with support intersecting this cell.

Parameters
cellthe cell to inquire

Global indices start at 0 and are consecutive - in the range returned here, an unordered subset is contained.

Definition at line 334 of file partitionedspace.hh.

◆ globalIndices() [2/2]

GlobalIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::globalIndices ( size_t  n) const
inlineinherited

Returns an immutable sequence containing the global indices of the ansatz functions with support intersecting this cell.

Parameters
nthe index of the cell in the grid view used.

Global indices start at 0 and are consecutive - in the range returned here, an unordered subset is contained.

Definition at line 348 of file partitionedspace.hh.

◆ gridView()

GridView const & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::gridView ( ) const
inlineinherited

Returns the grid view used.

Definition at line 300 of file partitionedspace.hh.

◆ initGlobalIndexRange()

GlobalIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::initGlobalIndexRange ( ) const
inlineinherited

Returns an empty range just for initialization purposes, since RangeView is not default constructible.

Definition at line 318 of file partitionedspace.hh.

◆ initSortedIndexRange()

static SortedIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::initSortedIndexRange ( )
inlinestaticinherited

Returns an empty range just for initialization, since RangeView is not default constructible.

Definition at line 357 of file partitionedspace.hh.

◆ lowerShapeFunctions()

ShapeFunctionSet const & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::lowerShapeFunctions ( Cell const &  cell) const
inlineinherited

Returns a shape function set for a lower ansatz order (or an empty shape function set if there is no lower order).

Definition at line 431 of file partitionedspace.hh.

◆ maxOrder()

int Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::maxOrder ( ) const
inlineinherited

Returns the maximal polynomial order of shape functions encountered in any cell.

Definition at line 308 of file partitionedspace.hh.

◆ shapefunctions() [1/2]

ShapeFunctionSet const & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::shapefunctions ( Cell const &  cell,
bool  contained = false 
) const
inlineinherited

Returns the set of shape functions defined on this cell.

Parameters
cellthe codim 0 entity of the grid for wich the shape functions are to be retrieved
containedif true, the method may assume that the cell is contained in the index set of the space. (The other case occurs during interpolation between different grids).

Definition at line 399 of file partitionedspace.hh.

◆ shapefunctions() [2/2]

ShapeFunctionSet const & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::shapefunctions ( size_t  n) const
inlineinherited

Returns the set of shape functions defined on this cell.

Definition at line 416 of file partitionedspace.hh.

◆ shapefunctions_non_const() [1/2]

ShapeFunctionSet & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::shapefunctions_non_const ( Cell const &  cell)
inlineinherited

Definition at line 408 of file partitionedspace.hh.

◆ shapefunctions_non_const() [2/2]

ShapeFunctionSet & Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::shapefunctions_non_const ( size_t  n)
inlineinherited

Definition at line 422 of file partitionedspace.hh.

◆ size()

size_t Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::size ( ) const
inlineinherited

Returns the number of global scalar degrees of freedom managed.

Note that this does not correspond directly to the number of scalar coefficients in a FE function (if the FE function has more than one component).

Definition at line 388 of file partitionedspace.hh.

◆ sortedIndices() [1/2]

SortedIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::sortedIndices ( Cell const &  cell) const
inlineinherited

Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order.

Definition at line 367 of file partitionedspace.hh.

◆ sortedIndices() [2/2]

SortedIndexRange Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::sortedIndices ( size_t  n) const
inlineinherited

Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index order.

Definition at line 376 of file partitionedspace.hh.

◆ unrestrictedToRestrictedIndex()

std::pair< bool, size_t > Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::unrestrictedToRestrictedIndex ( size_t  unrestrictedIndex)
inlineinherited

Definition at line 453 of file partitionedspace.hh.

◆ update()

(Re)computes the internal data.

This has to be called after grid modifications and on construction.

Definition at line 463 of file partitionedspace.hh.

Member Data Documentation

◆ continuity

template<class Tagger , class GV , class ScalarType = double>
int const Kaskade::PiecewiseContinuousLagrangeMapper< Tagger, GV, ScalarType >::continuity = -1
static

Definition at line 694 of file partitionedspace.hh.

◆ globalSupport

bool const Kaskade::UniformPartitionedMapper< ContinuousLagrangeMapperImplementation< double, GV > , Tagger, PartitionedSpaceDetail::Empty >::globalSupport
staticinherited

Whether the ansatz functions have global support (i.e. lead to dense matrices).

Definition at line 282 of file partitionedspace.hh.

◆ implementation

Definition at line 654 of file partitionedspace.hh.


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