LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Represents a reference element with all its properties. More...
#include <lf/base/ref_el.h>
Public Types | |
using | dim_t = unsigned int |
template<RefElType type> | |
using | NodeCoordVector |
Type of the node coordinate iterator that is returned from NodeCoords() | |
Public Member Functions | |
constexpr | RefEl (RefElType type) noexcept |
Create a RefEl from a lf::base::RefElType enum. | |
constexpr | RefEl (const RefEl &)=default |
Default copy constructor. | |
constexpr | RefEl (RefEl &&)=default |
Default move constructor. | |
constexpr RefEl & | operator= (const RefEl &rhs) |
Default copy assignment operator. | |
constexpr RefEl & | operator= (RefEl &&rhs) noexcept |
Default move assignment operator. | |
constexpr dim_t | Dimension () const |
Return the dimension of this reference element. | |
constexpr size_type | NumNodes () const |
The number of nodes of this reference element. | |
const Eigen::MatrixXd & | NodeCoords () const |
Get the coordinates of the nodes of this reference element. | |
constexpr size_type | NumSubEntities (dim_t sub_codim) const |
Get the number of sub-entities of this RefEl with the given codimension. | |
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 . | |
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. | |
std::string | ToString () const |
Return a string representation of this Reference element. | |
constexpr | operator RefElType () const |
Conversion operator, converts this RefEl to a lf::base::RefElType enum. | |
constexpr unsigned int | Id () const |
Return a unique id for this reference element. | |
~RefEl ()=default | |
Static Public Member Functions | |
static constexpr RefEl | kPoint () |
Returns the (0-dimensional) reference point. | |
static constexpr RefEl | kSegment () |
Returns the (1-dimensional) reference segment. | |
static constexpr RefEl | kTria () |
Returns the reference triangle. | |
static constexpr RefEl | kQuad () |
Returns the reference quadrilateral. | |
template<RefElType type> | |
static const std::vector< NodeCoordVector< type > > & | NodeCoords () |
Get the coordinates of the nodes of a reference element. | |
Private Attributes | |
RefElType | type_ |
Static Private Attributes | |
static const Eigen::MatrixXd | ncoords_point_dynamic_ = Eigen::VectorXd(0) |
static const Eigen::MatrixXd | ncoords_segment_dynamic_ |
static const Eigen::MatrixXd | ncoords_tria_dynamic_ |
static const Eigen::MatrixXd | ncoords_quad_dynamic_ |
static const std::vector< Eigen::Matrix< double, 0, 1 > > | ncoords_point_static_ |
static const std::vector< Eigen::Matrix< double, 1, 1 > > | ncoords_segment_static_ |
static const std::vector< Eigen::Vector2d > | ncoords_tria_static_ |
static const std::vector< Eigen::Vector2d > | ncoords_quad_static_ |
static constexpr std::array< std::array< base::dim_t, 2 >, 3 > | sub_sub_entity_index_tria_ = {{{0, 1}, {1, 2}, {2, 0}}} |
static constexpr std::array< std::array< base::dim_t, 2 >, 4 > | sub_sub_entity_index_quad_ = {{{0, 1}, {1, 2}, {2, 3}, {3, 0}}} |
Related Symbols | |
(Note that these are not member symbols.) | |
void | PrintInfo (std::ostream &stream, const RefEl &ref_el, int output_ctrl=0) |
Diagnostic output operator. Prints info about a reference element. | |
Represents a reference element with all its properties.
Every entity of a mesh in LehrFEM++ is the image of a reference element under an entity-specific (smooth) transformation (which is described by the lf::mesh::Geometry class). This transformation describes the shape of the actual entity, but also the algebraic relations between its sub-entities. For details consult Lecture Document Subsection 10.2.2
On a more fundamental level objects of this class distinguish the topological type of a mesh entity, see Lecture Document Paragraph 2.7.2.9.
Summing up, the reference element serves two main purposes:
There is a fixed number of reference elements. This class has a static member field for every type of reference element:
RefEl::kPoint()
is the Reference element of every point/node in a mesh. The point itself doesn't have any sub-entities.RefEl::kSegment()
is the reference element of every edge in a mesh. It connects two points with each other.RefEl::kTria()
is the reference element of every triangular element in the mesh. It has three segments (codim=1) and three points (codim=2) as sub-entities.RefEl::kQuad()
is the reference element of every quadrilateral element in the mesh. It has four segments (codim=1) and four points (codim=2) as sub-entities.sizeof(RefEl) == sizeof(char)
. It can be copied around as needed.using lf::base::RefEl::dim_t = unsigned int |
using lf::base::RefEl::NodeCoordVector |
Type of the node coordinate iterator that is returned from NodeCoords()
|
inlineexplicitconstexprnoexcept |
Create a RefEl from a lf::base::RefElType enum.
type | The type of the Reference Element |
Definition at line 175 of file ref_el.h.
Referenced by lf::assemble::DynamicFEDofHandler::GlobalDofIndices(), lf::assemble::DynamicFEDofHandler::InteriorGlobalDofIndices(), lf::mesh::test_utils::isWatertightMesh(), kPoint(), kQuad(), kSegment(), kTria(), lf::assemble::DynamicFEDofHandler::NumInteriorDofs(), lf::assemble::DynamicFEDofHandler::NumLocalDofs(), lf::mesh::Mesh::PrintInfo(), lf::mesh::Entity::PrintInfo(), lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >::RefEl(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
constexprdefault |
Default copy constructor.
|
constexprdefault |
Default move constructor.
|
default |
|
inlineconstexpr |
Return the dimension of this reference element.
Definition at line 204 of file ref_el.h.
References type_.
Referenced by lf::mesh::hybrid2d::MeshFactory::AddEntity(), lf::mesh::test_utils::checkGeometryOrientation(), lf::mesh::test_utils::checkLocalTopology(), lf::mesh::test_utils::checkRelCodim(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::fe::ScalarReferenceFiniteElement< SCALAR >::Dimension(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), NumSubEntities(), lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::operator()(), lf::mesh::utils::MeshFunctionGlobal< F >::operator()(), lf::mesh::Entity::PrintInfo(), PrintInfo(), lf::quad::QuadRule::QuadRule(), SubSubEntity2SubEntity(), and SubType().
|
inlineconstexpr |
Return a unique id for this reference element.
The id numbers are contiguous over the known entity types starting from 0. Thus, this identification number can be used as an array index, if one wants to store information for (all) type of entities.
Definition at line 494 of file ref_el.h.
References type_.
Referenced by lf::quad::QuadRuleCache::Get(), and lf::io::VtkWriter::WriteCellData().
|
inlinestaticconstexpr |
Returns the (0-dimensional) reference point.
Definition at line 144 of file ref_el.h.
References lf::base::kPoint, and RefEl().
Referenced by lf::mesh::hybrid2d::MeshFactory::AddPoint(), lf::geometry::test_utils::checkChildGeometry(), lf::geometry::test_utils::checkSubGeometry(), lf::geometry::Point::ChildGeometry(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace(), lf::uscalfe::UniformScalarFESpace< SCALAR >::InitDofHandler(), lf::io::GmshReader::InitGmshFile(), lf::io::GmshReader::InitGmshFile(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), lf::assemble::UniformFEDofHandler::NumCoveredDofs(), lf::assemble::UniformFEDofHandler::NumInteriorDofs(), lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::refinement::MeshHierarchy::PerformRefinement(), lf::mesh::hybrid2d::Point::Point(), lf::uscalfe::PrintInfo(), lf::mesh::Entity::PrintInfo(), PrintInfo(), lf::fe::FePoint< SCALAR >::RefEl(), lf::geometry::Point::RefEl(), lf::mesh::hybrid2d::Point::RefEl(), lf::io::RefElOf(), lf::io::RefElOf(), lf::uscalfe::UniformScalarFESpace< SCALAR >::ShapeFunctionLayout(), SubType(), lf::assemble::UniformFEDofHandler::UniformFEDofHandler(), lf::geometry::Volume(), lf::io::VtkWriter::WriteCellData(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
inlinestaticconstexpr |
Returns the reference quadrilateral.
Definition at line 169 of file ref_el.h.
References lf::base::kQuad, and RefEl().
Referenced by lf::mesh::utils::TorusMeshBuilder::Build(), lf::mesh::utils::TPQuadMeshBuilder::Build(), lf::geometry::test_utils::checkChildGeometry(), lf::geometry::QuadO1::ChildGeometry(), lf::geometry::Parallelogram::ChildGeometry(), lf::geometry::QuadO2::ChildGeometry(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::mesh::test_utils::GenerateHybrid2DTestMesh(), lf::fe::HierarchicScalarFESpace< SCALAR >::Init(), lf::uscalfe::UniformScalarFESpace< SCALAR >::InitDofHandler(), lf::io::GmshReader::InitGmshFile(), lf::io::GmshReader::InitGmshFile(), lf::assemble::UniformFEDofHandler::initIndexArrays(), lf::quad::make_QuadQR_EdgeMidpointRule(), lf::quad::make_QuadQR_MidpointRule(), lf::quad::make_QuadQR_P4O4(), lf::quad::make_QuadRule(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), lf::assemble::UniformFEDofHandler::NumCoveredDofs(), lf::assemble::UniformFEDofHandler::NumInteriorDofs(), lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::uscalfe::PrintInfo(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), lf::fe::FeHierarchicQuad< SCALAR >::RefEl(), lf::geometry::QuadO1::RefEl(), lf::geometry::Parallelogram::RefEl(), lf::geometry::QuadO2::RefEl(), lf::mesh::hybrid2d::Quadrilateral::RefEl(), lf::uscalfe::FeLagrangeO1Quad< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO2Quad< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO3Quad< SCALAR >::RefEl(), lf::io::RefElOf(), lf::io::RefElOf(), lf::refinement::MeshHierarchy::RefineMarked(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider(), lf::uscalfe::UniformScalarFESpace< SCALAR >::ShapeFunctionLayout(), lf::assemble::UniformFEDofHandler::UniformFEDofHandler(), lf::geometry::Volume(), lf::io::VtkWriter::WriteCellData(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
inlinestaticconstexpr |
Returns the (1-dimensional) reference segment.
Definition at line 153 of file ref_el.h.
References lf::base::kSegment, and RefEl().
Referenced by lf::mesh::hybrid2d::MeshFactory::AddEntity(), lf::mesh::utils::TPTriagMeshBuilder::Build(), lf::geometry::test_utils::checkChildGeometry(), lf::geometry::SegmentO1::ChildGeometry(), lf::geometry::SegmentO2::ChildGeometry(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::uscalfe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::FeHierarchicQuad< SCALAR >::FeHierarchicQuad(), lf::fe::FeHierarchicTria< SCALAR >::FeHierarchicTria(), lf::uscalfe::UniformScalarFESpace< SCALAR >::InitDofHandler(), lf::io::GmshReader::InitGmshFile(), lf::io::GmshReader::InitGmshFile(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::isActive(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::isActive(), lf::quad::make_QuadRule(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), lf::assemble::UniformFEDofHandler::NumCoveredDofs(), lf::assemble::UniformFEDofHandler::NumInteriorDofs(), lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::uscalfe::PrintInfo(), lf::fe::FeHierarchicSegment< SCALAR >::RefEl(), lf::geometry::SegmentO1::RefEl(), lf::geometry::SegmentO2::RefEl(), lf::mesh::hybrid2d::Segment::RefEl(), lf::uscalfe::FeLagrangeO1Segment< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO2Segment< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO3Segment< SCALAR >::RefEl(), lf::io::RefElOf(), lf::io::RefElOf(), lf::refinement::MeshHierarchy::RefineMarked(), lf::uscalfe::UniformScalarFESpace< SCALAR >::ShapeFunctionLayout(), SubType(), lf::assemble::UniformFEDofHandler::UniformFEDofHandler(), lf::geometry::Volume(), lf::io::VtkWriter::WriteCellData(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
inlinestaticconstexpr |
Returns the reference triangle.
Definition at line 161 of file ref_el.h.
References lf::base::kTria, and RefEl().
Referenced by lf::mesh::utils::TPTriagMeshBuilder::Build(), lf::geometry::test_utils::checkChildGeometry(), lf::geometry::TriaO1::ChildGeometry(), lf::geometry::TriaO2::ChildGeometry(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::mesh::test_utils::GenerateHybrid2DTestMesh(), lf::fe::HierarchicScalarFESpace< SCALAR >::Init(), lf::uscalfe::UniformScalarFESpace< SCALAR >::InitDofHandler(), lf::io::GmshReader::InitGmshFile(), lf::io::GmshReader::InitGmshFile(), lf::assemble::UniformFEDofHandler::initIndexArrays(), lf::quad::make_QuadRule(), lf::quad::make_QuadRuleNodal(), lf::quad::make_TriaQR_EdgeMidpointRule(), lf::quad::make_TriaQR_MidpointRule(), lf::quad::make_TriaQR_P6O4(), lf::refinement::MeshHierarchy::MeshHierarchy(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), lf::assemble::UniformFEDofHandler::NumCoveredDofs(), lf::assemble::UniformFEDofHandler::NumInteriorDofs(), lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::uscalfe::PrintInfo(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), lf::fe::FeHierarchicTria< SCALAR >::RefEl(), lf::geometry::TriaO1::RefEl(), lf::geometry::TriaO2::RefEl(), lf::mesh::hybrid2d::Triangle::RefEl(), lf::uscalfe::FeLagrangeO1Tria< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO2Tria< SCALAR >::RefEl(), lf::uscalfe::FeLagrangeO3Tria< SCALAR >::RefEl(), lf::io::RefElOf(), lf::io::RefElOf(), lf::refinement::MeshHierarchy::RefineMarked(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider(), lf::uscalfe::UniformScalarFESpace< SCALAR >::ShapeFunctionLayout(), lf::assemble::UniformFEDofHandler::UniformFEDofHandler(), lf::geometry::Volume(), lf::io::VtkWriter::WriteCellData(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), lf::io::VtkWriter::WritePointData(), and lf::io::writeTikZ().
|
inlinestatic |
Get the coordinates of the nodes of a reference element.
type | The RefElType of the reference element for which the node coordinates should be returned, see usage below. |
std::vector
with NumNodes()
elements. Every Element is a fixed-size column vector (e.g. Eigen::Matrix<double,1,1>
for a segment or Eigen::Matrix<double,2,1>
for a triangle/quadrilateral)Definition at line 274 of file ref_el.h.
References lf::base::kPoint, lf::base::kQuad, lf::base::kSegment, lf::base::kTria, ncoords_point_static_, ncoords_quad_static_, ncoords_segment_static_, and ncoords_tria_static_.
|
inline |
Get the coordinates of the nodes of this reference element.
Dimension() x NumNodes()
that contains the the coordinates of every node as a column vector.Definition at line 241 of file ref_el.h.
References lf::base::kPoint, lf::base::kQuad, lf::base::kSegment, lf::base::kTria, ncoords_point_dynamic_, ncoords_quad_dynamic_, ncoords_segment_dynamic_, ncoords_tria_dynamic_, and type_.
Referenced by lf::mesh::test_utils::checkGeometryOrientation(), lf::geometry::test_utils::checkSubGeometry(), lf::geometry::Corners(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::uscalfe::FeLagrangeO1Tria< SCALAR >::EvaluationNodes(), lf::uscalfe::FeLagrangeO1Quad< SCALAR >::EvaluationNodes(), lf::uscalfe::FeLagrangeO1Segment< SCALAR >::EvaluationNodes(), lf::mesh::test_utils::isWatertightMesh(), lf::quad::make_QuadRuleNodal(), lf::geometry::Geometry::PrintInfo(), lf::mesh::Mesh::PrintInfo(), lf::mesh::Entity::PrintInfo(), PrintInfo(), lf::io::VtkWriter::WriteCellData(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
inlineconstexpr |
The number of nodes of this reference element.
NumSubEntities(Dimension())
Definition at line 213 of file ref_el.h.
References lf::base::kPoint, lf::base::kQuad, lf::base::kSegment, lf::base::kTria, and type_.
Referenced by lf::mesh::hybrid2d::MeshFactory::AddEntity(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::mesh::test_utils::isWatertightMesh(), lf::quad::make_QuadRuleNodal(), lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumEvaluationNodes(), lf::uscalfe::FeLagrangeO1Quad< SCALAR >::NumEvaluationNodes(), lf::uscalfe::FeLagrangeO1Segment< SCALAR >::NumEvaluationNodes(), PrintInfo(), and lf::io::writeTikZ().
Get the number of sub-entities of this RefEl with the given codimension.
sub_codim | The codimension of the subEntities that should be counted. |
codim=1
sub-entities, therefore RefEl::kSegment().NumSubEntities(1) == 2
codim=1
(all Segments), therefore RefEl::kTria().NumSubEntities(1) == 3
codim=2
(all Points), therefore RefEl::kTria().NumSubEntities(2) == 3
Definition at line 308 of file ref_el.h.
References Dimension(), lf::base::kQuad, lf::base::kSegment, lf::base::kTria, and type_.
Referenced by lf::refinement::Hybrid2DRefinementPattern::Hybrid2DRefinementPattern(), lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::mesh::Entity::PrintInfo(), PrintInfo(), lf::refinement::Hybrid2DRefinementPattern::setAnchor(), SubSubEntity2SubEntity(), and SubType().
|
inlineconstexpr |
Conversion operator, converts this RefEl to a lf::base::RefElType enum.
Definition at line 481 of file ref_el.h.
References type_.
|
inlineconstexpr |
Identifies sub-entities of sub-entities (so-called sub-sub-entities) with sub-entities.
sub_codim | The codimension of the sub-entity. |
sub_index | The zero-based index of the sub-entity. |
sub_rel_codim | The codimension of the sub-sub-entity w.r.t. the sub-entity identified by sub_codim and sub_index . |
sub_rel_index | The index of the sub-sub-entity w.r.t. the sub-entity identified by sub_codim and sub_index . |
sub_rel_codim
is a relative co-dimension.sub_rel_index
is a local index in the sub-entity.RefEl::kTria()
with sub_codim=1
, sub_index=1
is a RefEl::kSegment()
that connects node 1 with node 2. The (sub-)sub-entity of this segment with codim sub_rel_codim=1
and sub-index sub_rel_index=0
(both w.r.t. to the segment) is the first point of the segment, i.e. node 1. Therefore SubSubEntity2SubEntity(1,1,1,0) == 1
sub_rel_index=1
: SubSubEntity2SubEntity(1,1,1,1) == 2
Definition at line 412 of file ref_el.h.
References Dimension(), lf::base::kPoint, lf::base::kQuad, lf::base::kTria, NumSubEntities(), sub_sub_entity_index_quad_, sub_sub_entity_index_tria_, SubType(), and type_.
Return the RefEl of the sub-entity with codim sub_codim
and index sub_index
.
sub_codim | The codimension of the sub-entity (w.r.t. Dimension() ). Should be <= Dimension() . |
sub_index | The zero-based index of the sub-entity. sub_index should be smaller than NumSumEntities(sub_codim) |
RefEl::kTria().SubType(2,i) == RefEl::kPoint()
for i=0,1,2.RefEl::kQuad().SubType(1,i) == RefEl::kSegment()
for i=0,1,2,3.RefEl::kTria().SubType(0,0) == RefEl::kTria()
.Definition at line 351 of file ref_el.h.
References Dimension(), kPoint(), kSegment(), and NumSubEntities().
Referenced by PrintInfo(), and SubSubEntity2SubEntity().
|
inline |
Return a string representation of this Reference element.
This string is supposed to described the topological type of the entity: NODE, EDGE (2 nodes), TRIA (3-node triangle), QUAD (4-node quadrilateral)
Definition at line 458 of file ref_el.h.
References lf::base::kPoint, lf::base::kQuad, lf::base::kSegment, lf::base::kTria, and type_.
Referenced by lf::mesh::hybrid2d::MeshFactory::AddEntity(), lf::geometry::SegmentO1::ChildGeometry(), lf::geometry::SegmentO2::ChildGeometry(), lf::geometry::QuadO1::ChildGeometry(), lf::geometry::Parallelogram::ChildGeometry(), lf::geometry::QuadO2::ChildGeometry(), lf::geometry::TriaO1::ChildGeometry(), lf::geometry::TriaO2::ChildGeometry(), lf::geometry::Point::ChildGeometry(), lf::refinement::Hybrid2DRefinementPattern::ChildPolygons(), lf::refinement::Hybrid2DRefinementPattern::Hybrid2DRefinementPattern(), lf::refinement::Hybrid2DRefinementPattern::NumChildren(), lf::base::operator<<(), PrintInfo(), lf::refinement::Hybrid2DRefinementPattern::setAnchor(), lf::geometry::Volume(), and lf::io::writeMatplotlib().
|
related |
Diagnostic output operator. Prints info about a reference element.
ref_el | The reference element to print info about |
stream | The stream to which this function should output |
output_ctrl | Determines the level of detail with which the output is printed to o |
output_ctrl
= 0: Type of reference element, dimension and number of nodesoutput_ctrl
> 0: The above and number of subentities and their types for each codimensionoutput_ctrl
> 10: The above and type of subentity for each subentities in each codimension.output_ctrl
> 20: The above and coordinates of the points of the reference element Definition at line 17 of file ref_el.cc.
References Dimension(), kPoint(), NodeCoords(), NumNodes(), NumSubEntities(), SubType(), and ToString().
|
staticprivate |
Definition at line 111 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 117 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 114 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 120 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 112 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 118 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 113 of file ref_el.h.
Referenced by NodeCoords().
|
staticprivate |
Definition at line 119 of file ref_el.h.
Referenced by NodeCoords().
|
staticconstexprprivate |
Definition at line 129 of file ref_el.h.
Referenced by SubSubEntity2SubEntity().
|
staticconstexprprivate |
Definition at line 127 of file ref_el.h.
Referenced by SubSubEntity2SubEntity().
|
private |
Definition at line 123 of file ref_el.h.
Referenced by Dimension(), Id(), NodeCoords(), NumNodes(), NumSubEntities(), operator RefElType(), operator=(), operator=(), SubSubEntity2SubEntity(), and ToString().