LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
check_local_topology.cc
1#include "check_local_topology.h"
2
3#include <gtest/gtest.h>
4
5namespace lf::mesh::test_utils {
6void checkLocalTopology(const Entity &e) {
7 using dim_t = lf::base::RefEl::dim_t;
8 using RefEl = lf::base::RefEl;
9 using size_type = lf::mesh::Mesh::size_type;
10
11 // Obtain basic information about current Entity
12 const RefEl ref_el = e.RefEl();
13 const dim_t dimension = ref_el.Dimension();
14
15 // What this function does in the case of a 2D cell entity (dimension = 2)
16 // It runs through all edges (co-dimension = 1), fetches their endnodes
17 // (co-dimension again 1 relative to an edge) and test whether they
18 // are also sub-entities of co-dimension 2 of the cell.
19
20 // Co-dimensions of sub-entities run from 1 to dimension
21 for (dim_t sub_codim = 1; sub_codim <= dimension; sub_codim++) {
22 // Number of sub-entities with relative co-dimensjon sub_codim
23 const size_type num_sub_entities = ref_el.NumSubEntities(sub_codim);
24 // Obtain sequence of sub-entities of co-dimensjon sub_codim
25 auto sub_ent_range = e.SubEntities(sub_codim);
26 // Run through sub-entities
27 size_type sub_ent_cnt{0};
28 for (const Entity *sub_ent : sub_ent_range) {
29 const RefEl sub_ref_el = sub_ent->RefEl();
30 const dim_t sub_dim = sub_ref_el.Dimension();
31 EXPECT_EQ(sub_dim, dimension - sub_codim) << "Dim/Co-dim mismatch";
32 // The sub-entity has further sub-entities of codimension 1 to sub_dim
33 for (dim_t sub_sub_codim = 1; sub_sub_codim <= sub_dim; sub_sub_codim++) {
34 auto sub_sub_ent_range = sub_ent->SubEntities(sub_sub_codim);
35 for (const Entity *sub_sub_ent : sub_sub_ent_range) {
36 // The entity referenced by sub_sub_ent has co-dimension
37 // sub_codim + sub_sub_codim w.r.t. the entity referenced by e
38 // Hence get all corresponding sub-entities of e
39 auto e_sub_sub_range = e.SubEntities(sub_codim + sub_sub_codim);
40 // See whether we find sub_sub_ent in this range
41 int found = 0; // Count how many times we find the sub-entity
42 for (const Entity *e_sub_sub : e_sub_sub_range) {
43 if (e_sub_sub == sub_sub_ent) {
44 found++;
45 }
46 }
47 EXPECT_EQ(found, 1) << "Sub-sub-entity hit " << found << " times";
48 }
49 }
50 sub_ent_cnt++;
51 } // end loop over sub-entities
52 EXPECT_EQ(num_sub_entities, sub_ent_cnt)
53 << "Subent cnt mismatch " << num_sub_entities << " <-> " << sub_ent_cnt;
54 }
55}
56
57/* SAM_LISTING_BEGIN_1 */
58void checkRelCodim(const Entity &e) {
59 using dim_t = lf::base::RefEl::dim_t;
60 using RefEl = lf::base::RefEl;
61 using size_type = lf::mesh::Mesh::size_type;
62
63 // Obtain basic information about current Entity
64 const RefEl ref_el = e.RefEl();
65 const dim_t dimension = ref_el.Dimension();
66 // Loop over all possible co-dimensions of sub-entities
67 for (dim_t sub_codim = 1; sub_codim <= dimension; ++sub_codim) {
68 // Obtain array of sub-entities of co-dimensjon sub\_codim
69 std::span<const lf::mesh::Entity *const> sub_ent_array{
70 e.SubEntities(sub_codim)};
71 // Query number of sub-entities
72 const size_type num_subent = ref_el.NumSubEntities(sub_codim);
73 // Index-based loop over sub-entities
74 for (int sub_ent_idx = 0; sub_ent_idx < num_subent; ++sub_ent_idx) {
75 // Test whether relative dimension matches absolute dimensions
76 EXPECT_EQ(sub_codim,
77 dimension - sub_ent_array[sub_ent_idx]->RefEl().Dimension())
78 << "Dimension mismatch: " << e << " <-> subent(" << sub_ent_idx << ")"
79 << std::endl;
80 }
81 }
82}
83/* SAM_LISTING_END_1 */
84
85} // namespace lf::mesh::test_utils
Represents a reference element with all its properties.
Definition ref_el.h:109
unsigned int dim_t
Definition ref_el.h:132
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition ref_el.h:204
Interface class representing a topological entity in a cellular complex
Definition entity.h:42
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 base::RefEl RefEl() const =0
Describes the reference element type of this entity.
lf::base::size_type size_type
Utilities for testing sanity of mesh data structures and tests involving meshes.
void checkRelCodim(const Entity &e)
consistence check for local (relative) co-dimensions
void checkLocalTopology(const Entity &e)
Function for testing consistency of subentities.