refine and fairing functions requires a triangle mesh

This commit is contained in:
Sébastien Loriot 2016-02-15 14:23:09 +01:00
parent 77b69cc881
commit 6c4f1274ef
3 changed files with 50 additions and 42 deletions

View File

@ -72,7 +72,7 @@ The visual results of aforementioned steps are depicted by \cgalFigureRef{Mech_s
\subsubsection Meshing
Refinement and fairing functions can be applied to an arbitrary region on a mesh, using :
Refinement and fairing functions can be applied to an arbitrary region on a triangle mesh, using :
- `CGAL::Polygon_mesh_processing::refine()` : given a set of facets on a mesh, refines the region.
- `CGAL::Polygon_mesh_processing::fair()` : given a set of vertices on a mesh, fairs the region.
@ -124,10 +124,10 @@ Isotropic remeshing. (a) Triangulated input surface mesh.
\subsection MeshingExamples Meshing Examples
\subsubsection MeshingExample_1 Refine and Fair a Region on a Polygon Mesh
\subsubsection MeshingExample_1 Refine and Fair a Region on a Triangle Mesh
The following example calls the functions `CGAL::Polygon_mesh_processing::refine()`
and `CGAL::Polygon_mesh_processing::fair()` for some selected regions on the input polygon mesh.
and `CGAL::Polygon_mesh_processing::fair()` for some selected regions on the input triangle mesh.
\cgalExample{Polygon_mesh_processing/refine_fair_example.cpp}
@ -477,8 +477,8 @@ For each requirement, the table answers a question :
<center>
Function or Class | Pure Triangle | Self-Intersection Free | Closed
:------------------------------------ | :------------: | :--------------------: | :-------:
`fair()` | no | no | no
`refine()` | no | no | no
`fair()` | yes | no | no
`refine()` | yes | no | no
`triangulate_faces()` | no | no | no
`isotropic_remeshing()` | no | no | no
`split_long_edges()` | no | no | no
@ -494,7 +494,7 @@ Function or Class | Pure Triangle | Self-Intersection Free
`stitch_borders()` | no | no | no
`polygon_soup_to_polygon_mesh()` | no | no | no
`orient_polygon_soup()` | no | no | no
`remove_isolated_vertices()` | no | no | no
`remove_isolated_vertices()` | no | no | no
`compute_face_normal()` | no | no | no
`compute_face_normals()` | no | no | no
`compute_vertex_normal()` | no | no | no

View File

@ -41,10 +41,10 @@ namespace internal {
// weight_calculator a function object to calculate weights, defaults to Cotangent weights and can be omitted
template<typename SparseLinearSolver,
typename WeightCalculator,
typename PolygonMesh,
typename TriangleMesh,
typename VertexRange,
typename VertexPointMap>
bool fair(PolygonMesh& pmesh,
bool fair(TriangleMesh& tmesh,
const VertexRange& vertices,
SparseLinearSolver solver,
WeightCalculator weight_calculator,
@ -52,8 +52,8 @@ namespace internal {
VertexPointMap vpmap)
{
CGAL::Polygon_mesh_processing::internal::Fair_Polyhedron_3
<PolygonMesh, SparseLinearSolver, WeightCalculator, VertexPointMap>
fair_functor(pmesh, vpmap, weight_calculator);
<TriangleMesh, SparseLinearSolver, WeightCalculator, VertexPointMap>
fair_functor(tmesh, vpmap, weight_calculator);
return fair_functor.fair(vertices, solver, continuity);
}
@ -61,7 +61,7 @@ namespace internal {
/*!
\ingroup PMP_meshing_grp
@brief fairs a region on a polygon mesh.
@brief fairs a region on a triangle mesh.
The points of the selected vertices are
relocated to yield an as-smooth-as-possible surface patch,
based on solving a linear bi-Laplacian system with boundary constraints,
@ -76,39 +76,43 @@ namespace internal {
Fairing might fail if fixed vertices, which are used as boundary conditions,
do not suffice to solve constructed linear system.
Note that if the vertex range to which fairing is applied contains all the vertices of the polygon mesh,
Note that if the vertex range to which fairing is applied contains all the vertices of the triangle mesh,
fairing does not fail, but the mesh gets shrinked to `CGAL::ORIGIN`.
@tparam PolygonMesh a model of `FaceGraph` and `MutableFaceGraph`
@tparam TriangleMesh a model of `FaceGraph` and `MutableFaceGraph`
that has an internal property map for `CGAL::vertex_point_t`
@tparam VertexRange a range of vertex descriptors of `PolygonMesh`, model of `Range`.
@tparam VertexRange a range of vertex descriptors of `TriangleMesh`, model of `Range`.
Its iterator type is `InputIterator`.
@tparam NamedParameters a sequence of \ref namedparameters
@param pmesh the polygon mesh with patches to be faired
@param tmesh the triangle mesh with patches to be faired
@param vertices the vertices of the patches to be faired (the positions of only those vertices will be changed)
@param np optional sequence of \ref namedparameters among the ones listed below
\cgalNamedParamsBegin
\cgalParamBegin{vertex_point_map} the property map with the points associated to the vertices of `pmesh` \cgalParamEnd
\cgalParamBegin{vertex_point_map} the property map with the points associated to the vertices of `tmesh` \cgalParamEnd
\cgalParamBegin{fairing_continuity} tangential continuity of the output surface patch. The larger `fairing_continuity` gets, the more fixed vertices are required \cgalParamEnd
\cgalParamBegin{sparse_linear_solver} an instance of the sparse linear solver used for fairing \cgalParamEnd
\cgalNamedParamsEnd
@return `true` if fairing is successful, otherwise no vertices are relocated
@pre `is_triangle_mesh(tmesh)`
@todo accuracy of solvers are not good, for example when there is no boundary condition pre_factor should fail, but it does not.
*/
template<typename PolygonMesh,
template<typename TriangleMesh,
typename VertexRange,
typename NamedParameters>
bool fair(PolygonMesh& pmesh,
bool fair(TriangleMesh& tmesh,
const VertexRange& vertices,
const NamedParameters& np)
{
using boost::get_param;
using boost::choose_param;
CGAL_precondition(is_triangle_mesh(tmesh));
#if defined(CGAL_EIGEN3_ENABLED)
#if EIGEN_VERSION_AT_LEAST(3,2,0)
typedef CGAL::Eigen_solver_traits<Eigen::SparseLU<
@ -133,27 +137,27 @@ namespace internal {
"The function `fair` requires Eigen3 version 3.2 or later.");
#endif
typedef typename GetVertexPointMap < PolygonMesh, NamedParameters>::type VPMap;
typedef CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<PolygonMesh, VPMap>
typedef typename GetVertexPointMap < TriangleMesh, NamedParameters>::type VPMap;
typedef CGAL::internal::Cotangent_weight_with_voronoi_area_fairing<TriangleMesh, VPMap>
Default_Weight_calculator;
VPMap vpmap_ = choose_param(get_param(np, vertex_point), get(CGAL::vertex_point, pmesh));
VPMap vpmap_ = choose_param(get_param(np, vertex_point), get(CGAL::vertex_point, tmesh));
return internal::fair(pmesh, vertices,
return internal::fair(tmesh, vertices,
choose_param(get_param(np, sparse_linear_solver), Default_solver()),
choose_param(get_param(np, weight_calculator), Default_Weight_calculator(pmesh, vpmap_)),
choose_param(get_param(np, weight_calculator), Default_Weight_calculator(tmesh, vpmap_)),
choose_param(get_param(np, fairing_continuity), 1),
vpmap_
);
}
template<typename PolygonMesh, typename VertexRange>
bool fair(PolygonMesh& pmesh, const VertexRange& vertices)
template<typename TriangleMesh, typename VertexRange>
bool fair(TriangleMesh& tmesh, const VertexRange& vertices)
{
return fair(pmesh, vertices,
return fair(tmesh, vertices,
CGAL::Polygon_mesh_processing::parameters::all_default());
}
} //end namespace Polygon_mesh_processing
} //end namespace CGAL

View File

@ -32,28 +32,28 @@ namespace Polygon_mesh_processing {
/*!
\ingroup PMP_meshing_grp
@brief refines a region of a polygon mesh
@brief refines a region of a triangle mesh
@tparam PolygonMesh model of `MutableFaceGraph`
@tparam TriangleMesh model of `MutableFaceGraph`
that has an internal property map for `CGAL::vertex_point_t`
@tparam FaceRange range of face descriptors, model of `Range`.
Its iterator type is `InputIterator`.
@tparam FaceOutputIterator model of `OutputIterator`
holding `boost::graph_traits<PolygonMesh>::%face_descriptor` for patch faces
holding `boost::graph_traits<TriangleMesh>::%face_descriptor` for patch faces
@tparam VertexOutputIterator model of `OutputIterator`
holding `boost::graph_traits<PolygonMesh>::%vertex_descriptor` for patch vertices
holding `boost::graph_traits<TriangleMesh>::%vertex_descriptor` for patch vertices
@tparam NamedParameters a sequence of \ref namedparameters
@param pmesh polygon mesh with patches to be refined
@param tmesh triangle mesh with patches to be refined
@param faces the range of faces defining the patches to refine
@param faces_out output iterator into which descriptors of new faces are recorded
@param vertices_out output iterator into which descriptors of new vertices are recorded
@param np optional sequence of \ref namedparameters among the ones listed below
\cgalNamedParamsBegin
\cgalParamBegin{vertex_point_map} the property map with the points associated to the vertices of `pmesh`
\cgalParamBegin{vertex_point_map} the property map with the points associated to the vertices of `tmesh`
Instance of a class model of `ReadWritePropertyMap` \cgalParamEnd
\cgalParamBegin{density_control_factor} factor to control density of the ouput mesh,
\cgalParamBegin{density_control_factor} factor to control density of the output mesh,
where larger values lead to denser refinements.
The density of vertices of `faces_out` is this factor times higher than the vertices of `faces.` \cgalParamEnd
\cgalNamedParamsEnd
@ -61,15 +61,17 @@ namespace Polygon_mesh_processing {
@return pair of `faces_out` and `vertices_out`
\pre `is_triangle_mesh(tmesh)`
@todo current algorithm iterates 10 times at most, since (I guess) there is no termination proof.
*/
template<typename PolygonMesh,
template<typename TriangleMesh,
typename FaceRange,
typename FaceOutputIterator,
typename VertexOutputIterator,
typename NamedParameters>
std::pair<FaceOutputIterator, VertexOutputIterator>
refine(PolygonMesh& pmesh,
refine(TriangleMesh& tmesh,
const FaceRange& faces,
FaceOutputIterator faces_out,
VertexOutputIterator vertices_out,
@ -79,12 +81,14 @@ namespace Polygon_mesh_processing {
using boost::choose_param;
using boost::get_param;
typedef typename GetVertexPointMap<PolygonMesh,NamedParameters>::type VPmap;
CGAL_precondition(is_triangle_mesh(tmesh) );
typedef typename GetVertexPointMap<TriangleMesh,NamedParameters>::type VPmap;
VPmap vpm = choose_pmap(get_param(np, boost::vertex_point),
pmesh,
tmesh,
boost::vertex_point);
internal::Refine_Polyhedron_3<PolygonMesh, VPmap> refine_functor(pmesh, vpm);
internal::Refine_Polyhedron_3<TriangleMesh, VPmap> refine_functor(tmesh, vpm);
refine_functor.refine(faces,
faces_out,
vertices_out,
@ -93,18 +97,18 @@ namespace Polygon_mesh_processing {
}
///\cond SKIP_IN_MANUAL
template<typename PolygonMesh,
template<typename TriangleMesh,
typename FaceRange,
typename FaceOutputIterator,
typename VertexOutputIterator>
std::pair<FaceOutputIterator, VertexOutputIterator>
refine(PolygonMesh& pmesh,
refine(TriangleMesh& tmesh,
const FaceRange& faces,
FaceOutputIterator faces_out,
VertexOutputIterator vertices_out)
{
return refine(pmesh, faces, faces_out, vertices_out,
return refine(tmesh, faces, faces_out, vertices_out,
CGAL::Polygon_mesh_processing::parameters::all_default());
}
///\endcond