Dune::Copasi
|
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. | |
This class describes a local basis cache.
This class cache results of the local basis for different local basis with same traits.
LocalBasisTraits | Local basis traits |
|
inline |
Constructs a new instance with empty values.
|
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.
dim | Dimension of the quadrature rule to apply |
finite_element | Finite element to apply cache to |
quad_rule | Quadrature rule to evaluate the finite element |
quad_proj | Quadrature positions projection into the finite element reference element |
force | Whether to force a recomputation of the results |
|
inline |
Binds a finite element to this cache.
|
inline |
Evaluate function with the bound finite element and quadrature rule.
[in] | position_id | The position id wrt quadrature rule |
|
inline |
Evaluate jacobian with the bound finite element and quadrature rule.
[in] | position_id | The position id wrt quadrature rule |
|
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:
dim | Dimension of the quadrature rule to apply |
finite_element | Finite element to apply cache to |
quad_rule | Quadrature rule to evaluate the finite element |
quad_proj | Quadrature positions projection into the finite element reference element |
|
inlinenoexcept |
Unbind finite element from this cache.