dolfin/refinement

Documentation for C++ code found in dolfin/refinement/*.h

Functions

p_refine

C++ documentation for p_refine from dolfin/refinement/refine.h:

void dolfin::p_refine(Mesh &refined_mesh, const Mesh &mesh)

Increase the polynomial order of the mesh from 1 to 2, i.e. add points at the Edge midpoints, to make a quadratic mesh.

Parameters:
  • refined_mesh – (Mesh ) The mesh that will be the quadratic mesh.
  • mesh – (Mesh ) The original linear mesh.

C++ documentation for p_refine from dolfin/refinement/refine.h:

Mesh dolfin::p_refine(const Mesh &mesh)

Return a p_refined mesh Increase the polynomial order of the mesh from 1 to 2, i.e. add points at the Edge midpoints, to make a quadratic mesh.

Parameters:mesh – (Mesh ) The original linear mesh.
Returns:Mesh

refine

C++ documentation for refine from dolfin/refinement/refine.h:

void dolfin::refine(Mesh &refined_mesh, const Mesh &mesh, bool redistribute = true)

Create uniformly refined mesh

Parameters:
  • refined_mesh – (Mesh ) The mesh that will be the refined mesh.
  • mesh – (Mesh ) The original mesh.
  • redistribute – (bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.

C++ documentation for refine from dolfin/refinement/refine.h:

void dolfin::refine(Mesh &refined_mesh, const Mesh &mesh, const MeshFunction<bool> &cell_markers, bool redistribute = true)

Create locally refined mesh

Parameters:
  • refined_mesh – (Mesh ) The mesh that will be the refined mesh.
  • mesh – (Mesh ) The original mesh.
  • cell_markers – (MeshFunction<bool>) A mesh function over booleans specifying which cells that should be refined (and which should not).
  • redistribute – (bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.

C++ documentation for refine from dolfin/refinement/refine.h:

Mesh dolfin::refine(const Mesh &mesh, bool redistribute = true)

Create uniformly refined mesh

mesh = refine(mesh);
Parameters:
  • mesh – (Mesh ) The mesh to refine.
  • redistribute – (bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
Returns:

Mesh The refined mesh.

C++ documentation for refine from dolfin/refinement/refine.h:

Mesh dolfin::refine(const Mesh &mesh, const MeshFunction<bool> &cell_markers, bool redistribute = true)

Create locally refined mesh

MeshFunction<bool> cell_markers(mesh, mesh->topology().dim());
cell_markers.set_all(false);
Point origin(0.0, 0.0, 0.0);
for (CellIterator cell(mesh); !cell.end(); ++cell)
{
Point p = cell->midpoint();
if (p.distance(origin) < 0.1)
cell_markers[*cell] = true;
}
mesh = refine(mesh, cell_markers);
Parameters:
  • mesh – (Mesh ) The mesh to refine.
  • cell_markers – (MeshFunction<bool>) A mesh function over booleans specifying which cells that should be refined (and which should not).
  • redistribute – (bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
Returns:

Mesh The locally refined mesh.

C++ documentation for refine from dolfin/refinement/refine.h:

std::shared_ptr<const MeshHierarchy> dolfin::refine(const MeshHierarchy &hierarchy, const MeshFunction<bool> &markers)

Refine a MeshHierarchy .

Parameters:
  • hierarchy
  • markers

Classes

BisectionRefinement1D

C++ documentation for BisectionRefinement1D from dolfin/refinement/BisectionRefinement1D.h:

class dolfin::BisectionRefinement1D

This class implements mesh refinement in 1D.

void dolfin::BisectionRefinement1D::refine(Mesh &refined_mesh, const Mesh &mesh, bool redistribute = false)

Refine mesh uniformly

Parameters:
  • refined_mesh – (Mesh )
  • mesh – (const Mesh )
  • redistribute – (bool)
void dolfin::BisectionRefinement1D::refine(Mesh &refined_mesh, const Mesh &mesh, const MeshFunction<bool> &cell_markers, bool redistribute = false)

Refine mesh based on cell markers

Parameters:
  • refined_mesh – (Mesh )
  • mesh – (const Mesh )
  • cell_markers – (const MeshFunction<bool>)
  • redistribute – (bool)

LocalMeshCoarsening

C++ documentation for LocalMeshCoarsening from dolfin/refinement/LocalMeshCoarsening.h:

class dolfin::LocalMeshCoarsening

This class implements local mesh coarsening for different mesh types.

bool dolfin::LocalMeshCoarsening::coarsen_cell(Mesh &mesh, Mesh &coarse_mesh, int cell_id, std::vector<int> &old2new_vertex, std::vector<int> &old2new_cell, bool coarsen_boundary = false)

Coarsen simplicial cell by edge collapse.

Parameters:
  • mesh
  • coarse_mesh
  • cell_id
  • old2new_vertex
  • old2new_cell
  • coarsen_boundary
void dolfin::LocalMeshCoarsening::coarsen_mesh_by_edge_collapse(Mesh &mesh, MeshFunction<bool> &cell_marker, bool coarsen_boundary = false)

Coarsen simplicial mesh locally by edge collapse.

Parameters:
  • mesh
  • cell_marker
  • coarsen_boundary
bool dolfin::LocalMeshCoarsening::coarsen_mesh_ok(Mesh &mesh, std::size_t edge_index, std::size_t *edge_vertex, MeshFunction<bool> &vertex_forbidden)

Check that edge collapse is ok.

Parameters:
  • mesh
  • edge_index
  • edge_vertex
  • vertex_forbidden
void dolfin::LocalMeshCoarsening::collapse_edge(Mesh &mesh, Edge &edge, Vertex &vertex_to_remove, MeshFunction<bool> &cell_to_remove, std::vector<int> &old2new_vertex, std::vector<int> &old2new_cell, MeshEditor &editor, std::size_t &current_cell)

Collapse edge by node deletion.

Parameters:
  • mesh
  • edge
  • vertex_to_remove
  • cell_to_remove
  • old2new_vertex
  • old2new_cell
  • editor
  • current_cell

ParallelRefinement

C++ documentation for ParallelRefinement from dolfin/refinement/ParallelRefinement.h:

class dolfin::ParallelRefinement

Data structure and methods for refining meshes in parallel. ParallelRefinement encapsulates two main features: a distributed MeshFunction defined over the mesh edes, which can be updated across processes, and storage for local mesh data, which can be used to construct the new Mesh

dolfin::ParallelRefinement::ParallelRefinement(const Mesh &mesh)

Constructor.

Parameters:mesh
void dolfin::ParallelRefinement::build_local(Mesh &new_mesh) const

Build local mesh from internal data when not running in parallel

Parameters:new_mesh – (Mesh )
void dolfin::ParallelRefinement::create_new_vertices()

Add new vertex for each marked edge, and create new_vertex_coordinates and global_edge->new_vertex mapping. Communicate new vertices with MPI to all affected processes.

std::shared_ptr<const std::map<std::size_t, std::size_t>> dolfin::ParallelRefinement::edge_to_new_vertex() const

Mapping of old edge (to be removed) to new global vertex number. Useful for forming new topology

bool dolfin::ParallelRefinement::is_marked(std::size_t edge_index) const

Return marked status of edge

Parameters:edge_index – (std::size_t)
std::shared_ptr<std::map<std::size_t, std::size_t>> dolfin::ParallelRefinement::local_edge_to_new_vertex
void dolfin::ParallelRefinement::mark(const MeshEntity &cell)

Mark all incident edges of an entity

Parameters:cell – (MeshEntity )
void dolfin::ParallelRefinement::mark(const MeshFunction<bool> &refinement_marker)

Mark all edges incident on entities indicated by refinement marker

Parameters:refinement_marker – (const MeshFunction<bool>)
void dolfin::ParallelRefinement::mark(std::size_t edge_index)

Mark edge by index

Parameters:edge_index – (std::size_t) Index of edge to mark
void dolfin::ParallelRefinement::mark_all()

Mark all edges in mesh.

std::vector<std::size_t> dolfin::ParallelRefinement::marked_edge_list(const MeshEntity &cell) const

Return list of marked edges incident on this MeshEntity - usually a cell

Parameters:cell – (const MeshEntity )
std::vector<bool> dolfin::ParallelRefinement::marked_edges
std::vector<std::vector<std::size_t>> dolfin::ParallelRefinement::marked_for_update
const Mesh &dolfin::ParallelRefinement::mesh() const

Original mesh associated with this refinement.

void dolfin::ParallelRefinement::new_cell(const Cell &cell)

Add a new cell to the list in 3D or 2D

Parameters:cell – (const Cell )
void dolfin::ParallelRefinement::new_cell(std::size_t i0, std::size_t i1, std::size_t i2)

Add a new cell with vertex indices

Parameters:
  • i0 – (std::size_t)
  • i1 – (std::size_t)
  • i2 – (std::size_t)
void dolfin::ParallelRefinement::new_cell(std::size_t i0, std::size_t i1, std::size_t i2, std::size_t i3)

Add a new cell with vertex indices

Parameters:
  • i0 – (std::size_t)
  • i1 – (std::size_t)
  • i2 – (std::size_t)
  • i3 – (std::size_t)
std::vector<std::size_t> dolfin::ParallelRefinement::new_cell_topology
void dolfin::ParallelRefinement::new_cells(const std::vector<std::size_t> &idx)

Add new cells with vertex indices

Parameters:idx – (const std::vector<std::size_t>)
std::vector<double> dolfin::ParallelRefinement::new_vertex_coordinates
void dolfin::ParallelRefinement::partition(Mesh &new_mesh, bool redistribute) const

Use vertex and topology data to partition new mesh across processes

Parameters:
  • new_mesh – (Mesh )
  • redistribute – (bool)
std::unordered_map<unsigned int, std::vector<std::pair<unsigned int, unsigned int>>> dolfin::ParallelRefinement::shared_edges
void dolfin::ParallelRefinement::update_logical_edgefunction()

Transfer marked edges between processes.

dolfin::ParallelRefinement::~ParallelRefinement()

Destructor.

PlazaRefinementND

C++ documentation for PlazaRefinementND from dolfin/refinement/PlazaRefinementND.h:

class dolfin::PlazaRefinementND

Implementation of the refinement method described in Plaza and Carey “Local refinement of simplicial grids based on the skeleton” (Applied Numerical Mathematics 32 (2000) 195-218)

void dolfin::PlazaRefinementND::do_refine(Mesh &new_mesh, const Mesh &mesh, ParallelRefinement &p_ref, const std::vector<unsigned int> &long_edge, const std::vector<bool> &edge_ratio_ok, bool redistribute, bool calculate_parent_facets, MeshRelation &mesh_relation)
Parameters:
  • new_mesh
  • mesh
  • p_ref
  • long_edge
  • edge_ratio_ok
  • redistribute
  • calculate_parent_facets
  • mesh_relation
void dolfin::PlazaRefinementND::enforce_rules(ParallelRefinement &p_ref, const Mesh &mesh, const std::vector<unsigned int> &long_edge)
Parameters:
  • p_ref
  • mesh
  • long_edge
void dolfin::PlazaRefinementND::face_long_edge(std::vector<unsigned int> &long_edge, std::vector<bool> &edge_ratio_ok, const Mesh &mesh)
Parameters:
  • long_edge
  • edge_ratio_ok
  • mesh
void dolfin::PlazaRefinementND::get_simplices(std::vector<std::size_t> &simplex_set, const std::vector<bool> &marked_edges, const std::vector<std::size_t> &longest_edge, std::size_t tdim, bool uniform)

Get the subdivision of an original simplex into smaller simplices, for a given set of marked edges, and the longest edge of each facet (cell local indexing). A flag indicates if a uniform subdivision is preferable in 2D.

Parameters:
  • simplex_set – Returned set of triangles/tets topological description
  • marked_edgesVector indicating which edges are to be split
  • longest_edgeVector indicating the longest edge for each triangle. For tdim=2, one entry, for tdim=3, four entries.
  • tdim – Topological dimension (2 or 3)
  • uniform – Make a “uniform” subdivision with all triangles being similar shape
void dolfin::PlazaRefinementND::get_tetrahedra(std::vector<std::size_t> &tet_set, const std::vector<bool> &marked_edges, const std::vector<std::size_t> &longest_edge)
Parameters:
  • tet_set
  • marked_edges
  • longest_edge
void dolfin::PlazaRefinementND::get_triangles(std::vector<std::size_t> &tri_set, const std::vector<bool> &marked_edges, const std::size_t longest_edge, bool uniform)
Parameters:
  • tri_set
  • marked_edges
  • longest_edge
  • uniform
void dolfin::PlazaRefinementND::refine(Mesh &new_mesh, const Mesh &mesh, bool redistribute, bool calculate_parent_facets)

Uniform refine, optionally redistributing and optionally calculating the parent-child relation for facets (in 2D)

Parameters:
  • new_mesh – New Mesh
  • mesh – Input mesh to be refined
  • redistribute – Flag to call the Mesh Partitioner to redistribute after refinement
  • calculate_parent_facets – Flag to build parent facet information, needed to propagate information on boundaries
void dolfin::PlazaRefinementND::refine(Mesh &new_mesh, const Mesh &mesh, const MeshFunction<bool> &refinement_marker, bool calculate_parent_facets, MeshRelation &mesh_relation)

Refine with markers, optionally calculating facet relations, and saving relation data in MeshRelation structure

Parameters:
  • new_mesh – New Mesh
  • mesh – Input mesh to be refined
  • refinement_markerMeshFunction listing MeshEntities which should be split by this refinement
  • calculate_parent_facets – Flag to build parent facet information, needed to propagate information on boundaries
  • mesh_relation – New relationship between the two meshes
void dolfin::PlazaRefinementND::refine(Mesh &new_mesh, const Mesh &mesh, const MeshFunction<bool> &refinement_marker, bool redistribute, bool calculate_parent_facets)

Refine with markers, optionally redistributing and optionally calculating the parent-child relation for facets (in 2D)

Parameters:
  • new_mesh – New Mesh
  • mesh – Input mesh to be refined
  • refinement_markerMeshFunction listing MeshEntities which should be split by this refinement
  • redistribute – Flag to call the Mesh Partitioner to redistribute after refinement
  • calculate_parent_facets – Flag to build parent facet information, needed to propagate information on boundaries
void dolfin::PlazaRefinementND::set_parent_facet_markers(const Mesh &mesh, Mesh &new_mesh, const std::map<std::size_t, std::size_t> &new_vertex_map)
Parameters:
  • mesh
  • new_mesh
  • new_vertex_map

RegularCutRefinement

C++ documentation for RegularCutRefinement from dolfin/refinement/RegularCutRefinement.h:

class dolfin::RegularCutRefinement

This class implements local mesh refinement by a regular cut of each cell marked for refinement in combination with propagation of cut edges to neighboring cells.

void dolfin::RegularCutRefinement::compute_markers(std::vector<int> &refinement_markers, IndexSet &marked_edges, const Mesh &mesh, const MeshFunction<bool> &cell_markers)
Parameters:
  • refinement_markers
  • marked_edges
  • mesh
  • cell_markers
std::size_t dolfin::RegularCutRefinement::count_markers(const std::vector<bool> &markers)
Parameters:markers
std::size_t dolfin::RegularCutRefinement::extract_edge(const std::vector<bool> &markers)
Parameters:markers
std::pair<std::size_t, std::size_t> dolfin::RegularCutRefinement::find_bisection_edges(const Cell &cell, const Mesh &mesh, std::size_t bisection_twin)
Parameters:
  • cell
  • mesh
  • bisection_twin
std::pair<std::size_t, std::size_t> dolfin::RegularCutRefinement::find_bisection_vertices(const Cell &cell, const Mesh &mesh, std::size_t bisection_twin, const std::pair<std::size_t, std::size_t> &bisection_edges)
Parameters:
  • cell
  • mesh
  • bisection_twin
  • bisection_edges
std::pair<std::size_t, std::size_t> dolfin::RegularCutRefinement::find_common_edges(const Cell &cell, const Mesh &mesh, std::size_t bisection_twin)
Parameters:
  • cell
  • mesh
  • bisection_twin
enum dolfin::RegularCutRefinement::marker_type
enumerator dolfin::RegularCutRefinement::marker_type::no_refinement = -1
enumerator dolfin::RegularCutRefinement::marker_type::regular_refinement = -2
enumerator dolfin::RegularCutRefinement::marker_type::backtrack_bisection = -3
enumerator dolfin::RegularCutRefinement::marker_type::backtrack_bisection_refine = -4
void dolfin::RegularCutRefinement::refine(Mesh &refined_mesh, const Mesh &mesh, const MeshFunction<bool> &cell_markers)

Refine mesh based on cell markers

Parameters:
  • refined_mesh – New Mesh
  • mesh – Input Mesh
  • cell_markers – Markers of MeshEntities to split (typically Cells)
void dolfin::RegularCutRefinement::refine_marked(Mesh &refined_mesh, const Mesh &mesh, const std::vector<int> &refinement_markers, const IndexSet &marked_edges)
Parameters:
  • refined_mesh
  • mesh
  • refinement_markers
  • marked_edges
bool dolfin::RegularCutRefinement::too_thin(const Cell &cell, const std::vector<bool> &edge_markers)
Parameters:
  • cell
  • edge_markers