CGAL 5.6 - 3D Simplicial Mesh Data Structures
|
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
A data structure to represent and maintain a 3D complex embedded in a 3D triangulation.
The class Mesh_complex_3_in_triangulation_3
implements a data structure to store the 3D restricted Delaunay triangulation used by a mesh generation process.
This class is a model of the concept MeshComplexWithFeatures_3InTriangulation_3
.
Tr | can be instantiated with any 3D triangulation of CGAL provided that its vertex and cell base class are models of the concepts SimplicialMeshVertexBase_3 and SimplicialMeshCellBase_3 , respectively. |
CornerIndex | Type of indices for corners (i.e. \( 0\)–dimensional features) of the discretized geometric domain. It must be a model of CopyConstructible , Assignable , DefaultConstructible and LessThanComparable . It must match the Corner_index of the model of the MeshDomainWithFeatures_3 concept when used for mesh generation. |
CurveIndex | Type of indices for curves (i.e. \( 1\)-dimensional features) of the discretized geometric domain. It must be a model of CopyConstructible , Assignable , DefaultConstructible and LessThanComparable . The default constructed value must be the value for an edge which does not approximate a 1-dimensional feature of the geometric domain. It must match the Curve_index types of the model of the MeshDomainWithFeatures_3 concept when used for mesh generation. |
Those two last template parameters default to int
, so that they can be ignored if the domain used for mesh generation does not include 0 and 1-dimensionnal features (i.e is only a model of the concept MeshDomain_3
).
Public Member Functions | |
void | remove_isolated_vertices () |
The tetrahedral mesh generation algorithm implemented in CGAL::make_mesh_3() and CGAL::refine_mesh_3() does not guarantee that all the points inserted by the algorithm are actually present in the final mesh. | |
template<typename OutputIterator > | |
OutputIterator | adjacent_vertices_in_complex (const Vertex_handle &v, OutputIterator out) const |
fills out with incident edges (1-dimensional features of v ). | |
bool | has_incident_facets_in_complex (const Vertex_handle &v) const |
Returns true if the vertex v has is incident to at least a facet of the complex. | |
Types | |
typedef Tr | Triangulation |
typedef Tr::size_type | size_type |
typedef Tr::Point | Point |
typedef Tr::Edge | Edge |
typedef Tr::Facet | Facet |
typedef Tr::Vertex_handle | Vertex_handle |
typedef Tr::Cell_handle | Cell_handle |
typedef Tr::Vertex::Index | Index |
Index type. | |
typedef Tr::Cell::Surface_patch_index | Surface_patch_index |
Surface index type. | |
typedef Tr::Cell::Subdomain_index | Subdomain_index |
Subdomain index type. | |
typedef CornerIndex | Corner_index |
Corner index type. | |
typedef CurveIndex | Curve_index |
Curve index type. | |
Traversal of the complex | |
typedef unspecified_type | Cells_in_complex_iterator |
Iterator type to visit the cells of the 3D complex. | |
typedef unspecified_type | Facets_in_complex_iterator |
Iterator type to visit the facets of the 2D complex. | |
typedef unspecified_type | Edges_in_complex_iterator |
Iterator type to visit the edges of the 1D complex. | |
typedef unspecified_type | Vertices_in_complex_iterator |
Iterator type to visit the vertices of the 0D complex. | |
typedef Iterator_range< unspecified_type > | Cells_in_complex |
Range type for iterating over all cells of the 3D complex, with a nested type iterator that has as value type Cell_handle . | |
typedef Iterator_range< unspecified_type > | Facets_in_complex |
Range type for iterating over all facets of the 2D complex, with a nested type iterator that has as value type Facet . | |
typedef Iterator_range< unspecified_type > | Edges_in_complex |
Range type for iterating over all cells of the 1D complex, with a nested type iterator that has as value type Edge . | |
typedef Iterator_range< unspecified_type > | Vertices_in_complex |
Range type for iterating over all vertices of the 0D complex, with a nested type iterator that has as value type Vertex_handle . | |
Creation | |
Mesh_complex_3_in_triangulation_3 () | |
Constructor builds an empty 3D complex. | |
Mesh_complex_3_in_triangulation_3 (const Self &rhs) | |
Copy constructor. | |
Mesh_complex_3_in_triangulation_3 (Self &&rhs) | |
Move constructor. | |
Self & | operator= (Self rhs) |
Assignment operator, also serves as move-assignment. | |
void | swap (Self &rhs) |
swaps this and rhs | |
Access Functions | |
const Triangulation & | triangulation () const |
returns a const reference to the triangulation | |
Non const access | |
Triangulation & | triangulation () |
returns a reference to the triangulation | |
Modifiers | |
void | clear () |
clears data of the complex | |
void | add_to_complex (const Cell_handle &cell, const Subdomain_index &index) |
adds cell cell to the 3D complex, with subdomain index index | |
void | add_to_complex (const Facet &facet, const Surface_patch_index &index) |
adds facet facet to the 2D complex, with surface index index | |
void | add_to_complex (const Cell_handle &cell, const int i, const Surface_patch_index &index) |
adds facet(cell , i ) to the 2D complex, with surface index index | |
void | add_to_complex (const Edge &e, const Curve_index &index) |
adds edge e to complex, with curve index index | |
void | add_to_complex (const Vertex_handle &v1, const Vertex_handle &v2, const Curve_index &index) |
adds edge (v1 , v2 ) to the 1D complex, with Curve_index index | |
void | add_to_complex (const Vertex_handle &v, const Corner_index &index) |
marks vertex v as a corner of the complex | |
void | remove_from_complex (const Cell_handle &cell) |
removes cell cell from the 3D complex | |
void | remove_from_complex (const Facet &facet) |
removes facet facet from the 2D complex | |
void | remove_from_complex (const Cell_handle &c, const int i) |
removes facet(cell , i ) from the 2D complex | |
void | remove_from_complex (const Edge &e) |
removes edge e from the 1D complex | |
void | remove_from_complex (const Vertex_handle &v1, const Vertex_handle &v2) |
removes edge (v1,v2) from the 1D complex | |
void | remove_from_complex (const Vertex_handle &v) |
removes vertex v from the complex | |
void | set_index (const Vertex_handle &vertex, const Index &index) const |
sets the index of vertex vertex to index | |
void | set_surface_patch_index (const Facet &f, const Surface_patch_index &index) |
sets the surface index of facet facet to index | |
void | set_surface_patch_index (const Cell_handle &cell, const int i, const Surface_patch_index &index) const |
sets the surface index of facet(cell , i ) to index | |
void | set_subdomain_index (const Cell_handle &cell, const Subdomain_index &index) const |
sets the subdomain index of cell cell to index | |
void | set_dimension (const Vertex_handle &vertex, int dimension) const |
sets the dimension of vertex vertex to dimension | |
Queries on the identifier of the face complex including triangulation cells, facets, and vertices. | |
Index | index (const Vertex_handle &v) const |
returns the index of vertex v | |
Subdomain_index | subdomain_index (const Cell_handle &cell) const |
returns the subdomain index of cell cell | |
Surface_patch_index | surface_patch_index (const Facet &f) const |
returns the surface index of facet f | |
Surface_patch_index | surface_patch_index (const Cell_handle &cell, const int i) const |
returns the surface index of facet(cell , i ) | |
int | in_dimension (const Vertex_handle &v) const |
returns the dimension of the lowest dimensional face of the input 3D complex that contains the vertex | |
Curve_index | curve_index (const Edge &e) const |
returns the curve index of edge e | |
Curve_index | curve_index (const Vertex_handle &v1, const Vertex_handle &v2) const |
returns the curve index of the edge formed by v1 and v2 | |
Corner_index | corner_index (const Vertex_handle &v) const |
returns the corner index of vertex v | |
I/O Functions | |
std::ostream & | output_boundary_to_off (std::ostream &out) const |
outputs the outer boundary of the entire domain, with facets oriented outward. | |
std::ostream & | output_boundary_to_off (std::ostream &out, Subdomain_index subdomain) const |
outputs the outer boundary of the selected subdomain, with facets oriented outward. | |
std::ostream & | output_facets_in_complex_to_off (std::ostream &out) const |
outputs the surface facets, with a consistent orientation at the interface of two subdomains. | |
void | output_to_medit (std::ostream &os) const |
outputs the mesh to os in Medit format. | |
void | output_to_maya (std::ostream &os, bool surfaceOnly=true) const |
Iterators | |
Cells_in_complex_iterator | cells_in_complex_begin () const |
returns a Cells_in_complex_iterator to the first cell of the 3D complex | |
Cells_in_complex_iterator | cells_in_complex_begin (const Subdomain_index &index) const |
returns a Cells_in_complex_iterator to the first cell of the 3D complex | |
Cells_in_complex_iterator | cells_in_complex_end () const |
returns the past-the-end iterator for the cells of the 3D complex | |
Facets_in_complex_iterator | facets_in_complex_begin () const |
returns a Facets_in_complex_iterator to the first facet of the 2D complex | |
Facets_in_complex_iterator | facets_in_complex_begin (const Surface_patch_index &index) const |
returns a Facets_in_complex_iterator to the first facet of the 2D complex | |
Facets_in_complex_iterator | facets_in_complex_end (const Surface_patch_index=Surface_patch_index()) const |
returns past-the-end iterator on facet of the 2D complex | |
Edges_in_complex_iterator | edges_in_complex_begin () const |
returns a Edges_in_complex_iterator to the first edge of the 1D complex | |
Edges_in_complex_iterator | edges_in_complex_begin (const Curve_index &index) const |
returns a Edges_in_complex_iterator to the first edge of the 1D complex | |
Edges_in_complex_iterator | edges_in_complex_end (const Curve_index &=Curve_index()) const |
returns past-the-end iterator on edges of the 1D complex | |
Vertices_in_complex_iterator | vertices_in_complex_begin () const |
returns a Vertices_in_complex_iterator to the first vertex of the 0D complex | |
Vertices_in_complex_iterator | vertices_in_complex_begin (const Corner_index &index) const |
returns a Vertices_in_complex_iterator to the first vertex of the 0D complex | |
Vertices_in_complex_iterator | vertices_in_complex_end () const |
returns past-the-end iterator on vertices of the 0D complex | |
Vertices_in_complex | vertices_in_complex () const |
returns a range of iterators over vertices of the 0D complex | |
Edges_in_complex | edges_in_complex () const |
returns a range of iterators over the edges of the 1D complex, starting at an arbitrary edge. | |
Facets_in_complex | facets_in_complex () const |
returns a range of iterators over the facets of the 2D complex, starting at an arbitrary facet. | |
Cells_in_complex | cells_in_complex () const |
returns a range of iterators over cells of the 3D complex. | |
OutputIterator CGAL::Mesh_complex_3_in_triangulation_3< Tr, CI_, CSI_ >::adjacent_vertices_in_complex | ( | const Vertex_handle & | v, |
OutputIterator | out | ||
) | const |
fills out
with incident edges (1-dimensional features of v
).
OutputIterator value type is std::pair<Vertex_handle,Curve_index>
Cells_in_complex CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex >::cells_in_complex | ( | ) | const |
returns a range of iterators over cells of the 3D complex.
Returns an empty range when triangulation().number_of_cells() == 0
or complex is empty.
Cells_in_complex::iterator
is Cell_handle
. Edges_in_complex CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex >::edges_in_complex | ( | ) | const |
returns a range of iterators over the edges of the 1D complex, starting at an arbitrary edge.
Returns an empty range when t.dimension() < 2
.
Facets_in_complex CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex >::facets_in_complex | ( | ) | const |
returns a range of iterators over the facets of the 2D complex, starting at an arbitrary facet.
Returns an empty range when t.dimension() < 2
.
void CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex >::remove_isolated_vertices | ( | ) |
The tetrahedral mesh generation algorithm implemented in CGAL::make_mesh_3()
and CGAL::refine_mesh_3()
does not guarantee that all the points inserted by the algorithm are actually present in the final mesh.
In most cases, all points are used, but if the geometry of the object has small features compared to the size of the simplices (triangles and tetrahedra), it might be that the Delaunay facets that are selected in the restricted Delaunay triangulation miss some vertices of the triangulation. The concurrent version of the tetrahedral mesh generation algorithm also inserts a small set of auxiliary vertices that belong to the triangulation but are isolated from the complex at the end of the meshing process.
This function removes these so-called isolated vertices, that belong to the triangulation but not to any simplex of the C3T3
, from the triangulation.
Vertices_in_complex CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex >::vertices_in_complex | ( | ) | const |
returns a range of iterators over vertices of the 0D complex
Vertices_in_complex::iterator
is Vertex_handle
.