LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
dofhandler.cc
1/***************************************************************************
2 * LehrFEM++ - A simple C++ finite element libray for teaching
3 * Developed from 2018 at the Seminar of Applied Mathematics of ETH Zurich,
4 * lead developers Dr. R. Casagrande and Prof. R. Hiptmair
5 ***************************************************************************/
6
14
15#include "dofhandler.h"
16
17#include "lf/mesh/entity.h"
18
19namespace lf::assemble {
20
21// Implementation of output operator for interface class
22std::ostream &operator<<(std::ostream &o, const DofHandler &dof_handler) {
23 PrintInfo(o, dof_handler, 0);
24 return o;
25}
26
27// NOLINTNEXTLINE(readability-function-cognitive-complexity)
28void PrintInfo(std::ostream &stream, const DofHandler &dof_handler,
29 unsigned ctrl) {
30 // Obtain pointer to the underlying mesh
31 auto mesh = dof_handler.Mesh();
32 // Number of degrees of freedom managed by the DofHandler object
33 const lf::assemble::size_type N_dofs(dof_handler.NumDofs());
34
35 stream << "DofHandler(" << dof_handler.NumDofs() << " dofs)";
36 if (ctrl > 0) {
37 // More detailed output
38 stream << '\n';
39 if (ctrl % 2 == 0) {
40 for (lf::base::dim_t codim = 0; codim <= mesh->DimMesh(); codim++) {
41 // Visit all entities of a specific codimension
42 for (const lf::mesh::Entity *e : mesh->Entities(codim)) {
43 // Fetch unique index of current entity supplied by mesh object
44 const lf::base::glb_idx_t e_idx = mesh->Index(*e);
45 // Number of shape functions covering current entity
46 const lf::assemble::size_type no_dofs(dof_handler.NumLocalDofs(*e));
47 // Obtain global indices of those shape functions ...
48 const std::span<const gdof_idx_t> doflist(
49 dof_handler.GlobalDofIndices(*e));
50 // and print them
51 stream << *e << ' ' << e_idx << ": " << no_dofs << " dofs = [";
52 for (const lf::assemble::gdof_idx_t &dof : doflist) {
53 stream << dof << ' ';
54 }
55 stream << ']';
56 if (ctrl % 5 == 0) {
57 // Also output indices of interior shape functions
58 const std::span<const gdof_idx_t> intdoflist(
59 dof_handler.InteriorGlobalDofIndices(*e));
60 stream << " int = [";
61 for (const lf::assemble::gdof_idx_t int_dof : intdoflist) {
62 stream << int_dof << ' ';
63 }
64 stream << ']';
65 }
66 stream << '\n';
67 }
68 }
69 }
70 if (ctrl % 3 == 0) {
71 // List entities associated wit the dofs managed by the current
72 // DofHandler object
73 for (lf::assemble::gdof_idx_t dof_idx = 0; dof_idx < N_dofs; dof_idx++) {
74 const lf::mesh::Entity &e(dof_handler.Entity(dof_idx));
75 stream << "dof " << dof_idx << " -> " << e << " " << mesh->Index(e)
76 << '\n';
77 }
78 }
79 }
80}
81
82// ----------------------------------------------------------------------
83// Implementation UniformFEDofHandler
84// ----------------------------------------------------------------------
85
87 std::shared_ptr<const lf::mesh::Mesh> mesh, dof_map_t dofmap,
88 bool check_edge_orientation)
89 : mesh_(std::move(mesh)),
90 num_dofs_(),
91 check_edge_orientation_(check_edge_orientation) {
92 LF_ASSERT_MSG((mesh_->DimMesh() == 2), "Can handle 2D meshes only");
93
94 // For checking whether a key was found
95 auto map_end = dofmap.end();
96
97 // Get no of interior dof specified for nodes
98 auto map_el_pt = dofmap.find(lf::base::RefEl::kPoint());
99 if (map_el_pt != map_end) {
100 num_loc_dof_point_ = map_el_pt->second;
101 }
102
103 // Get no of interior dof specified for edges
104 auto map_el_ed = dofmap.find(lf::base::RefEl::kSegment());
105 if (map_el_ed != map_end) {
106 num_loc_dof_segment_ = map_el_ed->second;
107 }
108
109 // Get no of interior dof specified for triangles
110 auto map_el_tr = dofmap.find(lf::base::RefEl::kTria());
111 if (map_el_tr != map_end) {
112 num_loc_dof_tria_ = map_el_tr->second;
113 }
114
115 // Get no of interior dof specified for quads
116 auto map_el_qd = dofmap.find(lf::base::RefEl::kQuad());
117 if (map_el_qd != map_end) {
118 num_loc_dof_quad_ = map_el_qd->second;
119 }
120
121 // If an entity type is not represented in the map, we assume that
122 // no shape functions are attached to those entities
123
124 // Initialize total number of shape functions covering an entity.
126
127 // Initializatin of dof index arrays
129}
130
140
141// NOLINTNEXTLINE(readability-function-cognitive-complexity)
143 // This method assumes a proper initialization
144 // of the data in no_loc_dof_* and num_dofs_, num_dof_tria, num_dofs_quad_
145 gdof_idx_t dof_idx = 0;
146
147 // Step I: Set indices for shape functions on nodes
148 // Total number of degrees of freedom on nodes (entities of co-dim = 2)
149 const size_type no_nodes = mesh_->NumEntities(2);
150 const size_type num_dofs_nodes = no_nodes * num_dofs_[kNodeOrd];
151 dofs_[kNodeOrd].resize(num_dofs_nodes);
152 // Run through nodes in the order given by their numbering
153 for (glb_idx_t node_idx = 0; node_idx < no_nodes; node_idx++) {
154 const mesh::Entity *node_p{mesh_->EntityByIndex(2, node_idx)};
155 LF_ASSERT_MSG(mesh_->Index(*node_p) == node_idx, "Node index mismatch");
156 // Beginning of section for concrete node in the dof index vector
157 // for entities of co-dimension 2
158 glb_idx_t node_dof_offset = node_idx * num_dofs_[kNodeOrd];
159 for (unsigned j = 0; j < num_loc_dof_point_; j++) {
160 dofs_[kNodeOrd][node_dof_offset++] = dof_idx;
161 dof_entities_.push_back(node_p); // Store entity for current dof
162 dof_idx++; // Move on to next index
163 }
164 }
165
166 // Step II: Set indices for shape functions on edges
167 // Total number of degrees of freedom belonging to edges (entities of co-dim =
168 // 1)
169 const size_type no_edges = mesh_->NumEntities(1);
170 const size_type num_dofs_edges = no_edges * num_dofs_[kEdgeOrd];
171 dofs_[kEdgeOrd].resize(num_dofs_edges);
172 // Visit all edges
173 // Old implementation, see remarks above
174 // for (const lf::mesh::Entity &edge : mesh_->Entities(1)) {
175 for (glb_idx_t edge_idx = 0; edge_idx < no_edges; edge_idx++) {
176 // Obtain pointer to edge entity
177 const mesh::Entity *edge_p{mesh_->EntityByIndex(1, edge_idx)};
178 LF_ASSERT_MSG(mesh_->Index(*edge_p) == edge_idx, "Edge index mismatch");
179 // Beginning of section for concrete edge in the dof index vector
180 // for entities of co-dimension 1
181 glb_idx_t edge_dof_offset = edge_idx * num_dofs_[kEdgeOrd];
182
183 // Obtain indices for basis functions sitting at endpoints
184 for (const lf::mesh::Entity *endpoint : edge_p->SubEntities(1)) {
185 const glb_idx_t ep_idx(mesh_->Index(*endpoint));
186 glb_idx_t ep_dof_offset = ep_idx * num_dofs_[kNodeOrd];
187 // Copy indices of shape functions from nodes to edge
188 for (unsigned j = 0; j < num_dofs_[kNodeOrd]; j++) {
189 dofs_[kEdgeOrd][edge_dof_offset++] = dofs_[kNodeOrd][ep_dof_offset++];
190 }
191 }
192 // Set indices for interior edge degrees of freedom
193 for (unsigned j = 0; j < num_loc_dof_segment_; j++) {
194 dofs_[kEdgeOrd][edge_dof_offset++] = dof_idx;
195 dof_entities_.push_back(edge_p);
196 dof_idx++;
197 }
198 }
199
200 // Step III: Set indices for shape functions on cells
201 const size_type no_cells = mesh_->NumEntities(0);
202 const size_type max_num_dof_cells = no_cells * num_dofs_[kCellOrd];
203 dofs_[kCellOrd].resize(max_num_dof_cells);
204
205 // Number of (non-)interior shape functins for edges
206 const size_type no_int_dof_edge = num_loc_dof_segment_;
207 const size_type num_ext_dof_edge = num_dofs_[kEdgeOrd] - no_int_dof_edge;
208
209 // Visit all cells
210 // Old implementation without strong link between cell
211 // indices and ordering of global shape functions
212 // for (const lf::mesh::Entity &cell : mesh_->Entities(0)) {
213 for (glb_idx_t cell_idx = 0; cell_idx < no_cells; cell_idx++) {
214 // Obtain pointer to current ell
215 const mesh::Entity *cell_p{mesh_->EntityByIndex(0, cell_idx)};
216 LF_ASSERT_MSG(cell_idx == mesh_->Index(*cell_p), "cell index mismatch");
217 // Offset for cell dof indices in large dof index vector
218 glb_idx_t cell_dof_offset = cell_idx * num_dofs_[kCellOrd];
219
220 // Obtain indices for basis functions in vertices
221 for (const lf::mesh::Entity *vertex : cell_p->SubEntities(2)) {
222 const glb_idx_t vt_idx(mesh_->Index(*vertex));
223 glb_idx_t vt_dof_offset = vt_idx * num_dofs_[kNodeOrd];
224 // Copy indices of shape functions from nodes to cell
225 for (unsigned j = 0; j < num_dofs_[kNodeOrd]; j++) {
226 dofs_[kCellOrd][cell_dof_offset++] = dofs_[kNodeOrd][vt_dof_offset++];
227 }
228 }
229
230 // Collect indices of interior shape functions of edges
231 // Internal ordering may depend on the orientation of the edge, if
232 // the check_edge_orientation_ flag is set
233 const std::span<const lf::mesh::Orientation> edge_orientations =
234 cell_p->RelativeOrientations();
235 auto edges = cell_p->SubEntities(1);
236 // Loop over edges
237 const size_type no_edges_cell = cell_p->RefEl().NumSubEntities(1);
238 for (int ed_sub_idx = 0; ed_sub_idx < no_edges_cell; ed_sub_idx++) {
239 const glb_idx_t edge_idx = mesh_->Index(*edges[ed_sub_idx]);
240 const glb_idx_t edge_int_dof_offset =
241 edge_idx * num_dofs_[kEdgeOrd] + num_ext_dof_edge;
242 // Copy indices of shape functions from edges to cell
243 // The order, in which they are copied depends on the relative orientation
244 // of the edge w.r.t. the cell, if the edge_orientation_flag is set
246 (edge_orientations[ed_sub_idx] == lf::mesh::Orientation::positive)) {
247 // Cell-internal and intrinsic orientation match, do not tinker with
248 // the numbering of local shape functions
249 for (int j = 0; j < no_int_dof_edge; j++) {
250 dofs_[kCellOrd][cell_dof_offset++] =
251 dofs_[kEdgeOrd][edge_int_dof_offset + j];
252 }
253 } else {
254 // lf::mesh::Orientation::negative: Mismatch of orientations
255 // reverse numbering of cell-internal d.o.f.s
256 for (int j = static_cast<int>(no_int_dof_edge - 1); j >= 0; j--) {
257 dofs_[kCellOrd][cell_dof_offset++] =
258 dofs_[kEdgeOrd][edge_int_dof_offset + j];
259 }
260 }
261 }
262
263 // Set indices for interior cell degrees of freedom depending on the type of
264 // cell. Here we add new degrees of freedom
265 size_type num_int_dofs_cell;
266 if (cell_p->RefEl() == lf::base::RefEl::kTria()) {
267 num_int_dofs_cell = num_loc_dof_tria_;
268 } else if (cell_p->RefEl() == lf::base::RefEl::kQuad()) {
269 num_int_dofs_cell = num_loc_dof_quad_;
270 } else {
271 LF_ASSERT_MSG(
272 false, "Illegal cell type; only triangles and quads are supported");
273 }
274
275 // enlist new interior cell-associated dofs
276 for (unsigned j = 0; j < num_int_dofs_cell; j++) {
277 dofs_[kCellOrd][cell_dof_offset++] = dof_idx;
278 dof_entities_.push_back(cell_p);
279 dof_idx++;
280 }
281 }
282 // Finally store total number of shape functions on the mesh.
283 num_dof_ = dof_idx;
284} // end constructor
285
286std::span<const gdof_idx_t> UniformFEDofHandler::GlobalDofIndices(
287 lf::base::RefEl ref_el_type, glb_idx_t entity_index) const {
288 // Co-dimension of entity in a 2D mesh
289 const dim_t codim = 2 - ref_el_type.Dimension();
290 const size_type no_covered_dofs = NumCoveredDofs(ref_el_type);
291
292 LF_ASSERT_MSG((mesh_->NumEntities(codim) > entity_index),
293 "Index " << entity_index << " out of range");
294 // Pointers to range of dof indices
295 const gdof_idx_t *begin =
296 dofs_[codim].data() +
297 (static_cast<std::size_t>(num_dofs_[codim]) * entity_index);
298 const gdof_idx_t *end = begin + no_covered_dofs;
299 return {begin, end};
300}
301
302std::span<const gdof_idx_t> UniformFEDofHandler::GlobalDofIndices(
303 const lf::mesh::Entity &entity) const {
304 return GlobalDofIndices(entity.RefEl(), mesh_->Index(entity));
305}
306
308 lf::base::RefEl ref_el_type, glb_idx_t entity_index) const {
309 // Co-dimension of entity in a 2D mesh
310 const dim_t codim = 2 - ref_el_type.Dimension();
311 const size_type no_covered_dofs = NumCoveredDofs(ref_el_type);
312 const size_type no_loc_dofs = NumInteriorDofs(ref_el_type);
313
314 LF_ASSERT_MSG((mesh_->NumEntities(codim) > entity_index),
315 "Index " << entity_index << " out of range");
316 // Pointers to range of dof indices
317 const gdof_idx_t *begin =
318 dofs_[codim].data() +
319 (static_cast<std::size_t>(num_dofs_[codim]) * entity_index);
320 const gdof_idx_t *end = begin + no_covered_dofs;
321 begin += (no_covered_dofs - no_loc_dofs);
322 return {begin, end};
323}
324
326 const lf::mesh::Entity &entity) const {
327 return InteriorGlobalDofIndices(entity.RefEl(), mesh_->Index(entity));
328}
329
331 const lf::mesh::Entity &entity) const {
332 return GetNumLocalDofs(entity.RefEl(), 0);
333}
334
336 const lf::mesh::Entity &entity) const {
337 return NumInteriorDofs(entity.RefEl());
338}
339
340// ----------------------------------------------------------------------
341// Implementation DynamicFEDofHandler
342// ----------------------------------------------------------------------
343
344std::span<const gdof_idx_t> DynamicFEDofHandler::GlobalDofIndices(
345 const lf::mesh::Entity &entity) const {
346 // Topological type
347 const lf::base::RefEl ref_el_type = entity.RefEl();
348 // Co-dimension of entity in a 2D mesh
349 const dim_t codim = 2 - ref_el_type.Dimension();
350 // Index of current mesh entity
351 const glb_idx_t entity_index = mesh_p_->Index(entity);
352 // Offset of dof indices for current entity
353 const size_type idx_offset = offsets_[codim][entity_index];
354 // Number of shape functions covering current entity
355 const size_type no_covered_dofs =
356 offsets_[codim][entity_index + 1] - idx_offset;
357 // Pointers to range of dof indices
358 const gdof_idx_t *begin = dofs_[codim].data() + idx_offset;
359 const gdof_idx_t *end = begin + no_covered_dofs;
360 return {begin, end};
361}
362
364 const lf::mesh::Entity &entity) const {
365 // Topological type
366 const lf::base::RefEl ref_el_type = entity.RefEl();
367 // Co-dimension of entity in a 2D mesh
368 const dim_t codim = 2 - ref_el_type.Dimension();
369 // Index of current mesh entity
370 const glb_idx_t entity_index = mesh_p_->Index(entity);
371 // Offset of indices for current entity
372 const size_type idx_offset = offsets_[codim][entity_index];
373 // Number of shape functions covering current entity
374 const size_type no_covered_dofs =
375 offsets_[codim][entity_index + 1] - idx_offset;
376 // Number of shape functions associated with the current entity
377 const size_type no_loc_dofs = num_int_dofs_[codim][entity_index];
378 // Pointers to range of dof indices
379 const gdof_idx_t *begin = dofs_[codim].data() + idx_offset;
380 const gdof_idx_t *end = begin + no_covered_dofs;
381 // Move pointer to location of indices for interior dofs
382 // This is the only difference to GlobalDofIndices()
383 begin += (no_covered_dofs - no_loc_dofs);
384 return {begin, end};
385}
386
388 const lf::mesh::Entity &entity) const {
389 const lf::base::RefEl ref_el_type = entity.RefEl();
390 const dim_t codim = 2 - ref_el_type.Dimension();
391 const glb_idx_t entity_index = mesh_p_->Index(entity);
392 return (offsets_[codim][entity_index + 1] - offsets_[codim][entity_index]);
393}
394
396 const lf::mesh::Entity &entity) const {
397 const lf::base::RefEl ref_el_type = entity.RefEl();
398 const dim_t codim = 2 - ref_el_type.Dimension();
399 const glb_idx_t entity_index = mesh_p_->Index(entity);
400 const size_type no_loc_dofs = num_int_dofs_[codim][entity_index];
401 return no_loc_dofs;
402}
403
404} // namespace lf::assemble
A general (interface) class for DOF handling, see Lecture Document Paragraph 2.7.4....
Definition dofhandler.h:112
virtual std::shared_ptr< const lf::mesh::Mesh > Mesh() const =0
Acess to underlying mesh object.
virtual const lf::mesh::Entity & Entity(gdof_idx_t dofnum) const =0
retrieve unique entity at which a basis function is located
virtual std::span< const gdof_idx_t > GlobalDofIndices(const lf::mesh::Entity &entity) const =0
access to indices of global dof's belonging to an entity
virtual std::span< const gdof_idx_t > InteriorGlobalDofIndices(const lf::mesh::Entity &entity) const =0
global indices of shape functions associated with an entity_
virtual size_type NumLocalDofs(const lf::mesh::Entity &entity) const =0
tells the number of degrees of freedom subordinate/_belonging_ to to an entity
virtual size_type NumDofs() const =0
total number of dof's handled by the object
size_type NumInteriorDofs(const lf::mesh::Entity &entity) const override
provides number of shape functions associated with an entity
std::shared_ptr< const lf::mesh::Mesh > mesh_p_
Definition dofhandler.h:765
std::span< const gdof_idx_t > InteriorGlobalDofIndices(const lf::mesh::Entity &entity) const override
global indices of shape functions associated with an entity_
std::span< const gdof_idx_t > GlobalDofIndices(const lf::mesh::Entity &entity) const override
access to indices of global dof's belonging to an entity
std::array< std::vector< size_type >, 3 > offsets_
Definition dofhandler.h:784
size_type NumLocalDofs(const lf::mesh::Entity &entity) const override
tells the number of degrees of freedom subordinate/_belonging_ to to an entity
std::array< std::vector< gdof_idx_t >, 3 > dofs_
Definition dofhandler.h:787
std::array< std::vector< size_type >, 3 > num_int_dofs_
Definition dofhandler.h:782
void initIndexArrays()
initialization of internal index arrays
size_type GetNumLocalDofs(lf::base::RefEl ref_el_type, glb_idx_t) const
Definition dofhandler.h:411
std::array< size_type, 3 > num_dofs_
Definition dofhandler.h:489
std::span< const gdof_idx_t > GlobalDofIndices(const lf::mesh::Entity &entity) const override
access to indices of global dof's belonging to an entity
std::array< std::vector< gdof_idx_t >, 3 > dofs_
Definition dofhandler.h:487
std::span< const gdof_idx_t > InteriorGlobalDofIndices(const lf::mesh::Entity &entity) const override
global indices of shape functions associated with an entity_
UniformFEDofHandler(std::shared_ptr< const lf::mesh::Mesh > mesh, dof_map_t dofmap, bool check_edge_orientation=true)
Construction from a map object.
Definition dofhandler.cc:86
size_type NumInteriorDofs(const lf::mesh::Entity &entity) const override
provides number of shape functions associated with an entity
size_type NumLocalDofs(const lf::mesh::Entity &entity) const override
tells the number of degrees of freedom subordinate/_belonging_ to to an entity
std::shared_ptr< const lf::mesh::Mesh > mesh_
Definition dofhandler.h:480
std::map< lf::base::RefEl, base::size_type > dof_map_t
Map data type for telling number of global shape functions associated with every topological kind of ...
Definition dofhandler.h:269
std::vector< const lf::mesh::Entity * > dof_entities_
Definition dofhandler.h:484
void InitTotalNumDofs()
compute number of shape functions covering an entity type
size_type NumCoveredDofs(lf::base::RefEl ref_el_type) const
Definition dofhandler.h:424
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition ref_el.h:204
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
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.
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition types.h:32
unsigned int glb_idx_t
type for global index of mesh entities (nodes, edges, cells)
Definition types.h:24
D.o.f. index mapping and assembly facilities.
Definition assemble.h:34
lf::base::glb_idx_t glb_idx_t
lf::base::size_type size_type
void PrintInfo(std::ostream &stream, const DofHandler &dof_handler, unsigned ctrl)
Definition dofhandler.cc:28
std::ostream & operator<<(std::ostream &o, const COOMatrix< SCALARTYPE > &mat)
Definition coomatrix.h:229
lf::base::dim_t dim_t
Eigen::Index gdof_idx_t
Defines a set of interface classes that define a mesh manager and provides mesh-related tools that bu...
Definition entity.cc:5