Subject: CGAL users discussion list
List archive
- From: Efi Fogel <>
- To:
- Subject: Re: [cgal-discuss] Obtain a specific part of the intersection of 2D polygons
- Date: Thu, 10 Feb 2022 02:16:44 +0200
- Authentication-results: mail2-smtp-roc.national.inria.fr; spf=None ; spf=Pass ; spf=None
- Ironport-data: A9a23:H4plHqNDbalXVSLvrR2OkcFynXyQoLVcMsEvi/4bfWQNrUpwhDIFm 2MWXmuFbquCa2b3eognb42w9h5VuZHcnNA2SnM5pCpnJ55ogZqcVI7Bdi8cHAvLc5adFBo/h yk6QoOdRCzhZiaE/n9BCpC48T8kk/jgqoPUUIYoAAgoLeNfYHpn2UILd9IR2NYy24DjW1LV4 rsenuWGULOb824sWo4rw/nbwP9flKyaVOQw4zTSzdgS1LPvvyF94KA3fcldHFOkKmVgJdNWc s6YpF2PEsw1yD92Yj+tuu6TnkTn2dc+NyDW4pZdc/DKbhSvOkXe345jXMfwZ3u7hB2Ovohfy 8hGq6CpSBobJaDQqqMWCjJxRnQW0a1uoNcrIFC6uM2XilLcKj7in6woA0YxMokVvO1wBAmi9 9RCcGFLPk3F3bjvhu7iIgVvrpxLwM3DMY0etHZvwDXxAvMvQJSFSKLPjTNd9G1o25AVRq2OD yYfQSJoZlPkQCVhAU4GMdFvkduJiSLYLyIN/Tp5ooJuuzSJpOBr65DmP9PRP9CLXs5IhV2wv XPD522/AxcANdXZxyDtz563rurGnCe+R5hLUbPkr7hlh1qcwmFVAxoTPbemnRWnokmfA/xcN kcKxjMvsJop+VCZCdKkdDTt9RZooSUgc9ZXFuQ77iSExazV/xuVCwA4othpOIxOWCgeFW1C6 7OZoz/6LWcw7+DNGBpx4p/R/GziY3FERYMXTXZcFVNt3jX1nG0kYvvyojtLFae0ipjqA2i1z WzT6ic5gLoXgIgA0KDTEbH7b9CE9sehou0dvF2/soeZAuVROtLNi2uAtwOz0Bq4BNzFJmRtR VBd8yRk0MgADIuWiAuGS/gXEbei6p6taWOA3QI0RsV5r2jyqxZPmLy8BhkudC+F1e5UKFfUj LP75Gu9GbcPZiLxNfYvC25PI55ykfO6fTgaahwkRoMWPsIZmP6v8yZpakqdt10BY2B9+ZzTz ayzKJ72ZV5DUfoP5GPvG481jOF2rghjmju7bc2ql3yPjOvFDFbIGOdtGAXVNYgRsvjUyDg5B v4FaKNmPT0EALOgCsQWmKZPRW03wY8TW8yo9ZwNK7Lrz8gPMDhJNsI9CIgJI+RN95m5XM+Rl p1kckMHmlf5m1PdLgCGNiJqZL/1DMRwqHs6OWonOlPxgyovZoOm7aE+cZorfOl/pLYzk6IsF /RVKd+dBvlvSyjc/2tPYJT4qrtkfkv5iA+LOR2jfzViLYVrQBbE+4O/cwa2rHsOAyO7uNEQu bql0g+HE5MPSx4zXsnTYfOriVi2uCFFyu51WkLJJPhVeVntoNA6cXyv0qdvLphVexvZxzac2 wKHOjsipLHA890v7d3EpaGYtIP2QeZzG0xtGWOEv7u7MC/t+HX6nd1NXeOOSjDqVG3u/ZKka +gIner3N+cKnQoTvodxT+RrwKY564e9rrNW1F4/TnDCblDuFak5Z3fahI9AsapCwrIfsgyzA xrd9t5fMLSPGcXkDF9Be1Z/P7rbjakZymvI8PA4AETm/ysrrrCJZkNfYkuXgytHIbopbY4on bU7tMgN51DtgxYmKIzd3CVd9mDJMWZZFqt779cVB4jkjgdtwVZHOMSOBijz6ZCJStNNLkh6f WPO1fSa3+xRlhjYbn4+NXnRxu4B154AjxZHkQ0ZLFOTl9uZ2/I60XW9K9jsoti5E/mG7w5yB oSvH0h8JKHL4C0xwcYfBiajHAZOABDf8Uv0o7fMeKs1UGHwPlEh7kVkUQpOwKzd221Zdzlfu rqfzQ4JlB70Kdrp0HJatVFN8pTeoB8YyuEGsM+iFsWBWZI9ZFIJR0NoiXUg83PaPC/6uKELS SSGMgq9hW0X+BP8e5EGNrQ=
- Ironport-hdrordr: A9a23:e8lMEK0EjZ4pO9RlbrlF8QqjBIskLtp133Aq2lEZdPU1SL3gqy nKpp4mPHDP+VMssR0b6LK90ey7MBDhHP1OgLX5X43SODUO0VHAROpfBMnZowEIcBeOkdK1u5 0QFZSWy+edMbG5t6vHCcWDfOrICePozJyV
- Ironport-phdr: A9a23:+yOxBBNZyQ5Aho4ztW0l6na/BxdPi9zP1u491JMrhvp0f7i5+Ny6Z QqDv68r0Q6CBNmDo9t/yMPu+5j6XmIB5ZvT+FsjS7drEyE/tMMNggY7C9SEA0CoZNTjbig9A dgQHAQ9pyLzPkdaAtvxaEPPqXOu8zESBg//NQ1oLejpB4Lelcu62/6s95HJYwhFgDWxba59I RmqsA7cqtQYjYx+J6gr1xDHuGFIe+NYxWNpIVKcgRPx7dqu8ZBg7ipdpesv+9ZPXqvmcas4S 6dYDCk9PGAu+MLrrxjDQhCR6XYaT24bjwBHAwnB7BH9Q5fxri73vfdz1SWGIcH7S60/VDK/5 KlpVRDokj8KODE38G7VisJ+gqFVrg+/qRNj2IPbep2ZOeBkc6/BYd8XR2xMVdtRWSxbBYO8a pMCAeUPPeZZsoLzp1wOrRSgCgmoGejizSFHhnH33a001OQhHh/J3Ag7EtIBtXTbttT1NKMIX e+py6nIyCzOYvVL0jjy9IbGaAouoe2QXb1ua8rRz1EiGgHbglmOqYHoIzyb2+sCvmaV8eZtV +yhhW4kpgx1rTWj28chh4jJiIwVy13K9Tl1zokoKdC4SEN1bsKpHYdfuiycKoB4TMQiQ2Ryt yY7zL0LoZ+7fC4QyJQm3RHTcfKHc5KO7xn+V+iROS91iGx5dL+7nRq/8kitxvfiWsWqzVpGt CVInsTKu3sQzRLc8NKHReF4/kq52TaAyQTT6uZcLEAxj6XbKpohzqc+l5oJrEjPByH2lUrrg KOMeUUk/e+o6+vjYrr4vJOTK4h0igTmPqQvnMywH/g4PxAQU2SH/emwzr7u8E3jTLlUkPE6j 7PVvZ/HKcgDo662GQ5V0oIt6xalCDem1cwVkmcJLFJEdhKHiIfpNE/KIP3jAve/hk6jkDZvx /zcIrLhBZDNImDFkLj6Zbl98VJTyBIvzdBD4JJZEqwNLOrpWkDtrNzYEgM5Mwuszun7B9Vyz IceVXuSDa+YK6PdrUKI5vk0I+SXf48UuDP9K+A/6PL0jH85n0Udfaiz0pcNZnC4BKcuHkOCf HC5gssdCXxY+U0lXenygRuDVyRSbjC8ReUn9zQjAcWnC4nEAYuiibjE0CagFYBNfTN7DEuRG 1f0coHRW+sQcDnAZYh6gzkcXP6gTZUg3Fegrkjh2r9/J63V/CMf8pns3dww6+zIngwp7m9JC d+A2V2AX30hnn8UXyRkm+dksEllwxGC17J5irpWD5tI9vZRW0A7M5DbiOd1AtS3VgPadcqSU wWaRYCtDjg1C94w2NQTeF1VGtO4jxmF0TD5LaUSkummCp183KXT0ny5c8tzynjB26Qlp1YjS 8pLc2ahg/gspEDoG4fVnhDBxO6RfqMG0XuVnI/i5W+HvUUCFRV1Tb2AR3cUIE3fsdX+4ErGC b6oE7UudAVbmoaZMqUfTNrvgB1dQev7fszEaje8lWa+AhmFwpuDaYPrfyMW2yCOQFMcnVUr9 G2dfRM7Gj/npmvfCDJ0Ele6bk3t/+5xpXeTQUo9zgXMZEpkhPKu4hBAo/uaRrsI264c/icsr zIhBFGmw9ffEMaNvSJkdaRYJMwnuRJJiTufuAt6MZitaatlgzbyaixRuEXjn1VyA4REy40xq W8yiRB1IuSe2U9AcDWR2dbxPKfWIy/85kLnbamewVzY3NuMn8VHoP0lt1Xuuh2oHUs+4j1m1 ddSyX6V+pTNCkIbT5vwVk898xUyqavdZ2Ex4Ibd1HskNqfR0HeK1twoCuwqxxKIcNJWMafCH wj3UoUbC8WoNO02ig2xdBtXdOtW9aMyI4anb67cgP/tbLsmxmj/yzges+UfmgqW+iFxS/DFx cMAyvCchE6cUivkyU2muYbxkJxFYjcbGiy+zzLlDchffP4XH85DBGGwLsmw3tg7iYTqXisS+ VCiCVQJ1cuBdh+bbli71gpVnxdywzTviW6jwjp4nit85KeR3SLJzOnmXBUCM29PAmJliB2/a ZjxhNccUk+yagEvnxbw/kf2yZ9Qo6FnJnXSS0NFF8TvB1lrSbD49r+LYsoVrYgtrT0SSuO3J 1aTVr/6pRIelSLlBWpXgj4hJXmmvZDwnhoyj2z4Tj47oXTYd8Z/yBP369nVRPoX1T0DDCV1k jjYAFGgMsLhp43F0cee9LrkDiT9D9VaamHzwJmFtTen6GEPY1X3hP21ltD9UEA73SL9y9h2R HDNpRf4bJPs0vfyOuZmc090QV7kvpAiS8cuz81q3shWhSFJ4/fdtWAKmmryL9hBjKf3bX5XA CUO38aQ+g/unktqMnOOwYv9EHSb2MpoIdegMQZ0kmow6d5HDKCM4flKhyxw9xC1owPfZvdwm h8SzPIv7DgRhORD629Phm2NR6sfG0VVJ3mmmhqJ4dexoaF/a2OmcLz2301714PpHPSJpQdSX 2z8c5EpEHpr78lxB1nL1WX69oDufNSDCLBb/g3RiRrLiPJZbY4gjvdfzzQyInrz5Dd2g/5+l xFl2ou2+ZSKO3k4trzsGQZWb1iXL4sS4m2/1vsYx5fOmdrzQdM5XW9XFJrwEaD2THRI7q+hb ljWVmV78ybTGKKDT1HBrh4+9TSXVcjsbSnyRjFRzM0+FkfDYhYD0UZEBHNi2cRhXgGymJ6+K gEgunZItwS+8lwVmqppL0WtDTuZ/V31LG9yEN/Gcn80pklD/xuHaJTOqLssQGcIuMXm9lXFK 3THNV0XVidQCxDCXxa7eeDwrdjYr7rCXrv4f6aIOO/e77QZDqjtp9rn05M6rWzVaIPSbj87X qd9gg0aDDh4A5iLwWxRDXZMxmSWNYjD407tsjt+qsT1mBjycCTo44bHS75bMNE1vguznb/GL embwiBwNTdf0JoIg37O0rkWmlAI2WlocHG2HLIMuDSoLuqYk7JLDxMddyJ4NddZp6M60A5XP MfHi9Tznrdmh/8xAl1BWBTvgMasLcANJmi8MhvACiPpfPyeIibXxsjsfa6mYbhZjeERqAfp/ DjHSgnsOTOMkzSvXBeqcKlNgCydIB1CqdS9fxJqWg2BBJrtbhy2NsMyjCVjm+VlwCOXczRGa H4gKhAozPXY9y5Tj/RhFnYU63NkKbPBgCOF96zCLZ1Qt/J3AyNynuYc4XIgyrIT4jsXIZ490 CbUsNNqpEmr1+eVzT8yGhdAqzhMi4+PlUpnMKTdsJJHXDyXmXBFpXXVEBkMq9Z/X5f3vLtMz 9HUiK/pADJL8taR4tFFQsaIcYSIN30uNRevEznRRlhgL3bjJSTUgEpTl+uX/3ueo80hq5Tir 5EJT6dSSF0/Ev5y4qVNGdUDJNJoRGphn+LCysEP4nW6oV/aQ8AI5vgvudqdBPzuLHCSirwWP nPgJJv3KI0SMsvw3EkwMjFH
- Ironport-sdr: XOUmd7RFAdiSHJYSF/S8yI2+zxqTXj3sr2p0pQa3UlIWWpI4xvBK6FBYQyINWSUkZocR3vV7/f PWyQkGAqELLmuxPso+91VuXQBOqUggpzNcop18l4EFk7OJ8P/9jfz9ZfoZp40dB6UgtzagDIcB C5ybqxjh63WVbgk/gL53pDoJirLcM0zUgewpd2KkH8fKN0DTVFzxwuqX4BfgR1S/jaqY0mqcAu ys0l8ZNFr/+P9/4hGWzVWffI5Wciz0iTP3mTbUqcdWobWoV61gyuJ+I3LBsFHymlJ5WXybpjj4 lgZDxWcvuge2QdD3NyuuTkKA
Hi Bob,
The Boolean operations are based on 2D arrangements (from the "2D Arrangements" package).
In general you have more options when you work directly with arrangements.
From your description I understand that you would like to obtain the extreme points on the boundary of the intersection that lie on the boundary of p2.
The code below does it.
Best,
Efi
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
#include <CGAL/Arr_segment_traits_2.h>
#include <CGAL/Arrangement_2.h>
#include <CGAL/Arr_extended_dcel.h>
#include <CGAL/Arr_observer.h>
#include <CGAL/Arr_overlay_2.h>
typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
typedef CGAL::Arr_segment_traits_2<Kernel> Traits;
typedef Traits::Point_2 Point;
typedef Traits::X_monotone_curve_2 X_monotone_curve;
typedef CGAL::Arr_extended_dcel<Traits, bool, bool, bool> Dcel;
typedef CGAL::Arrangement_2<Traits, Dcel> Arrangement;
template <typename Arrangement>
struct My_observer : CGAL::Arr_observer<Arrangement> {
typedef typename Arrangement::Vertex_handle V_handle;
typedef typename Arrangement::Halfedge_handle H_handle;
typedef typename Arrangement::Face_handle F_handle;
My_observer(Arrangement& arr, bool flag) :
CGAL::Arr_observer<Arrangement>(arr),
m_flag(flag)
{ arr.unbounded_face()->set_data(false); }
void after_create_vertex(V_handle v) { v->set_data(m_flag); }
void after_create_edge(H_handle h) {
h->set_data(m_flag);
h->twin()->set_data(m_flag);
}
void after_split_face(F_handle f1, F_handle f2, bool is_hole)
{ f2->set_data(true); }
private:
bool m_flag;
};
template <typename Arrangement>
struct Overlay_traits {
typedef typename Arrangement::Vertex_const_handle V_const_handle;
typedef typename Arrangement::Halfedge_const_handle H_const_handle;
typedef typename Arrangement::Face_const_handle F_const_handle;
typedef typename Arrangement::Vertex_handle V_handle;
typedef typename Arrangement::Halfedge_handle H_handle;
typedef typename Arrangement::Face_handle F_handle;
void create_face(F_const_handle f1, F_const_handle f2, F_handle f) const {
f->set_data(f1->data() && f2->data());
}
void create_vertex(H_const_handle h1, H_const_handle h2, V_handle v) const {
v->set_data(h2->data());
}
void create_vertex(V_const_handle v1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(V_const_handle v1, H_const_handle h2, V_handle v) const {
v->set_data(h2->data());
}
void create_vertex(H_const_handle h1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(F_const_handle f1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(V_const_handle v1, F_const_handle f2, V_handle v) const {}
void create_edge(H_const_handle h1, H_const_handle h2, H_handle h) const {}
void create_edge(H_const_handle h1, F_const_handle f2, H_handle h) const {}
void create_edge(F_const_handle f1, H_const_handle h2, H_handle h) const {}
};
int main(int argc, char* argv[]) {
const std::vector<Point> pts1{{-1, -1},{1, -1},{1, 1},{-1, 1}};
const std::vector<X_monotone_curve> cvs1 = {
{pts1[0], pts1[1]}, {pts1[1], pts1[2]},
{pts1[2], pts1[3]}, {pts1[3], pts1[0]},
};
Arrangement arr1;
My_observer<Arrangement> obs1(arr1, false);
CGAL::insert(arr1, cvs1.begin(), cvs1.end());
const std::vector<Point> pts2{{0, 1},{0.1, 0.2},{0.6, 0},{2, 0},{2, 1},{0, 2}};
const std::vector<X_monotone_curve> cvs2 = {
{pts2[0], pts2[1]}, {pts2[1], pts2[2]}, {pts2[2], pts2[3]},
{pts2[3], pts2[4]}, {pts2[4], pts2[5]}, {pts2[5], pts2[0]}
};
Arrangement arr2;
My_observer<Arrangement> obs2(arr2, true);
CGAL::insert(arr2, cvs2.begin(), cvs2.end());
Overlay_traits<Arrangement> ot;
Arrangement arr;
CGAL::overlay(arr1, arr2, arr, ot);
for (auto it = arr.faces_begin(); it != arr.faces_end(); ++it) {
if (! it->data()) continue;
auto curr = it->outer_ccb();
do {
if (curr->target()->data())
std::cout << curr->target()->point() << std::endl;
} while (++curr != it->outer_ccb());
}
return 0;
}
#include <CGAL/Arr_segment_traits_2.h>
#include <CGAL/Arrangement_2.h>
#include <CGAL/Arr_extended_dcel.h>
#include <CGAL/Arr_observer.h>
#include <CGAL/Arr_overlay_2.h>
typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
typedef CGAL::Arr_segment_traits_2<Kernel> Traits;
typedef Traits::Point_2 Point;
typedef Traits::X_monotone_curve_2 X_monotone_curve;
typedef CGAL::Arr_extended_dcel<Traits, bool, bool, bool> Dcel;
typedef CGAL::Arrangement_2<Traits, Dcel> Arrangement;
template <typename Arrangement>
struct My_observer : CGAL::Arr_observer<Arrangement> {
typedef typename Arrangement::Vertex_handle V_handle;
typedef typename Arrangement::Halfedge_handle H_handle;
typedef typename Arrangement::Face_handle F_handle;
My_observer(Arrangement& arr, bool flag) :
CGAL::Arr_observer<Arrangement>(arr),
m_flag(flag)
{ arr.unbounded_face()->set_data(false); }
void after_create_vertex(V_handle v) { v->set_data(m_flag); }
void after_create_edge(H_handle h) {
h->set_data(m_flag);
h->twin()->set_data(m_flag);
}
void after_split_face(F_handle f1, F_handle f2, bool is_hole)
{ f2->set_data(true); }
private:
bool m_flag;
};
template <typename Arrangement>
struct Overlay_traits {
typedef typename Arrangement::Vertex_const_handle V_const_handle;
typedef typename Arrangement::Halfedge_const_handle H_const_handle;
typedef typename Arrangement::Face_const_handle F_const_handle;
typedef typename Arrangement::Vertex_handle V_handle;
typedef typename Arrangement::Halfedge_handle H_handle;
typedef typename Arrangement::Face_handle F_handle;
void create_face(F_const_handle f1, F_const_handle f2, F_handle f) const {
f->set_data(f1->data() && f2->data());
}
void create_vertex(H_const_handle h1, H_const_handle h2, V_handle v) const {
v->set_data(h2->data());
}
void create_vertex(V_const_handle v1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(V_const_handle v1, H_const_handle h2, V_handle v) const {
v->set_data(h2->data());
}
void create_vertex(H_const_handle h1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(F_const_handle f1, V_const_handle v2, V_handle v) const {
v->set_data(v2->data());
}
void create_vertex(V_const_handle v1, F_const_handle f2, V_handle v) const {}
void create_edge(H_const_handle h1, H_const_handle h2, H_handle h) const {}
void create_edge(H_const_handle h1, F_const_handle f2, H_handle h) const {}
void create_edge(F_const_handle f1, H_const_handle h2, H_handle h) const {}
};
int main(int argc, char* argv[]) {
const std::vector<Point> pts1{{-1, -1},{1, -1},{1, 1},{-1, 1}};
const std::vector<X_monotone_curve> cvs1 = {
{pts1[0], pts1[1]}, {pts1[1], pts1[2]},
{pts1[2], pts1[3]}, {pts1[3], pts1[0]},
};
Arrangement arr1;
My_observer<Arrangement> obs1(arr1, false);
CGAL::insert(arr1, cvs1.begin(), cvs1.end());
const std::vector<Point> pts2{{0, 1},{0.1, 0.2},{0.6, 0},{2, 0},{2, 1},{0, 2}};
const std::vector<X_monotone_curve> cvs2 = {
{pts2[0], pts2[1]}, {pts2[1], pts2[2]}, {pts2[2], pts2[3]},
{pts2[3], pts2[4]}, {pts2[4], pts2[5]}, {pts2[5], pts2[0]}
};
Arrangement arr2;
My_observer<Arrangement> obs2(arr2, true);
CGAL::insert(arr2, cvs2.begin(), cvs2.end());
Overlay_traits<Arrangement> ot;
Arrangement arr;
CGAL::overlay(arr1, arr2, arr, ot);
for (auto it = arr.faces_begin(); it != arr.faces_end(); ++it) {
if (! it->data()) continue;
auto curr = it->outer_ccb();
do {
if (curr->target()->data())
std::cout << curr->target()->point() << std::endl;
} while (++curr != it->outer_ccb());
}
return 0;
}
____ _ ____ _
/_____/_) o /__________ __ //
(____ ( ( ( (_/ (_/-(-'_(/
_/
/_____/_) o /__________ __ //
(____ ( ( ( (_/ (_/-(-'_(/
_/
On Wed, 9 Feb 2022 at 19:11, Bob Bill <> wrote:
Dear all,let's say I have two polygons as in the attached image: P1 (black) and P2 (blue). After I obtained their intersection with CGAL::intersection, I can identify the outer boundary of the intersection.Unfortunately, I do not need the whole outer boundary, but only the vertices of the resulting intersected polygon that are coming from P2.In the picture they are precisely the ones described by the points (counterclockwise):(0 , 2)(0.1 , 0.2)(0.6 , 0)(1 , 0)How can I get those points only? The following snippet provides the whole intersection in the standard way, but then I don't know how to move from there.Best regards,Bob#include <iostream>#include <vector>#include <CGAL/Boolean_set_operations_2.h>#include <CGAL/Constrained_Delaunay_triangulation_2.h>#include <CGAL/Delaunay_mesh_face_base_2.h>#include <CGAL/Delaunay_mesh_size_criteria_2.h>#include <CGAL/Delaunay_mesher_2.h>#include <CGAL/Exact_predicates_exact_constructions_kernel.h>#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>#include <CGAL/IO/write_VTU.h>#include <CGAL/Polygon_2.h>#include <CGAL/Polygon_with_holes_2.h>#include <CGAL/Triangulation_2.h>#include <CGAL/draw_polygon_2.h>#include <CGAL/draw_polygon_with_holes_2.h>typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;typedef CGAL::Polygon_2<Kernel> Polygon;typedef CGAL::Polygon_with_holes_2<Kernel> Polygon_with_holes_2;typedef Polygon::Point_2 Point;typedef Polygon::Segment_2 Segment;typedef CGAL::Exact_predicates_inexact_constructions_kernel K;typedef CGAL::Triangulation_vertex_base_2<K> Vb;typedef CGAL::Delaunay_mesh_face_base_2<K> Fb;typedef CGAL::Triangulation_data_structure_2<Vb, Fb> Tds;typedef CGAL::Constrained_Delaunay_triangulation_2<K, Tds> CDT;typedef CGAL::Delaunay_mesh_size_criteria_2<CDT> Criteria;typedef CDT::Vertex_handle Vertex_handle;int main() {// define polygonconst std::vector<Point> pts{{-1, -1}, {1, -1}, {1, 1}, {-1, 1}};const Polygon p(pts.cbegin(), pts.cend());const std::vector<Point> pts2{{0, 1}, {0.1,0.2},{0.6, 0}, {2, 1}, {0, 2}};const Polygon p2(pts2.cbegin(), pts2.cend());std::vector<Polygon_with_holes_2> poly_list;auto oi = CGAL::intersection(p, p2, std::back_inserter(poly_list));Polygon intersect = poly_list[0].outer_boundary();return 0;}
--
You are currently subscribed to cgal-discuss.
To unsubscribe or access the archives, go to
https://sympa.inria.fr/sympa/info/cgal-discuss
- [cgal-discuss] Obtain a specific part of the intersection of 2D polygons, Bob Bill, 02/09/2022
- [cgal-discuss] having problems with Python-based installation on Windows in both vcpkg and CMake environments, Audette, Michel A., 02/09/2022
- Re: [cgal-discuss] Obtain a specific part of the intersection of 2D polygons, Efi Fogel, 02/10/2022
Archive powered by MHonArc 2.6.19+.