DOLFINx
0.5.1
DOLFINx C++ interface
|
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). More...
Enumerations | |
enum class | RefinementOptions : int { none = 0 , parent_cell = 1 , parent_facet = 2 , parent_cell_and_facet = 3 } |
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. | |
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 . More... | |
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 . More... | |
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. More... | |
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. More... | |
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).
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.
[in] | mesh | Input mesh to be refined |
[in] | edges | Indices of the edges that should be split by this refinement |
[in] | options | RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned. |
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.
[in] | mesh | Input mesh to be refined |
[in] | options | RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned. |
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
.
[in] | mesh | Input mesh to be refined |
[in] | redistribute | Flag to call the mesh partitioner to redistribute after refinement |
[in] | options | RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned. |
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
.
[in] | mesh | Input mesh to be refined |
[in] | edges | Indices of the edges that should be split by this refinement |
[in] | redistribute | Flag to call the Mesh Partitioner to redistribute after refinement |
[in] | options | RefinementOptions enum to choose the computation of parent facets, parent cells. If an option is unselected, an empty list is returned. |