CGAL / cgal

The public CGAL repository, see the README below
https://github.com/CGAL/cgal#readme
Other
4.88k stars 1.38k forks source link

2D Intersection of Curves #7137

Open Supranaturaler opened 1 year ago

Supranaturaler commented 1 year ago

Issue Details

CGAL provides the function "CGAL::compute_subcurves" to compute the non-intersecting sub-segments. However I want to get the composed sub-segments of each input segment. For example: the segment c1 is composed of 1, 2, 3; the segment c2 is composed of 4, 5, 6. Is there a fast way to get the result? image

Source Code

Environment

efifogel commented 1 year ago

You can use CGAL::Arr_curve_data_traits_2 to assign a distinct label to each curve. This label is carried over to the sub curves when the original curve is subdivided due to intersections. You can access the label using the data() member of the curve.

include

include <CGAL/Exact_predicates_exact_constructions_kernel.h>

include <CGAL/Arr_segment_traits_2.h>

include <CGAL/Surface_sweep_2_algorithms.h>

include <CGAL/Arr_curve_data_traits_2.h>

using Kernel = CGAL::Exact_predicates_exact_constructions_kernel; using Base_traits_2 = CGAL::Arr_segment_traits_2; using Base_x_monotone_curve_2 = Base_traits_2::X_monotone_curve_2; using Traits_2 = CGAL::Arr_curve_data_traits_2<Base_traits_2, std::size_t>; using Point_2 = Traits_2::Point_2; using X_monotone_curve_2 = Traits_2::X_monotone_curve_2;

int main() { Traits_2 traits; size_t id(0); Base_x_monotone_curve_2 base_segs[] = { Base_x_monotone_curve_2(Point_2(1, 5), Point_2(8, 5)), Base_x_monotone_curve_2(Point_2(1, 1), Point_2(8, 8)), Base_x_monotone_curve_2(Point_2(3, 1), Point_2(3, 8)), Base_x_monotone_curve_2(Point_2(8, 5), Point_2(8, 8)), }; auto size = sizeof(base_segs) / sizeof(Base_x_monotone_curve_2); X_monotone_curve_2 segs[size]; for (std::size_t i = 0; i < size; ++i) segs[i] = X_monotone_curve_2(base_segs[i], i); std::list sub_segs; CGAL::compute_subcurves(segs, segs + 4, std::back_inserter(sub_segs), false, traits); for (const auto& seg : sub_segs) std::cout << seg << "," << seg.data() << std::endl; return 0; }


//) o /__ // (__ ( ( ( (/ (/-(-'(/ /

On Thu, 22 Dec 2022 at 04:12, Supranaturaler @.***> wrote:

Issue Details

CGAL provides the function "CGAL::compute_subcurves" to compute the non-intersecting sub-segments. However I want to get the composed sub-segments of each input segment. For example: the segment c1 is composed of 1, 2, 3; the segment c2 is composed of 4, 5, 6. Is there a fast way to get the result? [image: image] https://user-images.githubusercontent.com/40469886/209039352-3d302a47-102b-4f3f-a713-24e20fc66fdc.png Source Code Environment

  • Operating system (Windows/Mac/Linux, 32/64 bits):
  • Compiler:
  • Release or debug mode:
  • Specific flags used (if any):
  • CGAL version:
  • Boost version:
  • Other libraries versions if used (Eigen, TBB, etc.):

— Reply to this email directly, view it on GitHub https://github.com/CGAL/cgal/issues/7137, or unsubscribe https://github.com/notifications/unsubscribe-auth/ABVBNOFNRBYRONRCYOILODLWOO2J7ANCNFSM6AAAAAATGHEI4Y . You are receiving this because you are subscribed to this thread.Message ID: @.***>