DofMap.h

Note

The documentation on this page was automatically extracted from the DOLFIN C++ code and may need to be edited or expanded.

class DofMap

Parent class(es)

Degree-of-freedom map This class handles the mapping of degrees of freedom. It builds a dof map based on a ufc::dofmap on a specific mesh. It will reorder the dofs when running in parallel. Sub-dofmaps, both views and copies, are supported.

DofMap(std::shared_ptr<const ufc::dofmap> ufc_dofmap, const Mesh &mesh)

Create dof map on mesh (mesh is not stored)

@param[in] ufc_dofmap (ufc::dofmap)
The ufc::dofmap.
@param[in] mesh (Mesh&)
The mesh.
DofMap(std::shared_ptr<const ufc::dofmap> ufc_dofmap, const Mesh &mesh, std::shared_ptr<const SubDomain> constrained_domain)

Create a periodic dof map on mesh (mesh is not stored)

@param[in] ufc_dofmap (ufc::dofmap)
The ufc::dofmap.
@param[in] mesh (Mesh)
The mesh.
@param[in] constrained_domain (SubDomain)
The subdomain marking the constrained (tied) boundaries.
bool is_view() const

True iff dof map is a view into another map

Returns
bool
True if the dof map is a sub-dof map (a view into another map).
std::size_t global_dimension() const

Return the dimension of the global finite element function space. Use index_map()->size() to get the local dimension.

Returns
std::size_t
The dimension of the global finite element function space.
std::size_t num_element_dofs(std::size_t cell_index) const

Return the dimension of the local finite element function space on a cell

@param cell_index (std::size_t)
Index of cell
@return std::size_t
Dimension of the local finite element function space.
std::size_t max_element_dofs() const

Return the maximum dimension of the local finite element function space

@return std::size_t
Maximum dimension of the local finite element function space.
std::size_t num_entity_dofs(std::size_t entity_dim) const

Return the number of dofs for a given entity dimension

@param entity_dim (std::size_t)
Entity dimension
@return std::size_t
Number of dofs associated with given entity dimension
std::size_t num_entity_closure_dofs(std::size_t entity_dim) const

Return the number of dofs for the closure of an entity of given dimension

Arguments
entity_dim (std::size_t)
Entity dimension
Returns
std::size_t
Number of dofs associated with closure of an entity of given dimension
std::size_t num_facet_dofs() const

Return number of facet dofs

@return std::size_t
The number of facet dofs.
std::pair<std::size_t, std::size_t> ownership_range() const

Return the ownership range (dofs in this range are owned by this process)

@return std::pair<std::size_t, std::size_t>
The ownership range.
const std::vector<int> &off_process_owner() const

Return map from nonlocal dofs that appear in local dof map to owning process

@return std::vector<unsigned int>
The map from non-local dofs.
const std::unordered_map<int, std::vector<int>> &shared_nodes() const

Return map from all shared nodes to the sharing processes (not including the current process) that share it.

@return std::unordered_map<std::size_t, std::vector<unsigned int>>
The map from dofs to list of processes
const std::set<int> &neighbours() const

Return set of processes that share dofs with this process

@return std::set<int>
The set of processes
void clear_sub_map_data()

Clear any data required to build sub-dofmaps (this is to reduce memory use)

Eigen::Map<const Eigen::Array<dolfin::la_index, Eigen::Dynamic, 1>> cell_dofs(std::size_t cell_index) const

Local-to-global mapping of dofs on a cell

@param cell_index (std::size_t)
The cell index.

@return ArrayView<const dolfin::la_index>

std::vector<dolfin::la_index> entity_dofs(const Mesh &mesh, std::size_t entity_dim, const std::vector<std::size_t> &entity_indices) const

Return the dof indices associated with entities of given dimension and entity indices

Arguments
entity_dim (std::size_t)
Entity dimension.
entity_indices (std::vector<dolfin::la_index>&)
Entity indices to get dofs for.
Returns
std::vector<dolfin::la_index>
Dof indices associated with selected entities.
std::vector<dolfin::la_index> entity_dofs(const Mesh &mesh, std::size_t entity_dim) const

Return the dof indices associated with all entities of given dimension

Arguments
entity_dim (std::size_t)
Entity dimension.
Returns
std::vector<dolfin::la_index>
Dof indices associated with selected entities.
std::vector<dolfin::la_index> entity_closure_dofs(const Mesh &mesh, std::size_t entity_dim, const std::vector<std::size_t> &entity_indices) const

Return the dof indices associated with the closure of entities of given dimension and entity indices

Arguments
entity_dim (std::size_t)
Entity dimension.
entity_indices (std::vector<dolfin::la_index>&)
Entity indices to get dofs for.
Returns
std::vector<dolfin::la_index>
Dof indices associated with selected entities and their closure.
std::vector<dolfin::la_index> entity_closure_dofs(const Mesh &mesh, std::size_t entity_dim) const

Return the dof indices associated with the closure of all entities of given dimension

@param mesh (Mesh)
Mesh
@param entity_dim (std::size_t)
Entity dimension.
@return std::vector<dolfin::la_index>
Dof indices associated with selected entities and their closure.
void tabulate_facet_dofs(std::vector<std::size_t> &element_dofs, std::size_t cell_facet_index) const

Tabulate local-local facet dofs

@param element_dofs (std::size_t)
Degrees of freedom on a single element.
@param cell_facet_index (std::size_t)
The local facet index on the cell.
void tabulate_entity_dofs(std::vector<std::size_t> &element_dofs, std::size_t entity_dim, std::size_t cell_entity_index) const

Tabulate local-local mapping of dofs on entity (dim, local_entity)

@param element_dofs (std::size_t)
Degrees of freedom on a single element.
@param entity_dim (std::size_t)
The entity dimension.
@param cell_entity_index (std::size_t)
The local entity index on the cell.
void tabulate_entity_closure_dofs(std::vector<std::size_t> &element_dofs, std::size_t entity_dim, std::size_t cell_entity_index) const

Tabulate local-local mapping of dofs on closure of entity (dim, local_entity)

@param element_dofs (std::size_t)
Degrees of freedom on a single element.
@param entity_dim (std::size_t)
The entity dimension.
@param cell_entity_index (std::size_t)
The local entity index on the cell.
void tabulate_global_dofs(std::vector<std::size_t> &element_dofs) const

Tabulate globally supported dofs

@param element_dofs (std::size_t)
Degrees of freedom.
std::shared_ptr<GenericDofMap> copy() const

Create a copy of the dof map

@return DofMap
The Dofmap copy.
std::shared_ptr<GenericDofMap> create(const Mesh &new_mesh) const

Create a copy of the dof map on a new mesh

@param new_mesh (Mesh)
The new mesh to create the dof map on.
@return DofMap
The new Dofmap copy.
std::shared_ptr<GenericDofMap> extract_sub_dofmap(const std::vector<std::size_t> &component, const Mesh &mesh) const

Extract subdofmap component

@param component (std::vector<std::size_t>)
The component.
@param mesh (Mesh)
The mesh.
@return DofMap
The subdofmap component.
std::shared_ptr<GenericDofMap> collapse(std::unordered_map<std::size_t, std::size_t> &collapsed_map, const Mesh &mesh) const

Create a “collapsed” dofmap (collapses a sub-dofmap)

@param collapsed_map (std::unordered_map<std::size_t, std::size_t>)
The “collapsed” map.
@param mesh (Mesh)
The mesh.
@return DofMap
The collapsed dofmap.
std::vector<dolfin::la_index> dofs(const Mesh &mesh, std::size_t dim) const

Return list of dof indices on this process that belong to mesh entities of dimension dim

void set(GenericVector &x, double value) const

Set dof entries in vector to a specified value. Parallel layout of vector must be consistent with dof map range. This function is typically used to construct the null space of a matrix operator.

@param x (GenericVector)
The vector to set.
@param value (double)
The value to set.
std::shared_ptr<const IndexMap> index_map() const

Return the map (const access)

int block_size() const

Return the block size for dof maps with components, typically used for vector valued functions.

void tabulate_local_to_global_dofs(std::vector<std::size_t> &local_to_global_map) const

Compute the map from local (this process) dof indices to global dof indices.

@param local_to_global_map (_std::vector<std::size_t>_)
The local-to-global map to fill.
std::size_t local_to_global_index(int local_index) const

Return global dof index for a given local (process) dof index

@param local_index (int)
The local local index.
@return std::size_t
The global dof index.
const std::vector<std::size_t> &local_to_global_unowned() const

Return indices of dofs which are owned by other processes

std::string str(bool verbose) const

Return informal string representation (pretty-print)

@param verbose (bool)
Flag to turn on additional output.
@return std::string
An informal representation of the function space.