LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
ref_el.h
1
2
3#ifndef INCG96e6ff0ee0034f4584fcdfc7e9c53f82
4#define INCG96e6ff0ee0034f4584fcdfc7e9c53f82
5
6#include <array>
7#include <vector>
8
9// #include <boost/range.hpp>
10
11#include <fmt/ostream.h>
12
13#include "eigen_tools.h"
14#include "lf_assert.h"
15#include "types.h"
16
17namespace lf::base {
18
31enum class RefElType : unsigned char {
32 kPoint = 1,
34 kSegment = 2,
36 kTria = 3,
38 kQuad = 4,
40};
41
42namespace internal {
43// Some utility methods that are needed to deduce compile time return types:
44constexpr dim_t DimensionImpl(RefElType type) {
45 switch (type) {
47 return 0;
49 return 1;
52 return 2;
53 default:
54 throw std::runtime_error(
55 "RefEl::Dimension() not implemented for this RefEl type.");
56 }
57}
58} // namespace internal
59
109class RefEl {
110 // Node coordinates as dynamic eigen matrices:
111 static const Eigen::MatrixXd ncoords_point_dynamic_;
112 static const Eigen::MatrixXd ncoords_segment_dynamic_;
113 static const Eigen::MatrixXd ncoords_tria_dynamic_;
114 static const Eigen::MatrixXd ncoords_quad_dynamic_;
115
116 // Node coordinates as fixed eigen matrices:
117 static const std::vector<Eigen::Matrix<double, 0, 1>> ncoords_point_static_;
118 static const std::vector<Eigen::Matrix<double, 1, 1>> ncoords_segment_static_;
119 static const std::vector<Eigen::Vector2d> ncoords_tria_static_;
120 static const std::vector<Eigen::Vector2d> ncoords_quad_static_;
121
122 // Member variable
124
125 // subSubEntities, used by SubSubEntity2SubEntity
126 static constexpr std::array<std::array<base::dim_t, 2>, 3>
127 sub_sub_entity_index_tria_ = {{{0, 1}, {1, 2}, {2, 0}}};
128 static constexpr std::array<std::array<base::dim_t, 2>, 4>
129 sub_sub_entity_index_quad_ = {{{0, 1}, {1, 2}, {2, 3}, {3, 0}}};
130
131 public:
132 using dim_t = unsigned int;
137 template <RefElType type>
139 Eigen::Matrix<double, internal::DimensionImpl(type), 1>;
140
144 static constexpr RefEl kPoint() { return RefEl(RefElType::kPoint); }
145
153 static constexpr RefEl kSegment() { return RefEl(RefElType::kSegment); }
154
161 static constexpr RefEl kTria() { return RefEl(RefElType::kTria); }
162
169 static constexpr RefEl kQuad() { return RefEl(RefElType::kQuad); }
170
175 explicit constexpr RefEl(RefElType type) noexcept : type_(type) {}
176
178 constexpr RefEl(const RefEl&) = default;
179
181 constexpr RefEl(RefEl&&) = default;
182
184 // NOLINTNEXTLINE(modernize-use-equals-default,hicpp-use-equals-default,cert-oop54-cpp)
185 constexpr RefEl& operator=(const RefEl& rhs) {
186 type_ = rhs.type_;
187 return *this;
188 }
189
191 constexpr RefEl& operator=(RefEl&& rhs) noexcept {
192 type_ = rhs.type_;
193 return *this;
194 }
195
204 [[nodiscard]] constexpr dim_t Dimension() const {
205 return internal::DimensionImpl(type_);
206 }
207
213 [[nodiscard]] constexpr size_type NumNodes() const {
214 switch (type_) {
216 return 1;
218 return 2;
219 case RefElType::kTria:
220 return 3;
221 case RefElType::kQuad:
222 return 4;
223 default:
224 throw std::runtime_error(
225 "RefEl::NumNodes() not implemented for this RefEl type.");
226 }
227 }
228
241 [[nodiscard]] const Eigen::MatrixXd& NodeCoords() const {
242 switch (type_) {
247 case RefElType::kTria:
249 case RefElType::kQuad:
251 default:
252 LF_VERIFY_MSG(
253 false, "RefEl::NodeCoords() not implemented for this RefEl type.");
254 }
255 }
256
273 template <RefElType type>
274 static const std::vector<NodeCoordVector<type>>& NodeCoords() {
275 // NOLINTNEXTLINE
276 if constexpr (type == RefElType::kPoint) {
278 }
279 // NOLINTNEXTLINE
280 if constexpr (type == RefElType::kSegment) {
282 }
283 // NOLINTNEXTLINE
284 if constexpr (type == RefElType::kTria) {
286 }
287 // NOLINTNEXTLINE
288 if constexpr (type == RefElType::kQuad) {
290 }
291 LF_VERIFY_MSG(false,
292 "RefEl::NodeCoords<>() not implemented for this RefEl type.");
293 }
294
308 [[nodiscard]] constexpr size_type NumSubEntities(dim_t sub_codim) const {
309 LF_ASSERT_MSG_CONSTEXPR(sub_codim >= 0, "sub_codim is negative");
310 LF_ASSERT_MSG_CONSTEXPR(sub_codim <= Dimension(),
311 "sub_codim > Dimension()");
312 if (sub_codim == 0) {
313 return 1;
314 }
315 switch (type_) {
317 return 2; // sub_codim=1
318 case RefElType::kTria:
319 return 3; // sub_codim=1,2
320 case RefElType::kQuad:
321 return 4; // sub_codim=1,2
322 default:
323 LF_ASSERT_MSG_CONSTEXPR(
324 false,
325 "RefEl::NumSubEntities() not implemented for this RefElType.");
326 }
327 return 0; // prevent warnings from compilers
328 }
329
351 [[nodiscard]] constexpr RefEl SubType(dim_t sub_codim,
352 dim_t sub_index) const {
353 LF_ASSERT_MSG_CONSTEXPR(sub_codim >= 0, "sub_codim is negative");
354 LF_ASSERT_MSG_CONSTEXPR(sub_codim <= Dimension(),
355 "sub_codim > Dimension()");
356 LF_ASSERT_MSG_CONSTEXPR(sub_index >= 0, "sub_index is negative");
357 LF_ASSERT_MSG_CONSTEXPR(sub_index < NumSubEntities(sub_codim),
358 "sub_index >= NumSubEntities");
359
360 if (sub_codim == 0) {
361 return *this;
362 }
363 if (sub_codim == Dimension()) {
364 return kPoint();
365 }
366 if (Dimension() - sub_codim == 1) {
367 return kSegment();
368 } else { // NOLINT(readability-else-after-return)
369 LF_ASSERT_MSG_CONSTEXPR(false, "This code should never be reached.");
370 }
371
372 return kPoint(); // prevent warnings from compiler
373 }
374
412 [[nodiscard]] constexpr sub_idx_t SubSubEntity2SubEntity(
413 dim_t sub_codim, sub_idx_t sub_index, dim_t sub_rel_codim,
414 sub_idx_t sub_rel_index) const {
415 LF_ASSERT_MSG_CONSTEXPR(sub_codim >= 0, "sub_codim negative");
416 LF_ASSERT_MSG_CONSTEXPR(sub_codim <= Dimension(), "sub_codim > Dimension");
417 LF_ASSERT_MSG_CONSTEXPR(sub_index >= 0, "sub_index negative");
418 LF_ASSERT_MSG_CONSTEXPR(sub_index <= NumSubEntities(sub_codim),
419 "sub_index >= NumSubEntities");
420 LF_ASSERT_MSG_CONSTEXPR(sub_rel_codim >= 0, "sub_rel_codim negative.");
421 LF_ASSERT_MSG_CONSTEXPR(sub_rel_codim <= Dimension() - sub_codim,
422 "subSubCodim out of bounds.");
423 LF_ASSERT_MSG_CONSTEXPR(sub_rel_index >= 0, "sub_rel_index negative.");
424 LF_ASSERT_MSG_CONSTEXPR(
425 sub_rel_index <
426 SubType(sub_codim, sub_index).NumSubEntities(sub_rel_codim),
427 "sub_sub_index out of bounds.");
428
429 if (type_ == RefElType::kPoint) {
430 return 0;
431 }
432 if (sub_codim == 0) {
433 return sub_rel_index;
434 }
435 if (sub_codim == Dimension()) {
436 return sub_index;
437 }
438
439 // from here on, it must be a segment
440 switch (type_) {
441 case RefElType::kTria:
442 return sub_sub_entity_index_tria_[sub_index][sub_rel_index];
443 case RefElType::kQuad:
444 return sub_sub_entity_index_quad_[sub_index][sub_rel_index];
445 default:
446 LF_ASSERT_MSG_CONSTEXPR(false, "This code should never be reached.");
447 }
448
449 return 0; // Prevent warnings from compiler...
450 }
451
458 [[nodiscard]] std::string ToString() const {
459 switch (type_) {
461 return "NODE";
463 return "EDGE";
464 case RefElType::kTria:
465 return "TRIA";
466 case RefElType::kQuad:
467 return "QUAD";
468 default:
469 LF_VERIFY_MSG(false, "ToString() not implemented for this RefElType");
470 }
471 }
472
480 // NOLINTNEXTLINE(google-explicit-constructor, hicpp-explicit-conversions)
481 [[nodiscard]] constexpr operator RefElType() const { return type_; }
482
493 // NOLINTNEXTLINE(google-explicit-constructor, hicpp-explicit-conversions)
494 [[nodiscard]] constexpr unsigned int Id() const {
495 return static_cast<unsigned int>(type_);
496 }
497
498 ~RefEl() = default;
499
500}; // class RefEl
501
502// Declare print function
522void PrintInfo(std::ostream& stream, const RefEl& ref_el, int output_ctrl = 0);
523
535inline std::ostream& operator<<(std::ostream& stream, const RefEl& ref_el) {
536 return stream << ref_el.ToString();
537}
538
539} // namespace lf::base
540
542
546template <>
547struct fmt::formatter<lf::base::RefEl> : ostream_formatter {};
549#endif // INCG96e6ff0ee0034f4584fcdfc7e9c53f82
Represents a reference element with all its properties.
Definition ref_el.h:109
static const Eigen::MatrixXd ncoords_quad_dynamic_
Definition ref_el.h:114
constexpr sub_idx_t SubSubEntity2SubEntity(dim_t sub_codim, sub_idx_t sub_index, dim_t sub_rel_codim, sub_idx_t sub_rel_index) const
Identifies sub-entities of sub-entities (so-called sub-sub-entities) with sub-entities.
Definition ref_el.h:412
static const std::vector< Eigen::Matrix< double, 0, 1 > > ncoords_point_static_
Definition ref_el.h:117
~RefEl()=default
RefElType type_
Definition ref_el.h:123
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
constexpr RefEl(RefElType type) noexcept
Create a RefEl from a lf::base::RefElType enum.
Definition ref_el.h:175
unsigned int dim_t
Definition ref_el.h:132
Eigen::Matrix< double, internal::DimensionImpl(type), 1 > NodeCoordVector
Type of the node coordinate iterator that is returned from NodeCoords()
Definition ref_el.h:138
static const std::vector< Eigen::Matrix< double, 1, 1 > > ncoords_segment_static_
Definition ref_el.h:118
static const std::vector< Eigen::Vector2d > ncoords_tria_static_
Definition ref_el.h:119
constexpr RefEl(RefEl &&)=default
Default move constructor.
constexpr RefEl(const RefEl &)=default
Default copy constructor.
static const Eigen::MatrixXd ncoords_segment_dynamic_
Definition ref_el.h:112
constexpr RefEl & operator=(const RefEl &rhs)
Default copy assignment operator.
Definition ref_el.h:185
constexpr unsigned int Id() const
Return a unique id for this reference element.
Definition ref_el.h:494
static constexpr std::array< std::array< base::dim_t, 2 >, 3 > sub_sub_entity_index_tria_
Definition ref_el.h:127
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
constexpr RefEl & operator=(RefEl &&rhs) noexcept
Default move assignment operator.
Definition ref_el.h:191
constexpr size_type NumSubEntities(dim_t sub_codim) const
Get the number of sub-entities of this RefEl with the given codimension.
Definition ref_el.h:308
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition ref_el.h:204
static const Eigen::MatrixXd ncoords_tria_dynamic_
Definition ref_el.h:113
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static const std::vector< NodeCoordVector< type > > & NodeCoords()
Get the coordinates of the nodes of a reference element.
Definition ref_el.h:274
static constexpr std::array< std::array< base::dim_t, 2 >, 4 > sub_sub_entity_index_quad_
Definition ref_el.h:129
static const Eigen::MatrixXd ncoords_point_dynamic_
Definition ref_el.h:111
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
constexpr size_type NumNodes() const
The number of nodes of this reference element.
Definition ref_el.h:213
static const std::vector< Eigen::Vector2d > ncoords_quad_static_
Definition ref_el.h:120
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
std::string ToString() const
Return a string representation of this Reference element.
Definition ref_el.h:458
constexpr RefEl SubType(dim_t sub_codim, dim_t sub_index) const
Return the RefEl of the sub-entity with codim sub_codim and index sub_index.
Definition ref_el.h:351
unsigned int size_type
general type for variables related to size of arrays
Definition types.h:20
unsigned int sub_idx_t
type for local indices of sub-entities
Definition types.h:28
Contains basic functionality that is used by other parts of LehrFEM++.
Definition base.h:15
RefElType
An enum that defines all possible RefEl types.
Definition ref_el.h:31
@ kSegment
Returns the (1-dimensional) reference segment.
@ kPoint
Returns the (0-dimensional) reference point.
@ kTria
Returns the reference triangle.
@ kQuad
Returns the reference quadrilateral.
void PrintInfo(std::ostream &stream, const RefEl &ref_el, int output_ctrl)
Definition ref_el.cc:17
std::ostream & operator<<(std::ostream &stream, const RefEl &ref_el)
Operator overload to print a RefEl to a stream, such as std::cout
Definition ref_el.h:535
Definition assemble.h:31