KASKADE 7 development version
Classes | Public Types | Static Public Attributes | Protected Attributes | Related Functions | List of all members
Kaskade::FunctionSpaceElement< FunctionSpace, m > Class Template Reference

A class for representing finite element functions. More...

#include <functionspace.hh>

Detailed Description

template<class FunctionSpace, int m>
class Kaskade::FunctionSpaceElement< FunctionSpace, m >

A class for representing finite element functions.

A class for representing finite element function space elements (also known as finite element functions). Finite element functions are piecewisely smooth functions (on each codim 0 entity, i.e. cell, of an associated grid) that may or may not be globally continuous or differentiable. FE functions support pointwise evaluation and differentiation, but on sets of measure zero (more precisely, on the codim 1 entities, i.e. faces, of the grid), these values are undefined, and calling the corresponding methods may yield arbitrary values.

Template Parameters
FunctionSpaceEach FunctionSpaceElement is associated with an FEFunctionSpace (it "lives" in this space). This function space is defined in the first template parameter
mFunctionSpaceElements may have multiple components. The number of components attached to each shape function is given by this template parameter. Note that even with m=1, the function can be vector-valued if the shape functions themselves are vector-valued.

Usually, FSElements are created as part of a VariableSet, which in turn is created via a VariableSetDescription, containing information about variables and their spaces (VariableDescription s). FSElements contain the data for the variables.

If a GridManager is used, FunctionSpaceElements are automatically prolongated, if the grid is refined, and approximated, if the grid is coarsened.

If a VariableSet vs is at hand, the FunctionSpaceElements can be accessed as public data-members in a boost::fusion::vector vs.data. This means: boost::fusion::at_c<N>(vs.data) gives a reference to the N'th FSElement in vs.

Simple assignments are performed by the assignment operator. For more complex specification of FE function values see spaceTransfer, interpolateGlobally, and interpolateGloballyWeak.

Definition at line 314 of file functionspace.hh.

Public Types

using Space = FunctionSpace
 FEFunctionSpace, this function is an element of. More...
 
using GridView = typename Space::GridView
 the grid view this function is defined on More...
 
using ctype = typename GridView::ctype
 the scalar type used for spatial coordinates More...
 
typedef FunctionSpaceElement< Space, m > Self
 own type More...
 
typedef Space::Scalar Scalar
 scalar field type More...
 
using field_type = typename Space::field_type
 scalar field type More...
 
typedef Dune::FieldVector< Scalar, m > StorageValueType
 type of the elements of the data vector More...
 
using SfValueType = typename Space::SfValueType
 type of shape function values More...
 
using block_type = StorageValueType
 type of coefficient vector entries More...
 
typedef Dune::BlockVector< StorageValueTypeStorageType
 type of the data vector More...
 
typedef Dune::FieldVector< Scalar, componentsValueType
 return type of the function value(...) More...
 
typedef Dune::FieldMatrix< Scalar, components, Space::dimworld > DerivativeType
 return type of the function derivative (...) More...
 
using HessianType = Tensor< Scalar, components, Space::dimworld, Space::dimworld >
 return type of the function hessian(...) More...
 

Public Member Functions

Constructors & Destructor
 FunctionSpaceElement (Space const &fs)
 Constructor. More...
 
 FunctionSpaceElement (Self const &fse)
 Copy constructor. More...
 
 ~FunctionSpaceElement ()
 
Assignment
Selfoperator= (Self const &fse)
 Copy assignment. More...
 
Selfoperator= (Self &&fse)
 Move assignment. More...
 
template<class OtherSpace >
Selfoperator= (FunctionSpaceElement< OtherSpace, m > const &f)
 Assignment from functions belonging to other function spaces. More...
 
template<template< class, class > class DomainMapper>
Selfoperator= (FunctionSpaceElement< FEFunctionSpace< BoundaryMapper< DomainMapper, Scalar, typename Space::GridView > >, m > const &bf)
 Assignment from an FE function which is restricted to the boundary to a suitable FE function. More...
 
template<class Mapper >
std::enable_if_t< std::is_same< Mapper, FunctionSpace_Detail::ChooseDomainMapper< typename Space::Mapper > >::value, Self & > operator= (FunctionSpaceElement< FEFunctionSpace< Mapper >, m > const &f)
 Assignment from a suitable FE function to an FE function which is restricted to the boundary. More...
 
template<class Functor >
Selfoperator= (WeakFunctionViewAdaptor< Functor > const &f)
 Assignment from a weak function view adaptor. More...
 
template<class Space , class Functor >
Selfoperator= (FunctionViewAdaptor< Space, Functor > const &f)
 Assignment from a function view adaptor. More...
 
Selfoperator= (StorageType const &v)
 Assignment from vector of coefficients. More...
 
Selfoperator= (StorageValueType const &a)
 Sets all coefficients to the given value. As StorageValueType is a vector one value for each component can be specified. More...
 
Selfoperator= (Scalar val)
 Sets all coefficients to the given value. More...
 
Basic linear algebra
Selfoperator+= (Self const &f)
 Addition of a function from the same (type of) space. Note that despite the same type the function spaces may actually be different (e.g. different polynomial order). In this case, interpolation is employed. More...
 
template<class OtherSpace >
Selfoperator+= (FunctionSpaceElement< OtherSpace, m > const &f)
 Addition of a function from a different function space, but with the same number of components. More...
 
Selfoperator+= (StorageType const &v)
 Addition of coefficient vectors. More...
 
Selfoperator+= (Scalar val)
 Addition of a constant value to the coefficient vector. Note that this does in general not add a constant function, except for the notable case of Lagrange finite elements. More...
 
Selfoperator-= (Self const &f)
 Subtraction of a function from the same (type of) space. More...
 
template<class OtherSpace >
Selfoperator-= (FunctionSpaceElement< OtherSpace, m > const &f)
 Subtraction of a function from a different function space, but with the same number of components. More...
 
Selfoperator-= (StorageType const &v)
 Subtraction of coefficient vectors. More...
 
Selfaxpy (Scalar a, Self const &f)
 Computes \( y \leftarrow ax + y \) for a function x from the same (type of) space. Note that despite the same type the function spaces may actually be different (e.g. different polynomial order). In this case, interpolation is employed. More...
 
template<class OtherSpace >
Selfaxpy (Scalar a, FunctionSpaceElement< OtherSpace, m > const &f)
 Computes \( y \leftarrow ax + y \) for a function x from a different function space, but with the same number of components. More...
 
Selfaxpy (Scalar a, StorageType const &v)
 Computes \( y \leftarrow av + y \) for coefficient vectors. More...
 
Selfoperator*= (Scalar a)
 Scaling. More...
 
Function and derivative evaluation
ValueType value (typename Space::Evaluator const &evaluator) const
 Evaluates the FE function at the position used to construct the evaluator.
More...
 
ValueType value (Cell< GridView > const &cell, LocalPosition< GridView > const &local) const
 Evaluates the FE function at the specified local position in the given cell. More...
 
ValueType value (GlobalPosition< GridView > const &x) const
 Evaluates the FE function at the specified global position. More...
 
DerivativeType derivative (typename FunctionSpace::Evaluator const &evaluator) const
 Evaluates the derivative of the FE function wrt global coordinates at the position used to construct the evaluator.
More...
 
DerivativeType derivative (Cell< GridView > const &cell, LocalPosition< GridView > const &local) const
 Evaluates the derivative of the FE function at the given local position inside the given cell.
More...
 
DerivativeType derivative (GlobalPosition< GridView > const &global) const
 Evaluates the FE function's gradient at the specified global position. More...
 
HessianType hessian (typename FunctionSpace::Evaluator const &evaluator) const
 Evaluates the Hessian of the FE function at the position used to construct the evaluator.
More...
 
HessianType hessian (Cell< GridView > const &cell, LocalPosition< GridView > const &local) const
 Evaluates the Hessian of the FE function at the given local position inside the given cell.
More...
 
HessianType hessian (GlobalPosition< GridView > const &global) const
 Evaluates the FE function's Hessian at the specified global position. More...
 
Access to coefficient vector
StorageType const & coefficients () const
 Provides read-only access to the coefficients of the finite element basis functions. More...
 
StorageTypecoefficients ()
 Provides access to the coefficients of the finite element basis functions. More...
 
StorageValueTypeoperator[] (size_t i)
 EXPERIMENTAL Provides access to the coefficients of the finite element basis functions. More...
 
StorageValueType const & operator[] (size_t i) const
 EXPERIMENTAL Provides read-only access to the coefficients of the finite element basis functions. More...
 
Miscellaneous
int dim () const
 Number of scalar degrees of freedom. More...
 
int size () const
 Number of coefficient blocks. More...
 
Space const & space () const
 Provides access to the FEFunctionSpace to which this function belongs. More...
 
int order (Cell< GridView > const &cell) const
 The polynomial order of the function on the provided cell. More...
 
template<class SFS >
int order (SFS const &sfs) const
 
void transfer (typename TransferData< Space >::TransferMatrix const &transferMatrix)
 

Static Public Attributes

static int const components = Space::sfComponents * m
 components at each point of evaluation More...
 

Protected Attributes

Dune::BlockVector< StorageValueTypedata
 
Space const * sp
 

Related Functions

(Note that these are not member functions.)

template<class Space , int m, class OutIter >
OutIter vectorToSequence (FunctionSpaceElement< Space, m > const &v, OutIter iter)
 Writes the coefficients into a flat sequence. <Space,m> More...
 

Member Typedef Documentation

◆ block_type

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::block_type = StorageValueType

type of coefficient vector entries

This is provided for compatibility with Dune::BlockVector.

Definition at line 351 of file functionspace.hh.

◆ ctype

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::ctype = typename GridView::ctype

the scalar type used for spatial coordinates

Definition at line 324 of file functionspace.hh.

◆ DerivativeType

template<class FunctionSpace , int m>
typedef Dune::FieldMatrix<Scalar,components,Space::dimworld> Kaskade::FunctionSpaceElement< FunctionSpace, m >::DerivativeType

return type of the function derivative (...)

The derivative is always computed wrt the global world coordinate system, even for functions defined on a lower dimensional manifold (e.g. in shell or plate problems). In that case, the function is extended to a neighborhood of the manifold in the canonical way - constant extension normal to the manifold. The derivative of that extension wrt the global world coordinate system is then returned. Hence, the number of entries in the derivative is always the world dimension.

Design rationale: The alternative would be to consider the derivative wrt a coordinate system of the tangent space. The choice of such a coordinate system is, however, not unique.

Definition at line 374 of file functionspace.hh.

◆ field_type

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::field_type = typename Space::field_type

scalar field type

This is provided for consistency with Dune::BlockVector.

Definition at line 336 of file functionspace.hh.

◆ GridView

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::GridView = typename Space::GridView

the grid view this function is defined on

Definition at line 321 of file functionspace.hh.

◆ HessianType

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::HessianType = Tensor<Scalar,components,Space::dimworld,Space::dimworld>

return type of the function hessian(...)

Definition at line 377 of file functionspace.hh.

◆ Scalar

template<class FunctionSpace , int m>
typedef Space::Scalar Kaskade::FunctionSpaceElement< FunctionSpace, m >::Scalar

scalar field type

Definition at line 329 of file functionspace.hh.

◆ Self

template<class FunctionSpace , int m>
typedef FunctionSpaceElement<Space, m> Kaskade::FunctionSpaceElement< FunctionSpace, m >::Self

own type

Definition at line 327 of file functionspace.hh.

◆ SfValueType

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::SfValueType = typename Space::SfValueType

type of shape function values

Definition at line 344 of file functionspace.hh.

◆ Space

template<class FunctionSpace , int m>
using Kaskade::FunctionSpaceElement< FunctionSpace, m >::Space = FunctionSpace

FEFunctionSpace, this function is an element of.

Definition at line 318 of file functionspace.hh.

◆ StorageType

template<class FunctionSpace , int m>
typedef Dune::BlockVector<StorageValueType> Kaskade::FunctionSpaceElement< FunctionSpace, m >::StorageType

type of the data vector

Definition at line 354 of file functionspace.hh.

◆ StorageValueType

template<class FunctionSpace , int m>
typedef Dune::FieldVector<Scalar,m> Kaskade::FunctionSpaceElement< FunctionSpace, m >::StorageValueType

type of the elements of the data vector

Definition at line 339 of file functionspace.hh.

◆ ValueType

template<class FunctionSpace , int m>
typedef Dune::FieldVector<Scalar,components> Kaskade::FunctionSpaceElement< FunctionSpace, m >::ValueType

return type of the function value(...)

Definition at line 360 of file functionspace.hh.

Constructor & Destructor Documentation

◆ FunctionSpaceElement() [1/2]

template<class FunctionSpace , int m>
Kaskade::FunctionSpaceElement< FunctionSpace, m >::FunctionSpaceElement ( Space const &  fs)
inline

Constructor.

The function is initialized to zero. Space is a FEFunctionSpace.

Definition at line 389 of file functionspace.hh.

◆ FunctionSpaceElement() [2/2]

template<class FunctionSpace , int m>
Kaskade::FunctionSpaceElement< FunctionSpace, m >::FunctionSpaceElement ( Self const &  fse)
inline

Copy constructor.

Definition at line 397 of file functionspace.hh.

◆ ~FunctionSpaceElement()

template<class FunctionSpace , int m>
Kaskade::FunctionSpaceElement< FunctionSpace, m >::~FunctionSpaceElement ( )
inline

Definition at line 403 of file functionspace.hh.

Member Function Documentation

◆ axpy() [1/3]

template<class FunctionSpace , int m>
template<class OtherSpace >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::axpy ( Scalar  a,
FunctionSpaceElement< OtherSpace, m > const &  f 
)
inline

Computes \( y \leftarrow ax + y \) for a function x from a different function space, but with the same number of components.

Here, interpolation is used.

Definition at line 720 of file functionspace.hh.

◆ axpy() [2/3]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::axpy ( Scalar  a,
Self const &  f 
)
inline

Computes \( y \leftarrow ax + y \) for a function x from the same (type of) space. Note that despite the same type the function spaces may actually be different (e.g. different polynomial order). In this case, interpolation is employed.

Definition at line 704 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::axpy().

◆ axpy() [3/3]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::axpy ( Scalar  a,
StorageType const &  v 
)
inline

Computes \( y \leftarrow av + y \) for coefficient vectors.

Definition at line 729 of file functionspace.hh.

◆ coefficients() [1/2]

template<class FunctionSpace , int m>
StorageType & Kaskade::FunctionSpaceElement< FunctionSpace, m >::coefficients ( )
inline

Provides access to the coefficients of the finite element basis functions.

Definition at line 924 of file functionspace.hh.

◆ coefficients() [2/2]

template<class FunctionSpace , int m>
StorageType const & Kaskade::FunctionSpaceElement< FunctionSpace, m >::coefficients ( ) const
inline

◆ derivative() [1/3]

template<class FunctionSpace , int m>
DerivativeType Kaskade::FunctionSpaceElement< FunctionSpace, m >::derivative ( Cell< GridView > const &  cell,
LocalPosition< GridView > const &  local 
) const
inline

Evaluates the derivative of the FE function at the given local position inside the given cell.

On codimension 1 entitites, the value is undefined if Space::continuity<1. The components are sorted as described for value().

Definition at line 830 of file functionspace.hh.

◆ derivative() [2/3]

template<class FunctionSpace , int m>
DerivativeType Kaskade::FunctionSpaceElement< FunctionSpace, m >::derivative ( GlobalPosition< GridView > const &  global) const
inline

Evaluates the FE function's gradient at the specified global position.

Warning
This method is comparatively DEAD SLOW (it has to locate the cell in which the position is).

Definition at line 845 of file functionspace.hh.

◆ derivative() [3/3]

template<class FunctionSpace , int m>
DerivativeType Kaskade::FunctionSpaceElement< FunctionSpace, m >::derivative ( typename FunctionSpace::Evaluator const &  evaluator) const
inline

Evaluates the derivative of the FE function wrt global coordinates at the position used to construct the evaluator.

On codimension 1 entitites, the value is undefined if Space::continuity<1. The components are sorted as described for value().

Definition at line 807 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::derivative().

◆ dim()

template<class FunctionSpace , int m>
int Kaskade::FunctionSpaceElement< FunctionSpace, m >::dim ( ) const
inline

Number of scalar degrees of freedom.

Definition at line 946 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::size().

◆ hessian() [1/3]

template<class FunctionSpace , int m>
HessianType Kaskade::FunctionSpaceElement< FunctionSpace, m >::hessian ( Cell< GridView > const &  cell,
LocalPosition< GridView > const &  local 
) const
inline

Evaluates the Hessian of the FE function at the given local position inside the given cell.

On codimension 1 entitites, the value is undefined if Space::continuity<2. The components are sorted as described for value().

Definition at line 886 of file functionspace.hh.

◆ hessian() [2/3]

template<class FunctionSpace , int m>
HessianType Kaskade::FunctionSpaceElement< FunctionSpace, m >::hessian ( GlobalPosition< GridView > const &  global) const
inline

Evaluates the FE function's Hessian at the specified global position.

Warning
This method is comparatively DEAD SLOW.

Definition at line 901 of file functionspace.hh.

◆ hessian() [3/3]

template<class FunctionSpace , int m>
HessianType Kaskade::FunctionSpaceElement< FunctionSpace, m >::hessian ( typename FunctionSpace::Evaluator const &  evaluator) const
inline

Evaluates the Hessian of the FE function at the position used to construct the evaluator.

On codimension 1 entitites, the value is undefined if Space::continuity<2. The components are sorted as described for value().

Definition at line 859 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::hessian().

◆ operator*=()

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator*= ( Scalar  a)
inline

Scaling.

Definition at line 732 of file functionspace.hh.

◆ operator+=() [1/4]

template<class FunctionSpace , int m>
template<class OtherSpace >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator+= ( FunctionSpaceElement< OtherSpace, m > const &  f)
inline

Addition of a function from a different function space, but with the same number of components.

Here, interpolation is used.

Definition at line 638 of file functionspace.hh.

◆ operator+=() [2/4]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator+= ( Scalar  val)
inline

Addition of a constant value to the coefficient vector. Note that this does in general not add a constant function, except for the notable case of Lagrange finite elements.

Definition at line 654 of file functionspace.hh.

◆ operator+=() [3/4]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator+= ( Self const &  f)
inline

Addition of a function from the same (type of) space. Note that despite the same type the function spaces may actually be different (e.g. different polynomial order). In this case, interpolation is employed.

Definition at line 622 of file functionspace.hh.

◆ operator+=() [4/4]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator+= ( StorageType const &  v)
inline

Addition of coefficient vectors.

Definition at line 647 of file functionspace.hh.

◆ operator-=() [1/3]

template<class FunctionSpace , int m>
template<class OtherSpace >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator-= ( FunctionSpaceElement< OtherSpace, m > const &  f)
inline

Subtraction of a function from a different function space, but with the same number of components.

Here, interpolation of the given function to our function space is used.

Definition at line 683 of file functionspace.hh.

◆ operator-=() [2/3]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator-= ( Self const &  f)
inline

Subtraction of a function from the same (type of) space.

Note that despite the same type the function spaces may actually be different (e.g. different polynomial order). In this case, interpolation is employed.

Definition at line 667 of file functionspace.hh.

◆ operator-=() [3/3]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator-= ( StorageType const &  v)
inline

Subtraction of coefficient vectors.

Definition at line 692 of file functionspace.hh.

◆ operator=() [1/10]

template<class FunctionSpace , int m>
template<template< class, class > class DomainMapper>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( FunctionSpaceElement< FEFunctionSpace< BoundaryMapper< DomainMapper, Scalar, typename Space::GridView > >, m > const &  bf)
inline

Assignment from an FE function which is restricted to the boundary to a suitable FE function.

If the underlying mappers are the same for both FE function, assignment can be done fast by simply copying the coefficients. Otherwise interpolation will be used. (TODO: Currently interpolateGlobally with a boundary FE function as second parameter does not work as it should. Therefore, before doing the interpolation the boundary FE function is assigned to a suitable FE function, which is then used for interpolation.)

Definition at line 484 of file functionspace.hh.

◆ operator=() [2/10]

template<class FunctionSpace , int m>
template<class Mapper >
std::enable_if_t< std::is_same< Mapper, FunctionSpace_Detail::ChooseDomainMapper< typename Space::Mapper > >::value, Self & > Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( FunctionSpaceElement< FEFunctionSpace< Mapper >, m > const &  f)
inline

Assignment from a suitable FE function to an FE function which is restricted to the boundary.

If the underlying mappers are the same for both FE function, assignment can be done fast by simply copying the coefficients. Otherwise interpolation will be used. This overload will only be available if this FE function belongs to a boundary space (otherwise ChooseDomainMapper yields void which is not the same as Mapper and type substitution will fail).

Definition at line 519 of file functionspace.hh.

◆ operator=() [3/10]

template<class FunctionSpace , int m>
template<class OtherSpace >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( FunctionSpaceElement< OtherSpace, m > const &  f)
inline

Assignment from functions belonging to other function spaces.

This assignment is defined by interpolation and hence results in pointwise equality only if the right hand side f lies inside the function space represented by our space. In cases where the coefficients of the involved spaces are just subsets (either way round), the interpolation approach is probably of suboptimal performance.

Definition at line 467 of file functionspace.hh.

◆ operator=() [4/10]

template<class FunctionSpace , int m>
template<class Space , class Functor >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( FunctionViewAdaptor< Space, Functor > const &  f)
inline

Assignment from a function view adaptor.

Template Parameters
Functorthe actual callable defining the function view

Usage example doubling an existing FE function:

f = makeFunctionView(space,[&f](auto const& eval)
{
return 2*f.value(eval);
});
Space const & space() const
Provides access to the FEFunctionSpace to which this function belongs.

Definition at line 575 of file functionspace.hh.

◆ operator=() [5/10]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( Scalar  val)
inline

Sets all coefficients to the given value.

Note that this does in general not result in a constant function (a notable exception being Lagrange elements, where this property indeed holds).

Definition at line 607 of file functionspace.hh.

◆ operator=() [6/10]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( Self &&  fse)
inline

Move assignment.

Note that even if the type of function (and hence of the function space) is the same, the actual function spaces may be different (e.g., by different polynomial order). If different spaces are involved, assignment is just interpolation.

Definition at line 444 of file functionspace.hh.

◆ operator=() [7/10]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( Self const &  fse)
inline

Copy assignment.

Note that even if the type of function (and hence of the function space) is the same, the actual function spaces may be different (e.g., by different polynomial order). If different spaces are involved, assignment is just interpolation.

Definition at line 425 of file functionspace.hh.

◆ operator=() [8/10]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( StorageType const &  v)
inline

Assignment from vector of coefficients.

The size of the provided coefficient vector has to match the size of the current coefficient vector.

Definition at line 586 of file functionspace.hh.

◆ operator=() [9/10]

template<class FunctionSpace , int m>
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( StorageValueType const &  a)
inline

Sets all coefficients to the given value. As StorageValueType is a vector one value for each component can be specified.

Note that this does in general not result in a constant function (a notable exception being Lagrange elements, where this property indeed holds).

Definition at line 599 of file functionspace.hh.

◆ operator=() [10/10]

template<class FunctionSpace , int m>
template<class Functor >
Self & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator= ( WeakFunctionViewAdaptor< Functor > const &  f)
inline

Assignment from a weak function view adaptor.

Template Parameters
Functorthe actual callable defining the function view

Usage example creating a conic function:

f = makeWeakFunctionView([&](auto const& cell, auto const& xi)
{
auto x = cell.geometry().global(xi);
return Dune::FieldVector<double,1>(x.two_norm());
});

Definition at line 555 of file functionspace.hh.

◆ operator[]() [1/2]

template<class FunctionSpace , int m>
StorageValueType & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator[] ( size_t  i)
inline

EXPERIMENTAL Provides access to the coefficients of the finite element basis functions.

Definition at line 929 of file functionspace.hh.

◆ operator[]() [2/2]

template<class FunctionSpace , int m>
StorageValueType const & Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator[] ( size_t  i) const
inline

EXPERIMENTAL Provides read-only access to the coefficients of the finite element basis functions.

Definition at line 934 of file functionspace.hh.

◆ order() [1/2]

template<class FunctionSpace , int m>
int Kaskade::FunctionSpaceElement< FunctionSpace, m >::order ( Cell< GridView > const &  cell) const
inline

The polynomial order of the function on the provided cell.

Definition at line 964 of file functionspace.hh.

◆ order() [2/2]

template<class FunctionSpace , int m>
template<class SFS >
int Kaskade::FunctionSpaceElement< FunctionSpace, m >::order ( SFS const &  sfs) const
inline

Definition at line 972 of file functionspace.hh.

◆ size()

template<class FunctionSpace , int m>
int Kaskade::FunctionSpaceElement< FunctionSpace, m >::size ( ) const
inline

Number of coefficient blocks.

This is dim()/m. Method is added for interface consistency with Dune::BlockVector.

Definition at line 954 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator=().

◆ space()

template<class FunctionSpace , int m>
Space const & Kaskade::FunctionSpaceElement< FunctionSpace, m >::space ( ) const
inline

◆ transfer()

template<class FunctionSpace , int m>
void Kaskade::FunctionSpaceElement< FunctionSpace, m >::transfer ( typename TransferData< Space >::TransferMatrix const &  transferMatrix)
inline

Group of members that are necessary for automatic grid transfer after refinement. The member transferMe is connected to a boost::signal in function space automatically at construction

Definition at line 982 of file functionspace.hh.

◆ value() [1/3]

template<class FunctionSpace , int m>
ValueType Kaskade::FunctionSpaceElement< FunctionSpace, m >::value ( Cell< GridView > const &  cell,
LocalPosition< GridView > const &  local 
) const
inline

Evaluates the FE function at the specified local position in the given cell.

If several values on the same cell are to be calculated and performance is important, consider creating an evaluator and compute the function's values using the evaluator.

Definition at line 776 of file functionspace.hh.

◆ value() [2/3]

template<class FunctionSpace , int m>
ValueType Kaskade::FunctionSpaceElement< FunctionSpace, m >::value ( GlobalPosition< GridView > const &  x) const
inline

Evaluates the FE function at the specified global position.

\WARNING This method is comparatively DEAD SLOW.
The components are sorted as described for value().

Definition at line 792 of file functionspace.hh.

◆ value() [3/3]

template<class FunctionSpace , int m>
ValueType Kaskade::FunctionSpaceElement< FunctionSpace, m >::value ( typename Space::Evaluator const &  evaluator) const
inline

Evaluates the FE function at the position used to construct the evaluator.

On codimension >0 entitites, the value is undefined if Space::continuity<0. For vectorial shape functions and m>1, the components are grouped by shape function component first. Example with shape functions with 3 components and m==2: [sf0m0,sf1m0,sf2m0,sf0m1,sf1m1,sf2m1]

Definition at line 752 of file functionspace.hh.

Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator=(), and Kaskade::FunctionSpaceElement< FunctionSpace, m >::value().

Member Data Documentation

◆ components

template<class FunctionSpace , int m>
int const Kaskade::FunctionSpaceElement< FunctionSpace, m >::components = Space::sfComponents * m
static

components at each point of evaluation

Definition at line 357 of file functionspace.hh.

◆ data

template<class FunctionSpace , int m>
Dune::BlockVector<StorageValueType> Kaskade::FunctionSpaceElement< FunctionSpace, m >::data
protected

◆ sp

template<class FunctionSpace , int m>
Space const* Kaskade::FunctionSpaceElement< FunctionSpace, m >::sp
protected

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