mirror of https://github.com/CGAL/cgal
Document PMP/curvature.h
This commit is contained in:
parent
db297d7eff
commit
c418aea7c2
|
|
@ -25,7 +25,7 @@
|
|||
/// \ingroup PkgPolygonMeshProcessingRef
|
||||
|
||||
/// \defgroup PMP_measure_grp Geometric Measure Functions
|
||||
/// Functions to compute lengths of edges and borders, areas of faces and patches, as well as volumes of closed meshes.
|
||||
/// Functions to compute discrete curvatures, lengths of edges and borders, areas of faces and patches, volumes of closed meshes.
|
||||
/// \ingroup PkgPolygonMeshProcessingRef
|
||||
|
||||
/// \defgroup PMP_orientation_grp Orientation Functions
|
||||
|
|
@ -239,6 +239,11 @@ The page \ref bgl_namedparameters "Named Parameters" describes their usage.
|
|||
- `CGAL::Polygon_mesh_processing::sample_triangle_mesh()`
|
||||
|
||||
\cgalCRPSection{Geometric Measure Functions}
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::angle_sum()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::discrete_Gaussian_curvatures()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::discrete_Gaussian_curvature()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::discrete_mean_curvature()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::discrete_mean_curvatures()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::edge_length()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::squared_edge_length()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::face_area()` \endlink
|
||||
|
|
@ -248,7 +253,6 @@ The page \ref bgl_namedparameters "Named Parameters" describes their usage.
|
|||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::face_border_length()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::longest_border()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::centroid()` \endlink
|
||||
- \link PMP_measure_grp `CGAL::Polygon_mesh_processing::match_faces()` \endlink
|
||||
|
||||
\cgalCRPSection{Feature Detection Functions}
|
||||
- `CGAL::Polygon_mesh_processing::sharp_edges_segmentation()`
|
||||
|
|
|
|||
|
|
@ -1225,6 +1225,17 @@ compute the curvatures on a specific vertex.
|
|||
|
||||
\cgalExample{Polygon_mesh_processing/interpolated_corrected_curvatures_vertex.cpp}
|
||||
|
||||
\subsection DCurvartures Discrete Curvatures
|
||||
|
||||
The package also provides methods to compute the standard, non-interpolated discrete mean and Gaussian
|
||||
curvatures on triangle meshes, based on the work of Meyer et al. \cgalCite{cgal:mdsb-ddgot-02}.
|
||||
These curvatures are computed at each vertex of the mesh, and are based on the angles of the incident
|
||||
triangles. The functions are:
|
||||
- `CGAL::Polygon_mesh_processing::discrete_mean_curvature()`
|
||||
- `CGAL::Polygon_mesh_processing::discrete_mean_curvatures()`
|
||||
- `CGAL::Polygon_mesh_processing::discrete_Gaussian_curvature()`
|
||||
- `CGAL::Polygon_mesh_processing::discrete_Gaussian_curvatures()`
|
||||
|
||||
****************************************
|
||||
\section PMPSlicer Slicer
|
||||
|
||||
|
|
|
|||
|
|
@ -8,10 +8,11 @@
|
|||
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
|
||||
//
|
||||
//
|
||||
// Author(s) : Mael Rouxel-Labbé
|
||||
// Author(s) : Andreas Fabri,
|
||||
// Mael Rouxel-Labbé
|
||||
|
||||
#ifndef CGAL_PMP_INTERNAL_CURVATURE_H
|
||||
#define CGAL_PMP_INTERNAL_CURVATURE_H
|
||||
#ifndef CGAL_PMP_CURVATURE_H
|
||||
#define CGAL_PMP_CURVATURE_H
|
||||
|
||||
#include <CGAL/license/Polygon_mesh_processing/measure.h>
|
||||
|
||||
|
|
@ -27,10 +28,12 @@ namespace CGAL {
|
|||
namespace Polygon_mesh_processing {
|
||||
|
||||
/**
|
||||
* \ingroup PMP_vertex_angle_grp
|
||||
* \ingroup PMP_measure_grp
|
||||
*
|
||||
* computes the sum of the angles around a vertex.
|
||||
*
|
||||
* The angle sum is given in degrees.
|
||||
*
|
||||
* @tparam PolygonMesh a model of `HalfedgeGraph`
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
|
|
@ -60,7 +63,6 @@ namespace Polygon_mesh_processing {
|
|||
* or the geometric traits class deduced from the point property map of `pmesh`.
|
||||
*
|
||||
* \warning This function involves trigonometry.
|
||||
*
|
||||
*/
|
||||
template<typename PolygonMesh,
|
||||
typename CGAL_NP_TEMPLATE_PARAMETERS>
|
||||
|
|
@ -105,6 +107,45 @@ angle_sum(typename boost::graph_traits<PolygonMesh>::vertex_descriptor v,
|
|||
|
||||
// Discrete Gaussian Curvature
|
||||
|
||||
/**
|
||||
* \ingroup PMP_measure_grp
|
||||
*
|
||||
* computes the discrete Gaussian curvature at a vertex.
|
||||
*
|
||||
* We refer to Meyer et al. \cgalCite{cgal:mdsb-ddgot-02} for the definition of <i>discrete Gaussian curvature</i>.
|
||||
*
|
||||
* @tparam TriangleMesh a model of `HalfedgeGraph`
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* @param v the vertex whose discrete Gaussian curvature is being computed
|
||||
* @param tmesh the triangle mesh to which `v` belongs
|
||||
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tmesh`}
|
||||
* \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
|
||||
* as key type and `%Point_3` as value type}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, tmesh)`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{geom_traits}
|
||||
* \cgalParamDescription{an instance of a geometric traits class}
|
||||
* \cgalParamType{The traits class must be a model of `Kernel`.}
|
||||
* \cgalParamDefault{a \cgal kernel deduced from the point type, using `CGAL::Kernel_traits`}
|
||||
* \cgalParamExtra{The geometric traits class must be compatible with the vertex point type.}
|
||||
* \cgalParamNEnd
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* @return the discrete Gaussian curvature at `v`. The return type `FT` is a number type either deduced
|
||||
* from the `geom_traits` \ref bgl_namedparameters "Named Parameters" if provided,
|
||||
* or the geometric traits class deduced from the point property map of `tmesh`.
|
||||
*
|
||||
* \warning This function involves trigonometry.
|
||||
* \warning The current formulation is not well defined for border vertices.
|
||||
*
|
||||
* \pre `tmesh` is a triangle mesh
|
||||
*/
|
||||
template <typename TriangleMesh,
|
||||
typename CGAL_NP_TEMPLATE_PARAMETERS>
|
||||
#ifdef DOXYGEN_RUNNING
|
||||
|
|
@ -113,7 +154,7 @@ FT
|
|||
typename GetGeomTraits<TriangleMesh, CGAL_NP_CLASS>::type::FT
|
||||
#endif
|
||||
discrete_Gaussian_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descriptor v,
|
||||
const TriangleMesh& tm,
|
||||
const TriangleMesh& tmesh,
|
||||
const CGAL_NP_CLASS& np = parameters::default_values())
|
||||
{
|
||||
using parameters::choose_parameter;
|
||||
|
|
@ -128,7 +169,7 @@ discrete_Gaussian_curvature(typename boost::graph_traits<TriangleMesh>::vertex_d
|
|||
|
||||
GeomTraits gt = choose_parameter<GeomTraits>(get_parameter(np, internal_np::geom_traits));
|
||||
VertexPointMap vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
|
||||
get_const_property_map(vertex_point, tm));
|
||||
get_const_property_map(vertex_point, tmesh));
|
||||
|
||||
typename GeomTraits::Construct_vector_3 vector =
|
||||
gt.construct_vector_3_object();
|
||||
|
|
@ -141,13 +182,13 @@ discrete_Gaussian_curvature(typename boost::graph_traits<TriangleMesh>::vertex_d
|
|||
|
||||
FT angle_sum = 0;
|
||||
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tm))
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tmesh))
|
||||
{
|
||||
if(is_border(h, tm))
|
||||
if(is_border(h, tmesh))
|
||||
continue;
|
||||
|
||||
const Vector_3 v0 = vector(get(vpm, v), get(vpm, target(next(h, tm), tm))); // p1p2
|
||||
const Vector_3 v1 = vector(get(vpm, v), get(vpm, source(h, tm))); // p1p0
|
||||
const Vector_3 v0 = vector(get(vpm, v), get(vpm, target(next(h, tmesh), tmesh))); // p1p2
|
||||
const Vector_3 v1 = vector(get(vpm, v), get(vpm, source(h, tmesh))); // p1p0
|
||||
|
||||
const FT dot = scalar_product(v0, v1);
|
||||
const Vector_3 cross = cross_product(v0, v1);
|
||||
|
|
@ -172,28 +213,108 @@ discrete_Gaussian_curvature(typename boost::graph_traits<TriangleMesh>::vertex_d
|
|||
}
|
||||
}
|
||||
|
||||
Weights::Secure_cotangent_weight_with_voronoi_area<TriangleMesh, VertexPointMap, GeomTraits> wc(tm, vpm, gt);
|
||||
Weights::Secure_cotangent_weight_with_voronoi_area<TriangleMesh, VertexPointMap, GeomTraits> wc(tmesh, vpm, gt);
|
||||
|
||||
const FT gaussian_curvature = (2 * CGAL_PI - angle_sum) / wc.voronoi(v);
|
||||
|
||||
return gaussian_curvature;
|
||||
}
|
||||
|
||||
/**
|
||||
* \ingroup PMP_measure_grp
|
||||
*
|
||||
* computes the discrete Gaussian curvatures at the vertices of a mesh.
|
||||
*
|
||||
* We refer to Meyer et al. \cgalCite{cgal:mdsb-ddgot-02} for the definition of <i>discrete Gaussian curvature</i>.
|
||||
*
|
||||
* @tparam TriangleMesh a model of `HalfedgeGraph`
|
||||
* @tparam VertexCurvatureMap must be a model of `WritablePropertyMap` with key type
|
||||
* `boost::graph_traits<TriangleMesh>::%vertex_descriptor` and value type `FT`,
|
||||
* which is either `geom_traits::FT` if this named parameter is provided,
|
||||
* or `kernel::FT` with the kernel deduced from from the point property map of `tmesh`.
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* @param tmesh the triangle mesh to which `v` belongs
|
||||
* @param vcm the property map that contains the computed discrete curvatures
|
||||
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tmesh`}
|
||||
* \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
|
||||
* as key type and `%Point_3` as value type}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, tmesh)`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{geom_traits}
|
||||
* \cgalParamDescription{an instance of a geometric traits class}
|
||||
* \cgalParamType{The traits class must be a model of `Kernel`.}
|
||||
* \cgalParamDefault{a \cgal kernel deduced from the point type, using `CGAL::Kernel_traits`}
|
||||
* \cgalParamExtra{The geometric traits class must be compatible with the vertex point type.}
|
||||
* \cgalParamNEnd
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* \warning This function involves trigonometry.
|
||||
* \warning The current formulation is not well defined for border vertices.
|
||||
*
|
||||
* \pre `tmesh` is a triangle mesh
|
||||
*/
|
||||
template <typename TriangleMesh,
|
||||
typename VertexCurvatureMap,
|
||||
typename CGAL_NP_TEMPLATE_PARAMETERS>
|
||||
void discrete_Gaussian_curvatures(const TriangleMesh& tm,
|
||||
void discrete_Gaussian_curvatures(const TriangleMesh& tmesh,
|
||||
VertexCurvatureMap vcm,
|
||||
const CGAL_NP_CLASS& np = parameters::default_values())
|
||||
{
|
||||
using vertex_descriptor = typename boost::graph_traits<TriangleMesh>::vertex_descriptor;
|
||||
|
||||
for(vertex_descriptor v : vertices(tm))
|
||||
put(vcm, v, discrete_Gaussian_curvature(v, tm, np));
|
||||
for(vertex_descriptor v : vertices(tmesh))
|
||||
{
|
||||
put(vcm, v, discrete_Gaussian_curvature(v, tmesh, np));
|
||||
// std::cout << "curvature: " << get(vcm, v) << std::endl;
|
||||
}
|
||||
}
|
||||
|
||||
// Discrete Mean Curvature
|
||||
|
||||
/**
|
||||
* \ingroup PMP_measure_grp
|
||||
*
|
||||
* computes the discrete mean curvature at a vertex.
|
||||
*
|
||||
* We refer to Meyer et al. \cgalCite{cgal:mdsb-ddgot-02} for the definition of <i>discrete mean curvature</i>.
|
||||
*
|
||||
* @tparam TriangleMesh a model of `HalfedgeGraph`
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* @param v the vertex whose discrete mean curvature is being computed
|
||||
* @param tmesh the triangle mesh to which `v` belongs
|
||||
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tmesh`}
|
||||
* \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
|
||||
* as key type and `%Point_3` as value type}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, tmesh)`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{geom_traits}
|
||||
* \cgalParamDescription{an instance of a geometric traits class}
|
||||
* \cgalParamType{The traits class must be a model of `Kernel`.}
|
||||
* \cgalParamDefault{a \cgal kernel deduced from the point type, using `CGAL::Kernel_traits`}
|
||||
* \cgalParamExtra{The geometric traits class must be compatible with the vertex point type.}
|
||||
* \cgalParamNEnd
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* @return the discrete mean curvature at `v`. The return type `FT` is a number type either deduced
|
||||
* from the `geom_traits` \ref bgl_namedparameters "Named Parameters" if provided,
|
||||
* or the geometric traits class deduced from the point property map of `tmesh`.
|
||||
*
|
||||
* \warning The current formulation is not well defined for border vertices.
|
||||
*
|
||||
* \pre `tmesh` is a triangle mesh
|
||||
*/
|
||||
template <typename TriangleMesh,
|
||||
typename CGAL_NP_TEMPLATE_PARAMETERS>
|
||||
#ifdef DOXYGEN_RUNNING
|
||||
|
|
@ -202,7 +323,7 @@ FT
|
|||
typename GetGeomTraits<TriangleMesh, CGAL_NP_CLASS>::type::FT
|
||||
#endif
|
||||
discrete_mean_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descriptor v,
|
||||
const TriangleMesh& tm,
|
||||
const TriangleMesh& tmesh,
|
||||
const CGAL_NP_CLASS& np = parameters::default_values())
|
||||
{
|
||||
using parameters::choose_parameter;
|
||||
|
|
@ -220,7 +341,7 @@ discrete_mean_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descr
|
|||
|
||||
GeomTraits gt = choose_parameter<GeomTraits>(get_parameter(np, internal_np::geom_traits));
|
||||
VertexPointMap vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
|
||||
get_const_property_map(vertex_point, tm));
|
||||
get_const_property_map(vertex_point, tmesh));
|
||||
|
||||
#if 0
|
||||
typename GeomTraits::Compute_squared_distance_3 squared_distance =
|
||||
|
|
@ -231,12 +352,12 @@ discrete_mean_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descr
|
|||
const FT two_pi = 2 * CGAL_PI;
|
||||
|
||||
FT hi = 0;
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tm))
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tmesh))
|
||||
{
|
||||
const Point_3& p = get(vpm, source(h, tm));
|
||||
const Point_3& q = get(vpm, target(h, tm));
|
||||
const Point_3& r = get(vpm, target(next(h, tm), tm));
|
||||
const Point_3& s = get(vpm, target(next(opposite(h, tm), tm), tm));
|
||||
const Point_3& p = get(vpm, source(h, tmesh));
|
||||
const Point_3& q = get(vpm, target(h, tmesh));
|
||||
const Point_3& r = get(vpm, target(next(h, tmesh), tmesh));
|
||||
const Point_3& s = get(vpm, target(next(opposite(h, tmesh), tmesh), tmesh));
|
||||
const FT l = squared_distance(p,q);
|
||||
|
||||
FT phi = CGAL_PI * dihedral_angle(p, q, r, s) / FT(180);
|
||||
|
|
@ -260,27 +381,27 @@ discrete_mean_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descr
|
|||
typename GeomTraits::Compute_squared_length_3 squared_length =
|
||||
gt.compute_squared_length_3_object();
|
||||
|
||||
Weights::Secure_cotangent_weight_with_voronoi_area<TriangleMesh, VertexPointMap, GeomTraits> wc(tm, vpm, gt);
|
||||
Weights::Secure_cotangent_weight_with_voronoi_area<TriangleMesh, VertexPointMap, GeomTraits> wc(tmesh, vpm, gt);
|
||||
|
||||
Vector_3 kh = vector(CGAL::NULL_VECTOR);
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tm))
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_target(v, tmesh))
|
||||
{
|
||||
const vertex_descriptor v1 = source(h, tm);
|
||||
const vertex_descriptor v1 = source(h, tmesh);
|
||||
|
||||
const Point_ref p0 = get(vpm, v);
|
||||
const Point_ref p1 = get(vpm, v1);
|
||||
|
||||
FT local_c = 0;
|
||||
if(!is_border(h, tm))
|
||||
if(!is_border(h, tmesh))
|
||||
{
|
||||
const vertex_descriptor v2 = target(next(h, tm), tm);
|
||||
const vertex_descriptor v2 = target(next(h, tmesh), tmesh);
|
||||
const Point_ref p2 = get(vpm, v2);
|
||||
local_c += Weights::cotangent_3_clamped(p0, p2, p1, gt);
|
||||
}
|
||||
|
||||
if(!is_border(opposite(h, tm), tm))
|
||||
if(!is_border(opposite(h, tmesh), tmesh))
|
||||
{
|
||||
const vertex_descriptor v3 = target(next(opposite(h, tm), tm), tm);
|
||||
const vertex_descriptor v3 = target(next(opposite(h, tmesh), tmesh), tmesh);
|
||||
const Point_ref p3 = get(vpm, v3);
|
||||
local_c += Weights::cotangent_3_clamped(p1, p3, p0, gt);
|
||||
}
|
||||
|
|
@ -297,20 +418,58 @@ discrete_mean_curvature(typename boost::graph_traits<TriangleMesh>::vertex_descr
|
|||
#endif
|
||||
}
|
||||
|
||||
/**
|
||||
* \ingroup PMP_measure_grp
|
||||
*
|
||||
* computes the discrete mean curvatures at the vertices of a mesh.
|
||||
*
|
||||
* We refer to Meyer et al. \cgalCite{cgal:mdsb-ddgot-02} for the definition of <i>discrete mean curvature</i>.
|
||||
*
|
||||
* @tparam TriangleMesh a model of `HalfedgeGraph`
|
||||
* @tparam VertexCurvatureMap must be a model of `WritablePropertyMap` with key type
|
||||
* `boost::graph_traits<TriangleMesh>::%vertex_descriptor` and value type `FT`,
|
||||
* which is either `geom_traits::FT` if this named parameter is provided,
|
||||
* or `kernel::FT` with the kernel deduced from from the point property map of `tmesh`.
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* @param tmesh the triangle mesh to which `v` belongs
|
||||
* @param vcm the property map that contains the computed discrete curvatures
|
||||
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tmesh`}
|
||||
* \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
|
||||
* as key type and `%Point_3` as value type}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, tmesh)`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{geom_traits}
|
||||
* \cgalParamDescription{an instance of a geometric traits class}
|
||||
* \cgalParamType{The traits class must be a model of `Kernel`.}
|
||||
* \cgalParamDefault{a \cgal kernel deduced from the point type, using `CGAL::Kernel_traits`}
|
||||
* \cgalParamExtra{The geometric traits class must be compatible with the vertex point type.}
|
||||
* \cgalParamNEnd
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* \warning The current formulation is not well defined for border vertices.
|
||||
*
|
||||
* \pre `tmesh` is a triangle mesh
|
||||
*/
|
||||
template <typename TriangleMesh,
|
||||
typename VertexCurvatureMap,
|
||||
typename CGAL_NP_TEMPLATE_PARAMETERS>
|
||||
void discrete_mean_curvatures(const TriangleMesh& tm,
|
||||
void discrete_mean_curvatures(const TriangleMesh& tmesh,
|
||||
VertexCurvatureMap vcm,
|
||||
const CGAL_NP_CLASS& np = parameters::default_values())
|
||||
{
|
||||
using vertex_descriptor = typename boost::graph_traits<TriangleMesh>::vertex_descriptor;
|
||||
|
||||
for(vertex_descriptor v : vertices(tm))
|
||||
put(vcm, v, discrete_mean_curvature(v, tm, np));
|
||||
for(vertex_descriptor v : vertices(tmesh))
|
||||
put(vcm, v, discrete_mean_curvature(v, tmesh, np));
|
||||
}
|
||||
|
||||
} // namespace Polygon_mesh_processing
|
||||
} // namespace CGAL
|
||||
|
||||
#endif //CGAL_PMP_INTERNAL_CURVATURE_H
|
||||
#endif //CGAL_PMP_CURVATURE_H
|
||||
|
|
|
|||
Loading…
Reference in New Issue