3#ifndef INCG96e6ff0ee0034f4584fcdfc7e9c53f82
4#define INCG96e6ff0ee0034f4584fcdfc7e9c53f82
11#include <fmt/ostream.h>
13#include "eigen_tools.h"
44constexpr dim_t DimensionImpl(
RefElType type) {
54 throw std::runtime_error(
55 "RefEl::Dimension() not implemented for this RefEl type.");
126 static constexpr std::array<std::array<base::dim_t, 2>, 3>
128 static constexpr std::array<std::array<base::dim_t, 2>, 4>
137 template <RefElType type>
139 Eigen::Matrix<double, internal::DimensionImpl(type), 1>;
205 return internal::DimensionImpl(
type_);
224 throw std::runtime_error(
225 "RefEl::NumNodes() not implemented for this RefEl type.");
253 false,
"RefEl::NodeCoords() not implemented for this RefEl type.");
273 template <RefElType type>
274 static const std::vector<NodeCoordVector<type>>&
NodeCoords() {
292 "RefEl::NodeCoords<>() not implemented for this RefEl type.");
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) {
323 LF_ASSERT_MSG_CONSTEXPR(
325 "RefEl::NumSubEntities() not implemented for this RefElType.");
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");
358 "sub_index >= NumSubEntities");
360 if (sub_codim == 0) {
369 LF_ASSERT_MSG_CONSTEXPR(
false,
"This code should never be reached.");
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");
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(
427 "sub_sub_index out of bounds.");
432 if (sub_codim == 0) {
433 return sub_rel_index;
446 LF_ASSERT_MSG_CONSTEXPR(
false,
"This code should never be reached.");
469 LF_VERIFY_MSG(
false,
"ToString() not implemented for this RefElType");
494 [[nodiscard]]
constexpr unsigned int Id()
const {
495 return static_cast<unsigned int>(
type_);
522void PrintInfo(std::ostream& stream,
const RefEl& ref_el,
int output_ctrl = 0);
547struct fmt::formatter<
lf::base::RefEl> : ostream_formatter {};
Represents a reference element with all its properties.
static const Eigen::MatrixXd ncoords_quad_dynamic_
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.
static const std::vector< Eigen::Matrix< double, 0, 1 > > ncoords_point_static_
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
constexpr RefEl(RefElType type) noexcept
Create a RefEl from a lf::base::RefElType enum.
Eigen::Matrix< double, internal::DimensionImpl(type), 1 > NodeCoordVector
Type of the node coordinate iterator that is returned from NodeCoords()
static const std::vector< Eigen::Matrix< double, 1, 1 > > ncoords_segment_static_
static const std::vector< Eigen::Vector2d > ncoords_tria_static_
constexpr RefEl(RefEl &&)=default
Default move constructor.
constexpr RefEl(const RefEl &)=default
Default copy constructor.
static const Eigen::MatrixXd ncoords_segment_dynamic_
constexpr RefEl & operator=(const RefEl &rhs)
Default copy assignment operator.
constexpr unsigned int Id() const
Return a unique id for this reference element.
static constexpr std::array< std::array< base::dim_t, 2 >, 3 > sub_sub_entity_index_tria_
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
constexpr RefEl & operator=(RefEl &&rhs) noexcept
Default move assignment operator.
constexpr size_type NumSubEntities(dim_t sub_codim) const
Get the number of sub-entities of this RefEl with the given codimension.
constexpr dim_t Dimension() const
Return the dimension of this reference element.
static const Eigen::MatrixXd ncoords_tria_dynamic_
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
static const std::vector< NodeCoordVector< type > > & NodeCoords()
Get the coordinates of the nodes of a reference element.
static constexpr std::array< std::array< base::dim_t, 2 >, 4 > sub_sub_entity_index_quad_
static const Eigen::MatrixXd ncoords_point_dynamic_
static constexpr RefEl kTria()
Returns the reference triangle.
constexpr size_type NumNodes() const
The number of nodes of this reference element.
static const std::vector< Eigen::Vector2d > ncoords_quad_static_
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
std::string ToString() const
Return a string representation of this Reference element.
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.
unsigned int size_type
general type for variables related to size of arrays
unsigned int sub_idx_t
type for local indices of sub-entities
Contains basic functionality that is used by other parts of LehrFEM++.
RefElType
An enum that defines all possible RefEl types.
@ 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)
std::ostream & operator<<(std::ostream &stream, const RefEl &ref_el)
Operator overload to print a RefEl to a stream, such as std::cout