dune-functions 2.9.0
Namespaces | Classes | Typedefs | Functions | Variables
Dune::Functions Namespace Reference

Namespaces

namespace  BasisBuilder
 
namespace  BasisFactory
 
namespace  Concept
 
namespace  Experimental
 
namespace  ImplDoc
 

Classes

class  AnalyticGridViewFunction
 
class  AnalyticGridViewFunction< Range(Domain), GV, F, DerivativeTraits >
 Class wrapping any differentiable function as grid function. More...
 
class  BasisNodeMixin
 
class  BrezziDouglasMariniNode
 
class  BrezziDouglasMariniPreBasis
 
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  BSplineLocalFiniteElement
 LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More...
 
class  BSplineLocalInterpolation
 Local interpolation in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More...
 
class  BSplineNode
 
class  BSplinePreBasis
 Pre-basis for B-spline basis. More...
 
class  CallableFunctionWrapper
 Wrap a Dune::VirtualFunction into a callable object. More...
 
class  ComposedGridFunction
 Composition of grid functions with another function. More...
 
class  CompositeBasisNode
 
class  CompositePreBasis
 A pre-basis for composite bases. More...
 
struct  DefaultDerivativeTraits
 Default implementation for derivative traits. More...
 
struct  DefaultDerivativeTraits< double(double) >
 Default implementation for derivative traits. More...
 
struct  DefaultDerivativeTraits< FieldMatrix< K, 1, m >(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< K(FieldVector< K, n >)>
 Default implementation for derivative traits. More...
 
class  DefaultGlobalBasis
 Global basis for given pre-basis. More...
 
class  DefaultLocalView
 The restriction of a finite element basis to a single element. More...
 
struct  DefaultNodeToRangeMap
 A simple node to range map using lexicographic ordering. 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  DiscreteGlobalBasisFunction
 A grid function induced by a global basis and a coefficient vector. More...
 
class  DiscreteGlobalBasisFunctionDerivative
 Derivative of a DiscreteGlobalBasisFunction More...
 
class  FaceNormalGridFunction
 Grid function implementing the piecewise element face normal. More...
 
class  FunctionFromCallable
 
class  FunctionFromCallable< Range(Domain), F, FunctionInterface >
 Wrap a callable object as Dune::Function or Dune::VirtualFunction. 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  HasStaticSize
 Check if type is a statically sized container. More...
 
class  HierarchicalLagrangeNode
 
class  HierarchicalLagrangePreBasis
 A pre-basis for a hierarchical basis. More...
 
struct  HierarchicNodeToRangeMap
 A simple node to range map using the nested tree indices. More...
 
class  HierarchicVectorWrapper
 A wrapper providing multiindex access to vector entries. More...
 
class  InvalidRange
 Dummy range class to be used if no proper type is available. More...
 
struct  IsCallable
 Helper class to check that F is callable. More...
 
class  LagrangeDGPreBasis
 
class  LagrangeNode
 
class  LagrangePreBasis
 A pre-basis for a PQ-lagrange bases with given order. More...
 
struct  LastType
 Get last entry of type list. More...
 
class  LeafBasisNode
 
struct  LocalDerivativeTraits
 Derivative traits for local functions. More...
 
class  LocalFunction
 
class  LocalFunction< Range(Domain), LocalContext, DerivativeTraits, bufferSize >
 Class storing local functions using type erasure. More...
 
class  NedelecNode
 
class  NedelecPreBasis
 
class  OverflowArray
 A dynamically sized array-like class with overflow. More...
 
class  PolymorphicSmallObject
 A wrapper providing small object optimization with polymorphic types. More...
 
class  PolymorphicType
 Base class with polymorphic type boiler plate code. More...
 
class  Polynomial
 A scalar polynomial implementation. More...
 
class  PowerBasisNode
 
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  ReservedDeque
 A double-ended queue (deque) class with statically reserved memory. More...
 
struct  RotateTuple
 Rotate type list by one, such that last entry is moved to first position. More...
 
struct  SignatureTag
 
struct  SignatureTag< Range(Domain), DerivativeTraitsT >
 Tag-class to encapsulate signature information. More...
 
struct  SignatureTraits
 Helper class to deduce the signature of a callable. More...
 
class  SizeInfo
 A class encapsulating size information. More...
 
class  StaticMultiIndex
 A statically sized MultiIndex type. More...
 
class  StaticMultiIndex< size_type, 1 >
 A statically sized MultiIndex type. More...
 
struct  StaticSize
 Obtain size of statically sized container. 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  TaylorHoodBasisTree
 
class  TaylorHoodPreBasis
 Pre-basis for lowest order Taylor-Hood basis. More...
 
class  TaylorHoodVelocityTree
 
class  TreeData
 Container allowing to attach data to each node of a tree. More...
 
class  TrigonometricFunction
 A linear combination of trigonomic functions. More...
 
class  TypeErasureBase
 Base class for type-erased interface wrapper. More...
 
struct  UniformNodeVisitor
 Mixin for visitors that should apply the same action on all nodes. More...
 

Typedefs

template<class T >
using ResolveRef_t = Dune::ResolveRef_t< T >
 This is an alias for Dune::ResolveRef_t. More...
 
template<class T , class... Args>
using enableIfConstructible = typename std::enable_if< std::is_constructible< T, Args... >::value, int >::type
 Helper to constrain forwarding constructors. More...
 
template<template< class... > class T, class ArgTuple >
using ExpandTuple = typename Imp::ExpandTupleHelper< T, ArgTuple >::Type
 Expand tuple arguments as template arguments. More...
 
template<template< class... > class F, class... Tuples>
using TransformTuple = typename Imp::TransformTupleHelper< F, Tuples... >::Type
 Transform tuple types argument using type-functor. More...
 
template<class IntegerSequence >
using IntegerSequenceTuple = typename Imp::IntegerSequenceTupleHelper< IntegerSequence >::Type
 Transform integer_sequence<I,k...> to tuple<integral_constant<I,k>...> More...
 
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. More...
 
template<typename GV >
using BSplineBasis = DefaultGlobalBasis< BSplinePreBasis< GV > >
 A global B-spline basis. More...
 
template<class size_type >
using FlatMultiIndex = StaticMultiIndex< size_type, 1 >
 A multi-index class with only one level. More...
 
template<typename GV , int k, typename R = double>
using HierarchicalLagrangeBasis = DefaultGlobalBasis< HierarchicalLagrangePreBasis< GV, k, R > >
 Basis of a scalar Hierarchical Lagrange finite element space. More...
 
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. More...
 
template<typename GV , int k>
using LagrangeDGNode = LagrangeNode< GV, k >
 
template<typename GV , int k>
using LagrangeDGBasis = DefaultGlobalBasis< LagrangeDGPreBasis< GV, k > >
 Basis of a scalar k-th-order Lagrangean-DG finite element space. More...
 
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. More...
 
template<typename GV >
using RannacherTurekBasis = DefaultGlobalBasis< RannacherTurekPreBasis< GV > >
 Rannacher-Turek basis. More...
 
template<typename GV , int k>
using RaviartThomasBasis = DefaultGlobalBasis< RaviartThomasPreBasis< GV, k > >
 Basis of a k-th-order Raviart Thomas finite element space. More...
 
template<typename GV >
using TaylorHoodBasis = DefaultGlobalBasis< TaylorHoodPreBasis< GV > >
 Nodal basis for a lowest order Taylor-Hood Lagrangean finite element space. More...
 

Functions

template<class K , int sinFactor, int cosFactor>
TrigonometricFunction< K, -cosFactor, sinFactor > derivative (const TrigonometricFunction< K, sinFactor, cosFactor > &f)
 Obtain derivative of TrigonometricFunction function. More...
 
template<class V >
constexpr auto fieldTypes ()
 Generate list of field types in container. More...
 
template<class V >
constexpr bool hasUniqueFieldType ()
 Check if container has a unique field type. More...
 
template<class Vector >
auto istlVectorBackend (Vector &v)
 Return a vector backend wrapping non-const ISTL like containers. More...
 
template<class Vector >
auto istlVectorBackend (const Vector &v)
 Return a vector backend wrapping const ISTL like containers. More...
 
template<class F >
CallableFunctionWrapper< F > callable (const F &f)
 Create a callable object from some Dune::VirtualFunction. More...
 
template<class F >
CallableFunctionWrapper< F > callable (const std::shared_ptr< F > &fp)
 Create a callable object from std::shared_ptr<F> More...
 
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. More...
 
template<class C , class I , class F , typename std::enable_if< Dune::models< Imp::Concept::HasDynamicIndexAccess< I >, C >(), int >::type = 0>
auto hybridIndexAccess (C &&c, const I &i, F &&f) -> decltype(f(c[i]))
 Provide operator[] index-access for containers. More...
 
template<class C , class I , class F , typename std::enable_if< not Dune::models< Imp::Concept::HasDynamicIndexAccess< I >, C >(), int >::type = 0>
decltype(auto) hybridIndexAccess (C &&c, const I &i, F &&f)
 Provide operator[] index-access for containers. More...
 
template<class Result , class C , class MultiIndex >
Result hybridMultiIndexAccess (C &&c, const MultiIndex &index)
 Provide multi-index access by chaining operator[]. More...
 
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[]. More...
 
template<class C , class MultiIndex >
constexpr decltype(auto) resolveDynamicMultiIndex (C &&c, const MultiIndex &multiIndex)
 Provide multi-index access by chaining operator[]. More...
 
template<class C , class MultiIndex >
constexpr decltype(auto) resolveStaticMultiIndex (C &&c, const MultiIndex &multiIndex)
 Provide multi-index access by chaining operator[]. More...
 
template<typename Stream , class size_type , std::size_t n>
Stream & operator<< (Stream &stream, const StaticMultiIndex< size_type, n > &c)
 
template<class T >
decltype(auto) resolveRef (T &&t)
 This is an alias for Dune::resolveRef. More...
 
template<class Range , class Domain , template< class > class DerivativeTraits>
auto derivativeSignatureTag (SignatureTag< Range(Domain), DerivativeTraits > tag)
 Construct SignatureTag for derivative. More...
 
template<std::size_t maxOrder, class Signature , template< class > class DerivativeTraits>
auto derivativeSignatureTags (Dune::Functions::SignatureTag< Signature, DerivativeTraits > tag)
 Construct SignatureTags for derivatives. More...
 
template<std::size_t begin_t, std::size_t end_t, class F , class... Args>
void staticFindInRange (F &&f, Args &&... args)
 Static find loop. More...
 
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. More...
 
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. More...
 
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. More...
 
template<class Expression >
auto callableCheck (Expression f)
 Create a predicate for checking validity of expressions. More...
 
template<class Check >
auto negatePredicate (Check check)
 Negate given predicate. More...
 
template<class T >
auto forwardCapture (T &&t)
 Create a capture object for perfect forwarding. More...
 
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. More...
 
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 Tree >
DefaultNodeToRangeMap< Tree > makeDefaultNodeToRangeMap (const Tree &tree)
 
template<class Basis , class TreePath >
auto makeDefaultNodeToRangeMap (const Basis &basis, TreePath &&treePath) -> decltype(makeDefaultNodeToRangeMap(TypeTree::child(basis.localView().tree(), treePath)))
 
template<class T >
auto flatVectorView (T &t)
 Create flat vector view of passed mutable container. More...
 
template<class T >
auto flatVectorView (const T &t)
 Create flat vector view of passed const container. More...
 
template<class T >
auto flatVectorView (T &&t)
 Create flat vector view of passed container temporary. More...
 
template<class V >
HierarchicVectorWrapper< V > hierarchicVector (V &v)
 
template<class MultiIndex , class V , typename std::enable_if< models< Concept::HasIndexAccess, V, MultiIndex >(), int >::type = 0>
V & makeHierarchicVectorForMultiIndex (V &v)
 
template<class MultiIndex , class V , typename std::enable_if< not models< Concept::HasIndexAccess, V, MultiIndex >(), int >::type = 0>
HierarchicVectorWrapper< V > makeHierarchicVectorForMultiIndex (V &v)
 
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. More...
 
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. More...
 
template<class B , class C , class F >
void interpolate (const B &basis, C &&coeff, const F &f)
 Interpolate given function in discrete function space. More...
 
template<typename Tree >
void clearSize (Tree &tree, std::size_t offset)
 
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 Basis >
SizeInfo< Basis > sizeInfo (const Basis &basis)
 
template<class T >
auto subEntityDOFs (const T &)
 Create SubEntityDOFs object. More...
 
template<class LocalView >
auto subEntityDOFs (const LocalView &localView, std::size_t subEntityIndex, std::size_t subEntityCodim)
 Create bound SubEntityDOFs object. More...
 
template<class LocalView , class Intersection >
auto subEntityDOFs (const LocalView &localView, const Intersection &intersection)
 Create bound SubEntityDOFs object. More...
 
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())>, Impl::JoinTreePath_t< InnerTP, OuterTP > >
 
template<class RootBasis , class... PrefixTreeIndices>
auto subspaceBasis (const RootBasis &rootBasis, const TypeTree::HybridTreePath< PrefixTreeIndices... > &prefixPath)
 Create SubspaceBasis from a root basis and a prefixPath. More...
 
template<class RootBasis , class... PrefixTreeIndices>
auto subspaceBasis (const RootBasis &rootBasis, const PrefixTreeIndices &... prefixTreeIndices)
 
template<class F , class GridView >
AnalyticGridViewFunction< typename std::invoke_result< F, typename GridView::template Codim< 0 >::Geometry::GlobalCoordinate >::type(typename GridView::template Codim< 0 >::Geometry::GlobalCoordinate), GridView, typename std::decay< F >::type > makeAnalyticGridViewFunction (F &&f, const GridView &gridView)
 Create an AnalyticGridViewFunction from a function and a grid view. More...
 
template<class OF , class... IF>
auto makeComposedGridFunction (OF &&outerFunction, IF &&... innerFunction)
 Create a ComposedGridFunction that composes grid-functions with another function. More...
 
template<typename R , typename B , typename V >
auto makeDiscreteGlobalBasisFunction (B &&basis, V &&vector)
 Generate a DiscreteGlobalBasisFunction. More...
 
template<class F , class GridView , typename std::enable_if< models< Imp::HasFreeLocalFunction, F >(), int >::type = 0>
std::decay< F >::type makeGridViewFunction (F &&f, const GridView &gridView)
 Construct a function modeling GridViewFunction from function and grid view. More...
 
template<class F , class GridView , typename std::enable_if< not(models< Imp::HasFreeLocalFunction, F >()), int >::type = 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. More...
 

Variables

template<class T >
constexpr bool IsReferenceWrapper_v = Dune::IsReferenceWrapper_v<T>
 This is an alias for Dune::IsReferenceWrapper_v. More...
 

Typedef Documentation

◆ BrezziDouglasMariniBasis

template<typename GV , int k>
using Dune::Functions::BrezziDouglasMariniBasis = typedef DefaultGlobalBasis<BrezziDouglasMariniPreBasis<GV, k> >

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 = typedef typename Imp::IntegerSequenceTupleHelper<IntegerSequence>::Type

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

◆ LagrangeDGNode

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

◆ NedelecBasis

template<typename GV , std::size_t kind, std::size_t order, typename Range = double>
using Dune::Functions::NedelecBasis = typedef 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 = typedef 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

◆ ResolveRef_t

template<class T >
using Dune::Functions::ResolveRef_t = typedef Dune::ResolveRef_t<T>

This is an alias for Dune::ResolveRef_t.

Deprecated:
Use Dune::resolveRef instead

Function Documentation

◆ bindTree()

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

◆ clearSize()

template<typename Tree >
void Dune::Functions::clearSize ( Tree &  tree,
std::size_t  offset 
)

◆ 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 >
constexpr 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.

◆ 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 >
constexpr bool Dune::Functions::hasUniqueFieldType ( )
constexpr

Check if container has a unique field type.

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

◆ hierarchicVector()

template<class V >
HierarchicVectorWrapper< V > Dune::Functions::hierarchicVector ( V &  v)

◆ hybridIndexAccess()

template<class C , class I , class F , typename std::enable_if< not Dune::models< Imp::Concept::HasDynamicIndexAccess< I >, C >(), int >::type = 0>
decltype(auto) Dune::Functions::hybridIndexAccess ( C &&  c,
const I &  i,
F &&  f 
)

Provide operator[] index-access for containers.

This is the overload for types providing a operator[] only for static arguments of type std::integral_constant<std::size_t,k>. This does a static linear search until a static index matching the given dynamic index is found. Since the result type will in general be different for different indices the method does not return the result directly but passes it to a given functor.

Parameters
cContainer to access
iThe index to use for accessing the container
fA functor to call with the result of operator[]

◆ 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
bitVectorA vector with flags marking all DOFs that should be interpolated
nodeToRangeEntryPolymorphic functor mapping local ansatz nodes to range-indices of given function

◆ makeDefaultNodeToRangeMap() [1/2]

template<class Basis , class TreePath >
auto Dune::Functions::makeDefaultNodeToRangeMap ( const Basis &  basis,
TreePath &&  treePath 
) -> decltype(makeDefaultNodeToRangeMap(TypeTree::child(basis.localView().tree(),treePath)))

◆ makeDefaultNodeToRangeMap() [2/2]

template<class Tree >
DefaultNodeToRangeMap< Tree > Dune::Functions::makeDefaultNodeToRangeMap ( const Tree &  tree)

◆ makeGridViewFunction() [1/2]

template<class F , class GridView , typename std::enable_if< models< Imp::HasFreeLocalFunction, F >(), int >::type = 0>
std::decay< F >::type 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 , typename std::enable_if< not(models< Imp::HasFreeLocalFunction, F >()), int >::type = 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.

◆ makeHierarchicVectorForMultiIndex() [1/2]

template<class MultiIndex , class V , typename std::enable_if< models< Concept::HasIndexAccess, V, MultiIndex >(), int >::type = 0>
V & Dune::Functions::makeHierarchicVectorForMultiIndex ( V &  v)

◆ makeHierarchicVectorForMultiIndex() [2/2]

template<class MultiIndex , class V , typename std::enable_if< not models< Concept::HasIndexAccess, V, MultiIndex >(), int >::type = 0>
HierarchicVectorWrapper< V > Dune::Functions::makeHierarchicVectorForMultiIndex ( V &  v)

◆ operator<<()

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

◆ resolveRef()

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

This is an alias for Dune::resolveRef.

Deprecated:
Use Dune::resolveRef instead

◆ sizeInfo()

template<class Basis >
SizeInfo< Basis > Dune::Functions::sizeInfo ( const Basis &  basis)

◆ SubspaceBasis() [1/2]

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

◆ 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::HybridTreePath< 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

◆ 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())>, Impl::JoinTreePath_t< InnerTP, OuterTP > >

Variable Documentation

◆ IsReferenceWrapper_v

template<class T >
constexpr bool Dune::Functions::IsReferenceWrapper_v = Dune::IsReferenceWrapper_v<T>
constexpr

This is an alias for Dune::IsReferenceWrapper_v.

Deprecated:
Use Dune::IsReferenceWrapper_v instead