Class

CGAL::Polyhedral_mesh_domain_with_features_3<IGT>

Definition

The class Polyhedral_mesh_domain_with_features_3<IGT> implements a domain whose boundary is a simplicial polyhedral surface. This surface must be closed and free of intersection. It is a model of the concept MeshDomainWithFeatures_3. It also provides a member function to automatically detect sharp features from the input polyhedral surface.

Parameters

The parameter IGT stands for a geometric traits class providing the types and functors required to implement the intersection tests and intersection computations for polyhedral boundary surfaces. This parameter has to be instantiated with a model of the concept IntersectionGeometricTraits_3.

#include <CGAL/Polyhedral_mesh_domain_with_features_3.h>

Is Model for the Concepts

MeshDomainWithFeatures_3

Inherits From

CGAL::Mesh_domain_with_polyline_features_3<CGAL::Polyhedral_mesh_domain_3<CGAL::Mesh_polyhedron_3<IGT>::type,IGT> >

Types

Polyhedral_mesh_domain_with_features_3<IGT>::FT
Numerical type.

Creation

template <typename Polyhedron>
Polyhedral_mesh_domain_with_features_3<IGT> md ( Polyhedron p);
Constructs a Polyhedral_mesh_domain_with_features_3 from a Polyhedron. The only requirement on type Polyhedron is that CGAL::Polyhedron_3 should be constructible from Polyhedron. No feature detection is done at this level. Note that a copy of p will be done.


Polyhedral_mesh_domain_with_features_3<IGT> md ( std::string filename);
Constructs a Polyhedral_mesh_domain_with_features_3 from an off file. No feature detection is done at this level.

Operations

void md.detect_features ( FT angle_bound=120)
Detects sharp features of the internal polyhedron and inserts them in as features of the domain. angle_bound gives the maximum dihedral angle (in degrees) between two triangles of the internal polyhedron which is used to consider that the edge between those triangles is a feature edge.

See Also

MeshDomainWithFeatures_3
CGAL::Mesh_domain_with_polyline_features_3<MeshDomain>
CGAL::Polyhedral_mesh_domain_3<Polyhedron,IGT,TriangleAccessor>
Mesh_polyhedron_3<Gt>