LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
uniform_scalar_fe_space.h
1#ifndef LF_FESPACE_H
2#define LF_FESPACE_H
3/***************************************************************************
4 * LehrFEM++ - A simple C++ finite element libray for teaching
5 * Developed from 2018 at the Seminar of Applied Mathematics of ETH Zurich,
6 * lead developers Dr. R. Casagrande and Prof. R. Hiptmair
7 ***************************************************************************/
8
17#include <lf/assemble/assemble.h>
18#include <lf/fe/scalar_fe_space.h>
19
20#include "lagr_fe.h"
21
22namespace lf::uscalfe {
23
49template <typename SCALAR>
51 public:
52 using Scalar = SCALAR;
53
56 UniformScalarFESpace &operator=(const UniformScalarFESpace &) = delete;
57 UniformScalarFESpace &operator=(UniformScalarFESpace &&) noexcept = default;
83 std::shared_ptr<const lf::mesh::Mesh> mesh_p,
84 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
85 rfs_tria_p,
86 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
87 rfs_quad_p,
88 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
89 rfs_edge_p = nullptr,
90 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
91 rfs_point_p = nullptr)
92 : lf::fe::ScalarFESpace<SCALAR>(),
93 rfs_tria_p_(std::move(rfs_tria_p)),
94 rfs_quad_p_(std::move(rfs_quad_p)),
95 rfs_edge_p_(std::move(rfs_edge_p)),
96 rfs_point_p_(std::move(rfs_point_p)),
97 dofh_(InitDofHandler(std::move(mesh_p))) {}
98
102 [[nodiscard]] std::shared_ptr<const lf::mesh::Mesh> Mesh() const override {
103 return dofh_.Mesh();
104 }
105
109 [[nodiscard]] const lf::assemble::DofHandler &LocGlobMap() const override {
110 LF_VERIFY_MSG(Mesh() != nullptr, "No valid FE space object: no mesh");
111 return dofh_;
112 }
113
118 ShapeFunctionLayout(const lf::mesh::Entity &entity) const override;
119
133 ShapeFunctionLayout(lf::base::RefEl ref_el_type) const;
134
138 [[nodiscard]] size_type NumRefShapeFunctions(
139 const lf::mesh::Entity &entity) const override;
140
144 [[nodiscard]] size_type NumRefShapeFunctions(
145 lf::base::RefEl ref_el_type) const;
146
148 ~UniformScalarFESpace() override = default;
149
150 private:
152 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
155 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
158 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
161 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
166
169
172 std::shared_ptr<const lf::mesh::Mesh> mesh_p);
174 [[nodiscard]] bool check_ptr() const {
175 LF_VERIFY_MSG(Mesh() != nullptr, "No valid FE space object: no mesh");
176 LF_VERIFY_MSG((rfs_quad_p_ != nullptr) && (rfs_quad_p_ != nullptr),
177 "No valid FE space object: no rsfs for cells");
178 return true;
179 }
180
181}; // end class definition UniformScalarFESpace
182
191const unsigned int kUniformScalarFESpaceOutMesh = 1;
192
199const unsigned int kUniformScalarFESpaceOutDofh = 2;
200
207const unsigned int kUniformScalarFESpaceOutRsfs = 4;
208
225template <class SCALAR>
226void PrintInfo(std::ostream &o, const UniformScalarFESpace<SCALAR> &fes,
227 unsigned int ctrl = 0);
228
234template <typename SCALAR>
235std::ostream &operator<<(std::ostream &o,
237
238// Initialization methods
239template <typename SCALAR>
240// NOLINTNEXTLINE(readability-function-cognitive-complexity)
242 std::shared_ptr<const lf::mesh::Mesh> mesh_p) {
243 // Check validity and consistency of mesh pointer
244 LF_VERIFY_MSG(mesh_p != nullptr, "Missing mesh!");
245 LF_VERIFY_MSG((rfs_quad_p_ != nullptr) || (rfs_tria_p_ != nullptr),
246 "Missing FE specification for cells");
247 LF_VERIFY_MSG((mesh_p->DimMesh() == 2), "Only for 2D meshes");
248
249 // Check whether all required finite element specifications are provided
250 LF_VERIFY_MSG(
251 (mesh_p->NumEntities(lf::base::RefEl::kTria()) == 0) ||
252 (rfs_tria_p_ != nullptr),
253 "Missing FE specification for triangles though mesh contains some");
254 LF_VERIFY_MSG((mesh_p->NumEntities(lf::base::RefEl::kQuad()) == 0) ||
255 (rfs_quad_p_ != nullptr),
256 "Missing FE specification for quads though mesh contains some");
257
258 // Compatibility checks and initialization of numbers of shape functions
259 // In particular only a single shape function may be associated to a node
260 // in the case of a SCALAR finite element space
261 if (rfs_tria_p_ != nullptr) {
262 // Probe local shape functions on a triangle
263 LF_VERIFY_MSG((*rfs_tria_p_).RefEl() == lf::base::RefEl::kTria(),
264 "Wrong type for triangle!");
265 LF_VERIFY_MSG((*rfs_tria_p_).NumRefShapeFunctions(2) <= 1,
266 "At most one shape function can be assigned to each vertex");
267 // Initialize numbers of shape functions associated to entities
268 num_rsf_node_ = (*rfs_tria_p_).NumRefShapeFunctions(2);
269 num_rsf_edge_ = (*rfs_tria_p_).NumRefShapeFunctions(1);
270 num_rsf_tria_ = (*rfs_tria_p_).NumRefShapeFunctions(0);
271 }
272 if (rfs_quad_p_ != nullptr) {
273 // Probe local shape functions for QUADs
274 LF_VERIFY_MSG((*rfs_quad_p_).RefEl() == lf::base::RefEl::kQuad(),
275 "Wrong type for quad!");
276 LF_VERIFY_MSG((*rfs_quad_p_).NumRefShapeFunctions(2) <= 1,
277 "At most one shape function can be assigned to each vertex");
278 // Initialize numbers of shape functions associated to entities
279 num_rsf_node_ = (*rfs_quad_p_).NumRefShapeFunctions(2);
280 num_rsf_edge_ = (*rfs_quad_p_).NumRefShapeFunctions(1);
281 num_rsf_quad_ = (*rfs_quad_p_).NumRefShapeFunctions(0);
282 }
283 if (rfs_edge_p_ != nullptr) {
284 LF_VERIFY_MSG((*rfs_edge_p_).RefEl() == lf::base::RefEl::kSegment(),
285 "Wrong type for edge!");
286 LF_VERIFY_MSG((*rfs_edge_p_).NumRefShapeFunctions(1) <= 1,
287 "At most one shape function can be assigned to each vertex");
288 num_rsf_node_ = (*rfs_edge_p_).NumRefShapeFunctions(1);
289 num_rsf_edge_ = (*rfs_edge_p_).NumRefShapeFunctions(0);
290 }
291
292 // Compatibility check for numbers of local shape functions associated with
293 // edges. Those must be the same for all reference shape function descriptions
294 // passed to the finite element space.
295 if ((rfs_tria_p_ != nullptr) && (rfs_quad_p_ != nullptr)) {
296 LF_ASSERT_MSG(((*rfs_tria_p_).NumRefShapeFunctions(2) ==
297 (*rfs_quad_p_).NumRefShapeFunctions(2)),
298 "#RSF mismatch on nodes "
299 << (*rfs_tria_p_).NumRefShapeFunctions(2) << " <-> "
300 << (*rfs_quad_p_).NumRefShapeFunctions(2));
301 LF_ASSERT_MSG(((*rfs_tria_p_).NumRefShapeFunctions(1) ==
302 (*rfs_quad_p_).NumRefShapeFunctions(1)),
303 "#RSF mismatch on edges "
304 << (*rfs_tria_p_).NumRefShapeFunctions(1) << " <-> "
305 << (*rfs_quad_p_).NumRefShapeFunctions(1));
306 }
307 if ((rfs_tria_p_ != nullptr) && (rfs_edge_p_ != nullptr)) {
308 LF_ASSERT_MSG(((*rfs_tria_p_).NumRefShapeFunctions(2) ==
309 (*rfs_edge_p_).NumRefShapeFunctions(1)),
310 "#RSF mismatch on nodes "
311 << (*rfs_tria_p_).NumRefShapeFunctions(2) << " <-> "
312 << (*rfs_edge_p_).NumRefShapeFunctions(1));
313 LF_ASSERT_MSG(((*rfs_tria_p_).NumRefShapeFunctions(1) ==
314 (*rfs_edge_p_).NumRefShapeFunctions(0)),
315 "#RSF mismatch on edges "
316 << (*rfs_tria_p_).NumRefShapeFunctions(1) << " <-> "
317 << (*rfs_edge_p_).NumRefShapeFunctions(0));
318 }
319 if ((rfs_quad_p_ != nullptr) && (rfs_edge_p_ != nullptr)) {
320 LF_ASSERT_MSG(((*rfs_quad_p_).NumRefShapeFunctions(2) ==
321 (*rfs_edge_p_).NumRefShapeFunctions(1)),
322 "#RSF mismatch on edges "
323 << (*rfs_quad_p_).NumRefShapeFunctions(2) << " <-> "
324 << (*rfs_edge_p_).NumRefShapeFunctions(1));
325 LF_ASSERT_MSG(((*rfs_quad_p_).NumRefShapeFunctions(1) ==
326 (*rfs_edge_p_).NumRefShapeFunctions(0)),
327 "#RSF mismatch on edges "
328 << (*rfs_quad_p_).NumRefShapeFunctions(1) << " <-> "
329 << (*rfs_edge_p_).NumRefShapeFunctions(0));
330 }
331
332 // Initialization of dof handler starting with collecting the number of
333 // interior reference shape functions
335 {lf::base::RefEl::kPoint(), num_rsf_node_},
336 {lf::base::RefEl::kSegment(), num_rsf_edge_},
337 {lf::base::RefEl::kTria(), num_rsf_tria_},
338 {lf::base::RefEl::kQuad(), num_rsf_quad_}};
339
340 // NOLINTNEXTLINE(modernize-return-braced-init-list)
341 return lf::assemble::UniformFEDofHandler(std::move(mesh_p), rsf_layout);
342}
343
344template <typename SCALAR>
347 const lf::mesh::Entity &entity) const {
348 return ShapeFunctionLayout(entity.RefEl());
349}
350
351template <typename SCALAR>
354 lf::base::RefEl ref_el_type) const {
355 // Retrieve specification of local shape functions
356 switch (ref_el_type) {
358 return rfs_point_p_.get();
359 }
361 // Null pointer valid return value: indicates that a shape function
362 // description for edges is missing.
363 // LF_ASSERT_MSG(rfs_edge_p_ != nullptr, "No RSF for edges!");
364 return rfs_edge_p_.get();
365 }
366 case lf::base::RefEl::kTria(): {
367 // Null pointer valid return value: indicates that a shape function
368 // description for triangular cells is missing.
369 // LF_ASSERT_MSG(rfs_tria_p_ != nullptr, "No RSF for triangles!");
370 return rfs_tria_p_.get();
371 }
372 case lf::base::RefEl::kQuad(): {
373 // Null pointer valid return value: indicates that a shape function
374 // description for quadrilaterals is missing.
375 // LF_ASSERT_MSG(rfs_quad_p_ != nullptr, "No RSF for quads!");
376 return rfs_quad_p_.get();
377 }
378 default: {
379 LF_VERIFY_MSG(false, "Illegal entity type");
380 }
381 }
382 return nullptr;
383}
384
385template <typename SCALAR>
387 const lf::mesh::Entity &entity) const {
388 return NumRefShapeFunctions(entity.RefEl());
389}
390
391/* number of _interior_ shape functions associated to entities of various types
392 */
393template <typename SCALAR>
395 lf::base::RefEl ref_el_type) const {
396 LF_ASSERT_MSG((rfs_quad_p_ != nullptr) && (rfs_quad_p_ != nullptr),
397 "No valid FE space object: no rsfs");
398 // Retrieve number of interior shape functions from rsf layouts
399 switch (ref_el_type) {
401 return num_rsf_node_;
402 }
404 return num_rsf_edge_;
405 }
406 case lf::base::RefEl::kTria(): {
407 return num_rsf_tria_;
408 }
409 case lf::base::RefEl::kQuad(): {
410 return num_rsf_quad_;
411 }
412 default: {
413 LF_VERIFY_MSG(false, "Illegal entity type");
414 }
415 }
416 return 0;
417}
418
419template <typename SCALAR>
420void PrintInfo(std::ostream &o, const UniformScalarFESpace<SCALAR> &fes,
421 unsigned ctrl) {
422 o << "Uniform scalar FE space, dim = " << fes.LocGlobMap().NumDofs()
423 << std::endl;
424
425 if ((ctrl & kUniformScalarFESpaceOutMesh) != 0) {
426 o << *fes.Mesh() << std::endl;
427 }
428 if ((ctrl & kUniformScalarFESpaceOutDofh) != 0) {
429 o << fes.LocGlobMap() << std::endl;
430 }
431 if ((ctrl & kUniformScalarFESpaceOutRsfs) != 0) {
432 o << fes.NumRefShapeFunctions(lf::base::RefEl::kPoint()) << " rsfs @ nodes"
433 << std::endl;
435 << " rsfs @ edges" << std::endl;
437 << " rsfs @ triangles" << std::endl;
438 o << fes.NumRefShapeFunctions(lf::base::RefEl::kQuad()) << " rsfs @ quads"
439 << std::endl;
440 }
441}
442
444template <typename SCALAR>
445std::ostream &operator<<(std::ostream &o,
446 const UniformScalarFESpace<SCALAR> &fes) {
447 fes.PrintInfo(o, 0);
448 return o;
449}
450
451} // namespace lf::uscalfe
452
454
458template <class SCALAR>
459struct fmt::formatter<lf::uscalfe::UniformScalarFESpace<SCALAR>>
460 : ostream_formatter {};
462
463#endif
A general (interface) class for DOF handling, see Lecture Document Paragraph 2.7.4....
Definition dofhandler.h:112
virtual size_type NumDofs() const =0
total number of dof's handled by the object
Dofhandler for uniform finite element spaces.
Definition dofhandler.h:260
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
Acess to underlying mesh object.
Definition dofhandler.h:385
std::map< lf::base::RefEl, base::size_type > dof_map_t
Map data type for telling number of global shape functions associated with every topological kind of ...
Definition dofhandler.h:269
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
Space of scalar valued finite element functions on a Mesh.
ScalarFESpace()=default
default constructor, needed by std::vector
Interface class for parametric scalar valued finite elements.
Interface class representing a topological entity in a cellular complex
Definition entity.h:42
virtual base::RefEl RefEl() const =0
Describes the reference element type of this entity.
Space of scalar valued finite element functions on a hybrid 2D mesh
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > rfs_quad_p_
void PrintInfo(std::ostream &o, const UniformScalarFESpace< SCALAR > &fes, unsigned int ctrl=0)
Print information about a UniformScalarFESpace to the given stream object.
const lf::assemble::DofHandler & LocGlobMap() const override
access to associated local-to-global map
assemble::UniformFEDofHandler dofh_
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > rfs_point_p_
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > rfs_tria_p_
UniformScalarFESpace(const UniformScalarFESpace &)=delete
~UniformScalarFESpace() override=default
No special destructor.
UniformScalarFESpace(UniformScalarFESpace &&) noexcept=default
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
access to underlying mesh
size_type NumRefShapeFunctions(const lf::mesh::Entity &entity) const override
number of interior shape functions associated to entities of various types
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > rfs_edge_p_
lf::fe::ScalarReferenceFiniteElement< SCALAR > const * ShapeFunctionLayout(const lf::mesh::Entity &entity) const override
access to shape function layout for cells
assemble::UniformFEDofHandler InitDofHandler(std::shared_ptr< const lf::mesh::Mesh > mesh_p)
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
const unsigned int kUniformScalarFESpaceOutDofh
information about the dof handler will be printed.
lf::assemble::size_type size_type
Definition lagr_fe.h:30
void PrintInfo(std::ostream &o, const UniformScalarFESpace< SCALAR > &fes, unsigned ctrl)
const unsigned int kUniformScalarFESpaceOutMesh
mesh information will be printed
const unsigned int kUniformScalarFESpaceOutRsfs
information about the reference shape functions will be printed
std::ostream & operator<<(std::ostream &o, const UniformScalarFESpace< SCALAR > &fes)
Definition assemble.h:31