Data structure and methods for refining meshes in parallel.
More...
#include <ParallelRefinement.h>
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
◆ build_local()
void ParallelRefinement::build_local |
( |
Mesh & |
new_mesh | ) |
const |
Build local mesh from internal data when not running in parallel
- Parameters
-
◆ create_new_vertices()
void 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.
◆ edge_to_new_vertex()
std::shared_ptr< const std::map< std::size_t, std::size_t > > ParallelRefinement::edge_to_new_vertex |
( |
| ) |
const |
Mapping of old edge (to be removed) to new global vertex number. Useful for forming new topology
◆ is_marked()
bool ParallelRefinement::is_marked |
( |
std::size_t |
edge_index | ) |
const |
Return marked status of edge
- Parameters
-
◆ mark() [1/3]
void ParallelRefinement::mark |
( |
std::size_t |
edge_index | ) |
|
Mark edge by index
- Parameters
-
edge_index | (std::size_t) Index of edge to mark |
◆ mark() [2/3]
void ParallelRefinement::mark |
( |
const MeshFunction< bool > & |
refinement_marker | ) |
|
Mark all edges incident on entities indicated by refinement marker
- Parameters
-
refinement_marker | (const MeshFunction<bool>) |
◆ mark() [3/3]
void ParallelRefinement::mark |
( |
const MeshEntity & |
cell | ) |
|
Mark all incident edges of an entity
- Parameters
-
◆ marked_edge_list()
std::vector< std::size_t > ParallelRefinement::marked_edge_list |
( |
const MeshEntity & |
cell | ) |
const |
Return list of marked edges incident on this MeshEntity - usually a cell
- Parameters
-
◆ new_cell() [1/3]
void ParallelRefinement::new_cell |
( |
const Cell & |
cell | ) |
|
Add a new cell to the list in 3D or 2D
- Parameters
-
◆ new_cell() [2/3]
void 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) |
◆ new_cell() [3/3]
void 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) |
◆ new_cells()
void ParallelRefinement::new_cells |
( |
const std::vector< std::size_t > & |
idx | ) |
|
Add new cells with vertex indices
- Parameters
-
idx | (const std::vector<std::size_t>) |
◆ partition()
void 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) |
The documentation for this class was generated from the following files:
- /home/fenics/shared/dolfin/dolfin/refinement/ParallelRefinement.h
- /home/fenics/shared/dolfin/dolfin/refinement/ParallelRefinement.cpp