LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
hierarchic_scalar_fe_space.h
1
9#ifndef LF_USCALFE_FE_SPACE_HP_H_
10#define LF_USCALFE_FE_SPACE_HP_H_
11
12#include <lf/mesh/utils/all_codim_mesh_data_set.h>
13
14#include "hierarchic_fe.h"
15#include "scalar_fe_space.h"
16
17namespace lf::fe {
18
42template <typename SCALAR>
43class HierarchicScalarFESpace : public ScalarFESpace<SCALAR> {
44 public:
45 using Scalar = SCALAR;
46
52 default;
53
60 const std::shared_ptr<const lf::mesh::Mesh> &mesh_p, unsigned degree)
62 [degree](const mesh::Entity &e) -> unsigned {
63 if (e.RefEl() == base::RefEl::kPoint()) {
64 return 1;
65 }
66 return degree;
67 }) {}
68
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 &&degree_functor)
90 : ScalarFESpace<SCALAR>(),
91 mesh_p_(mesh_p),
92 ref_el_(mesh_p),
93 dofh_(mesh_p,
94 [&degree_functor](const mesh::Entity &e) -> base::size_type {
95 if (e.RefEl() == base::RefEl::kPoint()) {
96 return 1;
97 }
98 auto degree = degree_functor(e);
99 switch (e.RefEl()) {
101 return degree - 1;
102 case base::RefEl::kTria():
103 return degree <= 2 ? 0 : (degree - 2) * (degree - 1) / 2;
104 case base::RefEl::kQuad():
105 return (degree - 1) * (degree - 1);
106 default:
107 LF_VERIFY_MSG(false, "Something went wrong.");
108 }
109 }) {
110 Init(std::forward<F>(degree_functor));
111 }
112
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");
118 return mesh_p_;
119 }
120
124 [[nodiscard]] const lf::assemble::DofHandler &LocGlobMap() const override {
125 LF_VERIFY_MSG(Mesh() != nullptr, "No valid FE space object: no mesh");
126 return dofh_;
127 }
128
133 const lf::mesh::Entity &entity) const override {
134 return std::visit(
135 [](const auto &fe) -> ScalarReferenceFiniteElement<SCALAR> const * {
136 if constexpr (std::is_same_v<decltype(fe),
137 const std::monostate &>) { // NOLINT
138 LF_VERIFY_MSG(false,
139 "Something is wrong, no "
140 "ScalarReferenceFiniteElement found for "
141 "this entity.");
142 return nullptr; // make compiler happy.
143 } else { // NOLINT
144 return &fe;
145 }
146 },
147 ref_el_(entity));
148 }
149
154 const lf::mesh::Entity &entity) const override {
155 return ShapeFunctionLayout(entity)->NumRefShapeFunctions();
156 }
157
158 ~HierarchicScalarFESpace() override = default;
159
160 private:
161 /* Underlying mesh */
162 std::shared_ptr<const lf::mesh::Mesh> mesh_p_;
164 // Store the ScalarReferenceFiniteElement for every entity.
166 std::variant<std::monostate, FePoint<SCALAR>, FeHierarchicSegment<SCALAR>,
170
171 // R.H. Detailed comment needed for this function
172 template <class F>
173 void Init(F &&degree_functor) {
174 // Initialize all shape function layouts for nodes
175 size_type num_rsf_node = 1;
176 for (auto entity : Mesh()->Entities(2)) {
177 ref_el_(*entity) = FePoint<SCALAR>(1);
178 }
179 // Initialize all shape function layouts for the edges
180 for (auto entity : Mesh()->Entities(1)) {
181 FeHierarchicSegment<SCALAR> fe(degree_functor(*entity), qr_cache_);
182
183 ref_el_(*entity) = std::move(fe);
184 }
185 // Initialize all shape function layouts for the cells
186 for (auto entity : Mesh()->Entities(0)) {
187 switch (entity->RefEl()) {
188 case lf::base::RefEl::kTria(): {
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])}};
193 FeHierarchicTria<SCALAR> fe(degree_functor(*entity), edge_degrees,
194 qr_cache_,
195 entity->RelativeOrientations());
196 ref_el_(*entity) = std::move(fe);
197 break;
198 }
199 case lf::base::RefEl::kQuad(): {
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])}};
205 FeHierarchicQuad<SCALAR> fe(degree_functor(*entity), edge_degrees,
206 qr_cache_,
207 entity->RelativeOrientations());
208 ref_el_(*entity) = std::move(fe);
209 break;
210 }
211 default:
212 LF_VERIFY_MSG(false, "Illegal entity type");
213 }
214 }
215 }
216};
217
218} // end namespace lf::fe
219
220#endif // LF_USCALFE_FE_SPACE_HP_H_
A general (interface) class for DOF handling, see Lecture Document Paragraph 2.7.4....
Definition: dofhandler.h:109
Dof handler allowing variable local dof layouts.
Definition: dofhandler.h:472
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition: ref_el.h:150
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition: ref_el.h:141
static constexpr RefEl kTria()
Returns the reference triangle.
Definition: ref_el.h:158
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition: ref_el.h:166
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.
Definition: fe_point.h:25
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
HierarchicScalarFESpace(const HierarchicScalarFESpace &)=delete
lf::assemble::DynamicFEDofHandler dofh_
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
access to underlying mesh
HierarchicScalarFESpace(const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, F &&degree_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
Definition: entity.h:39
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...
Definition: fe.h:47
lf::assemble::size_type size_type
Definition: fe.h:54
Definition: assemble.h:30