Dune::Copasi
Loading...
Searching...
No Matches
Classes | Public Member Functions | List of all members
Dune::Copasi::LocalBasisCache< LocalBasisTraits > Class Template Reference

This class describes a local basis cache. More...

#include <local_basis_cache.hh>

Classes

struct  Key
 Key for the local basis cache. More...
 

Public Member Functions

 LocalBasisCache ()
 Constructs a new instance with empty values.
 
template<int dim>
Key key (const auto &finite_element, const Dune::QuadratureRule< DomainField, dim > &quad_rule, auto &&quad_proj) const noexcept
 Key to match previews applications of this cache.
 
void bind (const auto &finite_element, const Dune::QuadratureRule< DomainField, dimDomain > &quad_rule, bool force=false)
 Binds a finite element to this cache.
 
template<int dim, typename ProjQuad >
requires std::is_invocable_r_v<Domain, ProjQuad, FieldVector<DomainField, dim>>
void bind (const auto &finite_element, const Dune::QuadratureRule< DomainField, dim > &quad_rule, ProjQuad &&quad_proj, bool force=false)
 Binds a finite element to this cache.
 
void unbind () noexcept
 Unbind finite element from this cache.
 
std::span< const Range > evaluateFunction (std::size_t position_id) const
 Evaluate function with the bound finite element and quadrature rule.
 
std::span< const Jacobian > evaluateJacobian (std::size_t position_id) const
 Evaluate jacobian with the bound finite element and quadrature rule.
 

Detailed Description

template<class LocalBasisTraits>
class Dune::Copasi::LocalBasisCache< LocalBasisTraits >

This class describes a local basis cache.

This class cache results of the local basis for different local basis with same traits.

Template Parameters
LocalBasisTraitsLocal basis traits

Constructor & Destructor Documentation

◆ LocalBasisCache()

template<class LocalBasisTraits >
Dune::Copasi::LocalBasisCache< LocalBasisTraits >::LocalBasisCache ( )
inline

Constructs a new instance with empty values.

Member Function Documentation

◆ bind() [1/2]

template<class LocalBasisTraits >
template<int dim, typename ProjQuad >
requires std::is_invocable_r_v<Domain, ProjQuad, FieldVector<DomainField, dim>>
void Dune::Copasi::LocalBasisCache< LocalBasisTraits >::bind ( const auto & finite_element,
const Dune::QuadratureRule< DomainField, dim > & quad_rule,
ProjQuad && quad_proj,
bool force = false )
inline

Binds a finite element to this cache.

The values of evaluateFunction and evaluateJacobian of the local basis in the finite element will be reused wrt a previous usage of this cache if the key has already been used. See key function for details. The force parameter will force recomputation of the values. This is useful if the key cannot uniquely the local basis evaluation. In either case the values resulting from evaluateFunction and evaluateJacobian will be stored within this cache. The quadrature projection is useful when the quadrature rule belongs to a sub-entity or an intersection of the finite element reference element.

Template Parameters
dimDimension of the quadrature rule to apply
Parameters
finite_elementFinite element to apply cache to
quad_ruleQuadrature rule to evaluate the finite element
quad_projQuadrature positions projection into the finite element reference element
forceWhether to force a recomputation of the results

◆ bind() [2/2]

template<class LocalBasisTraits >
void Dune::Copasi::LocalBasisCache< LocalBasisTraits >::bind ( const auto & finite_element,
const Dune::QuadratureRule< DomainField, dimDomain > & quad_rule,
bool force = false )
inline

Binds a finite element to this cache.

◆ evaluateFunction()

template<class LocalBasisTraits >
std::span< const Range > Dune::Copasi::LocalBasisCache< LocalBasisTraits >::evaluateFunction ( std::size_t position_id) const
inline

Evaluate function with the bound finite element and quadrature rule.

Parameters
[in]position_idThe position id wrt quadrature rule
Returns
The evaluation for each of the local basis at position_id

◆ evaluateJacobian()

template<class LocalBasisTraits >
std::span< const Jacobian > Dune::Copasi::LocalBasisCache< LocalBasisTraits >::evaluateJacobian ( std::size_t position_id) const
inline

Evaluate jacobian with the bound finite element and quadrature rule.

Parameters
[in]position_idThe position id wrt quadrature rule
Returns
The jacobian for each of the local basis at position_id

◆ key()

template<class LocalBasisTraits >
template<int dim>
Key Dune::Copasi::LocalBasisCache< LocalBasisTraits >::key ( const auto & finite_element,
const Dune::QuadratureRule< DomainField, dim > & quad_rule,
auto && quad_proj ) const
inlinenoexcept

Key to match previews applications of this cache.

The key is composed with the std::type_index of the finite element, the address of the quadrature rule, and the projected center of the quadrature geometry type. This means that:

  • Finite elements where shape functions change dynamically cannot be uniquely represented.
  • Quadrature rules on non-conforming intersections cannot be uniquely represented.
Template Parameters
dimDimension of the quadrature rule to apply
Parameters
finite_elementFinite element to apply cache to
quad_ruleQuadrature rule to evaluate the finite element
quad_projQuadrature positions projection into the finite element reference element
Returns
Key Key with the parameter ids

◆ unbind()

template<class LocalBasisTraits >
void Dune::Copasi::LocalBasisCache< LocalBasisTraits >::unbind ( )
inlinenoexcept

Unbind finite element from this cache.


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