LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Related Symbols | List of all members
lf::mesh::Entity Class Referenceabstract

Interface class representing a topological entity in a cellular complex More...

#include <lf/mesh/entity.h>

Inheritance diagram for lf::mesh::Entity:
lf::mesh::hybrid2d::Point lf::mesh::hybrid2d::Quadrilateral lf::mesh::hybrid2d::Segment lf::mesh::hybrid2d::Triangle

Public Member Functions

virtual unsigned Codim () const =0
 The codimension of this entity w.r.t. the Mesh.dimMesh() of the owning mesh manager.
 
virtual std::span< const Entity *const > SubEntities (unsigned rel_codim) const =0
 Return all sub entities of this entity that have the given codimension (w.r.t. this entity!)
 
virtual std::span< const OrientationRelativeOrientations () const =0
 return span of relative orientations of sub-entities of the next higher co-dimension.
 
virtual const geometry::GeometryGeometry () const =0
 Describes the geometry of this entity.
 
virtual base::RefEl RefEl () const =0
 Describes the reference element type of this entity.
 
virtual bool operator== (const Entity &rhs) const =0
 Check if two entities are the same.
 
bool operator!= (const Entity &rhs) const
 Check if two entities are different.
 
virtual ~Entity ()=default
 Virtual Destructor.
 

Protected Member Functions

 Entity ()=default
 
 Entity (const Entity &)=default
 
 Entity (Entity &&)=default
 
Entityoperator= (const Entity &)=default
 
Entityoperator= (Entity &&)=default
 

Related Symbols

(Note that these are not member symbols.)

void PrintInfo (std::ostream &stream, const lf::mesh::Entity &e, int output_ctrl=0)
 Prints info about an entity.
 

Detailed Description

Interface class representing a topological entity in a cellular complex

Example: a 2D hybrid mesh consists of cells (entities of co-dimension 0), edges (entities of co-dimension 1), and nodes (entities of co-dimension 2). The set of these entities endowed with incidence relations defines the topology of the mesh.

The core purpose of this class is to provide interface definitions for accessing incidence relations and geometry information. This interface applies to all (topological) types of entities.

Further information can be found in Subsection 2.7.2.2 of the Lecture Document

Definition at line 42 of file entity.h.

Constructor & Destructor Documentation

◆ Entity() [1/3]

lf::mesh::Entity::Entity ( )
protecteddefault

◆ Entity() [2/3]

lf::mesh::Entity::Entity ( const Entity & )
protecteddefault

◆ Entity() [3/3]

lf::mesh::Entity::Entity ( Entity && )
protecteddefault

◆ ~Entity()

virtual lf::mesh::Entity::~Entity ( )
virtualdefault

Virtual Destructor.

Member Function Documentation

◆ Codim()

virtual unsigned lf::mesh::Entity::Codim ( ) const
pure virtual

◆ Geometry()

virtual const geometry::Geometry * lf::mesh::Entity::Geometry ( ) const
pure virtual

Describes the geometry of this entity.

Returns
A pointer to a Geometry object that will remain valid for as long as the Mesh remains valid.

Why does this member function return a pointer instead of a reference? One reason is that entities without a geometric shape are conceivable as building block of a "mesh graph", where only topological information matters. Another reason is that during the construction of a mesh, it turns out to be convenient to build entities "without geometry" first and then endow them with geometric information. A nullptr is a good way to indicate missing geometric information.

Implemented in lf::mesh::hybrid2d::Point, lf::mesh::hybrid2d::Quadrilateral, lf::mesh::hybrid2d::Segment, and lf::mesh::hybrid2d::Triangle.

Referenced by lf::mesh::test_utils::checkGeometryOrientation(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::uscalfe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::operator()(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::mesh::utils::MeshFunctionGlobal< F >::operator()(), and PrintInfo().

◆ operator!=()

bool lf::mesh::Entity::operator!= ( const Entity & rhs) const
inline

Check if two entities are different.

See also
Entity::operator==

Definition at line 141 of file entity.h.

References operator==().

◆ operator=() [1/2]

Entity & lf::mesh::Entity::operator= ( const Entity & )
protecteddefault

◆ operator=() [2/2]

Entity & lf::mesh::Entity::operator= ( Entity && )
protecteddefault

◆ operator==()

virtual bool lf::mesh::Entity::operator== ( const Entity & rhs) const
pure virtual

Check if two entities are the same.

Parameters
rhsCheck if this entity is the same as the rhs entity.
Note
The behavior of this method is undefined if the rhs entity belongs to a different Mesh.

Implemented in lf::mesh::hybrid2d::Point, lf::mesh::hybrid2d::Quadrilateral, lf::mesh::hybrid2d::Segment, and lf::mesh::hybrid2d::Triangle.

Referenced by operator!=().

◆ RefEl()

virtual base::RefEl lf::mesh::Entity::RefEl ( ) const
pure virtual

Describes the reference element type of this entity.

Returns
An object of type lf::base::RefEl.

Implemented in lf::mesh::hybrid2d::Point, lf::mesh::hybrid2d::Quadrilateral, lf::mesh::hybrid2d::Segment, and lf::mesh::hybrid2d::Triangle.

Referenced by lf::mesh::test_utils::checkGeometryOrientation(), lf::mesh::test_utils::checkLocalTopology(), lf::mesh::test_utils::checkRelCodim(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::uscalfe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::isActive(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::isActive(), lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions(), lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::operator()(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::mesh::utils::MeshFunctionGlobal< F >::operator()(), lf::mesh::operator<<(), PrintInfo(), and lf::uscalfe::UniformScalarFESpace< SCALAR >::ShapeFunctionLayout().

◆ RelativeOrientations()

virtual std::span< const Orientation > lf::mesh::Entity::RelativeOrientations ( ) const
pure virtual

return span of relative orientations of sub-entities of the next higher co-dimension.

The corners of every entity are numbered and, thus, define a local orientation of the sub-entitites of co-dimension+1 as explained in Lecture Document Remark 2.7.2.15. This local orientation need not agree with the intrinsic orientation of that sub-entity as given by its own numbering of corners. Such a possible mismatch is detected by this function.

Implemented in lf::mesh::hybrid2d::Point, lf::mesh::hybrid2d::Quadrilateral, lf::mesh::hybrid2d::Segment, and lf::mesh::hybrid2d::Triangle.

◆ SubEntities()

virtual std::span< const Entity *const > lf::mesh::Entity::SubEntities ( unsigned rel_codim) const
pure virtual

Return all sub entities of this entity that have the given codimension (w.r.t. this entity!)

Parameters
rel_codimThe relative co-dimension w.r.t. this entity
Returns
A span of pointers to the sub-entities with the specified relative co-dimension

Implicitly this function defines the numbering of sub-entities, see lf::base::RefEl and Lecture Document Paragraph 2.7.2.14.

Note
For a mesh covering a manifold of dimension 2, we have the following cases
  • For a cell (co-dimension 0 entity), the cell itself is a subentity of relative co-dimension 0, the edges have relative co-dimension 1, and the vertices relative co-dimension 2: in this case the usual co-dimension agrees with the relative co-dimension.
  • For an edge (co-dimension 1 entity), the edge itself is the only sub-entity with relative co-dimension 0, and the endpoints are the sub-entitities of relative co-dimension 1.
The lifetime of the returned span equals the lifetime of the Parent Entity.

Demonstration of usage

bool checkLocalTopology(const Entity &e) {
// Obtain basic information about current Entity
const lf::base::RefEl ref_el = e.RefEl();
const lf::base::dim_t dimension = ref_el.Dimension();
// What this function does in the case of a 2D cell entity (dimension = 2):
// It runs through all edges (co-dimension = 1), fetches their endnodes
// (co-dimension again 1 relative to an edge) and test whether they
// are also sub-entities of co-dimension 2 of the cell.
// Co-dimensions of sub-entities run from 1 to dimension
for (lf::base::dim_t sub_codim = 1; sub_codim <= dimension; sub_codim++) {
// Number of sub-entities with relative co-dimensjon sub_codim
const lf::base::size_type num_sub_entities =
ref_el.NumSubEntities(sub_codim);
// Obtain sequence of pointers to sub-entities of co-dimension sub_codim
std::span<const Entity *const> sub_ent_range = e.SubEntities(sub_codim);
// Run through sub-entities
lf::base::size_type sub_ent_cnt{0};
for (const lf::mesh::Entity *sub_ent : sub_ent_range) {
const lf::base::RefEl sub_ref_el = sub_ent->RefEl();
const lf::base::dim_t sub_dim = sub_ref_el.Dimension();
if (sub_dim != dimension - sub_codim) return false;
// The sub-entity has further sub-entities of codimension 1 to sub_dim
for (lf::base::dim_t sub_sub_codim = 1; sub_sub_codim <= sub_dim;
sub_sub_codim++) {
std::span<const Entity *const> sub_sub_ent_range =
sub_ent->SubEntities(sub_sub_codim);
for (const lf::mesh::Entity *sub_sub_ent : sub_sub_ent_range) {
// The entity pointed to by sub_sub_ent has co-dimension
// sub_codim + sub_sub_codim w.r.t. the entity referenced by e
// Hence get all corresponding sub-entities of e
std::span<const Entity *const> e_sub_sub_range =
e.SubEntities(sub_codim + sub_sub_codim);
// See whether we find sub_sub_ent in this range
int found = 0; // Count how many times we find the sub-entity
for (const Entity *e_sub_sub : e_sub_sub_range) {
if (e_sub_sub == sub_sub_ent) found++;
}
// Any sub-entity should occur exactly once.
if (found != 1) return false;
}
}
sub_ent_cnt++;
} // end loop over sub-entities
if (num_sub_entities != sub_ent_cnt) return false;
}
return true;
} // end checklocalTopology

Use of this method is also shown in Lecture Document Example 2.7.2.12.

Implemented in lf::mesh::hybrid2d::Point, lf::mesh::hybrid2d::Quadrilateral, lf::mesh::hybrid2d::Segment, and lf::mesh::hybrid2d::Triangle.

Referenced by lf::mesh::test_utils::checkGeometryOrientation(), lf::mesh::test_utils::checkLocalTopology(), lf::mesh::test_utils::checkRelCodim(), and PrintInfo().

Friends And Related Symbol Documentation

◆ PrintInfo()

void PrintInfo ( std::ostream & stream,
const lf::mesh::Entity & e,
int output_ctrl = 0 )
related

Prints info about an entity.

Parameters
eThe entity to print info about
output_ctrlcontrols the level of output (see below)
streamThe stream to which this function should output

Output levels

  • output_ctrl == 0: Entity type is printed
  • output_ctrl > 10: The above and information of codimensions
  • output_ctrl > 50: The above and information about subentities
  • output_ctrl > 90: The above and coordinates

Definition at line 74 of file print_info.cc.

References lf::base::RefEl::Dimension(), Geometry(), lf::geometry::Geometry::Global(), lf::base::RefEl::kPoint(), lf::base::RefEl::NodeCoords(), lf::base::RefEl::NumSubEntities(), RefEl(), lf::base::RefEl::RefEl(), and SubEntities().


The documentation for this class was generated from the following files: