Refinement (dolfinx::refinement
)
-
namespace refinement
Mesh refinement algorithms.
Methods for refining meshes uniformly, or with markers, using edge bisection.
Functions
-
mesh::Mesh refine(const mesh::Mesh &mesh, bool redistribute = true)
Create a uniformly refined mesh.
- Parameters
mesh – [in] The mesh from which to build a refined Mesh
redistribute – [in] Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
- Returns
A refined mesh
-
mesh::Mesh refine(const mesh::Mesh &mesh, const std::span<const std::int32_t> &edges, bool redistribute = true)
Create a locally refined mesh.
- Parameters
mesh – [in] The mesh from which to build a refined Mesh
edges – [in] Indices of the edges that should be split by this refinement. mesh::compute_incident_entities can be used to compute the edges that are incident to other entities, e.g. incident to cells.
redistribute – [in] Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
- Returns
A locally refined mesh
-
void update_logical_edgefunction(MPI_Comm neighbor_comm, const std::vector<std::vector<std::int32_t>> &marked_for_update, std::vector<std::int8_t> &marked_edges, const common::IndexMap &map)
Communicate edge markers between processes that share edges.
- Parameters
neighbor_comm – [in] MPI Communicator for neighborhood
marked_for_update – [in] Lists of edges to be updated on each neighbor.
marked_for_update[r]
is the list of edge indices that are marked by the caller and are shared with local MPI rankr
.marked_edges – [inout] Marker for each edge on the calling process
map – [in] Index map for the mesh edges
-
std::tuple<std::map<std::int32_t, std::int64_t>, std::vector<double>, std::array<std::size_t, 2>> create_new_vertices(MPI_Comm neighbor_comm, const graph::AdjacencyList<int> &shared_edges, const mesh::Mesh &mesh, const std::vector<std::int8_t> &marked_edges)
Add new vertex for each marked edge, and create new_vertex_coordinates and global_edge->new_vertex map.
Communicate new vertices with MPI to all affected processes.
- Parameters
neighbor_comm – [in] MPI Communicator for neighborhood
shared_edges – [in]
mesh – [in] Existing mesh
marked_edges – [in]
- Returns
(0) map from local edge index to new vertex global index, (1) the coordinates of the new vertices (row-major storage) and (2) the shape of the new coordinates.
-
mesh::Mesh partition(const mesh::Mesh &old_mesh, const graph::AdjacencyList<std::int64_t> &cell_topology, std::span<const double> new_coords, std::array<std::size_t, 2> xshape, bool redistribute, mesh::GhostMode ghost_mode)
Use vertex and topology data to partition new mesh across processes.
- Parameters
old_mesh – [in]
cell_topology – [in] Topology of cells, (vertex indices)
new_coords – [in] New coordinates, row-major storage
xshape – [in] The shape of
new_coords
redistribute – [in] Call graph partitioner if true
ghost_mode – [in] None or shared_facet
- Returns
New mesh
-
std::vector<std::int64_t> adjust_indices(const common::IndexMap &map, std::int32_t n)
Add indices to account for extra n values on this process.
- Todo:
Fix docstring. It is unclear.
This is a utility to help add new topological vertices on each process into the space of the index map.
- Parameters
map – [in] Index map for the current mesh vertices
n – [in] Number of new entries to be accommodated on this process
- Returns
Global indices as if “n” extra values are appended on each process
Transfer facet MeshTags from coarse mesh to refined mesh.
Note
The refined mesh must not have been redistributed during refinement
Note
GhostMode must be GhostMode.none
- Parameters
meshtag – [in] Facet tags on parent mesh
refined_mesh – [in] Refined mesh based on parent mesh
cell – [in] Parent cell of each cell in refined mesh
facet – [in] Local facets of parent in each cell in refined mesh
- Returns
MeshTags on refined mesh
Transfer cell MeshTags from coarse mesh to refined mesh.
Note
The refined mesh must not have been redistributed during refinement
Note
GhostMode must be GhostMode.none
- Parameters
parent_meshtag – [in] Cell MeshTags on parent mesh
refined_mesh – [in] Refined mesh based on parent mesh
parent_cell – [in] Parent cell of each cell in refined mesh
- Returns
MeshTags on refined mesh, values copied over from coarse mesh
-
namespace plaza
Function in this namespace implement the refinement method described in Plaza and Carey “Local refinement of simplicial grids based on
the skeleton” (Applied Numerical Mathematics 32 (2000) 195-218).
Enums
-
enum class RefinementOptions : int
Selection of options when refining a Mesh.
parent_cell
will output a list containing the local parent cell index for each new cell,parent_facet
will output a list of the cell-local facet indices in the parent cell of each facet in each new cell (or -1 if no match).parent_cell_and_facet
will output both datasets.Values:
-
enumerator none
-
enumerator parent_cell
-
enumerator parent_facet
-
enumerator parent_cell_and_facet
-
enumerator none
Functions
-
std::tuple<mesh::Mesh, std::vector<std::int32_t>, std::vector<std::int8_t>> refine(const mesh::Mesh &mesh, bool redistribute, RefinementOptions options)
Uniform refine, optionally redistributing and optionally calculating the parent-child relationships, selected by
RefinementOptions
.- Parameters
mesh – [in] Input mesh to be refined
redistribute – [in] Flag to call the mesh partitioner to redistribute after refinement
options – [in] RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned.
- Returns
New Mesh and optional parent cell index, parent facet indices
-
std::tuple<mesh::Mesh, std::vector<std::int32_t>, std::vector<std::int8_t>> refine(const mesh::Mesh &mesh, const std::span<const std::int32_t> &edges, bool redistribute, RefinementOptions options)
Refine with markers, optionally redistributing, and optionally calculating the parent-child relationships, selected by
RefinementOptions
.- Parameters
mesh – [in] Input mesh to be refined
edges – [in] Indices of the edges that should be split by this refinement
redistribute – [in] Flag to call the Mesh Partitioner to redistribute after refinement
options – [in] RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned.
- Returns
New Mesh and optional parent cell index, parent facet indices
-
std::tuple<graph::AdjacencyList<std::int64_t>, std::vector<double>, std::array<std::size_t, 2>, std::vector<std::int32_t>, std::vector<std::int8_t>> compute_refinement_data(const mesh::Mesh &mesh, RefinementOptions options)
Refine mesh returning new mesh data.
- Parameters
mesh – [in] Input mesh to be refined
options – [in] RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned.
- Returns
New mesh data: cell topology, vertex coordinates, vertex coordinates shape, and optional parent cell index, and parent facet indices.
-
std::tuple<graph::AdjacencyList<std::int64_t>, std::vector<double>, std::array<std::size_t, 2>, std::vector<std::int32_t>, std::vector<std::int8_t>> compute_refinement_data(const mesh::Mesh &mesh, const std::span<const std::int32_t> &edges, RefinementOptions options)
Refine with markers returning new mesh data.
- Parameters
mesh – [in] Input mesh to be refined
edges – [in] Indices of the edges that should be split by this refinement
options – [in] RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned.
- Returns
New mesh data: cell topology, vertex coordinates and parent cell index, and stored parent facet indices (if requested).
-
enum class RefinementOptions : int
-
mesh::Mesh refine(const mesh::Mesh &mesh, bool redistribute = true)