mirror of https://github.com/CGAL/cgal
add cut_with_plane and new clip method
This commit is contained in:
parent
b19515e1f8
commit
5b4b19a1c8
|
|
@ -0,0 +1,649 @@
|
|||
// Copyright (c) 2024 GeometryFactory (France).
|
||||
// All rights reserved.
|
||||
//
|
||||
// This file is part of CGAL (www.cgal.org).
|
||||
//
|
||||
// $URL$
|
||||
// $Id$
|
||||
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
|
||||
//
|
||||
//
|
||||
// Author(s) : Sébastien Loriot
|
||||
|
||||
|
||||
#ifndef CGAL_POLYGON_MESH_PROCESSING_CUT_WITH_PLANE_H
|
||||
#define CGAL_POLYGON_MESH_PROCESSING_CUT_WITH_PLANE_H
|
||||
|
||||
#include <CGAL/license/Polygon_mesh_processing/corefinement.h>
|
||||
|
||||
#include <CGAL/Named_function_parameters.h>
|
||||
#include <CGAL/boost/graph/named_params_helper.h>
|
||||
#include <CGAL/Polygon_mesh_processing/connected_components.h>
|
||||
#include <CGAL/Polygon_mesh_processing/border.h>
|
||||
#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_TRIANGULATION
|
||||
#include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
|
||||
#endif
|
||||
#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_BOX_INTERSECTION_D
|
||||
#include <CGAL/Polygon_mesh_processing/self_intersections.h>
|
||||
#endif
|
||||
|
||||
namespace CGAL {
|
||||
namespace Polygon_mesh_processing {
|
||||
|
||||
template <class PolygonMesh>
|
||||
struct Default_cut_visitor
|
||||
{
|
||||
void before_edge_split(typename boost::graph_traits<PolygonMesh>::halfedge_descriptor, PolygonMesh&) {}
|
||||
void edge_split(typename boost::graph_traits<PolygonMesh>::halfedge_descriptor, PolygonMesh&) {}
|
||||
void vertices_on_cut(const std::vector<typename boost::graph_traits<PolygonMesh>::vertex_descriptor>&, PolygonMesh&){}
|
||||
};
|
||||
|
||||
// TODO: doc me or hide me in the np
|
||||
template <class Kernel>
|
||||
struct Orthogonal_cut_plane_traits
|
||||
{
|
||||
using FT = typename Kernel::FT;
|
||||
using Plane_3 = std::pair<int, FT>;
|
||||
using Point_3 = typename Kernel::Point_3;
|
||||
|
||||
struct Oriented_side_3
|
||||
{
|
||||
Oriented_side operator()(const Plane_3& plane, const Point_3& p) const
|
||||
{
|
||||
if (p[plane.first]==plane.second) return ON_ORIENTED_BOUNDARY;
|
||||
return p[plane.first]<plane.second ? ON_NEGATIVE_SIDE : ON_POSITIVE_SIDE;
|
||||
}
|
||||
};
|
||||
|
||||
// TODO here we should integrate the epsilon to reuse existing points (moving them on the grid as an option?)
|
||||
struct Construct_plane_line_intersection_point_3
|
||||
{
|
||||
Point_3 operator()(const Plane_3& plane, const Point_3& p, const Point_3& q)
|
||||
{
|
||||
//TODO: divide by largest value?
|
||||
FT alpha = (plane.second - q[plane.first]) / (p[plane.first] - q[plane.first]);
|
||||
std::array<FT,3> coords;
|
||||
for (int i=0;i<3;++i)
|
||||
{
|
||||
if (i==plane.first)
|
||||
coords[i] = plane.second;
|
||||
else
|
||||
coords[i] = (p[i]==q[i])?p[i]:p[i]*alpha +(1-alpha)*q[i];
|
||||
}
|
||||
return Point_3(coords[0], coords[1], coords[2]);
|
||||
}
|
||||
};
|
||||
|
||||
Oriented_side_3 oriented_side_3_object() const
|
||||
{
|
||||
return Oriented_side_3();
|
||||
}
|
||||
|
||||
Construct_plane_line_intersection_point_3 construct_plane_line_intersection_point_3_object() const
|
||||
{
|
||||
return Construct_plane_line_intersection_point_3();
|
||||
}
|
||||
|
||||
#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_BOX_INTERSECTION_D
|
||||
// for does self-intersect
|
||||
using Segment_3 = typename Kernel::Segment_3;
|
||||
using Triangle_3 = typename Kernel::Triangle_3;
|
||||
using Construct_segment_3 = typename Kernel::Construct_segment_3;
|
||||
using Construct_triangle_3 =typename Kernel::Construct_triangle_3;
|
||||
using Do_intersect_3 = typename Kernel::Do_intersect_3;
|
||||
Construct_segment_3 construct_segment_3_object() const { return Construct_segment_3(); }
|
||||
Construct_triangle_3 construct_triangle_3_object() const { return Construct_triangle_3(); }
|
||||
Do_intersect_3 do_intersect_3_object() const { return Do_intersect_3(); }
|
||||
#endif
|
||||
};
|
||||
|
||||
/*!
|
||||
* \ingroup PMP_corefinement_grp
|
||||
*
|
||||
* refines `pm` by inserting the intersection points of `plane` with its edges and faces.
|
||||
*
|
||||
* \tparam PolygonMesh a model of `HalfedgeListGraph`, `FaceListGraph`, and `MutableFaceGraph`
|
||||
* \tparam Plane_3 plane type, equal to `GeomTraits::Plane_3`, `GeomTraits` being the type of the parameter `geom_traits`.
|
||||
* \tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* \param pm input mesh to be refined
|
||||
* \param plane the plane used to refine the mesh
|
||||
* \param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below:
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
*
|
||||
* \cgalParamNBegin{concurrency_tag}
|
||||
* \cgalParamDescription{a tag indicating if the task should be done using one or several threads.}
|
||||
* \cgalParamType{Either `CGAL::Sequential_tag`, or `CGAL::Parallel_tag`, or `CGAL::Parallel_if_available_tag`}
|
||||
* \cgalParamDefault{`CGAL::Sequential_tag`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{edge_is_constrained_map}
|
||||
* \cgalParamDescription{a property map containing the constrained-or-not status of each edge of `pm`.
|
||||
* If an edge marked as constrained is split, the two resulting edges will be marked as constrained.}
|
||||
* \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<PolygonMesh>::%edge_descriptor`
|
||||
* as key type and `bool` as value type}
|
||||
* \cgalParamDefault{a constant property map returning `false` for any edge}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{edge_is_marked_map}
|
||||
* \cgalParamDescription{a property map filled by this function that puts `true` for all intersection edge of faces
|
||||
* of `pm` and `plane`, and `false` for all other edges.}
|
||||
* \cgalParamType{a class model of `WritablePropertyMap` with `boost::graph_traits<PolygonMesh>::%edge_descriptor`
|
||||
* as key type and `bool` as value type}
|
||||
* \cgalParamDefault{a constant property map returning `false` for any edge}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{vertex_oriented_side_map}
|
||||
* \cgalParamDescription{a property map filled by this function containing the position
|
||||
* of each vertex relative to the oriented plane `plane`.}
|
||||
* \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<PolygonMesh>::%vertex_descriptor`
|
||||
* as key type and `Oriented_side` as value type}
|
||||
* \cgalParamDefault{Dynamic vertex property map}
|
||||
* \cgalParamExtra{If the concurrenty tag is set to `Parallel_tag`, the property map might be filled by several thread at the same time.}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{visitor}
|
||||
* \cgalParamDescription{TODO add concept}
|
||||
* \cgalParamType{reference wrapper recommeded if it has state TODO}
|
||||
* \cgalParamDefault{None}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{do_not_triangulate_faces}
|
||||
* \cgalParamDescription{If the input mesh is triangulated and this parameter is set to `false`, the mesh will be kept triangulated.}
|
||||
* \cgalParamType{`bool`}
|
||||
* \cgalParamDefault{`true`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tm_in`}
|
||||
* \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<PolygonMesh>::%vertex_descriptor`
|
||||
* as key type and `GeomTraits::Point_3` as value type, `GeomTraits` being the type of the parameter `geom_traits`}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, pm)`}
|
||||
* \cgalParamExtra{If this parameter is omitted, an internal property map for `CGAL::vertex_point_t` must be available in `PolygonMesh`.}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{geom_traits}
|
||||
* \cgalParamDescription{an instance of a geometric traits class}
|
||||
* \cgalParamType{a class 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
|
||||
*/
|
||||
template <class PolygonMesh, class Plane_3, class NamedParameters = parameters::Default_named_parameters>
|
||||
void cut_with_plane(PolygonMesh& pm,
|
||||
const Plane_3& plane,
|
||||
const NamedParameters& np = parameters::default_values())
|
||||
{
|
||||
// TODO: concurrency tag
|
||||
// TODO: if you want to clip with many planes (**Kernel**),
|
||||
// it might be interesting to first classify all vertices with all planes
|
||||
// to limit the number of intersection points computed: several classify done lazily
|
||||
// on points not already eliminated (verices all out with adjacent vertices out too)
|
||||
// actually might be a global classifier to filter out edges, testing all planes at once per vertex and stop as soon as one is out
|
||||
using parameters::choose_parameter;
|
||||
using parameters::get_parameter;
|
||||
using parameters::get_parameter_reference;
|
||||
using parameters::is_default_parameter;
|
||||
|
||||
// graph typedefs
|
||||
using BGT = boost::graph_traits<PolygonMesh>;
|
||||
using face_descriptor = typename BGT::face_descriptor;
|
||||
using edge_descriptor = typename BGT::edge_descriptor;
|
||||
using halfedge_descriptor = typename BGT::halfedge_descriptor;
|
||||
using vertex_descriptor = typename BGT::vertex_descriptor;
|
||||
|
||||
// np typedefs
|
||||
using Default_ecm = Static_boolean_property_map<edge_descriptor, false>;
|
||||
using Default_visitor = Default_cut_visitor<PolygonMesh>;
|
||||
using Visitor_ref = typename internal_np::Lookup_named_param_def<internal_np::visitor_t, NamedParameters, Default_visitor>::reference;
|
||||
using GT = typename GetGeomTraits<PolygonMesh, NamedParameters>::type;
|
||||
GT traits = choose_parameter<GT>(get_parameter(np, internal_np::geom_traits));
|
||||
|
||||
static_assert(std::is_same_v<typename GT::Plane_3,Plane_3>);
|
||||
|
||||
auto ecm = choose_parameter<Default_ecm>(get_parameter(np, internal_np::edge_is_constrained));
|
||||
auto edge_is_marked = choose_parameter<Default_ecm>(get_parameter(np, internal_np::edge_is_marked_map));
|
||||
|
||||
Default_visitor default_visitor;
|
||||
Visitor_ref visitor = choose_parameter(get_parameter_reference(np, internal_np::visitor), default_visitor);
|
||||
auto vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
|
||||
get_property_map(vertex_point, pm));
|
||||
|
||||
bool triangulate = !choose_parameter(get_parameter(np, internal_np::do_not_triangulate_faces), true);
|
||||
if (triangulate && !is_triangle_mesh(pm))
|
||||
triangulate = false;
|
||||
|
||||
bool throw_on_self_intersection = choose_parameter(get_parameter(np, internal_np::use_compact_clipper), false);
|
||||
if (throw_on_self_intersection && !is_triangle_mesh(pm))
|
||||
throw_on_self_intersection = false;
|
||||
|
||||
typedef typename internal_np::Lookup_named_param_def <
|
||||
internal_np::concurrency_tag_t,
|
||||
NamedParameters,
|
||||
Sequential_tag
|
||||
> ::type Concurrency_tag;
|
||||
|
||||
// constexpr bool parallel_execution = std::is_same_v<Parallel_tag, Concurrency_tag>;
|
||||
|
||||
auto oriented_side = traits.oriented_side_3_object();
|
||||
auto intersection_point = traits.construct_plane_line_intersection_point_3_object();
|
||||
|
||||
// TODO: the default is not thread-safe for example for Polyhedron
|
||||
using V_os_tag = dynamic_vertex_property_t<Oriented_side>;
|
||||
static constexpr bool use_default_vosm =
|
||||
is_default_parameter<NamedParameters, internal_np::vertex_oriented_side_map_t>::value;
|
||||
|
||||
using Vertex_oriented_side_map =
|
||||
std::conditional_t<use_default_vosm,
|
||||
typename boost::property_map<PolygonMesh, V_os_tag>::type,
|
||||
typename internal_np::Get_param<typename NamedParameters::base,
|
||||
internal_np::vertex_oriented_side_map_t>::type>;
|
||||
|
||||
|
||||
Vertex_oriented_side_map vertex_os;
|
||||
if constexpr (use_default_vosm)
|
||||
vertex_os = get(V_os_tag(), pm);
|
||||
else
|
||||
vertex_os = get_parameter(np, internal_np::vertex_oriented_side_map);
|
||||
|
||||
std::vector<edge_descriptor> inters;
|
||||
|
||||
bool all_in = true;
|
||||
bool all_out = true;
|
||||
bool at_least_one_on = false;
|
||||
std::vector<vertex_descriptor> on_obnd;
|
||||
//TODO: parallel for
|
||||
for (vertex_descriptor v : vertices(pm))
|
||||
{
|
||||
Oriented_side os = oriented_side(plane, get(vpm, v));
|
||||
put(vertex_os,v,os);
|
||||
switch(os)
|
||||
{
|
||||
case ON_POSITIVE_SIDE:
|
||||
all_in = false;
|
||||
break;
|
||||
case ON_NEGATIVE_SIDE:
|
||||
all_out = false;
|
||||
break;
|
||||
case ON_ORIENTED_BOUNDARY:
|
||||
at_least_one_on=true;
|
||||
on_obnd.push_back(v);
|
||||
}
|
||||
}
|
||||
|
||||
if (at_least_one_on || (!all_in && !all_out))
|
||||
{
|
||||
//TODO: parallel for
|
||||
for(edge_descriptor e : edges(pm))
|
||||
{
|
||||
vertex_descriptor src = source(e,pm), tgt = target(e,pm);
|
||||
if (get(vertex_os, src)==CGAL::ON_ORIENTED_BOUNDARY)
|
||||
{
|
||||
if (get(vertex_os, tgt)==CGAL::ON_ORIENTED_BOUNDARY)
|
||||
{
|
||||
bool pure_coplanar=true;
|
||||
if (!is_border(e, pm))
|
||||
{
|
||||
halfedge_descriptor he=halfedge(e, pm);
|
||||
for (halfedge_descriptor h : halfedges_around_face(he,pm))
|
||||
if (get(vertex_os,target(h, pm))!=CGAL::ON_ORIENTED_BOUNDARY)
|
||||
{
|
||||
pure_coplanar=false;
|
||||
break;
|
||||
}
|
||||
if (pure_coplanar)
|
||||
{
|
||||
he=opposite(he, pm);
|
||||
for (halfedge_descriptor h : halfedges_around_face(he,pm))
|
||||
if (get(vertex_os, target(h, pm))!=CGAL::ON_ORIENTED_BOUNDARY)
|
||||
{
|
||||
pure_coplanar=false;
|
||||
break;
|
||||
}
|
||||
}
|
||||
}
|
||||
if (!pure_coplanar)
|
||||
put(edge_is_marked, e, true);
|
||||
}
|
||||
}
|
||||
else
|
||||
if (get(vertex_os, tgt)!=CGAL::ON_ORIENTED_BOUNDARY &&
|
||||
get(vertex_os, src)!=get(vertex_os, tgt))
|
||||
inters.push_back(e);
|
||||
}
|
||||
}
|
||||
|
||||
if (all_in || all_out)
|
||||
{
|
||||
visitor.vertices_on_cut(on_obnd, pm);
|
||||
return;
|
||||
}
|
||||
|
||||
std::unordered_map<face_descriptor, std::vector<halfedge_descriptor> > splitted_faces;
|
||||
|
||||
if (throw_on_self_intersection)
|
||||
{
|
||||
std::vector<face_descriptor> test_faces;
|
||||
for (edge_descriptor e : inters)
|
||||
{
|
||||
halfedge_descriptor h = halfedge(e, pm);
|
||||
if (!is_border(h,pm))
|
||||
test_faces.push_back(face(h,pm));
|
||||
h=opposite(h, pm);
|
||||
if (!is_border(h,pm))
|
||||
test_faces.push_back(face(h,pm));
|
||||
}
|
||||
std::sort(test_faces.begin(), test_faces.end());
|
||||
auto last = std::unique(test_faces.begin(), test_faces.end());
|
||||
test_faces.erase(last, test_faces.end());
|
||||
if (does_self_intersect<Concurrency_tag>(test_faces, pm, np))
|
||||
throw std::runtime_error("TODO Corefinement::Self_intersection_exception");
|
||||
}
|
||||
|
||||
//TODO: parallel for
|
||||
for (edge_descriptor e : inters)
|
||||
{
|
||||
halfedge_descriptor h = halfedge(e, pm);
|
||||
auto pts = CGAL::make_sorted_pair(get(vpm, source(h,pm)), get(vpm, target(h, pm)));
|
||||
typename GT::Point_3 ip = intersection_point(plane, pts.first, pts.second);
|
||||
|
||||
bool was_marked = get(ecm, edge(h, pm));
|
||||
visitor.before_edge_split(h, pm);
|
||||
h = CGAL::Euler::split_edge(h, pm);
|
||||
put(vpm, target(h, pm), ip);
|
||||
put(vertex_os, target(h, pm), ON_ORIENTED_BOUNDARY);
|
||||
visitor.edge_split(h, pm);
|
||||
if (was_marked)
|
||||
put(ecm, edge(h, pm), true);
|
||||
|
||||
if (!is_border(h, pm))
|
||||
splitted_faces[face(h, pm)].push_back(h);
|
||||
h=prev(opposite(h,pm),pm);
|
||||
if (!is_border(h, pm))
|
||||
splitted_faces[face(h, pm)].push_back(h);
|
||||
}
|
||||
|
||||
visitor.vertices_on_cut(on_obnd, pm);
|
||||
|
||||
// collect faces to be cut that have one vertex on the cut plane
|
||||
for (vertex_descriptor v : on_obnd)
|
||||
{
|
||||
halfedge_descriptor hv = halfedge(v, pm);
|
||||
for (halfedge_descriptor h : halfedges_around_target(hv, pm))
|
||||
{
|
||||
if (is_border(h, pm)) continue;
|
||||
Oriented_side prev_ori = get(vertex_os, source(h, pm)),
|
||||
next_ori = get(vertex_os, target(next(h, pm), pm));
|
||||
if ( prev_ori == ON_ORIENTED_BOUNDARY || next_ori == ON_ORIENTED_BOUNDARY) continue; // skip full edge
|
||||
if (prev_ori!=next_ori) splitted_faces[face(h, pm)].push_back(h); // skip tangency point
|
||||
}
|
||||
}
|
||||
|
||||
//TODO: parallel for
|
||||
for (std::pair<const face_descriptor, std::vector<halfedge_descriptor>>& f_and_hs : splitted_faces)
|
||||
{
|
||||
std::size_t nb_hedges = f_and_hs.second.size();
|
||||
|
||||
CGAL_assertion( nb_hedges%2 ==0 );
|
||||
|
||||
if (nb_hedges==2)
|
||||
{
|
||||
halfedge_descriptor h1=f_and_hs.second[0], h2=f_and_hs.second[1];
|
||||
CGAL_assertion(next(h1,pm)!=h2 && next(h2,pm)!=h1); // the edge does not already exist
|
||||
halfedge_descriptor res = CGAL::Euler::split_face(h1, h2, pm);
|
||||
put(edge_is_marked, edge(res, pm), true);
|
||||
|
||||
if (triangulate)
|
||||
{
|
||||
if (!is_triangle(res, pm))
|
||||
{
|
||||
// TODO: take the criteria in triangulate_faces ?
|
||||
halfedge_descriptor newh =
|
||||
CGAL::Euler::split_face(res, next(next(res, pm), pm), pm);
|
||||
put(edge_is_marked, edge(newh, pm), false);
|
||||
}
|
||||
else
|
||||
{
|
||||
res = opposite(res, pm);
|
||||
if (!is_triangle(res, pm))
|
||||
{
|
||||
// TODO: take the criteria in triangulate_faces ?
|
||||
halfedge_descriptor newh =
|
||||
CGAL::Euler::split_face(res, next(next(res, pm), pm), pm);
|
||||
put(edge_is_marked, edge(newh, pm), false);
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
// sort hedges to make them match
|
||||
CGAL_assertion(!triangulate);
|
||||
// TODO: need mechanism to make it robust even with EPICK
|
||||
auto less_hedge = [&pm, vpm](halfedge_descriptor h1, halfedge_descriptor h2)
|
||||
{
|
||||
return lexicographically_xyz_smaller(get(vpm,target(h1,pm)), get(vpm,target(h2,pm)));
|
||||
};
|
||||
std::sort(f_and_hs.second.begin(), f_and_hs.second.end(), less_hedge);
|
||||
|
||||
for (std::size_t i=0; i<nb_hedges; i+=2)
|
||||
{
|
||||
halfedge_descriptor h1=f_and_hs.second[i], h2=f_and_hs.second[i+1];
|
||||
CGAL_assertion(next(h1,pm)!=h2 && next(h2,pm)!=h1); // the edge does not already exist
|
||||
halfedge_descriptor res = CGAL::Euler::split_face(h1, h2, pm);
|
||||
put(edge_is_marked, edge(res, pm), true);
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
/**
|
||||
* \ingroup PMP_corefinement_grp
|
||||
*
|
||||
* \brief clips `tm` by keeping the part that is on the negative side of `plane` (side opposite to its normal vector).
|
||||
*
|
||||
* If `tm` is closed, the clipped part can be closed too if the named parameter `clip_volume` is set to `true`.
|
||||
* See Subsection \ref coref_clip for more details.
|
||||
*
|
||||
* \note `Plane_3` must be from the same %Kernel as the point of the internal vertex point map of `PolygonMesh`.
|
||||
* \note `Plane_3` must be from the same %Kernel as the point of the vertex point map of `tm`.
|
||||
*
|
||||
* \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(tm)` \endlink
|
||||
*
|
||||
* @tparam PolygonMesh a model of `MutableFaceGraph`, `HalfedgeListGraph` and `FaceListGraph`.
|
||||
* An internal property map for `CGAL::vertex_point_t` must be available.
|
||||
*
|
||||
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
|
||||
*
|
||||
* @param tm input triangulated surface mesh
|
||||
* @param plane plane whose negative side defines the half-space to intersect `tm` with.
|
||||
* `Plane_3` is the plane type for the same CGAL kernel as the point of the vertex point map of `tm`.
|
||||
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
*
|
||||
* \cgalParamNBegin{concurrency_tag}
|
||||
* \cgalParamDescription{a tag indicating if the task should be done using one or several threads.}
|
||||
* \cgalParamType{Either `CGAL::Sequential_tag`, or `CGAL::Parallel_tag`, or `CGAL::Parallel_if_available_tag`}
|
||||
* \cgalParamDefault{`CGAL::Sequential_tag`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{vertex_point_map}
|
||||
* \cgalParamDescription{a property map associating points to the vertices of `tm`}
|
||||
* \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<PolygonMesh>::%vertex_descriptor`
|
||||
* as key type and `%Point_3` as value type}
|
||||
* \cgalParamDefault{`boost::get(CGAL::vertex_point, tm)`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{visitor}
|
||||
* \cgalParamDescription{a visitor used to track the creation of new faces}
|
||||
* \cgalParamType{a class model of `PMPCorefinementVisitor`}
|
||||
* \cgalParamDefault{`Corefinement::Default_visitor<PolygonMesh>`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{throw_on_self_intersection}
|
||||
* \cgalParamDescription{If `true`, the set of triangles close to the intersection of `tm`
|
||||
* and `plane` will be checked for self-intersections
|
||||
* and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
|
||||
* will be thrown if at least one self-intersection is found. Always `false` if `pm` is not a triangle mesh.}
|
||||
* \cgalParamType{Boolean}
|
||||
* \cgalParamDefault{`false`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{clip_volume}
|
||||
* \cgalParamDescription{If `true`, and `tm` is closed, the clipping will be done on
|
||||
* the volume \link coref_def_subsec bounded \endlink by `tm`
|
||||
* rather than on its surface (i.e., `tm` will be kept closed).}
|
||||
* \cgalParamType{Boolean}
|
||||
* \cgalParamDefault{`false`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{use_compact_clipper}
|
||||
* \cgalParamDescription{if `false` the parts of `tm` coplanar with `plane` will not be part of the output}
|
||||
* \cgalParamType{Boolean}
|
||||
* \cgalParamDefault{`true`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{allow_self_intersections}
|
||||
* \cgalParamDescription{If `true`, self-intersections are accepted for `tm`.}
|
||||
* \cgalParamType{Boolean}
|
||||
* \cgalParamDefault{`false`}
|
||||
* \cgalParamExtra{If this option is set to `true`, `tm` is no longer required to be without self-intersection.
|
||||
* Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
|
||||
* and `clip_volume` to `false`.}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalParamNBegin{do_not_triangulate_faces}
|
||||
* \cgalParamDescription{If the input mesh is triangulated and this parameter is set to `false`, the mesh will be kept triangulated.
|
||||
* Always `true` if `pm` is not a triangle mesh.}
|
||||
* \cgalParamType{`bool`}
|
||||
* \cgalParamDefault{`true`}
|
||||
* \cgalParamNEnd
|
||||
*
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* @return `true` if the output surface mesh is manifold.
|
||||
* If `false` is returned `tm` is only refined by the intersection with `plane`.
|
||||
*
|
||||
* @see `split()`
|
||||
*/
|
||||
template <class PolygonMesh,
|
||||
class NamedParameters = parameters::Default_named_parameters>
|
||||
void new_clip(PolygonMesh& pm,
|
||||
#ifdef DOXYGEN_RUNNING
|
||||
const Plane_3& plane,
|
||||
#else
|
||||
const typename GetGeomTraits<PolygonMesh, NamedParameters>::type::Plane_3& plane,
|
||||
#endif
|
||||
const NamedParameters& np = parameters::default_values())
|
||||
{
|
||||
using parameters::choose_parameter;
|
||||
using parameters::get_parameter;
|
||||
|
||||
using face_descriptor = typename boost::graph_traits<PolygonMesh>::face_descriptor;
|
||||
using edge_descriptor = typename boost::graph_traits<PolygonMesh>::edge_descriptor;
|
||||
using halfedge_descriptor = typename boost::graph_traits<PolygonMesh>::halfedge_descriptor;
|
||||
using vertex_descriptor = typename boost::graph_traits<PolygonMesh>::vertex_descriptor;
|
||||
|
||||
using GT = typename GetGeomTraits<PolygonMesh, NamedParameters>::type;
|
||||
GT traits = choose_parameter<GT>(get_parameter(np, internal_np::geom_traits));
|
||||
auto vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
|
||||
get_property_map(vertex_point, pm));
|
||||
|
||||
typedef typename internal_np::Lookup_named_param_def <
|
||||
internal_np::concurrency_tag_t,
|
||||
NamedParameters,
|
||||
Sequential_tag
|
||||
> ::type Concurrency_tag;
|
||||
|
||||
// config flags
|
||||
bool clip_volume =
|
||||
parameters::choose_parameter(parameters::get_parameter(np, internal_np::clip_volume), false);
|
||||
bool use_compact_clipper =
|
||||
choose_parameter(get_parameter(np, internal_np::use_compact_clipper), true);
|
||||
const bool throw_on_self_intersection =
|
||||
choose_parameter(get_parameter(np, internal_np::use_compact_clipper), false);
|
||||
const bool allow_self_intersections =
|
||||
choose_parameter(get_parameter(np, internal_np::allow_self_intersections), false);
|
||||
bool triangulate = !choose_parameter(get_parameter(np, internal_np::do_not_triangulate_faces), true);
|
||||
|
||||
// TODO: generic versions
|
||||
auto vos = pm.template add_property_map<vertex_descriptor, CGAL::Oriented_side>("v:os").first;
|
||||
auto ecm = pm.template add_property_map<edge_descriptor, bool>("e:ecm", false).first;
|
||||
|
||||
if (triangulate && !is_triangle_mesh(pm))
|
||||
triangulate = false;
|
||||
|
||||
cut_with_plane(pm, plane, parameters::vertex_oriented_side_map(vos)
|
||||
.edge_is_marked_map(ecm)
|
||||
.vertex_point_map(vpm)
|
||||
.geom_traits(traits)
|
||||
.do_not_triangulate_faces(!triangulate)
|
||||
.throw_on_self_intersection(!allow_self_intersections &&
|
||||
throw_on_self_intersection)
|
||||
.concurrency_tag(Concurrency_tag()));
|
||||
|
||||
|
||||
if (allow_self_intersections)
|
||||
{
|
||||
clip_volume=false;
|
||||
triangulate=false;
|
||||
}
|
||||
|
||||
if (clip_volume && !is_closed(pm)) clip_volume=false;
|
||||
if (clip_volume && !use_compact_clipper) use_compact_clipper=true;
|
||||
|
||||
// TODO: dynamic map
|
||||
auto fcc = pm.template add_property_map<face_descriptor, std::size_t>("f:cc").first;
|
||||
|
||||
std::size_t nbcc = connected_components(pm, fcc, CGAL::parameters::edge_is_constrained_map(ecm));
|
||||
|
||||
std::vector<bool> classified(nbcc, false);
|
||||
std::vector<std::size_t> ccs_to_remove;
|
||||
|
||||
for (auto f : faces(pm))
|
||||
{
|
||||
std::size_t ccid = get(fcc, f);
|
||||
if (classified[ccid]) continue;
|
||||
halfedge_descriptor hf = halfedge(f, pm);
|
||||
for(halfedge_descriptor h : CGAL::halfedges_around_face(hf, pm))
|
||||
{
|
||||
CGAL::Oriented_side os = get(vos, target(h,pm));
|
||||
if (os==CGAL::ON_ORIENTED_BOUNDARY) continue;
|
||||
classified[ccid]=true;
|
||||
if (os==CGAL::ON_POSITIVE_SIDE) ccs_to_remove.push_back(ccid);
|
||||
}
|
||||
|
||||
if (!classified[ccid])
|
||||
{
|
||||
if (!use_compact_clipper) ccs_to_remove.push_back(ccid);
|
||||
classified[ccid]=true;
|
||||
}
|
||||
}
|
||||
|
||||
remove_connected_components(pm, ccs_to_remove, fcc);
|
||||
|
||||
if (clip_volume)
|
||||
{
|
||||
std::vector<halfedge_descriptor> borders;
|
||||
extract_boundary_cycles(pm, std::back_inserter(borders));
|
||||
|
||||
for (halfedge_descriptor h : borders)
|
||||
{
|
||||
Euler::fill_hole(h, pm);
|
||||
#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_TRIANGULATION
|
||||
if (triangulate)
|
||||
triangulate_face(face(h,pm), pm, parameters::vertex_point_map(vpm).geom_traits(traits));
|
||||
#endif
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
} } // CGAL::Polygon_mesh_processing
|
||||
|
||||
|
||||
#endif // CGAL_POLYGON_MESH_PROCESSING_CUT_WITH_PLANE_H
|
||||
|
|
@ -0,0 +1,11 @@
|
|||
OFF
|
||||
6 1 0
|
||||
|
||||
3 1 0
|
||||
3 0 0
|
||||
1 1 0
|
||||
0 0 0
|
||||
0 2 0
|
||||
2 2 0
|
||||
6 0 2 5 4 3 1
|
||||
|
||||
|
|
@ -0,0 +1,21 @@
|
|||
OFF
|
||||
16 1 0
|
||||
|
||||
1 5 0
|
||||
1 2 0
|
||||
4 1 0
|
||||
1 6 0
|
||||
1 1 0
|
||||
4 0 0
|
||||
4 4 0
|
||||
1 4 0
|
||||
1 3 0
|
||||
4 6 0
|
||||
4 7 0
|
||||
0 0 0
|
||||
4 2 0
|
||||
4 3 0
|
||||
0 7 0
|
||||
4 5 0
|
||||
16 8 7 6 15 0 3 9 10 14 11 5 2 4 1 12 13
|
||||
|
||||
|
|
@ -1,4 +1,5 @@
|
|||
#include <CGAL/Polygon_mesh_processing/clip.h>
|
||||
#include <CGAL/Polygon_mesh_processing/cut_with_plane.h>
|
||||
#include <CGAL/Polygon_mesh_processing/transform.h>
|
||||
|
||||
#include <CGAL/Surface_mesh.h>
|
||||
|
|
@ -856,6 +857,53 @@ void test_isocuboid()
|
|||
assert(vertices(meshes[0]).size() == 20);
|
||||
assert(vertices(meshes[1]).size() == 4);
|
||||
}
|
||||
|
||||
template <class TriangleMesh>
|
||||
void test_new_clip()
|
||||
{
|
||||
{
|
||||
TriangleMesh e;
|
||||
std::ifstream("data-clip/ee.off") >> e;
|
||||
PMP::cut_with_plane(e, K::Plane_3(1,0,0,-2));
|
||||
assert(faces(e).size()==5);
|
||||
assert(vertices(e).size()==24);
|
||||
}
|
||||
|
||||
{
|
||||
TriangleMesh c;
|
||||
std::ifstream("data-clip/c.off") >> c;
|
||||
PMP::cut_with_plane(c, K::Plane_3(1,0,0,-2));
|
||||
assert(faces(c).size()==2);
|
||||
assert(vertices(c).size()==8);
|
||||
}
|
||||
|
||||
{
|
||||
TriangleMesh e;
|
||||
std::ifstream("data-clip/ee.off") >> e;
|
||||
PMP::triangulate_faces(e);
|
||||
PMP::cut_with_plane(e, K::Plane_3(1,0,0,-2), CGAL::parameters::do_not_triangulate_faces(false));
|
||||
assert(faces(e).size()==30);
|
||||
assert(vertices(e).size()==28);
|
||||
}
|
||||
|
||||
{
|
||||
TriangleMesh c;
|
||||
std::ifstream("data-clip/c.off") >> c;
|
||||
PMP::triangulate_faces(c);
|
||||
PMP::cut_with_plane(c, K::Plane_3(1,0,0,-2), CGAL::parameters::do_not_triangulate_faces(false));
|
||||
assert(faces(c).size()==8);
|
||||
assert(vertices(c).size()==9);
|
||||
}
|
||||
|
||||
{
|
||||
TriangleMesh ele;
|
||||
std::ifstream(CGAL::data_file_path("meshes/elephant.off")) >> ele;
|
||||
PMP::new_clip(ele, K::Plane_3(1,0,0,0));
|
||||
PMP::new_clip(ele, K::Plane_3(0,1,0,0));
|
||||
PMP::new_clip(ele, K::Plane_3(0,0,1,0));
|
||||
}
|
||||
}
|
||||
|
||||
int main()
|
||||
{
|
||||
std::cout << "Surface Mesh" << std::endl;
|
||||
|
|
@ -878,5 +926,11 @@ int main()
|
|||
std::cout << "running test_iso_cuboid with Polyhedron\n";
|
||||
test_isocuboid<Polyhedron>();
|
||||
std::cout << "Done!" << std::endl;
|
||||
std::cout << "running test_new_clip with Surface_mesh\n";
|
||||
test_new_clip<Surface_mesh>();
|
||||
std::cout << "Done!" << std::endl;
|
||||
// std::cout << "running test_new_clip with Polyhedron\n";
|
||||
// test_new_clip<Polyhedron>();
|
||||
// std::cout << "Done!" << std::endl;
|
||||
return EXIT_SUCCESS;
|
||||
}
|
||||
|
|
|
|||
|
|
@ -20,6 +20,7 @@ CGAL_add_named_parameter(visitor_t, visitor, visitor)
|
|||
CGAL_add_named_parameter(point_t, point_map, point_map)
|
||||
|
||||
CGAL_add_named_parameter(edge_is_constrained_t, edge_is_constrained, edge_is_constrained_map)
|
||||
CGAL_add_named_parameter(edge_is_marked_map_t, edge_is_marked_map, edge_is_marked_map)
|
||||
CGAL_add_named_parameter(first_index_t, first_index, first_index)
|
||||
CGAL_add_named_parameter(number_of_iterations_t, number_of_iterations, number_of_iterations)
|
||||
CGAL_add_named_parameter(verbosity_level_t, verbosity_level, verbosity_level)
|
||||
|
|
@ -57,6 +58,8 @@ CGAL_add_named_parameter(face_color_output_iterator_t, face_color_output_iterato
|
|||
CGAL_add_named_parameter(vertex_normal_map_t, vertex_normal_map, vertex_normal_map)
|
||||
CGAL_add_named_parameter(vertex_color_map_t, vertex_color_map, vertex_color_map)
|
||||
CGAL_add_named_parameter(vertex_texture_map_t, vertex_texture_map, vertex_texture_map)
|
||||
CGAL_add_named_parameter(vertex_oriented_side_map_t, vertex_oriented_side_map, vertex_oriented_side_map)
|
||||
|
||||
CGAL_add_named_parameter(face_color_map_t, face_color_map, face_color_map)
|
||||
CGAL_add_named_parameter(repair_polygon_soup_t, repair_polygon_soup, repair_polygon_soup)
|
||||
CGAL_add_named_parameter(output_color_t, output_color, output_color)
|
||||
|
|
|
|||
Loading…
Reference in New Issue