Degree-of-freedom map.
More...
#include <DofMap.h>
|
template<typename E , typename U >
requires std::is_convertible_v<std::remove_cvref_t<E>, fem::ElementDofLayout> and std::is_convertible_v<std::remove_cvref_t<U>, std::vector<std::int32_t>> |
| DofMap (E &&element, std::shared_ptr< const common::IndexMap > index_map, int index_map_bs, U &&dofmap, int bs) |
| Create a DofMap from the layout of dofs on a reference element, an IndexMap defining the distribution of dofs across processes and a vector of indices.
|
|
| DofMap (const DofMap &dofmap)=delete |
|
| DofMap (DofMap &&dofmap)=default |
| Move constructor.
|
|
DofMap & | operator= (const DofMap &dofmap)=delete |
|
DofMap & | operator= (DofMap &&dofmap)=default |
| Move assignment.
|
|
bool | operator== (const DofMap &map) const |
| Equality operator.
|
|
std::span< const std::int32_t > | cell_dofs (std::int32_t c) const |
| Local-to-global mapping of dofs on a cell.
|
|
int | bs () const noexcept |
| Return the block size for the dofmap.
|
|
DofMap | extract_sub_dofmap (std::span< const int > component) const |
| Extract subdofmap component.
|
|
std::pair< DofMap, std::vector< std::int32_t > > | collapse (MPI_Comm comm, const mesh::Topology &topology, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &&reorder_fn=nullptr) const |
| Create a "collapsed" dofmap (collapses a sub-dofmap)
|
|
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > | map () const |
| Get dofmap data.
|
|
const ElementDofLayout & | element_dof_layout () const |
| Layout of dofs on an element.
|
|
int | index_map_bs () const |
| Block size associated with the index_map.
|
|
|
std::shared_ptr< const common::IndexMap > | index_map |
| Index map that describes the parallel distribution of the dofmap.
|
|
Degree-of-freedom map.
This class handles the mapping of degrees of freedom. It builds a dof map based on an ElementDofLayout on a specific mesh topology. It will reorder the dofs when running in parallel. Sub-dofmaps, both views and copies, are supported.
◆ DofMap()
template<typename E , typename U >
requires std::is_convertible_v<std::remove_cvref_t<E>,
fem::ElementDofLayout> and std::is_convertible_v<std::remove_cvref_t<U>, std::vector<std::int32_t>>
DofMap |
( |
E && | element, |
|
|
std::shared_ptr< const common::IndexMap > | index_map, |
|
|
int | index_map_bs, |
|
|
U && | dofmap, |
|
|
int | bs ) |
|
inline |
Create a DofMap from the layout of dofs on a reference element, an IndexMap defining the distribution of dofs across processes and a vector of indices.
- Parameters
-
[in] | element | The layout of the degrees of freedom on an element |
[in] | index_map | The map describing the parallel distribution of the degrees of freedom. |
[in] | index_map_bs | The block size associated with the index_map . |
[in] | dofmap | Adjacency list with the degrees-of-freedom for each cell. |
[in] | bs | The block size of the dofmap . |
◆ cell_dofs()
std::span< const std::int32_t > cell_dofs |
( |
std::int32_t | c | ) |
const |
|
inline |
Local-to-global mapping of dofs on a cell.
- Parameters
-
- Returns
- Local-global dof map for the cell (using process-local indices)
◆ collapse()
std::pair< DofMap, std::vector< std::int32_t > > collapse |
( |
MPI_Comm | comm, |
|
|
const mesh::Topology & | topology, |
|
|
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> && | reorder_fn = nullptr ) const |
Create a "collapsed" dofmap (collapses a sub-dofmap)
- Parameters
-
[in] | comm | MPI Communicator |
[in] | topology | Mesh topology that the dofmap is defined on |
[in] | reorder_fn | Graph re-ordering function to apply to the dof data |
- Returns
- The collapsed dofmap
◆ extract_sub_dofmap()
DofMap extract_sub_dofmap |
( |
std::span< const int > | component | ) |
const |
Extract subdofmap component.
- Parameters
-
[in] | component | The component indices |
- Returns
- The dofmap for the component
◆ map()
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > map |
( |
| ) |
const |
Get dofmap data.
- Returns
- The adjacency list with dof indices for each cell
◆ operator==()
bool operator== |
( |
const DofMap & | map | ) |
const |
Equality operator.
- Returns
- Returns true if the data for the two dofmaps are equal
The documentation for this class was generated from the following files:
- /__w/dolfinx/dolfinx/cpp/dolfinx/fem/DofMap.h
- /__w/dolfinx/dolfinx/cpp/dolfinx/fem/DofMap.cpp