LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
dofhandler.h
1#ifndef INCG_LF_DOFHD_H
2#define INCG_LF_DOFHD_H
3/***************************************************************************
4 * LehrFEM++ - A simple C++ finite element libray for teaching
5 * Developed from 2018 at the Seminar of Applied Mathematics of ETH Zurich,
6 * lead developers Dr. R. Casagrande and Prof. R. Hiptmair
7 ***************************************************************************/
8
17#include <lf/mesh/mesh.h>
18
19#include <map>
20
21#include "assembly_types.h"
22
23namespace lf::assemble {
113 protected:
117 DofHandler() = default;
118
123 DofHandler(const DofHandler &) = default;
124 DofHandler(DofHandler &&) = default;
125 DofHandler &operator=(const DofHandler &) = default;
128 public:
130 virtual ~DofHandler() = default;
131
135 [[nodiscard]] virtual size_type NumDofs() const = 0;
136
148 const lf::mesh::Entity &entity) const = 0;
149
161 const lf::mesh::Entity &entity) const = 0;
162
184 [[nodiscard]] virtual std::span<const gdof_idx_t> GlobalDofIndices(
185 const lf::mesh::Entity &entity) const = 0;
186
204 [[nodiscard]] virtual std::span<const gdof_idx_t> InteriorGlobalDofIndices(
205 const lf::mesh::Entity &entity) const = 0;
206
220 gdof_idx_t dofnum) const = 0;
221
227 [[nodiscard]] virtual std::shared_ptr<const lf::mesh::Mesh> Mesh() const = 0;
228};
229
231std::ostream &operator<<(std::ostream &o, const DofHandler &dof_handler);
232
248void PrintInfo(std::ostream &stream, const DofHandler &dof_handler,
249 unsigned int ctrl = 0);
250
251/* ====================================================================== */
252
261 public:
269 using dof_map_t = std::map<lf::base::RefEl, base::size_type>;
307 UniformFEDofHandler(std::shared_ptr<const lf::mesh::Mesh> mesh,
315
320
325
330
334 ~UniformFEDofHandler() override = default;
335
336 [[nodiscard]] size_type NumDofs() const override { return num_dof_; }
337
343 const lf::mesh::Entity &entity) const override;
344
350 const lf::mesh::Entity &entity) const override;
351
359 [[nodiscard]] std::span<const gdof_idx_t> GlobalDofIndices(
360 const lf::mesh::Entity &entity) const override;
361
369 [[nodiscard]] std::span<const gdof_idx_t> InteriorGlobalDofIndices(
370 const lf::mesh::Entity &entity) const override;
371
377 gdof_idx_t dofnum) const override {
378 LF_VERIFY_MSG(dofnum < dof_entities_.size(),
379 "Illegal dof index " << dofnum << ", max = " << num_dof_);
380 return *dof_entities_[dofnum];
381 }
382
385 [[nodiscard]] std::shared_ptr<const lf::mesh::Mesh> Mesh() const override {
386 return mesh_;
387 }
388
389 private:
393 void initIndexArrays();
402 void InitTotalNumDofs();
403
404 // Access method to numbers and values of indices of shape functions
405 [[nodiscard]] std::span<const gdof_idx_t> GlobalDofIndices(
407
408 [[nodiscard]] std::span<const gdof_idx_t> InteriorGlobalDofIndices(
410
415
426 switch (ref_el_type) {
429 break;
430 }
433 break;
434 }
435 case lf::base::RefEl::kTria(): {
437 break;
438 }
439 case lf::base::RefEl::kQuad(): {
441 break;
442 }
443 default: {
444 LF_VERIFY_MSG(false, "Illegal entity type");
445 break;
446 }
447 } // end switch
448 return no_covered_dofs;
449 }
450
454 switch (ref_el_type) {
457 break;
458 }
461 break;
462 }
463 case lf::base::RefEl::kTria(): {
465 break;
466 }
467 case lf::base::RefEl::kQuad(): {
469 break;
470 }
471 default: {
472 LF_VERIFY_MSG(false, "Illegal entity type");
473 break;
474 }
475 } // end switch
476 return no_loc_dofs;
477 }
478
480 std::shared_ptr<const lf::mesh::Mesh> mesh_;
484 std::vector<const lf::mesh::Entity *> dof_entities_;
487 std::array<std::vector<gdof_idx_t>, 3> dofs_;
489 std::array<size_type, 3> num_dofs_;
503};
504
505/* ====================================================================== */
506
515 public:
520
525
530
535
539 ~DynamicFEDofHandler() override = default;
540
564 template <typename LOCALDOFINFO>
565 // NOLINTNEXTLINE(readability-function-cognitive-complexity)
566 DynamicFEDofHandler(std::shared_ptr<const lf::mesh::Mesh> mesh_p,
568 : mesh_p_(std::move(mesh_p)) {
569 LF_ASSERT_MSG((mesh_p_->DimMesh() == 2), "Can handle 2D meshes only");
570
571 // Index counter for global shape functions = global dof
573
574 // Step I: Set indices for shape functions on nodes
575 // Run through node indices (entities of co-dimension 2)
576 const size_type no_nodes = mesh_p_->NumEntities(2);
577 num_int_dofs_[2].resize(no_nodes, 0);
578 offsets_[2].resize(no_nodes + 1, 0);
579 // Traverse nodes (co-dimension-2 entities) based on indices
581 // Obtain pointer to node entity
582 const mesh::Entity *node_p{mesh_p_->EntityByIndex(2, node_idx)};
583 LF_ASSERT_MSG(mesh_p_->Index(*node_p) == node_idx, "Node index mismatch");
584 // Offset for indices of node in index vector
587 // Request number of local shape functions associated with the node
590
591 // Store dof indices in array
592 for (unsigned j = 0; j < no_int_dof_node; j++) {
593 dofs_[2].push_back(dof_idx);
594 dof_entities_.push_back(node_p); // Store entity for current dof
595 dof_idx++; // Move on to next index
596 }
597 }
598 // Set sentinel
600
601 // Step II: Set indices for shape functions on edges (co-dimension = 1)
602 const size_type no_edges = mesh_p_->NumEntities(1);
603 // Set length of edge-related index vectors
604 num_int_dofs_[1].resize(no_edges, 0);
605 offsets_[1].resize(no_edges + 1, 0);
606 // points to the position of the current dof index
609 // Obtain pointer to edge entity
610 const mesh::Entity *edge{mesh_p_->EntityByIndex(1, edge_idx)};
611 LF_ASSERT_MSG(mesh_p_->Index(*edge) == edge_idx, "Edge index mismatch");
612 // Offset for indices of edge dof in index vector
616
617 // Obtain indices for basis functions sitting at endpoints
618 // Endpoints are mesh entities with co-dimension = 2
619 for (const lf::mesh::Entity *endpoint : edge->SubEntities(1)) {
620 const glb_idx_t ep_idx(mesh_p_->Index(*endpoint));
623 // Copy indices of shape functions from nodes to edge
624 for (unsigned j = 0; j < no_int_dofs_ep; j++) {
625 dofs_[1].push_back(dofs_[2][ep_dof_offset + j]);
627 }
628 }
629 // Set indices for interior edge degrees of freedom
630 for (unsigned j = 0; j < no_int_dof_edge; j++) {
631 dofs_[1].push_back(dof_idx);
633 dof_entities_.push_back(edge);
634 dof_idx++;
635 }
636 }
637 // Set sentinel
639
640 // Step III: Set indices for shape functions on cells (co-dimension = 0)
641 const size_type no_cells = mesh_p_->NumEntities(0);
642 // Set length of cell-related index vectors
643 num_int_dofs_[0].resize(no_cells, 0);
644 offsets_[0].resize(no_cells + 1, 0);
645 // points to the position of the current dof index
648 // Obtain pointer to current ell
649 const mesh::Entity *cell{mesh_p_->EntityByIndex(0, cell_idx)};
650 // Offset for indices of cell dof in index vector
654
655 // Obtain indices for basis functions in vertices
656 for (const lf::mesh::Entity *vertex : cell->SubEntities(2)) {
657 const glb_idx_t vt_idx(mesh_p_->Index(*vertex));
660 // Copy indices of shape functions from nodes to cell
661 for (unsigned j = 0; j < no_int_dofs_vt; j++) {
662 dofs_[0].push_back(dofs_[2][vt_dof_offset + j]);
664 }
665 }
666
667 // Collect indices of interior shape functions of edges
668 // Internal ordering will depend on the orientation of the edge
669 auto edge_orientations = cell->RelativeOrientations();
670 auto edges = cell->SubEntities(1);
671 // Loop over edges
672 for (int ed_sub_idx = 0; ed_sub_idx < cell->RefEl().NumSubEntities(1);
673 ed_sub_idx++) {
674 const glb_idx_t edge_idx = mesh_p_->Index(*edges[ed_sub_idx]);
678
679 // Copy indices of shape functions from edges to cell
680 // The order, in which they are copied depends on the relative
681 // orientation of the edge w.r.t. the cell
682 switch (edge_orientations[ed_sub_idx]) {
684 for (int j = 0; j < no_int_dof_edge; j++) {
685 dofs_[0].push_back(dofs_[1][edge_int_dof_offset + j]);
687 }
688 break;
689 }
691 for (int j = static_cast<int>(no_int_dof_edge - 1); j >= 0; j--) {
692 dofs_[0].push_back(dofs_[1][edge_int_dof_offset + j]);
694 }
695 break;
696 }
697 } // end switch
698 } // end loop over edges
699
700 // Set indices for interior shape functions of the cell
701 for (unsigned j = 0; j < no_int_dof_cell; j++) {
702 dofs_[0].push_back(dof_idx);
704 dof_entities_.push_back(cell);
705 dof_idx++;
706 } // end loop over interior dof of cell
707 } // end loop over cells
708 // Set sentinel
710
711 // Finally set number of global shape functions
713 } // end constructor
714
718 [[nodiscard]] size_type NumDofs() const override { return num_dof_; }
719
725 const lf::mesh::Entity &entity) const override;
726
732 const lf::mesh::Entity &entity) const override;
733
737 [[nodiscard]] std::span<const gdof_idx_t> GlobalDofIndices(
738 const lf::mesh::Entity &entity) const override;
739
743 [[nodiscard]] std::span<const gdof_idx_t> InteriorGlobalDofIndices(
744 const lf::mesh::Entity &entity) const override;
745
751 gdof_idx_t dofnum) const override {
752 LF_VERIFY_MSG(dofnum < dof_entities_.size(),
753 "Illegal dof index " << dofnum << ", max = " << num_dof_);
754 return *dof_entities_[dofnum];
755 }
756
759 [[nodiscard]] std::shared_ptr<const lf::mesh::Mesh> Mesh() const override {
760 return mesh_p_;
761 }
762
763 private:
765 std::shared_ptr<const lf::mesh::Mesh> mesh_p_;
769 std::vector<const lf::mesh::Entity *> dof_entities_;
782 std::array<std::vector<size_type>, 3> num_int_dofs_;
784 std::array<std::vector<size_type>, 3> offsets_;
787 std::array<std::vector<gdof_idx_t>, 3> dofs_;
789};
790
791} // namespace lf::assemble
792
794
798template <>
799struct fmt::formatter<lf::assemble::DofHandler> : ostream_formatter {};
800
802#endif
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.
void PrintInfo(std::ostream &stream, const DofHandler &dof_handler, unsigned int ctrl=0)
Output information about the given dof handler to the given stream object.
DofHandler(DofHandler &&)=default
DofHandler & operator=(DofHandler &&)=default
virtual size_type NumInteriorDofs(const lf::mesh::Entity &entity) const =0
provides number of shape functions associated with an entity
DofHandler()=default
Default constructor, can only be called from derived class.
DofHandler(const DofHandler &)=default
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
virtual ~DofHandler()=default
virtual Destructor
DofHandler & operator=(const DofHandler &)=default
Dof handler allowing variable local dof layouts.
Definition dofhandler.h:514
std::vector< const lf::mesh::Entity * > dof_entities_
Definition dofhandler.h:769
size_type NumInteriorDofs(const lf::mesh::Entity &entity) const override
provides number of shape functions associated with an entity
DynamicFEDofHandler & operator=(const DynamicFEDofHandler &)=delete
Copy assignment is forbidden.
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
DynamicFEDofHandler(DynamicFEDofHandler &&)=default
A DynamicFEDofHandler can be move constructed.
DynamicFEDofHandler & operator=(DynamicFEDofHandler &&)=default
A DynamicFEDofHandler can be moved into.
DynamicFEDofHandler(std::shared_ptr< const lf::mesh::Mesh > mesh_p, LOCALDOFINFO &&locdof)
Set-up of dof handler.
Definition dofhandler.h:566
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
const lf::mesh::Entity & Entity(gdof_idx_t dofnum) const override
retrieve unique entity at which a basis function is located
Definition dofhandler.h:750
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
Acess to underlying mesh object.
Definition dofhandler.h:759
~DynamicFEDofHandler() override=default
Virtual destructor.
std::array< std::vector< size_type >, 3 > num_int_dofs_
Definition dofhandler.h:782
DynamicFEDofHandler(const DynamicFEDofHandler &)=delete
It doesn't make much sense to copy construct a DynamicFEDofHandler.
size_type NumDofs() const override
total number of dof's handled by the object
Definition dofhandler.h:718
Dofhandler for uniform finite element spaces.
Definition dofhandler.h:260
~UniformFEDofHandler() override=default
Virtual Destructor.
void initIndexArrays()
initialization of internal index arrays
UniformFEDofHandler(UniformFEDofHandler &&)=default
UniformFEDofHandler can be moved.
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
UniformFEDofHandler(const UniformFEDofHandler &)=delete
Copy Construction doesn't make much sense for UniformFEDofHandler.
std::shared_ptr< const lf::mesh::Mesh > Mesh() const override
Acess to underlying mesh object.
Definition dofhandler.h:385
const lf::mesh::Entity & Entity(gdof_idx_t dofnum) const override
retrieve unique entity at which a basis function is located
Definition dofhandler.h:376
std::shared_ptr< const lf::mesh::Mesh > mesh_
Definition dofhandler.h:480
size_type NumInteriorDofs(lf::base::RefEl ref_el_type) const
Definition dofhandler.h:452
UniformFEDofHandler & operator=(UniformFEDofHandler &&)=default
Uniform FEDofHandler can be move assigned to.
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
size_type NumDofs() const override
total number of dof's handled by the object
Definition dofhandler.h:336
void InitTotalNumDofs()
compute number of shape functions covering an entity type
UniformFEDofHandler & operator=(const UniformFEDofHandler &)=delete
Copy assigning a UniformFEDofHandler doesn't make much sense.
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
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
void AssembleMatrixLocally(dim_t codim, const DofHandler &dof_handler_trial, const DofHandler &dof_handler_test, ENTITY_MATRIX_PROVIDER &entity_matrix_provider, TMPMATRIX &matrix)
Assembly function for standard assembly of finite element matrices.
Definition assembler.h:115
D.o.f. index mapping and assembly facilities.
Definition assemble.h:31
lf::base::glb_idx_t glb_idx_t
lf::base::size_type size_type
std::ostream & operator<<(std::ostream &o, const COOMatrix< SCALARTYPE > &mat)
Definition coomatrix.h:229
Eigen::Index gdof_idx_t
Definition assemble.h:31