9#ifndef LF_USCALFE_FE_SPACE_HP_H_
10#define LF_USCALFE_FE_SPACE_HP_H_
12#include <lf/mesh/utils/all_codim_mesh_data_set.h>
14#include "hierarchic_fe.h"
15#include "scalar_fe_space.h"
42template <
typename SCALAR>
60 const std::shared_ptr<const
lf::mesh::
Mesh> &mesh_p,
unsigned degree)
62 [degree](const mesh::Entity &e) ->
unsigned {
86 template <
class F,
class = std::enable_if_t<
87 std::is_invocable_v<F, const mesh::Entity &>>>
89 const std::shared_ptr<const lf::mesh::Mesh> &mesh_p, F &°ree_functor)
94 [°ree_functor](const mesh::Entity &e) -> base::
size_type {
95 if (e.RefEl() == base::RefEl::kPoint()) {
98 auto degree = degree_functor(e);
103 return degree <= 2 ? 0 : (degree - 2) * (degree - 1) / 2;
105 return (degree - 1) * (degree - 1);
107 LF_VERIFY_MSG(
false,
"Something went wrong.");
110 Init(std::forward<F>(degree_functor));
116 [[nodiscard]] std::shared_ptr<const lf::mesh::Mesh>
Mesh()
const override {
117 LF_VERIFY_MSG(mesh_p_ !=
nullptr,
"No valid FE space object: no mesh");
125 LF_VERIFY_MSG(Mesh() !=
nullptr,
"No valid FE space object: no mesh");
136 if constexpr (std::is_same_v<
decltype(fe),
137 const std::monostate &>) {
139 "Something is wrong, no "
140 "ScalarReferenceFiniteElement found for "
155 return ShapeFunctionLayout(entity)->NumRefShapeFunctions();
162 std::shared_ptr<const lf::mesh::Mesh>
mesh_p_;
173 void Init(F &°ree_functor) {
176 for (
auto entity : Mesh()->Entities(2)) {
180 for (
auto entity : Mesh()->Entities(1)) {
183 ref_el_(*entity) = std::move(fe);
186 for (
auto entity : Mesh()->Entities(0)) {
187 switch (entity->RefEl()) {
189 std::array<unsigned, 3> edge_degrees{
190 {degree_functor(*entity->SubEntities(1)[0]),
191 degree_functor(*entity->SubEntities(1)[1]),
192 degree_functor(*entity->SubEntities(1)[2])}};
195 entity->RelativeOrientations());
196 ref_el_(*entity) = std::move(fe);
200 std::array<unsigned, 4> edge_degrees{
201 {degree_functor(*entity->SubEntities(1)[0]),
202 degree_functor(*entity->SubEntities(1)[1]),
203 degree_functor(*entity->SubEntities(1)[2]),
204 degree_functor(*entity->SubEntities(1)[3])}};
207 entity->RelativeOrientations());
208 ref_el_(*entity) = std::move(fe);
212 LF_VERIFY_MSG(
false,
"Illegal entity type");
A general (interface) class for DOF handling, see Lecture Document Paragraph 2.7.4....
Dof handler allowing variable local dof layouts.
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
static constexpr RefEl kTria()
Returns the reference triangle.
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Hierarchic Finite Elements of arbitrary degree on quadrilaterals.
Hierarchic Finite Elements of arbitrary degree on segments.
Hierarchic Finite Elements of arbitrary degree on triangles.
Linear finite element on a point.
Finite Element Space that supports arbitrary, local degrees.
const lf::assemble::DofHandler & LocGlobMap() const override
access to associated local-to-global map
size_type NumRefShapeFunctions(const lf::mesh::Entity &entity) const override
number of interior shape functions associated to entities of various types
lf::mesh::utils::AllCodimMeshDataSet< std::variant< std::monostate, FePoint< SCALAR >, FeHierarchicSegment< SCALAR >, FeHierarchicTria< SCALAR >, FeHierarchicQuad< SCALAR > > > ref_el_
std::shared_ptr< const lf::mesh::Mesh > mesh_p_
HierarchicScalarFESpace(HierarchicScalarFESpace &&) noexcept=default
ScalarReferenceFiniteElement< SCALAR > const * ShapeFunctionLayout(const lf::mesh::Entity &entity) const override
access to shape function layout for cells
quad::QuadRuleCache qr_cache_
HierarchicScalarFESpace()=delete
HierarchicScalarFESpace(const HierarchicScalarFESpace &)=delete
lf::assemble::DynamicFEDofHandler dofh_
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
access to underlying mesh
void Init(F &°ree_functor)
HierarchicScalarFESpace(const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, F &°ree_functor)
Construct a new Hierarchic FESpace with (possibly) varying polynomial degrees.
~HierarchicScalarFESpace() override=default
Space of scalar valued finite element functions on a hybrid 2D mesh
Interface class for parametric scalar valued finite elements.
Interface class representing a topological entity in a cellular complex
Assigns to every entity(all codims) in a mesh a value of type T
A cache for make_QuadRule()
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
lf::assemble::size_type size_type