33#include < CGAL/Polygon_2.h>
44#include < CGAL/create_straight_skeleton_2.h>
55#include < CGAL/create_straight_skeleton_from_polygon_with_holes_2.h>
6+ #include < CGAL/create_offset_polygons_2.h>
7+ #include < CGAL/create_weighted_offset_polygons_from_polygon_with_holes_2.h>
8+ #include < CGAL/create_weighted_straight_skeleton_2.h>
69
710typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
811typedef K::Point_2 Point;
@@ -12,6 +15,8 @@ typedef CGAL::Straight_skeleton_2<K> Ss;
1215typedef boost::shared_ptr<Ss> SsPtr;
1316typedef CGAL::Straight_skeleton_2<K>::Halfedge_const_handle Halfedge_const_handle;
1417typedef CGAL::Straight_skeleton_2<K>::Vertex_const_handle Vertex_const_handle;
18+ typedef boost::shared_ptr<Polygon_2> PolygonPtr;
19+ typedef std::vector<PolygonPtr> PolygonPtrVector;
1520
1621compas::Edges pmp_create_interior_straight_skeleton (
1722 Eigen::Ref<const compas::RowMatrixXd> &V)
@@ -86,6 +91,108 @@ compas::Edges pmp_create_interior_straight_skeleton_with_holes(
8691
8792}
8893
94+ std::vector<compas::RowMatrixXd> pmp_create_offset_polygons_2_inner (Eigen::Ref<const compas::RowMatrixXd> &V, double &offset){
95+ Polygon_2 poly;
96+ for (int i = 0 ; i < V.rows (); i++){
97+ poly.push_back (Point (V (i, 0 ), V (i, 1 )));
98+ }
99+ PolygonPtrVector offset_polygons = CGAL::create_interior_skeleton_and_offset_polygons_2 (offset, poly);
100+
101+ std::vector<compas::RowMatrixXd> result;
102+ for (auto pi = offset_polygons.begin (); pi != offset_polygons.end (); ++pi){
103+ std::size_t n = (*pi)->size ();
104+ compas::RowMatrixXd points (n, 3 );
105+ int j = 0 ;
106+ for (auto vi = (*pi)->vertices_begin (); vi != (*pi)->vertices_end (); ++vi){
107+ points (j, 0 ) = (double )(*vi).x ();
108+ points (j, 1 ) = (double )(*vi).y ();
109+ points (j, 2 ) = 0 ;
110+ j++;
111+ }
112+ result.push_back (points);
113+ }
114+ return result;
115+ }
116+
117+ std::vector<compas::RowMatrixXd> pmp_create_offset_polygons_2_outer (Eigen::Ref<const compas::RowMatrixXd> &V, double &offset){
118+ Polygon_2 poly;
119+ for (int i = 0 ; i < V.rows (); i++){
120+ poly.push_back (Point (V (i, 0 ), V (i, 1 )));
121+ }
122+ PolygonPtrVector offset_polygons = CGAL::create_exterior_skeleton_and_offset_polygons_2 (offset, poly);
123+
124+ std::vector<compas::RowMatrixXd> result;
125+ for (auto pi = offset_polygons.begin (); pi != offset_polygons.end (); ++pi){
126+ std::size_t n = (*pi)->size ();
127+ compas::RowMatrixXd points (n, 3 );
128+ int j = 0 ;
129+ for (auto vi = (*pi)->vertices_begin (); vi != (*pi)->vertices_end (); ++vi){
130+ points (j, 0 ) = (double )(*vi).x ();
131+ points (j, 1 ) = (double )(*vi).y ();
132+ points (j, 2 ) = 0 ;
133+ j++;
134+ }
135+ result.push_back (points);
136+ }
137+ return result;
138+ }
139+
140+ std::vector<compas::RowMatrixXd> pmp_create_weighted_offset_polygons_2_inner (Eigen::Ref<const compas::RowMatrixXd> &V, double &offset, Eigen::Ref<const compas::RowMatrixXd> &weights){
141+ Polygon_2 poly;
142+ for (int i = 0 ; i < V.rows (); i++){
143+ poly.push_back (Point (V (i, 0 ), V (i, 1 )));
144+ }
145+ std::vector<double > weights_vec;
146+ for (int i = 0 ; i < weights.rows (); i++){
147+ weights_vec.push_back (weights (i, 0 ));
148+ }
149+ SsPtr iss = CGAL::create_interior_weighted_straight_skeleton_2 (poly, weights_vec);
150+ PolygonPtrVector offset_polygons = CGAL::create_offset_polygons_2<Polygon_2>(offset, *iss);
151+
152+ std::vector<compas::RowMatrixXd> result;
153+ for (auto pi = offset_polygons.begin (); pi != offset_polygons.end (); ++pi){
154+ std::size_t n = (*pi)->size ();
155+ compas::RowMatrixXd points (n, 3 );
156+ int j = 0 ;
157+ for (auto vi = (*pi)->vertices_begin (); vi != (*pi)->vertices_end (); ++vi){
158+ points (j, 0 ) = (double )(*vi).x ();
159+ points (j, 1 ) = (double )(*vi).y ();
160+ points (j, 2 ) = 0 ;
161+ j++;
162+ }
163+ result.push_back (points);
164+ }
165+ return result;
166+ }
167+
168+ std::vector<compas::RowMatrixXd> pmp_create_weighted_offset_polygons_2_outer (Eigen::Ref<const compas::RowMatrixXd> &V, double &offset, Eigen::Ref<const compas::RowMatrixXd> &weights){
169+ Polygon_2 poly;
170+ for (int i = 0 ; i < V.rows (); i++){
171+ poly.push_back (Point (V (i, 0 ), V (i, 1 )));
172+ }
173+ std::vector<double > weights_vec;
174+ for (int i = 0 ; i < weights.rows (); i++){
175+ weights_vec.push_back (weights (i, 0 ));
176+ }
177+ SsPtr iss = CGAL::create_exterior_weighted_straight_skeleton_2 (offset, weights_vec, poly);
178+ PolygonPtrVector offset_polygons = CGAL::create_offset_polygons_2<Polygon_2>(offset, *iss);
179+
180+ std::vector<compas::RowMatrixXd> result;
181+ for (auto pi = offset_polygons.begin (); pi != offset_polygons.end (); ++pi){
182+ std::size_t n = (*pi)->size ();
183+ compas::RowMatrixXd points (n, 3 );
184+ int j = 0 ;
185+ for (auto vi = (*pi)->vertices_begin (); vi != (*pi)->vertices_end (); ++vi){
186+ points (j, 0 ) = (double )(*vi).x ();
187+ points (j, 1 ) = (double )(*vi).y ();
188+ points (j, 2 ) = 0 ;
189+ j++;
190+ }
191+ result.push_back (points);
192+ }
193+ return result;
194+ }
195+
89196// ===========================================================================
90197// PyBind11
91198// ===========================================================================
@@ -104,4 +211,31 @@ void init_straight_skeleton_2(pybind11::module &m)
104211 &pmp_create_interior_straight_skeleton_with_holes,
105212 pybind11::arg (" V" ).noconvert (),
106213 pybind11::arg (" holes" ).noconvert ());
214+
215+ submodule.def (
216+ " create_offset_polygons_2_inner" ,
217+ &pmp_create_offset_polygons_2_inner,
218+ pybind11::arg (" V" ).noconvert (),
219+ pybind11::arg (" offset" ).noconvert ());
220+
221+ submodule.def (
222+ " create_offset_polygons_2_outer" ,
223+ &pmp_create_offset_polygons_2_outer,
224+ pybind11::arg (" V" ).noconvert (),
225+ pybind11::arg (" offset" ).noconvert ());
226+
227+ submodule.def (
228+ " create_weighted_offset_polygons_2_inner" ,
229+ &pmp_create_weighted_offset_polygons_2_inner,
230+ pybind11::arg (" V" ).noconvert (),
231+ pybind11::arg (" offset" ).noconvert (),
232+ pybind11::arg (" weights" ).noconvert ());
233+
234+ submodule.def (
235+ " create_weighted_offset_polygons_2_outer" ,
236+ &pmp_create_weighted_offset_polygons_2_outer,
237+ pybind11::arg (" V" ).noconvert (),
238+ pybind11::arg (" offset" ).noconvert (),
239+ pybind11::arg (" weights" ).noconvert ());
240+
107241};
0 commit comments