DOLFIN
DOLFIN C++ interface
Public Member Functions | Friends | List of all members
dolfin::Mesh Class Reference

#include <Mesh.h>

Inheritance diagram for dolfin::Mesh:
Inheritance graph
[legend]
Collaboration diagram for dolfin::Mesh:
Collaboration graph
[legend]

Public Member Functions

 Mesh ()
 Create empty mesh.
 
 Mesh (MPI_Comm comm)
 Create empty mesh.
 
 Mesh (const Mesh &mesh)
 
 Mesh (std::string filename)
 
 Mesh (MPI_Comm comm, std::string filename)
 
 Mesh (MPI_Comm comm, LocalMeshData &local_mesh_data)
 
 ~Mesh ()
 Destructor.
 
const Meshoperator= (const Mesh &mesh)
 
std::size_t num_vertices () const
 
std::size_t num_edges () const
 
std::size_t num_faces () const
 
std::size_t num_facets () const
 
std::size_t num_cells () const
 
std::size_t num_entities (std::size_t d) const
 
std::vector< double > & coordinates ()
 
const std::vector< double > & coordinates () const
 
const std::vector< unsigned int > & cells () const
 
std::size_t num_entities_global (std::size_t dim) const
 
MeshTopologytopology ()
 
const MeshTopologytopology () const
 
MeshGeometrygeometry ()
 
const MeshGeometrygeometry () const
 
MeshDomainsdomains ()
 
const MeshDomainsdomains () const
 
std::shared_ptr< BoundingBoxTreebounding_box_tree () const
 
MeshDatadata ()
 
const MeshDatadata () const
 
CellTypetype ()
 
const CellTypetype () const
 Get mesh cell type (const version).
 
std::size_t init (std::size_t dim) const
 
void init (std::size_t d0, std::size_t d1) const
 
void init () const
 Compute all entities and connectivity.
 
void init_global (std::size_t dim) const
 Compute global indices for entity dimension dim.
 
void clean ()
 
void order ()
 
bool ordered () const
 
Mesh renumber_by_color () const
 
void scale (double factor)
 
void translate (const Point &point)
 
void rotate (double angle, std::size_t axis=2)
 
void rotate (double angle, std::size_t axis, const Point &point)
 
void smooth (std::size_t num_iterations=1)
 
void smooth_boundary (std::size_t num_iterations=1, bool harmonic_smoothing=true)
 
void snap_boundary (const SubDomain &sub_domain, bool harmonic_smoothing=true)
 
const std::vector< std::size_t > & color (std::string coloring_type) const
 
const std::vector< std::size_t > & color (std::vector< std::size_t > coloring_type) const
 
double hmin () const
 
double hmax () const
 
double rmin () const
 
double rmax () const
 
std::size_t hash () const
 
std::string str (bool verbose) const
 
const std::vector< int > & cell_orientations () const
 
void init_cell_orientations (const Expression &global_normal)
 
MPI_Comm mpi_comm () const
 
std::string ghost_mode () const
 
- Public Member Functions inherited from dolfin::Variable
 Variable ()
 Create unnamed variable.
 
 Variable (const std::string name, const std::string label)
 Create variable with given name and label.
 
 Variable (const Variable &variable)
 Copy constructor.
 
virtual ~Variable ()
 Destructor.
 
const Variableoperator= (const Variable &variable)
 Assignment operator.
 
void rename (const std::string name, const std::string label)
 Rename variable.
 
std::string name () const
 Return name.
 
std::string label () const
 Return label (description)
 
std::size_t id () const
 
- Public Member Functions inherited from dolfin::Hierarchical< Mesh >
 Hierarchical (Mesh &self)
 Constructor.
 
virtual ~Hierarchical ()
 Destructor.
 
std::size_t depth () const
 
bool has_parent () const
 
bool has_child () const
 
Meshparent ()
 
const Meshparent () const
 Return parent in hierarchy (const version).
 
std::shared_ptr< Meshparent_shared_ptr ()
 
std::shared_ptr< const Meshparent_shared_ptr () const
 Return shared pointer to parent (const version).
 
Meshchild ()
 
const Meshchild () const
 Return child in hierarchy (const version).
 
std::shared_ptr< Meshchild_shared_ptr ()
 
std::shared_ptr< const Meshchild_shared_ptr () const
 Return shared pointer to child (const version).
 
Meshroot_node ()
 
const Meshroot_node () const
 Return root node object in hierarchy (const version).
 
std::shared_ptr< Meshroot_node_shared_ptr ()
 
std::shared_ptr< const Meshroot_node_shared_ptr () const
 Return shared pointer to root node object in hierarchy (const version).
 
Meshleaf_node ()
 
const Meshleaf_node () const
 Return leaf node object in hierarchy (const version).
 
std::shared_ptr< Meshleaf_node_shared_ptr ()
 
std::shared_ptr< const Meshleaf_node_shared_ptr () const
 Return shared pointer to leaf node object in hierarchy (const version).
 
void set_parent (std::shared_ptr< Mesh > parent)
 Set parent.
 
void clear_child ()
 Clear child.
 
void set_child (std::shared_ptr< Mesh > child)
 Set child.
 
const Hierarchicaloperator= (const Hierarchical &hierarchical)
 Assignment operator.
 
void _debug () const
 Function useful for debugging the hierarchy.
 

Friends

class MeshEditor
 
class TopologyComputation
 
class MeshPartitioning
 
Mesh create_mesh (Function &coordinates)
 

Additional Inherited Members

- Public Attributes inherited from dolfin::Variable
Parameters parameters
 Parameters.
 

Detailed Description

A Mesh consists of a set of connected and numbered mesh entities.

Both the representation and the interface are dimension-independent, but a concrete interface is also provided for standard named mesh entities:

Entity Dimension Codimension
Vertex 0
Edge 1
Face 2
Facet 1
Cell 0

When working with mesh iterators, all entities and connectivity are precomputed automatically the first time an iterator is created over any given topological dimension or connectivity.

Note that for efficiency, only entities of dimension zero (vertices) and entities of the maximal dimension (cells) exist when creating a Mesh. Other entities must be explicitly created by calling init(). For example, all edges in a mesh may be created by a call to mesh.init(1). Similarly, connectivities such as all edges connected to a given vertex must also be explicitly created (in this case by a call to mesh.init(0, 1)).

Constructor & Destructor Documentation

◆ Mesh() [1/4]

Mesh::Mesh ( const Mesh mesh)

Copy constructor.

Parameters
mesh(Mesh) Object to be copied.

◆ Mesh() [2/4]

Mesh::Mesh ( std::string  filename)
explicit

Create mesh from data file.

Parameters
filename(std::string) Name of file to load.

◆ Mesh() [3/4]

Mesh::Mesh ( MPI_Comm  comm,
std::string  filename 
)

Create mesh from data file.

Parameters
comm(MPI_Comm) The MPI communicator
filename(std::string) Name of file to load.

◆ Mesh() [4/4]

Mesh::Mesh ( MPI_Comm  comm,
LocalMeshData local_mesh_data 
)

Create a distributed mesh from local (per process) data.

Parameters
comm(MPI_Comm) MPI communicator for the mesh.
local_mesh_data(LocalMeshData) Data from which to build the mesh.

Member Function Documentation

◆ bounding_box_tree()

std::shared_ptr< BoundingBoxTree > Mesh::bounding_box_tree ( ) const

Get bounding box tree for mesh. The bounding box tree is initialized and built upon the first call to this function. The bounding box tree can be used to compute collisions between the mesh and other objects. It is the responsibility of the caller to use (and possibly rebuild) the tree. It is stored as a (mutable) member of the mesh to enable sharing of the bounding box tree data structure.

Returns
std::shared_ptr<BoundingBoxTree>

◆ cell_orientations()

const std::vector< int > & Mesh::cell_orientations ( ) const

Return cell_orientations (const version)

Returns
std::vector<int>
    Map from cell index to orientation of cell. Is empty
    if cell orientations have not been computed.  

◆ cells()

const std::vector<unsigned int>& dolfin::Mesh::cells ( ) const
inline

Get cell connectivity.

Returns
std::vector<unsigned int>& Connectivity for all cells.

◆ clean()

void Mesh::clean ( )

Clean out all auxiliary topology data. This clears all topological data, except the connectivity between cells and vertices.

◆ color() [1/2]

const std::vector< std::size_t > & Mesh::color ( std::string  coloring_type) const

Color the cells of the mesh such that no two neighboring cells share the same color. A colored mesh keeps a MeshFunction<std::size_t> named "cell colors" as mesh data which holds the colors of the mesh.

Parameters
coloring_type(std::string) Coloring type, specifying what relation makes two cells neighbors, can be one of "vertex", "edge" or "facet".
Returns
std::vector<std::size_t>& The colors as a mesh function over the cells of the mesh.

◆ color() [2/2]

const std::vector< std::size_t > & Mesh::color ( std::vector< std::size_t >  coloring_type) const

Color the cells of the mesh such that no two neighboring cells share the same color. A colored mesh keeps a MeshFunction<std::size_t> named "cell colors" as mesh data which holds the colors of the mesh.

Parameters
coloring_type(std::vector<std::size_t>&) Coloring type given as list of topological dimensions, specifying what relation makes two mesh entities neighbors.
Returns
std::vector<std::size_t>& The colors as a mesh function over entities of the mesh.

◆ coordinates() [1/2]

std::vector<double>& dolfin::Mesh::coordinates ( )
inline

Get vertex coordinates.

Returns
std::vector<double>& Coordinates of all vertices.

◆ coordinates() [2/2]

const std::vector<double>& dolfin::Mesh::coordinates ( ) const
inline

Return coordinates of all vertices (const version).

Returns
std::vector<double>& Coordinates of all vertices.

◆ data() [1/2]

MeshData & Mesh::data ( )

Get mesh data.

Returns
MeshData& The mesh data object associated with the mesh.

◆ data() [2/2]

const MeshData & Mesh::data ( ) const

Get mesh data (const version).

Returns
MeshData& The mesh data object associated with the mesh.

◆ domains() [1/2]

MeshDomains& dolfin::Mesh::domains ( )
inline

Get mesh (sub)domains.

Returns
MeshDomains The (sub)domains associated with the mesh.

◆ domains() [2/2]

const MeshDomains& dolfin::Mesh::domains ( ) const
inline

Get mesh (sub)domains (const).

Returns
MeshDomains The (sub)domains associated with the mesh.

◆ geometry() [1/2]

MeshGeometry& dolfin::Mesh::geometry ( )
inline

Get mesh geometry.

Returns
MeshGeometry The geometry object associated with the mesh.

◆ geometry() [2/2]

const MeshGeometry& dolfin::Mesh::geometry ( ) const
inline

Get mesh geometry (const version).

Returns
MeshGeometry The geometry object associated with the mesh.

◆ ghost_mode()

std::string Mesh::ghost_mode ( ) const

Ghost mode used for partitioning. Possible values are same as parameters["ghost_mode"]. WARNING: the interface may change in future without deprecation; the method is now intended for internal library use.

◆ hash()

std::size_t Mesh::hash ( ) const

Compute hash of mesh, currently based on the has of the mesh geometry and mesh topology.

Returns
std::size_t A tree-hashed value of the coordinates over all MPI processes

◆ hmax()

double Mesh::hmax ( ) const

Compute maximum cell size in mesh, measured greatest distance between any two vertices of a cell.

Returns
double The maximum cell size. The size is computed using Cell::h()

◆ hmin()

double Mesh::hmin ( ) const

Compute minimum cell size in mesh, measured greatest distance between any two vertices of a cell.

Returns
double The minimum cell size. The size is computed using Cell::h()

◆ init() [1/2]

std::size_t Mesh::init ( std::size_t  dim) const

Compute entities of given topological dimension.

Parameters
dim(std::size_t) Topological dimension.
Returns
std::size_t Number of created entities.

◆ init() [2/2]

void Mesh::init ( std::size_t  d0,
std::size_t  d1 
) const

Compute connectivity between given pair of dimensions.

Parameters
d0(std::size_t) Topological dimension.
d1(std::size_t) Topological dimension.

◆ init_cell_orientations()

void Mesh::init_cell_orientations ( const Expression global_normal)

Compute and initialize cell_orientations relative to a given global outward direction/normal/orientation. Only defined if mesh is orientable.

Parameters
global_normal(Expression) A global normal direction to the mesh

◆ mpi_comm()

MPI_Comm dolfin::Mesh::mpi_comm ( ) const
inline

Mesh MPI communicator

Returns
MPI_Comm

◆ num_cells()

std::size_t dolfin::Mesh::num_cells ( ) const
inline

Get number of cells in mesh.

Returns
std::size_t Number of cells.

◆ num_edges()

std::size_t dolfin::Mesh::num_edges ( ) const
inline

Get number of edges in mesh.

Returns
std::size_t Number of edges.

◆ num_entities()

std::size_t dolfin::Mesh::num_entities ( std::size_t  d) const
inline

Get number of entities of given topological dimension.

Parameters
d(std::size_t) Topological dimension.
Returns
std::size_t Number of entities of topological dimension d.

◆ num_entities_global()

std::size_t dolfin::Mesh::num_entities_global ( std::size_t  dim) const
inline

Get global number of entities of given topological dimension.

Parameters
dim(std::size_t) Topological dimension.
Returns
std::size_t Global number of entities of topological dimension d.

◆ num_faces()

std::size_t dolfin::Mesh::num_faces ( ) const
inline

Get number of faces in mesh.

Returns
std::size_t Number of faces.

◆ num_facets()

std::size_t dolfin::Mesh::num_facets ( ) const
inline

Get number of facets in mesh.

Returns
std::size_t Number of facets.

◆ num_vertices()

std::size_t dolfin::Mesh::num_vertices ( ) const
inline

Get number of vertices in mesh.

Returns
std::size_t Number of vertices.

◆ operator=()

const Mesh & Mesh::operator= ( const Mesh mesh)

Assignment operator

Parameters
mesh(Mesh) Another Mesh object.

◆ order()

void Mesh::order ( )

Order all mesh entities.

See also: UFC documentation (put link here!)

◆ ordered()

bool Mesh::ordered ( ) const

Check if mesh is ordered according to the UFC numbering convention.

Returns
bool The return values is true iff the mesh is ordered.

◆ renumber_by_color()

dolfin::Mesh Mesh::renumber_by_color ( ) const

Renumber mesh entities by coloring. This function is currently restricted to renumbering by cell coloring. The cells (cell-vertex connectivity) and the coordinates of the mesh are renumbered to improve the locality within each color. It is assumed that the mesh has already been colored and that only cell-vertex connectivity exists as part of the mesh.

Returns
Mesh

◆ rmax()

double Mesh::rmax ( ) const

Compute maximum cell inradius.

Returns
double The maximum of cells' inscribed sphere radii

◆ rmin()

double Mesh::rmin ( ) const

Compute minimum cell inradius.

Returns
double The minimum of cells' inscribed sphere radii

◆ rotate() [1/2]

void Mesh::rotate ( double  angle,
std::size_t  axis = 2 
)

Rotate mesh around a coordinate axis through center of mass of all mesh vertices

Parameters
angle(double) The number of degrees (0-360) of rotation.
axis(std::size_t) The coordinate axis around which to rotate the mesh.

◆ rotate() [2/2]

void Mesh::rotate ( double  angle,
std::size_t  axis,
const Point point 
)

Rotate mesh around a coordinate axis through a given point

Parameters
angle(double) The number of degrees (0-360) of rotation.
axis(std::size_t) The coordinate axis around which to rotate the mesh.
point(Point) The point around which to rotate the mesh.

◆ scale()

void Mesh::scale ( double  factor)

Scale mesh coordinates with given factor.

Arguments factor (double) The factor defining the scaling.

◆ smooth()

void Mesh::smooth ( std::size_t  num_iterations = 1)

Smooth internal vertices of mesh by local averaging.

Parameters
num_iterations(std::size_t) Number of iterations to perform smoothing, default value is 1.

◆ smooth_boundary()

void Mesh::smooth_boundary ( std::size_t  num_iterations = 1,
bool  harmonic_smoothing = true 
)

Smooth boundary vertices of mesh by local averaging.

Parameters
num_iterations(std::size_t) Number of iterations to perform smoothing, default value is 1.
harmonic_smoothing(bool) Flag to turn on harmonics smoothing, default value is true.

◆ snap_boundary()

void Mesh::snap_boundary ( const SubDomain sub_domain,
bool  harmonic_smoothing = true 
)

Snap boundary vertices of mesh to match given sub domain.

Parameters
sub_domain(SubDomain) A SubDomain object.
harmonic_smoothing(bool) Flag to turn on harmonics smoothing, default value is true.

◆ str()

std::string Mesh::str ( bool  verbose) const
virtual

Informal string representation.

Parameters
verbose(bool) Flag to turn on additional output.
Returns
std::string An informal representation of the mesh.

Reimplemented from dolfin::Variable.

◆ topology() [1/2]

MeshTopology& dolfin::Mesh::topology ( )
inline

Get mesh topology.

Returns
MeshTopology The topology object associated with the mesh.

◆ topology() [2/2]

const MeshTopology& dolfin::Mesh::topology ( ) const
inline

Get mesh topology (const version).

Returns
MeshTopology The topology object associated with the mesh.

◆ translate()

void Mesh::translate ( const Point point)

Translate mesh according to a given vector.

Parameters
point(Point) The vector defining the translation.

◆ type()

CellType& dolfin::Mesh::type ( )
inline

Get mesh cell type.

Returns
CellType& The cell type object associated with the mesh.

Friends And Related Function Documentation

◆ create_mesh

Mesh create_mesh ( Function coordinates)
friend

Creates mesh from coordinate function

Topology is given by underlying mesh of the function space and geometry is given by function values. Hence resulting mesh geometry has a degree of the function space degree. Geometry of function mesh is ignored.

Mesh connectivities d-0, d-1, ..., d-r are built on function mesh (where d is topological dimension of the mesh and r is maximal dimension of entity associated with any coordinate node). Consider clearing unneeded connectivities when finished.

Parameters
coordinates(Function) Vector Lagrange function of any degree
Returns
Mesh The mesh

The documentation for this class was generated from the following files: