Dune Core Modules (2.8.0)

Dune::LocalFiniteElementVariant< Implementations > Class Template Reference

Type erasure class for wrapping LocalFiniteElement classes. More...

#include <dune/localfunctions/common/localfiniteelementvariant.hh>

Public Types

using Traits = typename Dune::LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation >
 Export LocalFiniteElementTraits.
 

Public Member Functions

 LocalFiniteElementVariant ()=default
 Construct empty LocalFiniteElementVariant.
 
 LocalFiniteElementVariant (const std::monostate &monostate)
 Construct empty LocalFiniteElementVariant.
 
template<class Implementation , std::enable_if_t< std::disjunction< std::is_same< std::decay_t< Implementation >, Implementations >... >::value, int > = 0>
 LocalFiniteElementVariant (Implementation &&impl)
 Construct LocalFiniteElementVariant. More...
 
 LocalFiniteElementVariant (const LocalFiniteElementVariant &other)
 Copy constructor.
 
 LocalFiniteElementVariant (LocalFiniteElementVariant &&other)
 Move constructor.
 
LocalFiniteElementVariantoperator= (const LocalFiniteElementVariant &other)
 Copy assignment.
 
LocalFiniteElementVariantoperator= (LocalFiniteElementVariant &&other)
 Move assignment.
 
template<class Implementation , std::enable_if_t< std::disjunction< std::is_same< std::decay_t< Implementation >, Implementations >... >::value, int > = 0>
LocalFiniteElementVariantoperator= (Implementation &&impl)
 Assignment from implementation.
 
const Traits::LocalBasisType & localBasis () const
 Provide access to LocalBasis implementation of this LocalFiniteElement.
 
const Traits::LocalCoefficientsType & localCoefficients () const
 Provide access to LocalCoefficients implementation of this LocalFiniteElement.
 
const Traits::LocalInterpolationType & localInterpolation () const
 Provide access to LocalInterpolation implementation of this LocalFiniteElement.
 
unsigned int size () const
 Number of shape functions.
 
constexpr GeometryType type () const
 Number of shape functions.
 
const auto & variant () const
 Provide access to underlying std::variant. More...
 
 operator bool () const
 Check if LocalFiniteElementVariant stores an implementation. More...
 

Detailed Description

template<class... Implementations>
class Dune::LocalFiniteElementVariant< Implementations >

Type erasure class for wrapping LocalFiniteElement classes.

This is a type erasure wrapper class for types implementing the LocalFiniteElement interface. The types of the LocalFiniteElement implementations that this class can hold have to be provided as template parameter.

The implementation is based on std::variant. Notice that this prepends std::monostate to the Implementations list for the internally stored std::variant such that LocalFiniteElementVariant can be empty and is default-constructible. As a consequence providing std::monostate manually to LocalFiniteElementVariant is neither necessary nor allowed. Access to the stored implementation is internally implemented using std::visit(). To avoid multiple trivial std::visit() calls, the results of size(), order(), and type() are cached on creation and assignment.

In empty state accessing any method beyond operator bool(), variant(), or assignment leads to undefined behavior.

The LocalBasisTraits are extracted from the implementation provided as first template parameter. The other implementations are required to be compatible with this one.

Template Parameters
ImplementationsList of supported LocalFiniteElement implementations

Constructor & Destructor Documentation

◆ LocalFiniteElementVariant()

template<class... Implementations>
template<class Implementation , std::enable_if_t< std::disjunction< std::is_same< std::decay_t< Implementation >, Implementations >... >::value, int > = 0>
Dune::LocalFiniteElementVariant< Implementations >::LocalFiniteElementVariant ( Implementation &&  impl)
inline

Construct LocalFiniteElementVariant.

The created LocalFiniteElementVariant will store a copy of the provided implementation.

Member Function Documentation

◆ operator bool()

template<class... Implementations>
Dune::LocalFiniteElementVariant< Implementations >::operator bool ( ) const
inline

Check if LocalFiniteElementVariant stores an implementation.

This returns true iff variant() does not store a monostate.

References Dune::LocalFiniteElementVariant< Implementations >::variant().

◆ variant()

template<class... Implementations>
const auto & Dune::LocalFiniteElementVariant< Implementations >::variant ( ) const
inline

Provide access to underlying std::variant.

This allows to use std::visit on a higher level which allows to avoid the indirection of the std::variant - polymorphism inside the visitor code. Notice that the provided std::variant contains std::monostate in its type list. Hence any visitor used to access the variant has to be std::monostate-aware.

Referenced by Dune::LocalFiniteElementVariant< Implementations >::operator bool().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)