Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/dfe/namespacedolfinx_1_1refinement_1_1plaza.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Enumerations | Functions
dolfinx::refinement::plaza Namespace Reference

Plaza mesh refinement. More...

Enumerations

enum class  Option : int { none = 0 , parent_cell = 1 , parent_facet = 2 , parent_cell_and_facet = 3 }
 Options for data to compute during mesh refinement. More...
 

Functions

template<std::floating_point T>
std::tuple< mesh::Mesh< T >, std::vector< std::int32_t >, std::vector< std::int8_t > > refine (const mesh::Mesh< T > &mesh, bool redistribute, Option option)
 Uniform refine, optionally redistributing and optionally calculating the parent-child relationships`.
 
template<std::floating_point T>
std::tuple< mesh::Mesh< T >, std::vector< std::int32_t >, std::vector< std::int8_t > > refine (const mesh::Mesh< T > &mesh, std::span< const std::int32_t > edges, bool redistribute, Option option)
 Refine with markers, optionally redistributing, and optionally calculating the parent-child relationships.
 
template<std::floating_point T>
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< T >, std::array< std::size_t, 2 >, std::vector< std::int32_t >, std::vector< std::int8_t > > compute_refinement_data (const mesh::Mesh< T > &mesh, Option option)
 Refine mesh returning new mesh data.
 
template<std::floating_point T>
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< T >, std::array< std::size_t, 2 >, std::vector< std::int32_t >, std::vector< std::int8_t > > compute_refinement_data (const mesh::Mesh< T > &mesh, std::span< const std::int32_t > edges, Option option)
 

Detailed Description

Plaza mesh refinement.

Functions for the refinement method described in Plaza and Carey "Local refinement of simplicial grids based on the skeleton", Applied Numerical Mathematics 32 (2000), 195-218.

Enumeration Type Documentation

◆ Option

enum class Option : int
strong

Options for data to compute during mesh refinement.

Enumerator
none 

No extra data

parent_cell 

Compute list with the parent cell index for each new cell

parent_facet 

Compute 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 

Both cell and facet parent data

Function Documentation

◆ compute_refinement_data() [1/2]

template<std::floating_point T>
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< T >, std::array< std::size_t, 2 >, std::vector< std::int32_t >, std::vector< std::int8_t > > compute_refinement_data ( const mesh::Mesh< T > & mesh,
Option option )

Refine mesh returning new mesh data.

Parameters
[in]meshInput mesh to be refined
[in]optionControl computation of parent facets and 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.

◆ compute_refinement_data() [2/2]

template<std::floating_point T>
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< T >, std::array< std::size_t, 2 >, std::vector< std::int32_t >, std::vector< std::int8_t > > compute_refinement_data ( const mesh::Mesh< T > & mesh,
std::span< const std::int32_t > edges,
Option option )

Refine with markers returning new mesh data.

Parameters
[in]meshInput mesh to be refined
[in]edgesIndices of the edges that should be split by this refinement
[in]optionControl 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).

◆ refine() [1/2]

template<std::floating_point T>
std::tuple< mesh::Mesh< T >, std::vector< std::int32_t >, std::vector< std::int8_t > > refine ( const mesh::Mesh< T > & mesh,
bool redistribute,
Option option )

Uniform refine, optionally redistributing and optionally calculating the parent-child relationships`.

Parameters
[in]meshInput mesh to be refined
[in]redistributeFlag to call the mesh partitioner to redistribute after refinement
[in]optionControl the computation of parent facets, parent cells. If an option is unselected, an empty list is returned.
Returns
Refined mesh and optional parent cell index, parent facet indices

◆ refine() [2/2]

template<std::floating_point T>
std::tuple< mesh::Mesh< T >, std::vector< std::int32_t >, std::vector< std::int8_t > > refine ( const mesh::Mesh< T > & mesh,
std::span< const std::int32_t > edges,
bool redistribute,
Option option )

Refine with markers, optionally redistributing, and optionally calculating the parent-child relationships.

Parameters
[in]meshInput mesh to be refined
[in]edgesIndices of the edges that should be split by this refinement
[in]redistributeFlag to call the Mesh Partitioner to redistribute after refinement
[in]optionControl 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