Reference documentation for deal.II version 9.6.0
|
Namespaces | |
namespace | internal |
Classes | |
struct | AdditionalData |
struct | AdditionalData< 2 > |
struct | FaceInfo2 |
Typedefs | |
using | ConcurrencyTag = CGAL::Sequential_tag |
using | K = CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt |
using | K_exact = CGAL::Exact_predicates_exact_constructions_kernel |
using | CGALPolygon = CGAL::Polygon_2<K> |
using | Polygon_with_holes_2 = CGAL::Polygon_with_holes_2<K> |
using | CGALTriangle2 = K::Triangle_2 |
using | CGALTriangle3 = K::Triangle_3 |
using | CGALTriangle3_exact = K_exact::Triangle_3 |
using | CGALPoint2 = K::Point_2 |
using | CGALPoint3 = K::Point_3 |
using | CGALPoint3_exact = K_exact::Point_3 |
using | CGALSegment2 = K::Segment_2 |
using | Surface_mesh = CGAL::Surface_mesh<K_exact::Point_3> |
using | CGALSegment3 = K::Segment_3 |
using | CGALSegment3_exact = K_exact::Segment_3 |
using | CGALTetra = K::Tetrahedron_3 |
using | CGALTetra_exact = K_exact::Tetrahedron_3 |
using | Triangulation2 = CGAL::Triangulation_2<K> |
using | Triangulation3 = CGAL::Triangulation_3<K> |
using | Triangulation3_exact = CGAL::Triangulation_3<K_exact> |
using | Vb = CGAL::Triangulation_vertex_base_2<K> |
using | Fbb = CGAL::Triangulation_face_base_with_info_2<FaceInfo2, K> |
using | CFb = CGAL::Constrained_triangulation_face_base_2<K, Fbb> |
using | Fb = CGAL::Delaunay_mesh_face_base_2<K, CFb> |
using | Tds = CGAL::Triangulation_data_structure_2<Vb, Fb> |
using | Itag = CGAL::Exact_predicates_tag |
using | CDT = CGAL::Constrained_Delaunay_triangulation_2<K, Tds, Itag> |
using | Criteria = CGAL::Delaunay_mesh_size_criteria_2<CDT> |
using | Vertex_handle = CDT::Vertex_handle |
using | Face_handle = CDT::Face_handle |
Enumerations | |
enum class | FacetTopology { facet_vertices_on_surface = CGAL::FACET_VERTICES_ON_SURFACE , facet_vertices_on_same_surface_patch , facet_vertices_on_same_surface_patch_with_adjacency_check } |
enum class | BooleanOperation { compute_corefinement = 1 << 0 , compute_difference = 1 << 1 , compute_intersection = 1 << 2 , compute_union = 1 << 3 } |
Functions | |
template<int structdim0, int structdim1, int spacedim> | |
std::vector< std::array< Point< spacedim >, structdim1+1 > > | compute_intersection_of_cells (const typename Triangulation< structdim0, spacedim >::cell_iterator &cell0, const typename Triangulation< structdim1, spacedim >::cell_iterator &cell1, const Mapping< structdim0, spacedim > &mapping0, const Mapping< structdim1, spacedim > &mapping1, const double tol=1e-9) |
template<int structdim0, int structdim1, int spacedim> | |
std::vector< std::array< Point< spacedim >, structdim1+1 > > | compute_intersection_of_cells (const ArrayView< const Point< spacedim > > &vertices0, const ArrayView< const Point< spacedim > > &vertices1, const double tol=1e-9) |
template<typename CGALPointType , int dim, int spacedim> | |
void | dealii_cell_to_cgal_surface_mesh (const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const ::Mapping< dim, spacedim > &mapping, CGAL::Surface_mesh< CGALPointType > &mesh) |
template<typename CGALPointType , int dim, int spacedim> | |
void | dealii_tria_to_cgal_surface_mesh (const ::Triangulation< dim, spacedim > &triangulation, CGAL::Surface_mesh< CGALPointType > &mesh) |
template<typename C3t3 > | |
void | cgal_surface_mesh_to_cgal_triangulation (CGAL::Surface_mesh< typename C3t3::Point::Point > &surface_mesh, C3t3 &triangulation, const AdditionalData< 3 > &data=AdditionalData< 3 >{}) |
template<typename CGALPointType > | |
void | compute_boolean_operation (const CGAL::Surface_mesh< CGALPointType > &surface_mesh_1, const CGAL::Surface_mesh< CGALPointType > &surface_mesh_2, const BooleanOperation &boolean_operation, CGAL::Surface_mesh< CGALPointType > &output_surface_mesh) |
template<typename CGALTriangulationType > | |
::Quadrature< CGALTriangulationType::Point::Ambient_dimension::value > | compute_quadrature (const CGALTriangulationType &tria, const unsigned int degree) |
template<int dim0, int dim1, int spacedim> | |
::Quadrature< spacedim > | compute_quadrature_on_boolean_operation (const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const BooleanOperation &bool_op, const Mapping< dim0, spacedim > &mapping0=(ReferenceCells::get_hypercube< dim0 >() .template get_default_linear_mapping< dim0, spacedim >()), const Mapping< dim1, spacedim > &mapping1=(ReferenceCells::get_hypercube< dim1 >() .template get_default_linear_mapping< dim1, spacedim >())) |
template<int dim0, int dim1, int spacedim> | |
::Quadrature< spacedim > | compute_quadrature_on_intersection (const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const Mapping< dim0, spacedim > &mapping0, const Mapping< dim1, spacedim > &mapping1) |
template<typename CGALPointType > | |
void | remesh_surface (CGAL::Surface_mesh< CGALPointType > &surface_mesh, const AdditionalData< 3 > &data=AdditionalData< 3 >{}) |
template<class T , class... Types> | |
const T * | get_if_ (const std::variant< Types... > *v) |
template<class T , class... Types> | |
const T * | get_if_ (const boost::variant< Types... > *v) |
Interface to the Computational Geometry Algorithm Library (CGAL).
CGAL is a software project that provides easy access to efficient and reliable geometric algorithms. The library offers data structures and algorithms like triangulations, Voronoi diagrams, Boolean operations on polygons and polyhedra, point set processing, arrangements of curves, surface and volume mesh generation, geometry processing, alpha shapes, convex hull algorithms, shape reconstruction, AABB and KD trees...
You can learn more about the CGAL library at https://www.cgal.org/
using CGALWrappers::ConcurrencyTag = CGAL::Sequential_tag |
Definition at line 82 of file utilities.h.
using CGALWrappers::K = CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt |
Definition at line 61 of file intersections.cc.
using CGALWrappers::K_exact = CGAL::Exact_predicates_exact_constructions_kernel |
Definition at line 62 of file intersections.cc.
using CGALWrappers::CGALPolygon = CGAL::Polygon_2<K> |
Definition at line 63 of file intersections.cc.
using CGALWrappers::Polygon_with_holes_2 = CGAL::Polygon_with_holes_2<K> |
Definition at line 64 of file intersections.cc.
using CGALWrappers::CGALTriangle2 = K::Triangle_2 |
Definition at line 65 of file intersections.cc.
using CGALWrappers::CGALTriangle3 = K::Triangle_3 |
Definition at line 66 of file intersections.cc.
using CGALWrappers::CGALTriangle3_exact = K_exact::Triangle_3 |
Definition at line 67 of file intersections.cc.
using CGALWrappers::CGALPoint2 = K::Point_2 |
Definition at line 68 of file intersections.cc.
using CGALWrappers::CGALPoint3 = K::Point_3 |
Definition at line 69 of file intersections.cc.
using CGALWrappers::CGALPoint3_exact = K_exact::Point_3 |
Definition at line 70 of file intersections.cc.
using CGALWrappers::CGALSegment2 = K::Segment_2 |
Definition at line 71 of file intersections.cc.
using CGALWrappers::Surface_mesh = CGAL::Surface_mesh<K_exact::Point_3> |
Definition at line 72 of file intersections.cc.
using CGALWrappers::CGALSegment3 = K::Segment_3 |
Definition at line 73 of file intersections.cc.
using CGALWrappers::CGALSegment3_exact = K_exact::Segment_3 |
Definition at line 74 of file intersections.cc.
using CGALWrappers::CGALTetra = K::Tetrahedron_3 |
Definition at line 75 of file intersections.cc.
using CGALWrappers::CGALTetra_exact = K_exact::Tetrahedron_3 |
Definition at line 76 of file intersections.cc.
using CGALWrappers::Triangulation2 = CGAL::Triangulation_2<K> |
Definition at line 77 of file intersections.cc.
using CGALWrappers::Triangulation3 = CGAL::Triangulation_3<K> |
Definition at line 78 of file intersections.cc.
using CGALWrappers::Triangulation3_exact = CGAL::Triangulation_3<K_exact> |
Definition at line 79 of file intersections.cc.
using CGALWrappers::Vb = CGAL::Triangulation_vertex_base_2<K> |
Definition at line 92 of file intersections.cc.
using CGALWrappers::Fbb = CGAL::Triangulation_face_base_with_info_2<FaceInfo2, K> |
Definition at line 93 of file intersections.cc.
using CGALWrappers::CFb = CGAL::Constrained_triangulation_face_base_2<K, Fbb> |
Definition at line 94 of file intersections.cc.
using CGALWrappers::Fb = CGAL::Delaunay_mesh_face_base_2<K, CFb> |
Definition at line 95 of file intersections.cc.
using CGALWrappers::Tds = CGAL::Triangulation_data_structure_2<Vb, Fb> |
Definition at line 96 of file intersections.cc.
using CGALWrappers::Itag = CGAL::Exact_predicates_tag |
Definition at line 97 of file intersections.cc.
using CGALWrappers::CDT = CGAL::Constrained_Delaunay_triangulation_2<K, Tds, Itag> |
Definition at line 98 of file intersections.cc.
using CGALWrappers::Criteria = CGAL::Delaunay_mesh_size_criteria_2<CDT> |
Definition at line 99 of file intersections.cc.
using CGALWrappers::Vertex_handle = CDT::Vertex_handle |
Definition at line 100 of file intersections.cc.
using CGALWrappers::Face_handle = CDT::Face_handle |
Definition at line 101 of file intersections.cc.
|
strong |
Definition at line 31 of file additional_data.h.
|
strong |
An enum type given to functions that compute boolean operations between geometrical objects, defined by triangulated surface grids.
As an example, we show all supported boolean operations applied to a hedra (the union of two pyramids) and a cube.
Definition at line 94 of file utilities.h.
std::vector< std::array< Point< spacedim >, structdim1+1 > > CGALWrappers::compute_intersection_of_cells | ( | const typename Triangulation< structdim0, spacedim >::cell_iterator & | cell0, |
const typename Triangulation< structdim1, spacedim >::cell_iterator & | cell1, | ||
const Mapping< structdim0, spacedim > & | mapping0, | ||
const Mapping< structdim1, spacedim > & | mapping1, | ||
const double | tol = 1e-9 ) |
Given two deal.II affine cells, compute their intersection and return a vector of simplices, each one identified by an array of deal.II Points. All the simplices together are a subdivision of the intersection. If cells are non-affine, a geometrical error is introduced. If the measure of one of the simplices is below a certain threshold which defaults to 1e-9, then it is discarded. In case the two cells are disjoint, an empty array is returned.
cell0 | Iterator to the first cell. |
cell1 | Iterator to the second cell. |
mapping0 | Mapping for the first cell. |
mapping1 | Mapping for the second cell. |
tol | Threshold to decide whether or not a simplex is included. |
Definition at line 842 of file intersections.cc.
std::vector< std::array< Point< spacedim >, structdim1+1 > > CGALWrappers::compute_intersection_of_cells | ( | const ArrayView< const Point< spacedim > > & | vertices0, |
const ArrayView< const Point< spacedim > > & | vertices1, | ||
const double | tol = 1e-9 ) |
Same function as above, but working directly with vertices. It's assumed that the first vertices are coming from a Triangulation<dim0,spacedim>, while the other vertices are from a Triangulation<dim1,spacedim>, with dim0
> dim1
.
Definition at line 772 of file intersections.cc.
void CGALWrappers::dealii_cell_to_cgal_surface_mesh | ( | const typename ::Triangulation< dim, spacedim >::cell_iterator & | cell, |
const ::Mapping< dim, spacedim > & | mapping, | ||
CGAL::Surface_mesh< CGALPointType > & | mesh ) |
Build a CGAL::Surface_mesh from a deal.II cell.
The class Surface_mesh implements a halfedge data structure and can be used to represent polyhedral surfaces. It is an edge-centered data structure capable of maintaining incidence information of vertices, edges, and faces. Each edge is represented by two halfedges with opposite orientation. The orientation of a face is chosen so that the halfedges around a face are oriented counterclockwise.
More information on this class is available at https://doc.cgal.org/latest/Surface_mesh/index.html
The function will throw an exception in dimension one. In dimension two, it generates a surface mesh of the quadrilateral cell or of the triangle cell, while in dimension three it will generate the surface mesh of the cell, i.e., a polyhedral mesh containing the faces of the input cell.
The generated mesh is useful when performing geometric operations using CGAL::Polygon_mesh_processing, i.e., to compute boolean operations on cells, splitting, cutting, slicing, etc.
For examples on how to use the resulting CGAL::Surface_mesh see https://doc.cgal.org/latest/Polygon_mesh_processing/
[in] | cell | The input deal.II cell iterator |
[in] | mapping | The mapping used to map the vertices of the cell |
[out] | mesh | The output CGAL::Surface_mesh |
void CGALWrappers::dealii_tria_to_cgal_surface_mesh | ( | const ::Triangulation< dim, spacedim > & | triangulation, |
CGAL::Surface_mesh< CGALPointType > & | mesh ) |
Convert a deal.II triangulation to a CGAL::Surface_mesh. The output depends on the intrinsic dimension of the input deal.II triangulation.
In 2d, i.e. with a Triangulation<2> or a Triangulation<2,3>, the output is the CGAL::Surface_mesh describing the whole triangulation.
In 3d, the boundary the of the deal.II Triangulation is converted to a CGAL::Surface_mesh by looping over all the boundary faces.
[in] | triangulation | The input deal.II triangulation. |
[out] | mesh | The output CGAL::Surface_mesh. |
void CGALWrappers::cgal_surface_mesh_to_cgal_triangulation | ( | CGAL::Surface_mesh< typename C3t3::Point::Point > & | surface_mesh, |
C3t3 & | triangulation, | ||
const AdditionalData< 3 > & | data = AdditionalData< 3 >{} ) |
Given a closed CGAL::Surface_mesh, this function fills the region bounded by the surface with tets.
[in] | surface_mesh | The (closed) surface mesh bounding the volume that has to be filled. |
[out] | triangulation | The output triangulation filled with tetrahedra. |
[in] | data | AdditionalData object to pass to the CGAL::make_mesh_3 function. See the documentation of that struct for a description of those parameters. |
void CGALWrappers::compute_boolean_operation | ( | const CGAL::Surface_mesh< CGALPointType > & | surface_mesh_1, |
const CGAL::Surface_mesh< CGALPointType > & | surface_mesh_2, | ||
const BooleanOperation & | boolean_operation, | ||
CGAL::Surface_mesh< CGALPointType > & | output_surface_mesh ) |
Given two triangulated surface meshes that bound two volumes, execute a boolean operation on them, and store the result in a third surface mesh.
Quoting from CGAL documentation (https://doc.cgal.org/latest/Polygon_mesh_processing/index.html#title14):
Given a closed triangulated surface mesh, each connected component splits the 3d space into two subspaces. The vertex sequence of each face of a component is seen either clockwise or counterclockwise from these two subspaces. The subspace that sees the sequence clockwise (resp. counterclockwise) is on the negative (resp. positive) side of the component.
Given a closed triangulated surface mesh
surface
with no self-intersections, we say thatsurface
bounds a volume if each subspace lies exclusively on the positive (or negative) side of all the incident connected components ofsurface
. The volume bounded bysurface
is the union of all subspaces that are on negative sides of their incident connected components ofsurface
.
There is no restriction on the topology of the input volumes. However, there are some requirements on the input to guarantee that the operation is possible. First, the input meshes must not self-intersect. Second, the operation is possible only if the output can be bounded by a manifold triangulated surface mesh. In particular this means that the output volume has no part with zero thickness. Mathematically speaking, the intersection with an infinitesimally small ball centered in the output volume is a topological ball. At the surface level this means that no non-manifold vertex or edge is allowed in the output. For example, it is not possible to compute the union of two cubes that are disjoint but sharing an edge.
See BooleanOperation for a list of available operations, with the corresponding examples.
[in] | surface_mesh_1 | The first surface mesh. |
[in] | surface_mesh_2 | The second surface mesh. |
[in] | boolean_operation | See BooleanOperation for the list of the allowed operations. |
[out] | output_surface_mesh | The surface mesh with the result of the boolean operation. |
::Quadrature< CGALTriangulationType::Point::Ambient_dimension::value > CGALWrappers::compute_quadrature | ( | const CGALTriangulationType & | tria, |
const unsigned int | degree ) |
Given a CGAL Triangulation describing a polyhedral region, create a Quadrature rule to integrate over the polygon by looping through all the vertices and exploiting QGaussSimplex.
[in] | tria | The CGAL triangulation object describing the polyhedral region. |
[in] | degree | Desired degree of the Quadrature rule on each element of the polyhedral. |
::Quadrature< spacedim > CGALWrappers::compute_quadrature_on_boolean_operation | ( | const typename ::Triangulation< dim0, spacedim >::cell_iterator & | cell0, |
const typename ::Triangulation< dim1, spacedim >::cell_iterator & | cell1, | ||
const unsigned int | degree, | ||
const BooleanOperation & | bool_op, | ||
const Mapping< dim0, spacedim > & | mapping0 = (ReferenceCells::get_hypercube< dim0 >() .template get_default_linear_mapping< dim0, spacedim >()), | ||
const Mapping< dim1, spacedim > & | mapping1 = (ReferenceCells::get_hypercube< dim1 >() .template get_default_linear_mapping< dim1, spacedim >()) ) |
Compute a Quadrature formula over the polygonal/polyhedral region described by a BooleanOperation between two deal.II cell_iterator objects. The degree of the Quadrature rule is specified by the argument degree
. This function splits the polygon/polyhedron into simplices and collects QGaussSimplex quadrature rules.
[in] | cell0 | A cell_iterator to the first deal.II cell. |
[in] | cell1 | A cell_iterator to the second deal.II cell. |
[in] | degree | The degree of accuracy you wish to get for the global quadrature formula. |
[in] | bool_op | The BooleanOperation to be performed. |
[in] | mapping0 | Mapping object for the first cell. |
[in] | mapping1 | Mapping object for the first cell. |
::Quadrature< spacedim > CGALWrappers::compute_quadrature_on_intersection | ( | const typename ::Triangulation< dim0, spacedim >::cell_iterator & | cell0, |
const typename ::Triangulation< dim1, spacedim >::cell_iterator & | cell1, | ||
const unsigned int | degree, | ||
const Mapping< dim0, spacedim > & | mapping0, | ||
const Mapping< dim1, spacedim > & | mapping1 ) |
A specialization of the function above when the BooleanOperation is an intersection. The rationale behind this specialization is that deal.II affine cells are convex sets, and as the intersection of convex sets is itself convex, this function internally exploits this to use a cheaper way to mesh the inside.
[in] | cell0 | A cell_iterator to the first deal.II cell. |
[in] | cell1 | A cell_iterator to the second deal.II cell. |
[in] | mapping0 | Mapping object for the first cell. |
[in] | mapping1 | Mapping object for the first cell. |
[in] | degree | The degree of accuracy you wish to get for the global quadrature formula. |
void CGALWrappers::remesh_surface | ( | CGAL::Surface_mesh< CGALPointType > & | surface_mesh, |
const AdditionalData< 3 > & | data = AdditionalData< 3 >{} ) |
Remesh a CGAL::Surface_mesh.
If the domain has 1-dimensional exposed features, the criteria includes a sizing field to guide the sampling of 1-dimensional features with protecting balls centers.
This utility should be exploited to improve the quality of a mesh coming from boolean operations. As an example, consider two CGAL::Surface_mesh objects describing two hyper_spheres that intersect non-trivially. After computing their corefinement and union with compute_boolean_operation(), the resulting mesh is the following:
Clearly, elements (triangles) like the ones arising along the intersection of the two spheres should be avoided as they're quite degenerate and would decrease the accuracy of the results. A call to the present function will try to remesh the surface, according to the given named parameters, to improve its quality. In the present case, the result is the following:
surface_mesh | The input CGAL::Surface_mesh. | |
[in] | data | AdditionalData object to pass to the CGAL::make_mesh_3 function. See the documentation of that struct for a description of those parameters. |
const T * CGALWrappers::get_if_ | ( | const std::variant< Types... > * | v | ) |
Definition at line 105 of file intersections.cc.
const T * CGALWrappers::get_if_ | ( | const boost::variant< Types... > * | v | ) |
Definition at line 112 of file intersections.cc.