DUNE-FEM (unstable)
interface for local functions More...
#include <dune/fem/function/localfunction/localfunction.hh>
Public Types | |
| typedef BasisFunctionSet | BasisFunctionSetType |
| type of basis function set | |
| typedef LocalDofVector | LocalDofVectorType |
| type of local Dof Vector | |
| typedef LocalDofVectorType::value_type | DofType |
| type of DoF use with the discrete function | |
| typedef LocalDofVectorType::size_type | SizeType |
| type of index | |
| typedef BasisFunctionSetType::EntityType | EntityType |
| type of the entity, the local function lives on is given by the space | |
| typedef EntityType::Geometry | Geometry |
| type of the geometry, the local function lives on is given by the space | |
| typedef BasisFunctionSetType::FunctionSpaceType | FunctionSpaceType |
| type of functionspace | |
| typedef FunctionSpaceType::DomainFieldType | DomainFieldType |
| field type of the domain | |
| typedef FunctionSpaceType::RangeFieldType | RangeFieldType |
| field type of the range | |
| typedef FunctionSpaceType::DomainType | DomainType |
| type of domain vectors, i.e., type of coordinates | |
| typedef FunctionSpaceType::RangeType | RangeType |
| type of range vectors, i.e., type of function values | |
| typedef FunctionSpaceType::JacobianRangeType | JacobianRangeType |
| type of the Jacobian, i.e., type of evaluated Jacobian matrix | |
| typedef FunctionSpaceType::HessianRangeType | HessianRangeType |
| type of the Hessian | |
| typedef Geometry::LocalCoordinate | LocalCoordinateType |
| type of local coordinates | |
| typedef Traits::derived_type | derived_type |
| type of derived vector class | |
| typedef Traits::value_type | value_type |
| export the type representing the field | |
| typedef FieldTraits< value_type >::field_type | field_type |
| export the type representing the field | |
| typedef Traits::value_type | block_type |
| export the type representing the components | |
| typedef Traits::size_type | size_type |
| The type used for the index access and size operation. | |
| typedef DenseIterator< DenseVector, value_type > | Iterator |
| Iterator class for sequential access. | |
| typedef Iterator | iterator |
| typedef for stl compliant access | |
| typedef DenseIterator< const DenseVector, const value_type > | ConstIterator |
| ConstIterator class for sequential access. | |
| typedef ConstIterator | const_iterator |
| typedef for stl compliant access | |
Public Member Functions | |
| LocalFunction () | |
| default constructor, calls default ctor of BasisFunctionSetType and LocalDofVectorType | |
| LocalFunction (const BasisFunctionSetType &basisFunctionSet) | |
| ctor taking a basisFunctionSet, calling default ctor for LocalDofVectorType, and resize | |
| LocalFunction (const LocalDofVectorType &localDofVector) | |
| ctor taking a localDofVector, calling default ctor for BasisFunctionSetType | |
| LocalFunction (const BasisFunctionSetType &basisFunctionSet, const LocalDofVector &localDofVector) | |
| copy given agruments | |
| LocalFunction (LocalDofVectorType &&localDofVector) | |
| half move ctor | |
| LocalFunction (const BasisFunctionSetType &basisFunctionSet, LocalDofVector &&localDofVector) | |
| half move ctor | |
| LocalFunction (ThisType &&other) | |
| move constructor | |
| LocalFunction (const ThisType &other) | |
| copy constructor | |
| const DofType & | operator[] (SizeType num) const |
| access to local dofs (read-only) More... | |
| DofType & | operator[] (SizeType num) |
| access to local dofs (read-write) More... | |
| template<class T > | |
| void | assign (const LocalFunction< BasisFunctionSet, T > &other) |
| assign all DoFs of this local function More... | |
| void | clear () |
| set all DoFs to zero | |
| template<class PointType > | |
| void | axpy (const PointType &x, const RangeType &factor) |
| axpy operation for local function More... | |
| template<class PointType > | |
| void | axpy (const PointType &x, const JacobianRangeType &factor) |
| axpy operation for local function More... | |
| template<class PointType > | |
| void | axpy (const PointType &x, const RangeType &factor1, const JacobianRangeType &factor2) |
| axpy operation for local function More... | |
| int | order () const |
| obtain the order of this local function More... | |
| const BasisFunctionSetType & | basisFunctionSet () const |
| obtain the basis function set for this local function More... | |
| const EntityType & | entity () const |
| obtain the entity, this local function lives on More... | |
| const Geometry & | geometry () const |
| obtain the geometry, this local function lives on More... | |
| template<class PointType > | |
| void | evaluate (const PointType &x, RangeType &ret) const |
| evaluate the local function More... | |
| template<class PointType > | |
| void | jacobian (const PointType &x, JacobianRangeType &ret) const |
| evaluate Jacobian of the local function More... | |
| template<class PointType > | |
| void | hessian (const PointType &x, HessianRangeType &ret) const |
| evaluate Hessian of the local function More... | |
| int | numDofs () const |
| obtain the number of local DoFs More... | |
| SizeType | size () const |
| obtain the number of local DoFs More... | |
| template<class QuadratureType , class ... Vectors> | |
| void | axpyQuadrature (const QuadratureType &quad, const Vectors &... values) |
| evaluate all basisfunctions for all quadrature points, multiply with the given factor and add the result to the local coefficients | |
| template<class QuadratureType , class RangeVectorType , class JacobianRangeVectorType > | |
| void | axpyQuadrature (const QuadratureType &quad, const RangeVectorType &rangeVector, const JacobianRangeVectorType &jacobianVector) |
| evaluate all basisfunctions for all quadrature points, multiply with the given factor and add the result to the local coefficients | |
| template<class QuadratureType , class ... Vectors> | |
| void | evaluateQuadrature (const QuadratureType &quad, Vectors &... vec) const |
| evaluate all basisfunctions for all quadrature points and store the results in the result vector | |
| template<class QuadratureType , class ... Vectors> | |
| void | jacobianQuadrature (const QuadratureType &quad, Vectors &... vec) const |
| evaluate all Jacobians for all basis functions for all quadrature points and store the results in the result vector | |
| template<class QuadratureType , class ... Vectors> | |
| void | hessianQuadrature (const QuadratureType &quad, Vectors &... vec) const |
| evaluate all hessians of all basis functions for all quadrature points and store the results in the result vector | |
| const LocalDofVectorType & | localDofVector () const |
| return const reference to local Dof Vector | |
| LocalDofVectorType & | localDofVector () |
| return mutable reference to local Dof Vector | |
| bool | valid () const |
| Returns true if local function if bind or init was previously called. | |
| constexpr value_type & | operator[] (size_type i) |
| random access | |
| constexpr value_type & | front () |
| return reference to first element | |
| constexpr const value_type & | front () const |
| return reference to first element | |
| constexpr value_type & | back () |
| return reference to last element | |
| constexpr const value_type & | back () const |
| return reference to last element | |
| constexpr bool | empty () const |
| checks whether the container is empty | |
| constexpr Iterator | begin () |
| begin iterator | |
| constexpr ConstIterator | begin () const |
| begin ConstIterator | |
| constexpr Iterator | end () |
| end iterator | |
| constexpr ConstIterator | end () const |
| end ConstIterator | |
| constexpr Iterator | beforeEnd () |
| constexpr ConstIterator | beforeEnd () const |
| constexpr Iterator | beforeBegin () |
| constexpr ConstIterator | beforeBegin () const |
| constexpr Iterator | find (size_type i) |
| return iterator to given element or end() | |
| constexpr ConstIterator | find (size_type i) const |
| return iterator to given element or end() | |
| constexpr derived_type & | operator+= (const DenseVector< Other > &x) |
| vector space addition | |
| constexpr std::enable_if< std::is_convertible< ValueType, value_type >::value, derived_type >::type & | operator+= (const ValueType &kk) |
| vector space add scalar to all comps More... | |
| constexpr derived_type & | operator-= (const DenseVector< Other > &x) |
| vector space subtraction | |
| constexpr std::enable_if< std::is_convertible< ValueType, value_type >::value, derived_type >::type & | operator-= (const ValueType &kk) |
| vector space subtract scalar from all comps More... | |
| constexpr derived_type | operator+ (const DenseVector< Other > &b) const |
| Binary vector addition. | |
| constexpr derived_type | operator- (const DenseVector< Other > &b) const |
| Binary vector subtraction. | |
| constexpr derived_type | operator- () const |
| Vector negation. | |
| constexpr std::enable_if< std::is_convertible< FieldType, field_type >::value, derived_type >::type & | operator*= (const FieldType &kk) |
| vector space multiplication with scalar More... | |
| constexpr std::enable_if< std::is_convertible< FieldType, field_type >::value, derived_type >::type & | operator/= (const FieldType &kk) |
| vector space division by scalar More... | |
| constexpr bool | operator== (const DenseVector< Other > &x) const |
| Binary vector comparison. | |
| constexpr bool | operator!= (const DenseVector< Other > &x) const |
| Binary vector incomparison. | |
| constexpr derived_type & | axpy (const field_type &a, const DenseVector< Other > &x) |
| vector space axpy operation ( *this += a x ) | |
| constexpr PromotionTraits< field_type, typenameDenseVector< Other >::field_type >::PromotedType | operator* (const DenseVector< Other > &x) const |
| indefinite vector dot product \(\left (x^T \cdot y \right)\) which corresponds to Petsc's VecTDot More... | |
| constexpr PromotionTraits< field_type, typenameDenseVector< Other >::field_type >::PromotedType | dot (const DenseVector< Other > &x) const |
| vector dot product \(\left (x^H \cdot y \right)\) which corresponds to Petsc's VecDot More... | |
| constexpr FieldTraits< value_type >::real_type | one_norm () const |
| one norm (sum over absolute values of entries) | |
| constexpr FieldTraits< value_type >::real_type | one_norm_real () const |
| simplified one norm (uses Manhattan norm for complex values) | |
| constexpr FieldTraits< value_type >::real_type | two_norm () const |
| two norm sqrt(sum over squared values of entries) | |
| constexpr FieldTraits< value_type >::real_type | two_norm2 () const |
| square of two norm (sum over squared values of entries), need for block recursion | |
| constexpr FieldTraits< vt >::real_type | infinity_norm () const |
| infinity norm (maximum of absolute values of entries) | |
| constexpr FieldTraits< vt >::real_type | infinity_norm () const |
| infinity norm (maximum of absolute values of entries) | |
| constexpr FieldTraits< vt >::real_type | infinity_norm_real () const |
| simplified infinity norm (uses Manhattan norm for complex values) | |
| constexpr FieldTraits< vt >::real_type | infinity_norm_real () const |
| simplified infinity norm (uses Manhattan norm for complex values) | |
| constexpr size_type | N () const |
| number of blocks in the vector (are of size 1 here) | |
| constexpr size_type | dim () const |
| dimension of the vector space | |
Protected Member Functions | |
| void | init (const EntityType &entity) |
| initialize the local function for an entity More... | |
| void | bind (const EntityType &entity) |
| initialize the local function for an entity More... | |
| void | unbind () |
| clears the local function by removing the basisFunctionSet | |
| void | init (const BasisFunctionSetType &basisFunctionSet) |
| initialize the local function for an basisFunctionSet More... | |
Related Functions | |
(Note that these are not member functions.) | |
| std::ostream & | operator<< (std::ostream &s, const DenseVector< LocalFunction< BasisFunctionSet, LocalDofVector > > &v) |
| Write a DenseVector to an output stream. More... | |
Detailed Description
class Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >
interface for local functions
Local functions are used to represend a discrete function on one entity. The LocalFunctionInterface defines the functionality that can be expected from such a local function.
Member Function Documentation
◆ assign()
|
inline |
assign all DoFs of this local function
- Parameters
-
[in] lf local function to assign DoFs from
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
Referenced by Dune::Fem::ProlongFunction< LRP >::operator()(), Dune::Fem::RestrictFunction< LRP >::operator()(), and Dune::Fem::hpDG::DefaultDataProjection< DiscreteFunction >::operator()().
◆ axpy() [1/3]
|
inline |
axpy operation for local function
Denoting the DoFs of the local function by \(u_i\) and the basis functions by \(\varphi_i\), this function performs the following operation:
\[ u_i = u_i + factor \cdot \nabla\varphi_i( x ) \]
- Parameters
-
[in] x point to evaluate jacobian of basis functions in [in] factor axpy factor
References Dune::Fem::BasisFunctionSet< Entity, Range >::axpy(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
◆ axpy() [2/3]
|
inline |
axpy operation for local function
Denoting the DoFs of the local function by \(u_i\) and the basis functions by \(\varphi_i\), this function performs the following operation:
\[ u_i = u_i + factor \cdot \varphi_i( x ) \]
- Parameters
-
[in] x point to evaluate basis functions in [in] factor axpy factor
References Dune::Fem::BasisFunctionSet< Entity, Range >::axpy(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
Referenced by Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::axpyQuadrature().
◆ axpy() [3/3]
|
inline |
axpy operation for local function
Denoting the DoFs of the local function by \(u_i\) and the basis functions by \(\varphi_i\), this function performs the following operation:
\[ u_i = u_i + factor1 \cdot \varphi_i( x ) + factor2 \cdot \nabla\varphi_i( x ) \]
- Parameters
-
[in] x point to evaluate basis functions in [in] factor1 axpy factor for \(\varphi( x )\) [in] factor2 axpy factor for \(\nabla\varphi( x )\)
References Dune::Fem::BasisFunctionSet< Entity, Range >::axpy(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
◆ basisFunctionSet()
|
inline |
obtain the basis function set for this local function
- Returns
- reference to the basis function set
Referenced by Dune::Fem::LocalMassMatrixImplementationDgOrthoNormal< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse(), Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::axpy(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::axpyQuadrature(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::entity(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::evaluate(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::geometry(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::hessian(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::init(), Dune::Fem::ConstLocalDiscreteFunction< DiscreteFunction >::init(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::jacobian(), Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::LocalFunction(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::order().
◆ beforeBegin() [1/2]
|
inlineconstexprinherited |
- Returns
- an iterator that is positioned before the first entry of the vector.
◆ beforeBegin() [2/2]
|
inlineconstexprinherited |
- Returns
- an iterator that is positioned before the first entry of the vector.
◆ beforeEnd() [1/2]
|
inlineconstexprinherited |
- Returns
- an iterator that is positioned before the end iterator of the vector, i.e. at the last entry.
◆ beforeEnd() [2/2]
|
inlineconstexprinherited |
- Returns
- an iterator that is positioned before the end iterator of the vector. i.e. at the last element
◆ bind()
|
inlineprotected |
initialize the local function for an entity
Binds the local function to an basisFunctionSet and entity.
- Note
- Must be overloaded on the derived implementation class.
- Parameters
-
[in] entity to bind the local function to
References DUNE_THROW.
◆ dot()
|
inlineconstexprinherited |
vector dot product \(\left (x^H \cdot y \right)\) which corresponds to Petsc's VecDot
http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecDot.html
- Parameters
-
x other vector
- Returns
◆ entity()
|
inline |
obtain the entity, this local function lives on
- Returns
- reference to the entity
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::BasisFunctionSet< Entity, Range >::entity().
Referenced by Dune::Fem::LocalOrthonormalL2Projection< GridPart, BasisFunctionSet, Quadrature >::apply(), Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse(), Dune::Fem::LocalMassMatrixImplementationDgOrthoNormal< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverse(), Dune::Fem::ConstLocalDiscreteFunction< DiscreteFunction >::ConstLocalDiscreteFunction(), Dune::Fem::ConstLocalDiscreteFunction< DiscreteFunction >::init(), Dune::Fem::RestrictProlongDefault< DiscreteFunction >::restrictFinalize(), and Dune::Fem::RestrictProlongDefault< DiscreteFunction >::restrictLocal().
◆ evaluate()
|
inline |
evaluate the local function
- Parameters
-
[in] x evaluation point in local coordinates [out] ret value of the function in the given point
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), Dune::Fem::BasisFunctionSet< Entity, Range >::evaluateAll(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
Referenced by Dune::Fem::HdivProjection< DiscreteFunctionType >::normalJump().
◆ geometry()
|
inline |
obtain the geometry, this local function lives on
- Returns
- reference to the geometry
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet().
◆ hessian()
|
inline |
evaluate Hessian of the local function
- Note
- Though the Hessian is evaluated on the reference element, the return value is the Hessian with respect to the actual entity.
- Parameters
-
[in] x evaluation point in local coordinates [out] ret Hessian of the function in the evaluation point
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), Dune::Fem::BasisFunctionSet< Entity, Range >::hessianAll(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
◆ init() [1/2]
|
inlineprotected |
initialize the local function for an basisFunctionSet
Binds the local function to an basisFunctionSet and entity.
- Note
- A local function must be initialized to an entity before it can be used.
- This function can be called multiple times to use the local function for more than one entity.
- Parameters
-
[in] basisFunctionSet to bind the local function to
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::BasisFunctionSet< Entity, Range >::size().
◆ init() [2/2]
|
inlineprotected |
initialize the local function for an entity
Binds the local function to an basisFunctionSet and entity.
- Note
- Must be overloaded on the derived implementation class.
- Parameters
-
[in] entity to bind the local function to
References DUNE_THROW.
Referenced by Dune::Fem::BasicTemporaryLocalFunction< DiscreteFunctionSpace, DofVector >::init().
◆ jacobian()
|
inline |
evaluate Jacobian of the local function
- Note
- Though the Jacobian is evaluated on the reference element, the return value is the Jacobian with respect to the actual entity.
- Parameters
-
[in] x evaluation point in local coordinates [out] ret Jacobian of the function in the evaluation point
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), Dune::Fem::BasisFunctionSet< Entity, Range >::jacobianAll(), and Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
◆ numDofs()
|
inline |
obtain the number of local DoFs
Obtain the number of local DoFs of this local function. The value is identical to the number of basis functons on the entity.
- Returns
- number of local DoFs
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
◆ operator*()
|
inlineconstexprinherited |
indefinite vector dot product \(\left (x^T \cdot y \right)\) which corresponds to Petsc's VecTDot
http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecTDot.html
- Parameters
-
x other vector
- Returns
◆ operator*=()
|
inlineconstexprinherited |
vector space multiplication with scalar
we use enable_if to avoid an ambiguity, if the function parameter can be converted to field_type implicitly. (see FS#1457)
The function is only enabled, if the parameter is directly convertible to field_type.
◆ operator+=()
|
inlineconstexprinherited |
vector space add scalar to all comps
we use enable_if to avoid an ambiguity, if the function parameter can be converted to value_type implicitly. (see FS#1457)
The function is only enabled, if the parameter is directly convertible to value_type.
◆ operator-=()
|
inlineconstexprinherited |
vector space subtract scalar from all comps
we use enable_if to avoid an ambiguity, if the function parameter can be converted to value_type implicitly. (see FS#1457)
The function is only enabled, if the parameter is directly convertible to value_type.
◆ operator/=()
|
inlineconstexprinherited |
vector space division by scalar
we use enable_if to avoid an ambiguity, if the function parameter can be converted to field_type implicitly. (see FS#1457)
The function is only enabled, if the parameter is directly convertible to field_type.
◆ operator[]() [1/2]
|
inline |
access to local dofs (read-write)
- Parameters
-
[in] num local DoF number
- Returns
- reference to DoF
◆ operator[]() [2/2]
|
inline |
access to local dofs (read-only)
- Parameters
-
[in] num local dof number
- Returns
- reference to dof
◆ order()
|
inline |
obtain the order of this local function
The order of a local function refers to the polynomial order required to integrate it exactly.
- Note
- It is not completely clear what this value should be, e.g., for bilinear basis functions.
- Returns
- order of the local function
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::basisFunctionSet(), and Dune::Fem::BasisFunctionSet< Entity, Range >::order().
Referenced by Dune::Fem::LocalOrthonormalL2Projection< GridPart, BasisFunctionSet, Quadrature >::apply().
◆ size()
|
inline |
obtain the number of local DoFs
Obtain the number of local DoFs of this local function. The value is identical to the number of basis functons on the entity.
- Returns
- number of local DoFs
References Dune::Fem::LocalFunction< BasisFunctionSet, LocalDofVector >::localDofVector().
Referenced by Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverseDefault(), Dune::Fem::LocalMassMatrixImplementation< DiscreteFunctionSpace, VolumeQuadrature, refElemScaling >::applyInverseLocally(), and Dune::Fem::hpDG::DefaultDataProjection< DiscreteFunction >::operator()().
Friends And Related Function Documentation
◆ operator<<()
|
related |
Write a DenseVector to an output stream.
- Parameters
-
[in] s std :: ostream to write to [in] v DenseVector to write
- Returns
- the output stream (s)
The documentation for this class was generated from the following file:
- dune/fem/function/localfunction/localfunction.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden & Uni Heidelberg |
generated with Hugo v0.111.3
(Jan 9, 23:34, 2026)