Dune-Functions 2.11
Loading...
Searching...
No Matches
Dune::Functions Namespace Reference

Namespaces

namespace  Concept
namespace  ContainerDescriptors
namespace  Experimental
namespace  BasisFactory
namespace  ImplDoc

Classes

struct  MonomialSet
 Function, which evaluates all monomials up to degree maxDegree in a given coordinate. More...
class  Polynomial
 A univariate polynomial implementation. More...
class  TrigonometricFunction
 A linear combination of trigonomic functions. More...
class  InvalidRange
 Dummy range class to be used if no proper type is available. More...
struct  DefaultDerivativeTraits
 Default implementation for derivative traits. More...
struct  DefaultDerivativeTraits< double(double) >
 Default implementation for derivative traits. More...
struct  DefaultDerivativeTraits< K(FieldVector< K, n >)>
 Default implementation for derivative traits. More...
struct  DefaultDerivativeTraits< FieldVector< K, m >(FieldVector< K, n >)>
 Default implementation for derivative traits. More...
struct  DefaultDerivativeTraits< FieldMatrix< K, 1, m >(FieldVector< K, n >)>
 Default implementation for derivative traits. More...
class  DifferentiableFunction
class  DifferentiableFunction< Range(Domain), DerivativeTraits, bufferSize >
 Class storing differentiable functions using type erasure. More...
class  DifferentiableFunctionFromCallables
class  DifferentiableFunctionFromCallables< Range(Domain), DerivativeTraits, F >
 Wrap a list of callable objects as derivative sequence modelling Concept::DifferentiableFunction<Range(Domain), DerivativeTraits>. More...
class  DifferentiableFunctionFromCallables< Range(Domain), DerivativeTraits, F, DF, Derivatives... >
 Wrap a list of callable objects as derivative sequence modelling Concept::DifferentiableFunction<Range(Domain), DerivativeTraits>. More...
class  GeometryInAncestor
 A geometry embedding a descendent element into an ancestor. More...
class  PolymorphicType
 Base class with polymorphic type boiler plate code. More...
class  LocalFunction
class  LocalFunction< Range(Domain), LocalContext, DerivativeTraits, bufferSize >
 Class storing local functions using type erasure. More...
class  StaticMultiIndex
 A statically sized MultiIndex type. More...
class  StaticMultiIndex< size_type, 1 >
 A statically sized MultiIndex type. More...
class  OverflowArray
 A dynamically sized array-like class with overflow. More...
class  PolymorphicSmallObject
 A wrapper providing small object optimization with polymorphic types. More...
class  ReservedDeque
 A double-ended queue (deque) class with statically reserved memory. More...
struct  IsCallable
 Helper class to check that F is callable. More...
struct  SignatureTraits
 Helper class to deduce the signature of a callable. More...
struct  SignatureTag
struct  SignatureTag< Range(Domain), DerivativeTraitsT >
 Tag-class to encapsulate signature information. More...
struct  HasStaticSize
 Check if type is a statically sized container. More...
struct  StaticSizeOrZero
 Obtain size of statically sized container, or 0 if dynamic size. More...
class  TypeErasureBase
 Base class for type-erased interface wrapper. More...
struct  LastType
 Get last entry of type list. More...
struct  RotateTuple
 Rotate type list by one, such that last entry is moved to first position. More...
class  ArgyrisPreBasis
 A pre-basis for a Argyrisbasis. More...
class  ArgyrisNode
class  BrezziDouglasMariniNode
class  BrezziDouglasMariniPreBasis
class  BSplineLocalFiniteElement
 LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More...
class  BSplinePreBasis
 Pre-basis for B-spline basis. More...
class  BSplineLocalBasis
 LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch to a knot span. More...
class  BSplineLocalCoefficients
 Attaches a shape function to an entity. More...
class  BSplineLocalInterpolation
 Local interpolation in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More...
class  BSplineNode
class  CompositePreBasis
 A pre-basis for composite bases. More...
class  CubicHermitePreBasis
 A pre-basis for a Hermitebasis. More...
struct  H2LocalBasisTraits
class  CubicHermiteNode
class  DefaultGlobalBasis
 Global basis for given pre-basis. More...
class  DefaultLocalView
 The restriction of a finite element basis to a single element. More...
class  DynamicPowerPreBasis
 A pre-basis for dynamic power bases. More...
class  HierarchicalLagrangePreBasis
 A pre-basis for a hierarchical Lagrange basis. More...
class  HierarchicalLagrangePreBasis< GV, 1, R >
class  HierarchicalLagrangePreBasis< GV, 2, R >
class  HierarchicalLagrangeWithElementBubblePreBasis
 A pre-basis for a hierarchical basis. More...
class  HierarchicalLagrangeWithElementBubblePreBasis< GV, 1, R >
class  HierarchicalLagrangeWithElementBubblePreBasis< GV, 2, R >
struct  HierarchicNodeToRangeMap
 A simple node to range map using the nested tree indices. More...
class  LagrangeNode
class  LagrangePreBasis
 A pre-basis for a PQ-lagrange bases with given order. More...
class  LagrangeDGPreBasis
 PreBasis implementation for a Lagrangean-DG finite element space. More...
class  LeafPreBasisMapperMixin
 A generic MixIn class for PreBasis with flat indices computed from a mapper. More...
class  LeafPreBasisMixin
 A generic MixIn class for PreBasis. More...
class  LFEPreBasisMixin
 A pre-basis mixin class parametrized with a local finite-element and a DOF layout. More...
class  MorleyPreBasis
 A pre-basis for a Morleybasis. More...
class  MorleyNode
class  NedelecNode
class  NedelecPreBasis
class  BasisNodeMixin
class  LeafBasisNode
class  InnerBasisNodeMixin
class  PowerBasisNode
class  DynamicPowerBasisNode
class  CompositeBasisNode
class  PowerPreBasis
 A pre-basis for power bases. More...
class  RannacherTurekNode
class  RannacherTurekPreBasis
 Pre-basis for a Rannacher-Turek basis. More...
class  RaviartThomasNode
class  RaviartThomasPreBasis
class  RefinedLagrangeNode
class  RefinedLagrangePreBasis
 A pre-basis for a refined Lagrange bases. More...
class  SubEntityDOFs
 Range of DOFs associated to sub-entity. More...
class  SubspaceBasis
class  SubspaceLocalView
 The restriction of a finite element basis to a single element. More...
class  TaylorHoodVelocityTree
class  TaylorHoodBasisTree
class  TaylorHoodPreBasis
 Pre-basis for lowest order Taylor-Hood basis. More...
class  AnalyticGridViewFunction
class  AnalyticGridViewFunction< Range(Domain), GV, F, DerivativeTraits >
 Class wrapping any differentiable function as grid function. More...
class  CoarseFunctionOnFineGridView
 A wrapper representing a coarse grid function on a fine gridview. More...
class  ComposedGridFunction
 Composition of grid functions with another function. More...
class  DiscreteGlobalBasisFunctionDerivative
 Derivative of a DiscreteGlobalBasisFunction. More...
class  DiscreteGlobalBasisFunction
 A grid function induced by a global basis and a coefficient vector. More...
class  FaceNormalGridFunction
 Grid function implementing the piecewise element face normal. More...
class  FineFunctionOnCoarseGridView
 A wrapper representing a fine grid function on a gridview. More...
class  GridFunction
class  GridFunction< Range(Domain), ES, DerivativeTraits, bufferSize >
 Wrapper class for functions defined on a Grid. More...
class  GridViewEntitySet
 An entity set for all entities of given codim in a grid view. More...
class  GridViewFunction
class  GridViewFunction< Range(Domain), GV, DerivativeTraits, bufferSize >
 Wrapper class for functions defined on a GridView. More...
struct  LocalDerivativeTraits
 Derivative traits for local functions. More...

Typedefs

template<class T, class... Args>
using enableIfConstructible
 Helper to constrain forwarding constructors.
template<class T>
using StaticSize
 Obtain size of statically sized container as integral_constant, or fail.
template<template< class... > class T, class ArgTuple>
using ExpandTuple = typename Imp::ExpandTupleHelper<T, ArgTuple>::Type
 Expand tuple arguments as template arguments.
template<template< class... > class F, class... Tuples>
using TransformTuple = typename Imp::TransformTupleHelper<F, Tuples...>::Type
 Transform tuple types argument using type-functor.
template<class IntegerSequence>
using IntegerSequenceTuple = typename Imp::IntegerSequenceTupleHelper<IntegerSequence>::Type
 Transform integer_sequence<I,k...> to tuple<integral_constant<I,k>...>.
template<class GV, class R = double>
using ArgyrisBasis = DefaultGlobalBasis<ArgyrisPreBasis<GV, R> >
 Nodal basis of a scalar cubic Argyris finite element space.
template<typename GV, int k>
using BrezziDouglasMariniBasis = DefaultGlobalBasis<BrezziDouglasMariniPreBasis<GV, k> >
 Basis of a scalar k-th-order BDM finite element space on simplex and cube grids.
template<typename GV>
using BSplineBasis = DefaultGlobalBasis<BSplinePreBasis<GV> >
 A global B-spline basis.
template<class GV, bool reduced = false, class R = double>
using CubicHermiteBasis = DefaultGlobalBasis<CubicHermitePreBasis<GV, R, reduced> >
 Nodal basis of a scalar cubic Hermite finite element space.
template<class size_type>
using FlatMultiIndex = StaticMultiIndex<size_type, 1>
 A multi-index class with only one level.
template<typename GV, int k, typename R = double>
using HierarchicalLagrangeBasis = DefaultGlobalBasis<HierarchicalLagrangePreBasis<GV, k, R> >
 Basis of a scalar Hierarchical Lagrange finite element space.
template<typename GV, int k, typename R = double>
using HierarchicalLagrangeWithElementBubbleBasis = DefaultGlobalBasis<HierarchicalLagrangeWithElementBubblePreBasis<GV, k, R> >
 Basis of a Hierarchical Lagrange finite element space with element bubble functions.
template<typename GV, int k = -1, typename R = double>
using LagrangeBasis = DefaultGlobalBasis<LagrangePreBasis<GV, k, R> >
 Nodal basis of a scalar k-th-order Lagrangean finite element space.
template<typename GV, int k, typename R = double>
using LagrangeDGNode = LagrangeNode<GV, k, R>
template<typename GV, int k = -1, typename R = double>
using LagrangeDGBasis = DefaultGlobalBasis<LagrangeDGPreBasis<GV, k, R> >
 Basis of a scalar k-th-order Lagrangean-DG finite element space.
template<class GV, class R = double>
using MorleyBasis = DefaultGlobalBasis<MorleyPreBasis<GV, R> >
 Nodal basis of a scalar quadratic Morley finite element space.
template<typename GV, std::size_t kind, std::size_t order, typename Range = double>
using NedelecBasis = DefaultGlobalBasis<NedelecPreBasis<GV, Range, kind, order > >
 Basis of a k-th-order Nédélec finite element space.
template<typename GV>
using RannacherTurekBasis = DefaultGlobalBasis<RannacherTurekPreBasis<GV> >
 Rannacher-Turek basis.
template<typename GV, int k>
using RaviartThomasBasis = DefaultGlobalBasis<RaviartThomasPreBasis<GV, k> >
 Basis of a k-th-order Raviart Thomas finite element space.
template<typename GV, int k, typename R = double>
using RefinedLagrangeBasis = DefaultGlobalBasis<RefinedLagrangePreBasis<GV,k,R> >
 Nodal basis of a continuous Lagrange finite-element space on a uniformly refined simplex element.
template<typename GV>
using TaylorHoodBasis = DefaultGlobalBasis<TaylorHoodPreBasis<GV> >
 Nodal basis for a lowest order Taylor-Hood Lagrangean finite element space.

Functions

template<class K>
 Polynomial (std::vector< K >) -> Polynomial< K, std::vector< K > >
template<class K, std::size_t n>
 Polynomial (std::array< K, n >) -> Polynomial< K, std::array< K, n > >
template<class K, K... ci>
 Polynomial (std::integer_sequence< K, ci... >) -> Polynomial< K, std::integer_sequence< K, ci... > >
template<class K>
 Polynomial (std::initializer_list< K >) -> Polynomial< K, std::vector< K > >
template<class K, class Coefficients>
auto makePolynomial (Coefficients coefficients)
 Create Polynomial.
template<class K, class C>
auto makePolynomial (std::initializer_list< C > coefficients)
 Create Polynomial.
template<class K, int sinFactor, int cosFactor>
TrigonometricFunction< K, -cosFactor, sinFactor > derivative (const TrigonometricFunction< K, sinFactor, cosFactor > &f)
 Obtain derivative of TrigonometricFunction function.
template<class T, class ContainerDescriptor>
auto makeContainer (const ContainerDescriptor &descriptor, const T &defaultValue)
 Construct a nested random access container compatible with the container descriptor.
template<class T, class ContainerDescriptor>
auto makeContainer (const ContainerDescriptor &descriptor)
 Construct a nested random access container compatible with the container descriptor.
template<class V>
constexpr auto fieldTypes ()
 Generate list of field types in container.
template<class V>
constexpr bool hasUniqueFieldType ()
 Check if container has a unique field type.
template<class Vector>
auto istlVectorBackend (Vector &v)
 Return a vector backend wrapping non-const ISTL like containers.
template<class Vector>
auto istlVectorBackend (const Vector &v)
 Return a vector backend wrapping const ISTL like containers.
template<class T = double, class ContainerDescriptor>
auto makeISTLVector (const ContainerDescriptor &tree)
 Construct an istl vector type compatible with the container descriptor.
template<class Signature, template< class > class DerivativeTraits, class... F>
DifferentiableFunctionFromCallables< Signature, DerivativeTraits, F... > makeDifferentiableFunctionFromCallables (const SignatureTag< Signature, DerivativeTraits > &signatureTag, F &&... f)
 Create a DifferentiableFunction from callables.
template<class C, class I, class F, std::enable_if_t< Dune::models< Imp::Concept::HasDynamicIndexAccess< I >, C >(), int > = 0>
auto hybridIndexAccess (C &&c, const I &i, F &&f) -> decltype(f(c[i]))
 Provide operator[] index-access for containers.
template<class C, class I, class F, std::enable_if_t< not Dune::models< Imp::Concept::HasDynamicIndexAccess< I >, C >(), int > = 0>
decltype(auto) hybridIndexAccess (C &&c, const I &i, F &&f)
 Provide operator[] index-access for containers.
template<class Result, class C, class MultiIndex>
Result hybridMultiIndexAccess (C &&c, const MultiIndex &index)
 Provide multi-index access by chaining operator[].
template<class C, class MultiIndex, class IsFinal>
constexpr decltype(auto) resolveDynamicMultiIndex (C &&c, const MultiIndex &multiIndex, const IsFinal &isFinal)
 Provide multi-index access by chaining operator[].
template<class C, class MultiIndex>
constexpr decltype(auto) resolveDynamicMultiIndex (C &&c, const MultiIndex &multiIndex)
 Provide multi-index access by chaining operator[].
template<class C, class MultiIndex>
constexpr decltype(auto) resolveStaticMultiIndex (C &&c, const MultiIndex &multiIndex)
 Provide multi-index access by chaining operator[].
template<class size_type, std::size_t n>
std::ostream & operator<< (std::ostream &stream, const StaticMultiIndex< size_type, n > &c)
template<class Range, class Domain, template< class > class DerivativeTraits>
auto derivativeSignatureTag (SignatureTag< Range(Domain), DerivativeTraits > tag)
 Construct SignatureTag for derivative.
template<std::size_t maxOrder, class Signature, template< class > class DerivativeTraits>
auto derivativeSignatureTags (Dune::Functions::SignatureTag< Signature, DerivativeTraits > tag)
 Construct SignatureTags for derivatives.
template<std::size_t begin_t, std::size_t end_t, class F, class... Args>
void staticFindInRange (F &&f, Args &&... args)
 Static find loop.
template<class F, class size_type, size_type firstValue, class... Args>
auto forwardAsStaticInteger (std::integer_sequence< size_type, firstValue > values, const size_type &i, F &&f, Args &&... args) -> decltype(f(std::integral_constant< size_type, firstValue >(), std::forward< Args >(args)...))
template<class F, class size_type, size_type firstValue, size_type secondValue, size_type... otherValues, class... Args>
auto forwardAsStaticInteger (std::integer_sequence< size_type, firstValue, secondValue, otherValues... > values, const size_type i, F &&f, Args &&... args) -> decltype(f(std::integral_constant< size_type, firstValue >(), std::forward< Args >(args)...))
template<std::size_t end, class F, class size_type, class... Args>
auto forwardAsStaticIndex (const size_type &i, F &&f, Args &&... args) -> decltype(f(Dune::Indices::_0, std::forward< Args >(args)...))
 Transform dynamic index to static index_constant.
template<class F, class... T>
auto transformTuple (F &&f, const std::tuple< T... > &tuple) -> decltype(Imp::transformTupleHelper(std::forward< F >(f), tuple, std::index_sequence_for< T... >{}))
 Transform tuple value using a functor.
template<class F, class... T1, class... T2>
auto transformTuple (F &&f, const std::tuple< T1... > &tuple1, const std::tuple< T2... > &tuple2) -> decltype(Imp::transformTupleHelper(std::forward< F >(f), tuple1, tuple2, std::index_sequence_for< T1... >{}))
 Transform tuple value using a binary functor.
template<class Expression>
auto callableCheck (Expression f)
 Create a predicate for checking validity of expressions.
template<class Check>
auto negatePredicate (Check check)
 Negate given predicate.
template<class T>
auto forwardCapture (T &&t)
 Create a capture object for perfect forwarding.
template<class Basis, class F, decltype(std::declval< std::decay_t< F > >()(0, std::declval< typename Basis::LocalView >(), std::declval< typename Basis::GridView::Intersection >()), 0) = 0>
void forEachBoundaryDOF (const Basis &basis, F &&f)
 Loop over all DOFs on the boundary.
template<class Basis, class F, decltype(std::declval< std::decay_t< F > >()(0, std::declval< typename Basis::LocalView >()), 0) = 0>
void forEachBoundaryDOF (const Basis &basis, F &&f)
 Loop over all DOFs on the boundary.
template<class Basis, class F, decltype(std::declval< std::decay_t< F > >()(std::declval< typename Basis::MultiIndex >()), 0) = 0>
void forEachBoundaryDOF (const Basis &basis, F &&f)
 Loop over all DOFs on the boundary.
template<class PreBasis>
auto containerDescriptor (const PreBasis &preBasis)
 Return the container descriptor of the pre-basis, if defined, otherwise ContainerDescriptor::Unknown.
template<class PreBasis>
 DefaultGlobalBasis (PreBasis &&) -> DefaultGlobalBasis< std::decay_t< PreBasis > >
template<class GridView, class PreBasisFactory>
 DefaultGlobalBasis (const GridView &gv, PreBasisFactory &&f) -> DefaultGlobalBasis< std::decay_t< decltype(f(gv))> >
template<class T>
auto flatVectorView (T &t)
 Create flat vector view of passed mutable container.
template<class T>
auto flatVectorView (const T &t)
 Create flat vector view of passed const container.
template<class T>
auto flatVectorView (T &&t)
 Create flat vector view of passed container temporary.
template<class B, class C, class F, class BV, class NTRE>
void interpolate (const B &basis, C &&coeff, const F &f, const BV &bv, const NTRE &nodeToRangeEntry)
 Interpolate given function in discrete function space.
template<class B, class C, class F, class BV>
void interpolate (const B &basis, C &&coeff, const F &f, const BV &bitVector)
 Interpolate given function in discrete function space.
template<class B, class C, class F>
void interpolate (const B &basis, C &&coeff, const F &f)
 Interpolate given function in discrete function space.
template<class GridView, class LocalCoefficients>
auto subIndexRange (const Dune::MultipleCodimMultipleGeomTypeMapper< GridView > &mapper, const typename GridView::template Codim< 0 >::Entity &element, const LocalCoefficients &localCoefficients)
template<class GV, class LFE>
 LFEPreBasisMixin (const GV &, const LFE &, MCMGLayout) -> LFEPreBasisMixin< GV, LFE >
template<typename Tree, typename Entity>
void bindTree (Tree &tree, const Entity &entity, std::size_t offset=0)
template<typename Tree>
void initializeTree (Tree &tree, std::size_t treeIndexOffset=0)
template<class T>
auto subEntityDOFs (const T &)
 Create SubEntityDOFs object.
template<class LocalView>
auto subEntityDOFs (const LocalView &localView, std::size_t subEntityIndex, std::size_t subEntityCodim)
 Create bound SubEntityDOFs object.
template<class LocalView, class Intersection>
auto subEntityDOFs (const LocalView &localView, const Intersection &intersection)
 Create bound SubEntityDOFs object.
template<class RB, class TP>
 SubspaceBasis (const RB &, const TP) -> SubspaceBasis< RB, TP >
template<class RootRootBasis, class InnerTP, class OuterTP>
 SubspaceBasis (const SubspaceBasis< RootRootBasis, InnerTP > &rootBasis, const OuterTP &prefixPath) -> SubspaceBasis< std::decay_t< decltype(rootBasis.rootBasis())>, decltype(Dune::TypeTree::join(rootBasis.prefixPath(), prefixPath))>
template<class RootBasis, class... PrefixTreeIndices>
auto subspaceBasis (const RootBasis &rootBasis, const TypeTree::TreePath< PrefixTreeIndices... > &prefixPath)
 Create SubspaceBasis from a root basis and a prefixPath.
template<class RootBasis, class... PrefixTreeIndices>
auto subspaceBasis (const RootBasis &rootBasis, const PrefixTreeIndices &... prefixTreeIndices)
template<class F, class GridView, class Domain = typename GridView::template Codim<0>::Geometry::GlobalCoordinate, class Range = std::invoke_result_t<F,Domain>>
 AnalyticGridViewFunction (const F &, const GridView &) -> AnalyticGridViewFunction< Range(Domain), GridView, F >
template<class F, class GridView>
auto makeAnalyticGridViewFunction (F &&f, const GridView &gridView)
 Create an AnalyticGridViewFunction from a function and a grid view.
template<class OF, class... IF>
 ComposedGridFunction (const OF &, const IF &...) -> ComposedGridFunction< OF, IF... >
template<class OF, class... IF>
auto makeComposedGridFunction (OF &&outerFunction, IF &&... innerFunction)
 Create a ComposedGridFunction that composes grid-functions with another function.
template<typename R, typename B, typename V>
auto makeDiscreteGlobalBasisFunction (B &&basis, V &&vector)
 Generate a DiscreteGlobalBasisFunction.
template<class F, class GridView, std::enable_if_t< models< Imp::HasFreeLocalFunction, F >(), int > = 0>
std::decay_t< F > makeGridViewFunction (F &&f, const GridView &gridView)
 Construct a function modeling GridViewFunction from function and grid view.
template<class F, class GridView, std::enable_if_t< not(models< Imp::HasFreeLocalFunction, F >()), int > = 0>
auto makeGridViewFunction (F &&f, const GridView &gridView) -> decltype(makeAnalyticGridViewFunction(std::forward< F >(f), gridView))
 Construct a function modeling GridViewFunction from function and grid view.
template<class F, class GridView>
auto makeAnalyticGridViewFunction (F &&f, const GridView &gridView)
 Create an AnalyticGridViewFunction from a function and a grid view.
template<class OF, class... IF>
auto makeComposedGridFunction (OF &&outerFunction, IF &&... innerFunction)
 Create a ComposedGridFunction that composes grid-functions with another function.
template<typename R, typename B, typename V>
auto makeDiscreteGlobalBasisFunction (B &&basis, V &&vector)
 Generate a DiscreteGlobalBasisFunction.

Variables

template<class T>
constexpr bool HasStaticSize_v = HasStaticSize<T>::value
 A variable template representing the value of HasStaticSize.

Typedef Documentation

◆ BrezziDouglasMariniBasis

Basis of a scalar k-th-order BDM finite element space on simplex and cube grids.

TODO: Fix this for grids with more than one element type

Template Parameters
GVThe GridView that the space is defined on
kThe order of the basis

◆ IntegerSequenceTuple

template<class IntegerSequence>
using Dune::Functions::IntegerSequenceTuple = typename Imp::IntegerSequenceTupleHelper<IntegerSequence>::Type

Transform integer_sequence<I,k...> to tuple<integral_constant<I,k>...>.

◆ LagrangeDGNode

template<typename GV, int k, typename R = double>
using Dune::Functions::LagrangeDGNode = LagrangeNode<GV, k, R>

◆ NedelecBasis

template<typename GV, std::size_t kind, std::size_t order, typename Range = double>
using Dune::Functions::NedelecBasis = DefaultGlobalBasis<NedelecPreBasis<GV, Range, kind, order > >

Basis of a k-th-order Nédélec finite element space.

Template Parameters
GVThe GridView that the space is defined on
RangeNumber type used for shape function values
kindKind of the basis: 1 (for Nédélec element of the first kind) or 2
orderThe order of the basis (lowest order is '1')

◆ RaviartThomasBasis

template<typename GV, int k>
using Dune::Functions::RaviartThomasBasis = DefaultGlobalBasis<RaviartThomasPreBasis<GV, k> >

Basis of a k-th-order Raviart Thomas finite element space.

TODO: Fix this for grids with more than one element type

Template Parameters
GVThe GridView that the space is defined on
kThe order of the basis

Function Documentation

◆ AnalyticGridViewFunction()

template<class F, class GridView, class Domain = typename GridView::template Codim<0>::Geometry::GlobalCoordinate, class Range = std::invoke_result_t<F,Domain>>
Dune::Functions::AnalyticGridViewFunction ( const F & ,
const GridView &  )->AnalyticGridViewFunction< Range(Domain), GridView, F >

◆ bindTree()

template<typename Tree, typename Entity>
void Dune::Functions::bindTree ( Tree & tree,
const Entity & entity,
std::size_t offset = 0 )

◆ ComposedGridFunction()

template<class OF, class... IF>
Dune::Functions::ComposedGridFunction ( const OF & ,
const IF & ... )->ComposedGridFunction< OF, IF... >

◆ containerDescriptor()

template<class PreBasis>
auto Dune::Functions::containerDescriptor ( const PreBasis & preBasis)

Return the container descriptor of the pre-basis, if defined, otherwise ContainerDescriptor::Unknown.

◆ DefaultGlobalBasis() [1/2]

template<class GridView, class PreBasisFactory>
Dune::Functions::DefaultGlobalBasis ( const GridView & gv,
PreBasisFactory && f )->DefaultGlobalBasis< std::decay_t< decltype(f(gv))> >

◆ DefaultGlobalBasis() [2/2]

template<class PreBasis>
Dune::Functions::DefaultGlobalBasis ( PreBasis && ) ->DefaultGlobalBasis< std::decay_t< PreBasis > >

◆ fieldTypes()

template<class V>
auto Dune::Functions::fieldTypes ( )
constexpr

Generate list of field types in container.

This generates a Dune::TypeList of the field types in the given container type. To determine the field types, operator[] is called as often as passible with std::size_t or Dune::index_constant arguments. The return types obtained if no more operator[] call is available are returned as Dune::TypeList. Notice that possible duplicate entries are removed. However, const and reference qualifiers are deliberately preserved.

◆ flatVectorView() [1/3]

template<class T>
auto Dune::Functions::flatVectorView ( const T & t)

Create flat vector view of passed const container.

When passed a nested container, the resulting value is a flat-vector-like view object. It provides an operator[] method to access all entries of the underlying nested container using flat consecutive indices and a size() method to compute the corresponding total size.

This method will create a view object storing a pointer to the passed const container.

◆ flatVectorView() [2/3]

template<class T>
auto Dune::Functions::flatVectorView ( T && t)

Create flat vector view of passed container temporary.

When passed a nested container, the resulting value is a flat-vector-like view object. It provides an operator[] method to access all entries of the underlying nested container using flat consecutive indices and a size() method to compute the corresponding total size.

This method will create a 'view' object storing the provided temporary container by value.

◆ flatVectorView() [3/3]

template<class T>
auto Dune::Functions::flatVectorView ( T & t)

Create flat vector view of passed mutable container.

When passed a nested container, the resulting value is a flat-vector-like view object. It provides an operator[] method to access all entries of the underlying nested container using flat consecutive indices and a size() method to compute the corresponding total size.

This method will create a view object storing a pointer to the passed mutable container.

◆ forEachBoundaryDOF() [1/2]

template<class Basis, class F, decltype(std::declval< std::decay_t< F > >()(std::declval< typename Basis::MultiIndex >()), 0) = 0>
void Dune::Functions::forEachBoundaryDOF ( const Basis & basis,
F && f )

Loop over all DOFs on the boundary.

This loops over all DOFs of a basis associated to sub-entities on the boundary. This overload will pass a single argument to the given loop callback: The global (multi-)index of the boundary DOF. Notice that this may visit the same DOF multiple times.

If this callback signature is not suitable you can use one of the another variants of forEachBoundaryDOF.

Parameters
basisA function space basis
fA callback that will be called with the global index of the visited boundary DOF

◆ forEachBoundaryDOF() [2/2]

template<class Basis, class F, decltype(std::declval< std::decay_t< F > >()(0, std::declval< typename Basis::LocalView >()), 0) = 0>
void Dune::Functions::forEachBoundaryDOF ( const Basis & basis,
F && f )

Loop over all DOFs on the boundary.

This loops over all DOFs of a basis associated to sub-entities on the boundary. This overload will pass two arguments to the given loop callback: The local index of the boundary DOF and a bound local view this local index belongs to. Notice that this may visit the same DOF multiple times.

If this callback signature is not suitable you can use one of the another variants of forEachBoundaryDOF.

Parameters
basisA function space basis
fA callback that will be called with a local index and a bound local view of the visited boundary DOF

◆ forwardAsStaticInteger() [1/2]

template<class F, class size_type, size_type firstValue, class... Args>
auto Dune::Functions::forwardAsStaticInteger ( std::integer_sequence< size_type, firstValue > values,
const size_type & i,
F && f,
Args &&... args )->decltype(f(std::integral_constant< size_type, firstValue >(), std::forward< Args >(args)...))

◆ forwardAsStaticInteger() [2/2]

template<class F, class size_type, size_type firstValue, size_type secondValue, size_type... otherValues, class... Args>
auto Dune::Functions::forwardAsStaticInteger ( std::integer_sequence< size_type, firstValue, secondValue, otherValues... > values,
const size_type i,
F && f,
Args &&... args )->decltype(f(std::integral_constant< size_type, firstValue >(), std::forward< Args >(args)...))

◆ forwardCapture()

template<class T>
auto Dune::Functions::forwardCapture ( T && t)

Create a capture object for perfect forwarding.

The returned object will capture the passed argument t. If t is passed as r-value, then it is captured by value, otherwise by reference. The captured value is accessible once using the forward() method which either returns the catured reference or moves the captured value.

This allows to capture values for perfect forwarding in lambda functions using [t=forwardCapture(std::forward<T>(t))]() -> decltype(auto) { return t.forward(); }

◆ hasUniqueFieldType()

template<class V>
bool Dune::Functions::hasUniqueFieldType ( )
constexpr

Check if container has a unique field type.

This returns if fieldTypes<V>() has exactly one entry.

◆ initializeTree()

template<typename Tree>
void Dune::Functions::initializeTree ( Tree & tree,
std::size_t treeIndexOffset = 0 )

◆ interpolate() [1/3]

template<class B, class C, class F>
void Dune::Functions::interpolate ( const B & basis,
C && coeff,
const F & f )

Interpolate given function in discrete function space.

Notice that this will only work if the range type of f and the block type of coeff are compatible and supported by flatVectorView.

This function will only work, if the local ansatz tree of the basis is trivial, i.e., a single leaf node.

Parameters
basisGlobal function space basis of discrete function space
coeffCoefficient vector to represent the interpolation
fFunction to interpolate

◆ interpolate() [2/3]

template<class B, class C, class F, class BV>
void Dune::Functions::interpolate ( const B & basis,
C && coeff,
const F & f,
const BV & bitVector )

Interpolate given function in discrete function space.

Only vector coefficients marked as 'true' in the bitVector argument are interpolated. Use this, e.g., to interpolate Dirichlet boundary values.

Notice that this will only work if the range type of f and the block type of coeff are compatible and supported by flatVectorView.

Parameters
basisGlobal function space basis of discrete function space
coeffCoefficient vector to represent the interpolation
fFunction to interpolate
bitVectorA vector with flags marking all DOFs that should be interpolated

◆ interpolate() [3/3]

template<class B, class C, class F, class BV, class NTRE>
void Dune::Functions::interpolate ( const B & basis,
C && coeff,
const F & f,
const BV & bv,
const NTRE & nodeToRangeEntry )

Interpolate given function in discrete function space.

Only vector coefficients marked as 'true' in the bitVector argument are interpolated. Use this, e.g., to interpolate Dirichlet boundary values.

Notice that this will only work if the range type of f and the block type of coeff are compatible and supported by flatVectorView.

Parameters
basisGlobal function space basis of discrete function space
coeffCoefficient vector to represent the interpolation
fFunction to interpolate
bvA bit vector with flags marking all DOFs that should be interpolated
nodeToRangeEntryPolymorphic functor mapping local ansatz nodes to range-indices of given function

◆ LFEPreBasisMixin()

template<class GV, class LFE>
Dune::Functions::LFEPreBasisMixin ( const GV & ,
const LFE & ,
MCMGLayout  )->LFEPreBasisMixin< GV, LFE >

◆ makeContainer() [1/2]

template<class T, class ContainerDescriptor>
auto Dune::Functions::makeContainer ( const ContainerDescriptor & descriptor)

Construct a nested random access container compatible with the container descriptor.

Parameters
descriptorA ContainerDescriptor provided by a global basis

The constructed container mimics the nested structure of the container descriptor, but uses data structures like std::vector, std::array, and Dune::TupleVector to represent the block levels. The entries in the vector are of type T and default initialized.

◆ makeContainer() [2/2]

template<class T, class ContainerDescriptor>
auto Dune::Functions::makeContainer ( const ContainerDescriptor & descriptor,
const T & defaultValue )

Construct a nested random access container compatible with the container descriptor.

Parameters
descriptorA ContainerDescriptor provided by a global basis
defaultValueThe default value to initialize the container entries.

The constructed container mimics the nested structure of the container descriptor, but uses data structures like std::vector, std::array, and Dune::TupleVector to represent the block levels. The entries in the vector are of type T and initialized with the provided default value.

◆ makeGridViewFunction() [1/2]

template<class F, class GridView, std::enable_if_t< models< Imp::HasFreeLocalFunction, F >(), int > = 0>
std::decay_t< F > Dune::Functions::makeGridViewFunction ( F && f,
const GridView & gridView )

Construct a function modeling GridViewFunction from function and grid view.

This specialization is used for functions that already support localFunction(). It will simply return a copy of f.

Parameters
fA function object supporting argument compatible with global coordinates
gridViewThe GridView the function should act on.
Returns
A function that models the GridViewFunction interface.

◆ makeGridViewFunction() [2/2]

template<class F, class GridView, std::enable_if_t< not(models< Imp::HasFreeLocalFunction, F >()), int > = 0>
auto Dune::Functions::makeGridViewFunction ( F && f,
const GridView & gridView )->decltype(makeAnalyticGridViewFunction(std::forward< F >(f), gridView))

Construct a function modeling GridViewFunction from function and grid view.

This specialization is used for functions that do not support localFunction() themselves. It will forward to makeAnalyticGridViewFunction().

Notice that the returned function will store a copy of the original function and the GridView.

Parameters
fA function object supporting argument compatible with global coordinates
gridViewThe GridView the function should act on.
Returns
A function that models the GridFunction interface.

◆ makeISTLVector()

template<class T = double, class ContainerDescriptor>
auto Dune::Functions::makeISTLVector ( const ContainerDescriptor & tree)

Construct an istl vector type compatible with the container descriptor.

The constructed vector mimics the nested structure of the container descriptor, but uses data structures like BlockVector and FieldVector to represent the block levels. The entries in the vector are of type T and initialized with the default value 0.

◆ makePolynomial() [1/2]

template<class K, class Coefficients>
auto Dune::Functions::makePolynomial ( Coefficients coefficients)

Create Polynomial.

Template Parameters
KScalar type. The polynomial will map K to K
CCoefficient container type

This helper function allows to specify K, but lets C be deduced from the argument. If the scalar type K is the same as the coefficient type (i.e. the entry type of the coefficient container C), then you can also rely on class template argument deduction and call Polynomial(coefficients) directly.

◆ makePolynomial() [2/2]

template<class K, class C>
auto Dune::Functions::makePolynomial ( std::initializer_list< C > coefficients)

Create Polynomial.

Template Parameters
KScalar type. The polynomial will map K to K
CCoefficient type

The initializer list will be stored as std::vector in the created object of type Polynomial<K,std::vector<C>>.

◆ operator<<()

template<class size_type, std::size_t n>
std::ostream & Dune::Functions::operator<< ( std::ostream & stream,
const StaticMultiIndex< size_type, n > & c )
inline

◆ Polynomial() [1/4]

template<class K, std::size_t n>
Dune::Functions::Polynomial ( std::array< K, n > ) ->Polynomial< K, std::array< K, n > >

◆ Polynomial() [2/4]

template<class K>
Dune::Functions::Polynomial ( std::initializer_list< K > ) ->Polynomial< K, std::vector< K > >

◆ Polynomial() [3/4]

template<class K, K... ci>
Dune::Functions::Polynomial ( std::integer_sequence< K, ci... > ) ->Polynomial< K, std::integer_sequence< K, ci... > >

◆ Polynomial() [4/4]

template<class K>
Dune::Functions::Polynomial ( std::vector< K > ) ->Polynomial< K, std::vector< K > >

◆ subIndexRange()

template<class GridView, class LocalCoefficients>
auto Dune::Functions::subIndexRange ( const Dune::MultipleCodimMultipleGeomTypeMapper< GridView > & mapper,
const typename GridView::template Codim< 0 >::Entity & element,
const LocalCoefficients & localCoefficients )

◆ SubspaceBasis() [1/2]

template<class RB, class TP>
Dune::Functions::SubspaceBasis ( const RB & ,
const TP  )->SubspaceBasis< RB, TP >

◆ SubspaceBasis() [2/2]

template<class RootRootBasis, class InnerTP, class OuterTP>
Dune::Functions::SubspaceBasis ( const SubspaceBasis< RootRootBasis, InnerTP > & rootBasis,
const OuterTP & prefixPath )->SubspaceBasis< std::decay_t< decltype(rootBasis.rootBasis())>, decltype(Dune::TypeTree::join(rootBasis.prefixPath(), prefixPath))>

◆ subspaceBasis() [1/2]

template<class RootBasis, class... PrefixTreeIndices>
auto Dune::Functions::subspaceBasis ( const RootBasis & rootBasis,
const PrefixTreeIndices &... prefixTreeIndices )

◆ subspaceBasis() [2/2]

template<class RootBasis, class... PrefixTreeIndices>
auto Dune::Functions::subspaceBasis ( const RootBasis & rootBasis,
const TypeTree::TreePath< PrefixTreeIndices... > & prefixPath )

Create SubspaceBasis from a root basis and a prefixPath.

This will not return a nested SubspaceBasis if rootBasis is already a SubspaceBasis. Instead it will join the tree paths and use them to construct a non-nested SubspaceBasis relative to rootBasis.rootBasis().

Parameters
rootBasisCreate a subspace basis relative to this basis
prefixPathA prefix path of the subspace within the root basis

Variable Documentation

◆ HasStaticSize_v

template<class T>
bool Dune::Functions::HasStaticSize_v = HasStaticSize<T>::value
inlineconstexpr

A variable template representing the value of HasStaticSize.