From 47c7cbda17182565a2e98c3132af9969d8e5d287 Mon Sep 17 00:00:00 2001 From: lbartoletti Date: Tue, 10 Oct 2023 11:44:26 +0000 Subject: [PATCH] Extrude straight skeleton --- CMakeLists.txt | 2 +- NEWS | 1 + src/PolyhedralSurface.cpp | 18 + src/PolyhedralSurface.h | 8 + src/algorithm/extrude.cpp | 28 +- src/algorithm/extrude.h | 2 + src/algorithm/straightSkeleton.cpp | 49 ++ src/algorithm/straightSkeleton.h | 22 +- src/capi/sfcgal_c.cpp | 57 ++- src/capi/sfcgal_c.h | 20 +- .../SFCGAL/algorithm/StraightSkeletonTest.cpp | 476 +++++++++++++----- 11 files changed, 537 insertions(+), 146 deletions(-) diff --git a/CMakeLists.txt b/CMakeLists.txt index de737bc0..10e2e9db 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -68,7 +68,7 @@ if( ${CGAL_USE_AUTOLINK} ) add_definitions( "-DCGAL_NO_AUTOLINK" ) endif() -find_package( CGAL 5.3 COMPONENTS Core REQUIRED ) +find_package( CGAL 5.6 COMPONENTS Core REQUIRED ) message( STATUS "CGAL ${CGAL_VERSION} found" ) include_directories( ${CMAKE_BINARY_DIR}/include ) diff --git a/NEWS b/NEWS index c7615f7f..41bd566f 100644 --- a/NEWS +++ b/NEWS @@ -1,6 +1,7 @@ 1.5.0 (2023-xx-xx): * Add polygon partition (Raphaël Delhome, Loïc Bartoletti) * WKT: Fix triangle code (Loïc Bartoletti) + * Straight Skeleton: Add a version with extrusion 1.4.1 (2022-01-27): * Add alpha-shapes algorithm (Loïc Bartoletti, Hugo Mercier) * Fix build and tests for MSYS2/MinGW (Loïc Bartoletti) diff --git a/src/PolyhedralSurface.cpp b/src/PolyhedralSurface.cpp index da291bf3..f12bc148 100644 --- a/src/PolyhedralSurface.cpp +++ b/src/PolyhedralSurface.cpp @@ -54,6 +54,24 @@ PolyhedralSurface::PolyhedralSurface(const MarkedPolyhedron &poly) : Surface() } } +/// +/// +/// +PolyhedralSurface::PolyhedralSurface(const Mesh &sm) : Surface() +{ + + using vertex_descriptor = Mesh::Vertex_index; + for (auto face : sm.faces()) { + auto *new_face = new LineString(); + for (vertex_descriptor vd : vertices_around_face(sm.halfedge(face), sm)) { + new_face->addPoint(Point(sm.point(vd))); + } + + new_face->addPoint(new_face->startPoint().clone()); + _polygons.push_back(new Polygon(new_face)); + } +} + /// /// /// diff --git a/src/PolyhedralSurface.h b/src/PolyhedralSurface.h index 69b77c41..db731399 100644 --- a/src/PolyhedralSurface.h +++ b/src/PolyhedralSurface.h @@ -11,12 +11,16 @@ #include #include +#include #include #include #include #include #include +#include + +using Mesh = CGAL::Surface_mesh; namespace SFCGAL { @@ -42,6 +46,10 @@ class SFCGAL_API PolyhedralSurface : public Surface { * Constructor from a CGAL::Polyhedron_3 */ PolyhedralSurface(const detail::MarkedPolyhedron &poly); + /** + * Constructor from a CGAL::Surface_mesh + */ + PolyhedralSurface(const Mesh &sm); /** * Copy constructor */ diff --git a/src/algorithm/extrude.cpp b/src/algorithm/extrude.cpp index cc0da855..6a305f4a 100644 --- a/src/algorithm/extrude.cpp +++ b/src/algorithm/extrude.cpp @@ -38,7 +38,8 @@ extrude(const Point &g, const Kernel::Vector_3 &v) -> LineString *; auto extrude(const LineString &g, const Kernel::Vector_3 &v) -> PolyhedralSurface *; auto -extrude(const Polygon &g, const Kernel::Vector_3 &v) -> Solid *; +extrude(const Polygon &g, const Kernel::Vector_3 &v, bool addTop = true) + -> Solid *; auto extrude(const Triangle &g, const Kernel::Vector_3 &v) -> Solid *; @@ -120,7 +121,7 @@ extrude(const LineString &g, const Kernel::Vector_3 &v) -> PolyhedralSurface * /// /// auto -extrude(const Polygon &g, const Kernel::Vector_3 &v) -> Solid * +extrude(const Polygon &g, const Kernel::Vector_3 &v, bool addTop) -> Solid * { if (g.isEmpty()) { return new Solid(); @@ -142,11 +143,12 @@ extrude(const Polygon &g, const Kernel::Vector_3 &v) -> Solid * polyhedralSurface.addPolygon(bottom); // "top" - Polygon top(bottom); - top.reverse(); - translate(top, v); - polyhedralSurface.addPolygon(top); - + if (addTop) { + Polygon top(bottom); + top.reverse(); + translate(top, v); + polyhedralSurface.addPolygon(top); + } // exterior ring and interior rings extruded for (size_t i = 0; i < bottom.numRings(); i++) { std::unique_ptr boundaryExtruded( @@ -391,5 +393,17 @@ extrude(const Geometry &g, const double &dx, const double &dy, const double &dz) return extrude(g, Kernel::FT(dx), Kernel::FT(dy), Kernel::FT(dz)); } +SFCGAL_API auto +extrude(const Polygon &g, const double &height) -> std::unique_ptr +{ + + if (!std::isfinite(height)) { + BOOST_THROW_EXCEPTION(NonFiniteValueException( + "trying to extrude with non finite value in direction")); + } + + return std::unique_ptr( + extrude(g, Kernel::Vector_3(0.0, 0.0, height), false)); +} } // namespace algorithm } // namespace SFCGAL diff --git a/src/algorithm/extrude.h b/src/algorithm/extrude.h index 79264963..e426e47b 100644 --- a/src/algorithm/extrude.h +++ b/src/algorithm/extrude.h @@ -111,6 +111,8 @@ SFCGAL_API std::unique_ptr SFCGAL_API std::unique_ptr extrude(const Geometry &g, const Kernel::Vector_3 &v); +SFCGAL_API std::unique_ptr + extrude(const Polygon &g, const double &height); } // namespace algorithm } // namespace SFCGAL diff --git a/src/algorithm/straightSkeleton.cpp b/src/algorithm/straightSkeleton.cpp index e10f460e..1db7b996 100644 --- a/src/algorithm/straightSkeleton.cpp +++ b/src/algorithm/straightSkeleton.cpp @@ -10,6 +10,8 @@ #include #include #include +#include +#include #include #include @@ -17,10 +19,13 @@ #include #include #include +#include #include #include +#include #include +#include #include @@ -28,9 +33,11 @@ namespace SFCGAL { namespace algorithm { using Point_2 = Kernel::Point_2; +using Point_3 = Kernel::Point_3; using Polygon_2 = CGAL::Polygon_2; using Polygon_with_holes_2 = CGAL::Polygon_with_holes_2; using Straight_skeleton_2 = CGAL::Straight_skeleton_2; +using Mesh = CGAL::Surface_mesh; namespace { // anonymous @@ -379,5 +386,47 @@ approximateMedialAxis(const Geometry &g) -> std::unique_ptr return mx; } +auto +extrudeStraightSkeleton(const Polygon &g, double height) + -> std::unique_ptr +{ + + Mesh sm; + CGAL::extrude_skeleton(g.toPolygon_with_holes_2(), sm, + CGAL::parameters::maximum_height(height)); + std::unique_ptr polys(new PolyhedralSurface(sm)); + + return polys; +} + +auto +extrudeStraightSkeleton(const Geometry &g, double height) + -> std::unique_ptr +{ + SFCGAL_ASSERT_GEOMETRY_VALIDITY_2D(g); + + if (g.geometryTypeId() != TYPE_POLYGON) { + BOOST_THROW_EXCEPTION(Exception("Geometry must be a Polygon")); + } + std::unique_ptr result( + extrudeStraightSkeleton(g.as(), height)); + propagateValidityFlag(*result, true); + return result; +} + +auto +extrudeStraightSkeleton(const Geometry &g, double building_height, + double roof_height) + -> std::unique_ptr +{ + std::unique_ptr roof{ + extrudeStraightSkeleton(g, roof_height)}; + translate(*roof, 0.0, 0.0, building_height); + std::unique_ptr building(extrude(g.as(), building_height)); + std::unique_ptr result{ + new PolyhedralSurface(building->as().exteriorShell())}; + result->addPolygons(*roof); + return result; +}; } // namespace algorithm } // namespace SFCGAL diff --git a/src/algorithm/straightSkeleton.h b/src/algorithm/straightSkeleton.h index 0181c93c..49f29087 100644 --- a/src/algorithm/straightSkeleton.h +++ b/src/algorithm/straightSkeleton.h @@ -5,8 +5,10 @@ #ifndef _SFCGAL_ALGORITHM_STRAIGHTSKELETON_H_ #define _SFCGAL_ALGORITHM_STRAIGHTSKELETON_H_ +#include +#include +#include #include - #include namespace SFCGAL { @@ -84,6 +86,24 @@ SFCGAL_API std::unique_ptr bool innerOnly = false, bool outputDistanceInM = false, const double &toleranceAbs = 1e-8); +/** + * @brief build a 3D straight skeleton extruded for a Polygon + * @ingroup detail + * @throws NotImplementedException If g is a Polygon with point touching rings. + */ +SFCGAL_API auto +extrudedStraightSkeleton(const Polygon &g, double height) + -> std::unique_ptr; + +SFCGAL_API auto +extrudeStraightSkeleton(const Geometry &g, double height) + -> std::unique_ptr; + +SFCGAL_API auto +extrudeStraightSkeleton(const Geometry &g, double building_height, + double roof_height) + -> std::unique_ptr; + } // namespace algorithm } // namespace SFCGAL diff --git a/src/capi/sfcgal_c.cpp b/src/capi/sfcgal_c.cpp index 350bbf4c..df83687c 100644 --- a/src/capi/sfcgal_c.cpp +++ b/src/capi/sfcgal_c.cpp @@ -1170,6 +1170,49 @@ sfcgal_geometry_has_validity_flag(const sfcgal_geometry_t *geom) -> int return g1->hasValidityFlag() ? 1 : 0; } +extern "C" auto +sfcgal_geometry_extrude_straight_skeleton(const sfcgal_geometry_t *geom, + double height) -> sfcgal_geometry_t * +{ + const auto *g1 = reinterpret_cast(geom); + std::unique_ptr polys; + + try { + polys = SFCGAL::algorithm::extrudeStraightSkeleton(*g1, height); + } catch (std::exception &e) { + SFCGAL_WARNING("During straight_extrude_skeleton_distance(A):"); + SFCGAL_WARNING(" with A: %s", + ((const SFCGAL::Geometry *)(geom))->asText().c_str()); + SFCGAL_ERROR("%s", e.what()); + return nullptr; + } + + return polys.release(); +} + +extern "C" auto +sfcgal_geometry_extrude_polygon_straight_skeleton(const sfcgal_geometry_t *geom, + double building_height, + double roof_height) + -> sfcgal_geometry_t * +{ + const auto *g1 = reinterpret_cast(geom); + std::unique_ptr polys; + + try { + polys = SFCGAL::algorithm::extrudeStraightSkeleton(*g1, building_height, + roof_height); + } catch (std::exception &e) { + SFCGAL_WARNING("During straight_extrude_skeleton_distance(A):"); + SFCGAL_WARNING(" with A: %s", + ((const SFCGAL::Geometry *)(geom))->asText().c_str()); + SFCGAL_ERROR("%s", e.what()); + return nullptr; + } + + return polys.release(); +} + extern "C" auto sfcgal_geometry_straight_skeleton_distance_in_m(const sfcgal_geometry_t *geom) -> sfcgal_geometry_t * @@ -1259,7 +1302,6 @@ sfcgal_geometry_optimal_alpha_shapes(const sfcgal_geometry_t *geom, return result.release(); } - extern "C" sfcgal_geometry_t * sfcgal_y_monotone_partition_2(const sfcgal_geometry_t *geom) { @@ -1267,7 +1309,8 @@ sfcgal_y_monotone_partition_2(const sfcgal_geometry_t *geom) std::unique_ptr result; try { - result = SFCGAL::algorithm::partition_2(g1->as(), SFCGAL::algorithm::y_monotone); + result = SFCGAL::algorithm::partition_2(g1->as(), + SFCGAL::algorithm::y_monotone); } catch (std::exception &e) { SFCGAL_WARNING("During y_monotone_partition_2(A):"); SFCGAL_WARNING(" with A: %s", @@ -1286,7 +1329,8 @@ sfcgal_approx_convex_partition_2(const sfcgal_geometry_t *geom) std::unique_ptr result; try { - result = SFCGAL::algorithm::partition_2(g1->as(), SFCGAL::algorithm::approx_convex); + result = SFCGAL::algorithm::partition_2(g1->as(), + SFCGAL::algorithm::approx_convex); } catch (std::exception &e) { SFCGAL_WARNING("During approx_convex_partition_2(A):"); SFCGAL_WARNING(" with A: %s", @@ -1305,7 +1349,9 @@ sfcgal_greene_approx_convex_partition_2(const sfcgal_geometry_t *geom) std::unique_ptr result; try { - result = SFCGAL::algorithm::partition_2(g1->as(), SFCGAL::algorithm::greene_approx_convex); + result = + SFCGAL::algorithm::partition_2(g1->as(), + SFCGAL::algorithm::greene_approx_convex); } catch (std::exception &e) { SFCGAL_WARNING("During greene_approx_convex_partition_2(A):"); SFCGAL_WARNING(" with A: %s", @@ -1323,7 +1369,8 @@ sfcgal_optimal_convex_partition_2(const sfcgal_geometry_t *geom) std::unique_ptr result; try { - result = SFCGAL::algorithm::partition_2(g1->as(), SFCGAL::algorithm::optimal_convex); + result = SFCGAL::algorithm::partition_2(g1->as(), + SFCGAL::algorithm::optimal_convex); } catch (std::exception &e) { SFCGAL_WARNING("During optimal_convex_partition_2(A):"); SFCGAL_WARNING(" with A: %s", diff --git a/src/capi/sfcgal_c.h b/src/capi/sfcgal_c.h index 5d386425..2752365d 100644 --- a/src/capi/sfcgal_c.h +++ b/src/capi/sfcgal_c.h @@ -57,7 +57,7 @@ typedef enum { // TYPE_SURFACE = 14, //abstract SFCGAL_TYPE_POLYHEDRALSURFACE = 15, SFCGAL_TYPE_TRIANGULATEDSURFACE = 16, - SFCGAL_TYPE_TRIANGLE = 17, + SFCGAL_TYPE_TRIANGLE = 17, //-- not official codes SFCGAL_TYPE_SOLID = 101, @@ -980,6 +980,21 @@ sfcgal_geometry_straight_skeleton(const sfcgal_geometry_t *geom); SFCGAL_API sfcgal_geometry_t * sfcgal_geometry_straight_skeleton_distance_in_m(const sfcgal_geometry_t *geom); +/** + * Returns the extrude straight skeleton of the given Polygon + * @pre geom must be a Polygon + * @pre isValid(geom) == true + * @post isValid(return) == true + * @ingroup capi + */ +SFCGAL_API sfcgal_geometry_t * +sfcgal_geometry_extrude_straight_skeleton(const sfcgal_geometry_t *geom, + double height); +SFCGAL_API sfcgal_geometry_t * +sfcgal_geometry_extrude_polygon_straight_skeleton(const sfcgal_geometry_t *geom, + double building_height, + double roof_height); + /** * Returns the approximate medial axis for the given Polygon * Approximate medial axis is based on straight skeleton @@ -1064,7 +1079,8 @@ SFCGAL_API sfcgal_geometry_t * sfcgal_approx_convex_partition_2(const sfcgal_geometry_t *geom); /** - * Returns the greene approximal convex partition of a geometry (polygon without hole) + * Returns the greene approximal convex partition of a geometry (polygon without + * hole) * @pre isValid(geom) == true * @post isValid(return) == true * @ingroup capi diff --git a/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp b/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp index ce745317..db2515d8 100644 --- a/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp +++ b/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp @@ -15,190 +15,406 @@ * Library General Public License for more details. * You should have received a copy of the GNU Library General Public - * License along with this library; if not, see . + * License along with this library; if not, see + . */ #include +#include #include -#include #include -#include -#include -#include -#include -#include -#include -#include #include +#include #include #include -#include -#include +#include +#include +#include +#include +#include +#include #include #include +#include +#include #include -using namespace boost::unit_test ; -using namespace SFCGAL ; +using namespace boost::unit_test; +using namespace SFCGAL; -BOOST_AUTO_TEST_SUITE( SFCGAL_algorithm_StraightSkeletonTest ) +BOOST_AUTO_TEST_SUITE(SFCGAL_algorithm_StraightSkeletonTest) - -BOOST_AUTO_TEST_CASE( testTriangle ) +BOOST_AUTO_TEST_CASE(testTriangle) { - std::unique_ptr< Geometry > g( io::readWkt( "TRIANGLE((1 1,2 1,2 2,1 1))" ) ); - - std::string expectedWKT( "MULTILINESTRING((1.0 1.0,1.7 1.3),(2.0 1.0,1.7 1.3),(2.0 2.0,1.7 1.3))" ); - { - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 3U ); - BOOST_CHECK_EQUAL( result->asText( 1 ), expectedWKT ); - } - // check orientation insensitive - { - g->as< Triangle >().reverse() ; - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 3U ); - BOOST_CHECK_EQUAL( result->asText( 1 ), expectedWKT ); - } -} + std::unique_ptr g(io::readWkt("TRIANGLE((1 1,2 1,2 2,1 1))")); + std::string expectedWKT( + "MULTILINESTRING((1.0 1.0,1.7 1.3),(2.0 1.0,1.7 1.3),(2.0 2.0,1.7 1.3))"); + { + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 3U); + BOOST_CHECK_EQUAL(result->asText(1), expectedWKT); + } + // check orientation insensitive + { + g->as().reverse(); + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 3U); + BOOST_CHECK_EQUAL(result->asText(1), expectedWKT); + } +} -BOOST_AUTO_TEST_CASE( testPolygon ) +BOOST_AUTO_TEST_CASE(testPolygon) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON((1 1,11 1,11 11,1 11,1 1))" ) ); - - std::string expectedWKT( "MULTILINESTRING((1 1,6 6),(11 1,6 6),(11 11,6 6),(1 11,6 6))" ); - { - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 4U ); - BOOST_CHECK_EQUAL( result->asText( 0 ), expectedWKT ); - } - // check orientation insensitive - { - g->as< Polygon >().exteriorRing().reverse() ; - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 4U ); - BOOST_CHECK_EQUAL( result->asText( 0 ), expectedWKT ); - } + std::unique_ptr g( + io::readWkt("POLYGON((1 1,11 1,11 11,1 11,1 1))")); + + std::string expectedWKT( + "MULTILINESTRING((1 1,6 6),(11 1,6 6),(11 11,6 6),(1 11,6 6))"); + { + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 4U); + BOOST_CHECK_EQUAL(result->asText(0), expectedWKT); + } + // check orientation insensitive + { + g->as().exteriorRing().reverse(); + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 4U); + BOOST_CHECK_EQUAL(result->asText(0), expectedWKT); + } } -BOOST_AUTO_TEST_CASE( testPolygonWithHole ) +BOOST_AUTO_TEST_CASE(testPolygonWithHole) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON( (-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 -1.0)" - ", (-0.5 -0.5,-0.5 0.5,0.5 0.5,-0.5 -0.5)" - ")" - ) - ); - - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 13 ); - - std::unique_ptr< Geometry > expected( io::readWkt( "MULTILINESTRING( (-1/1 -1/1,-3/4 -3/4)" - ", (1/1 -1/1,1865452045155277/4503599627370496 -3730904090310553/9007199254740992)" - ", (1/1 1/1,3/4 3/4)" - ", (-1/1 1/1,-3/4 3/4)" - ", (-1/2 -1/2,-5822673418478105/9007199254740992 -7688125463633383/9007199254740992)" - ", (-1/2 1/2,-3/4 3/4)" - ", (1/2 1/2,3844062731816691/4503599627370496 727834177309763/1125899906842624)" - ", (-5822673418478105/9007199254740992 -7688125463633383/9007199254740992,1865452045155277/4503599627370496 -3730904090310553/9007199254740992)" - ", (-5822673418478105/9007199254740992 -7688125463633383/9007199254740992,-3/4 -3/4)" - ", (3844062731816691/4503599627370496 727834177309763/1125899906842624,3/4 3/4)" - ", (3844062731816691/4503599627370496 727834177309763/1125899906842624,1865452045155277/4503599627370496 -3730904090310553/9007199254740992)" - ", (-3/4 3/4,3/4 3/4)" - ", (-3/4 3/4,-3/4 -3/4)" - ")" - ) - ); - // results for gcc and clang differ due to rounding - // To avoid rounding errors in the results (-3730904090310553/9007199254740992 vs -466363011288819/1125899906842624), a text comparison is used. This is not optimal. - std::unique_ptr< Geometry > r( io::readWkt( result->asText( 10 ) ) ); - std::unique_ptr< Geometry > e( io::readWkt( expected->asText( 10 ) ) ); - BOOST_CHECK( algorithm::covers( *r, *e ) ); + std::unique_ptr g( + io::readWkt("POLYGON( (-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 -1.0)" + ", (-0.5 -0.5,-0.5 0.5,0.5 0.5,-0.5 -0.5)" + ")")); + + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 13); + std::unique_ptr expected( + io::readWkt("MULTILINESTRING( (-1/1 -1/1,-3/4 -3/4)" + ", (1/1 -1/1,1865452045155277/4503599627370496 " + "-3730904090310553/9007199254740992)" + ", (1/1 1/1,3/4 3/4)" + ", (-1/1 1/1,-3/4 3/4)" + ", (-1/2 -1/2,-5822673418478105/9007199254740992 " + "-7688125463633383/9007199254740992)" + ", (-1/2 1/2,-3/4 3/4)" + ", (1/2 1/2,3844062731816691/4503599627370496 " + "727834177309763/1125899906842624)" + ", (-5822673418478105/9007199254740992 " + "-7688125463633383/9007199254740992,1865452045155277/" + "4503599627370496 -3730904090310553/9007199254740992)" + ", (-5822673418478105/9007199254740992 " + "-7688125463633383/9007199254740992,-3/4 -3/4)" + ", (3844062731816691/4503599627370496 " + "727834177309763/1125899906842624,3/4 3/4)" + ", (3844062731816691/4503599627370496 " + "727834177309763/1125899906842624,1865452045155277/" + "4503599627370496 -3730904090310553/9007199254740992)" + ", (-3/4 3/4,3/4 3/4)" + ", (-3/4 3/4,-3/4 -3/4)" + ")")); + // results for gcc and clang differ due to rounding + // To avoid rounding errors in the results (-3730904090310553/9007199254740992 + // vs -466363011288819/1125899906842624), a text comparison is used. This is + // not optimal. + std::unique_ptr r(io::readWkt(result->asText(10))); + std::unique_ptr e(io::readWkt(expected->asText(10))); + BOOST_CHECK(algorithm::covers(*r, *e)); } -BOOST_AUTO_TEST_CASE( testPolygonWithHoleTouchingShell ) +BOOST_AUTO_TEST_CASE(testPolygonWithHoleTouchingShell) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON((-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 -1.0),(-0.5 -0.5,-0.5 0.5,0.5 0.5,1.0 -0.5,-0.5 -0.5))" ) ); - BOOST_CHECK_THROW( algorithm::straightSkeleton( *g ), NotImplementedException ); + std::unique_ptr g( + io::readWkt("POLYGON((-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 " + "-1.0),(-0.5 -0.5,-0.5 0.5,0.5 0.5,1.0 -0.5,-0.5 -0.5))")); + BOOST_CHECK_THROW(algorithm::straightSkeleton(*g), NotImplementedException); } -BOOST_AUTO_TEST_CASE( testPolygonWithTouchingHoles ) +BOOST_AUTO_TEST_CASE(testPolygonWithTouchingHoles) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON((-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 -1.0),(-0.5 -0.5,-0.5 0.5,-0.1 0.5,0.1 -0.5,-0.5 -0.5),(0.1 -0.5,0.1 0.5,0.5 0.5,0.5 -0.5,0.1 -0.5))" ) ); - - BOOST_CHECK_THROW( algorithm::straightSkeleton( *g ), NotImplementedException ); + std::unique_ptr g( + io::readWkt("POLYGON((-1.0 -1.0,1.0 -1.0,1.0 1.0,-1.0 1.0,-1.0 " + "-1.0),(-0.5 -0.5,-0.5 0.5,-0.1 0.5,0.1 -0.5,-0.5 -0.5),(0.1 " + "-0.5,0.1 0.5,0.5 0.5,0.5 -0.5,0.1 -0.5))")); + BOOST_CHECK_THROW(algorithm::straightSkeleton(*g), NotImplementedException); } -BOOST_AUTO_TEST_CASE( testMultiPolygon ) +BOOST_AUTO_TEST_CASE(testMultiPolygon) { - std::unique_ptr< Geometry > g( io::readWkt( "MULTIPOLYGON(((3.000000 0.000000,2.875000 0.484123,2.750000 0.661438,2.625000 0.780625,2.500000 0.866025,2.375000 0.927025,2.250000 0.968246,2.125000 0.992157,2.000000 1.000000,1.875000 1.484123,1.750000 1.661438,1.625000 1.780625,1.500000 1.866025,1.375000 1.927025,1.250000 1.968246,1.125000 1.992157,1.000000 2.000000,0.750000 2.661438,0.500000 2.866025,0.250000 2.968246,0.000000 3.000000,-0.250000 2.968246,-0.500000 2.866025,-0.750000 2.661438,-1.000000 2.000000,-1.125000 1.992157,-1.250000 1.968246,-1.375000 1.927025,-1.500000 1.866025,-1.625000 1.780625,-1.750000 1.661438,-1.875000 1.484123,-2.000000 1.000000,-2.125000 0.992157,-2.250000 0.968246,-2.375000 0.927025,-2.500000 0.866025,-2.625000 0.780625,-2.750000 0.661438,-2.875000 0.484123,-3.000000 0.000000,-2.875000 -0.484123,-2.750000 -0.661438,-2.625000 -0.780625,-2.500000 -0.866025,-2.375000 -0.927025,-2.250000 -0.968246,-2.125000 -0.992157,-2.000000 -1.000000,-1.875000 -1.484123,-1.750000 -1.661438,-1.625000 -1.780625,-1.500000 -1.866025,-1.375000 -1.927025,-1.250000 -1.968246,-1.125000 -1.992157,-1.000000 -2.000000,-0.750000 -2.661438,-0.500000 -2.866025,-0.250000 -2.968246,0.000000 -3.000000,0.250000 -2.968246,0.500000 -2.866025,0.750000 -2.661438,1.000000 -2.000000,1.125000 -1.992157,1.250000 -1.968246,1.375000 -1.927025,1.500000 -1.866025,1.625000 -1.780625,1.750000 -1.661438,1.875000 -1.484123,2.000000 -1.000000,2.125000 -0.992157,2.250000 -0.968246,2.375000 -0.927025,2.500000 -0.866025,2.625000 -0.780625,2.750000 -0.661438,2.875000 -0.484123,3.000000 0.000000),(0.000000 1.000000,0.125000 0.515877,0.250000 0.338562,0.375000 0.219375,0.500000 0.133975,0.625000 0.072975,0.750000 0.031754,0.875000 0.007843,1.000000 0.000000,0.875000 -0.007843,0.750000 -0.031754,0.625000 -0.072975,0.500000 -0.133975,0.375000 -0.219375,0.250000 -0.338562,0.125000 -0.515877,0.000000 -1.000000,-0.125000 -0.515877,-0.250000 -0.338562,-0.375000 -0.219375,-0.500000 -0.133975,-0.625000 -0.072975,-0.750000 -0.031754,-0.875000 -0.007843,-1.000000 0.000000,-0.875000 0.007843,-0.750000 0.031754,-0.625000 0.072975,-0.500000 0.133975,-0.375000 0.219375,-0.250000 0.338562,-0.125000 0.515877,0.000000 1.000000)))" ) ); - std::unique_ptr< MultiLineString > result( algorithm::straightSkeleton( *g ) ) ; - BOOST_CHECK_EQUAL( result->numGeometries(), 220U ); + std::unique_ptr g(io::readWkt( + "MULTIPOLYGON(((3.000000 0.000000,2.875000 0.484123,2.750000 " + "0.661438,2.625000 0.780625,2.500000 0.866025,2.375000 0.927025,2.250000 " + "0.968246,2.125000 0.992157,2.000000 1.000000,1.875000 1.484123,1.750000 " + "1.661438,1.625000 1.780625,1.500000 1.866025,1.375000 1.927025,1.250000 " + "1.968246,1.125000 1.992157,1.000000 2.000000,0.750000 2.661438,0.500000 " + "2.866025,0.250000 2.968246,0.000000 3.000000,-0.250000 " + "2.968246,-0.500000 2.866025,-0.750000 2.661438,-1.000000 " + "2.000000,-1.125000 1.992157,-1.250000 1.968246,-1.375000 " + "1.927025,-1.500000 1.866025,-1.625000 1.780625,-1.750000 " + "1.661438,-1.875000 1.484123,-2.000000 1.000000,-2.125000 " + "0.992157,-2.250000 0.968246,-2.375000 0.927025,-2.500000 " + "0.866025,-2.625000 0.780625,-2.750000 0.661438,-2.875000 " + "0.484123,-3.000000 0.000000,-2.875000 -0.484123,-2.750000 " + "-0.661438,-2.625000 -0.780625,-2.500000 -0.866025,-2.375000 " + "-0.927025,-2.250000 -0.968246,-2.125000 -0.992157,-2.000000 " + "-1.000000,-1.875000 -1.484123,-1.750000 -1.661438,-1.625000 " + "-1.780625,-1.500000 -1.866025,-1.375000 -1.927025,-1.250000 " + "-1.968246,-1.125000 -1.992157,-1.000000 -2.000000,-0.750000 " + "-2.661438,-0.500000 -2.866025,-0.250000 -2.968246,0.000000 " + "-3.000000,0.250000 -2.968246,0.500000 -2.866025,0.750000 " + "-2.661438,1.000000 -2.000000,1.125000 -1.992157,1.250000 " + "-1.968246,1.375000 -1.927025,1.500000 -1.866025,1.625000 " + "-1.780625,1.750000 -1.661438,1.875000 -1.484123,2.000000 " + "-1.000000,2.125000 -0.992157,2.250000 -0.968246,2.375000 " + "-0.927025,2.500000 -0.866025,2.625000 -0.780625,2.750000 " + "-0.661438,2.875000 -0.484123,3.000000 0.000000),(0.000000 " + "1.000000,0.125000 0.515877,0.250000 0.338562,0.375000 0.219375,0.500000 " + "0.133975,0.625000 0.072975,0.750000 0.031754,0.875000 0.007843,1.000000 " + "0.000000,0.875000 -0.007843,0.750000 -0.031754,0.625000 " + "-0.072975,0.500000 -0.133975,0.375000 -0.219375,0.250000 " + "-0.338562,0.125000 -0.515877,0.000000 -1.000000,-0.125000 " + "-0.515877,-0.250000 -0.338562,-0.375000 -0.219375,-0.500000 " + "-0.133975,-0.625000 -0.072975,-0.750000 -0.031754,-0.875000 " + "-0.007843,-1.000000 0.000000,-0.875000 0.007843,-0.750000 " + "0.031754,-0.625000 0.072975,-0.500000 0.133975,-0.375000 " + "0.219375,-0.250000 0.338562,-0.125000 0.515877,0.000000 1.000000)))")); + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 220U); } - -BOOST_AUTO_TEST_CASE( testInvalidTypes ) +BOOST_AUTO_TEST_CASE(testInvalidTypes) { - std::vector< std::string > wkt; - wkt.push_back( "POINT(1 2)" ); - wkt.push_back( "LINESTRING(0 0,1 1)" ); - - for ( std::vector< std::string >::const_iterator it=wkt.begin(), - itE=wkt.end(); - it != itE; ++it ) - { - std::unique_ptr< Geometry > g( io::readWkt( *it ) ); - std::unique_ptr< MultiLineString > result( - algorithm::straightSkeleton( *g ) - ); - BOOST_CHECK_EQUAL( result->numGeometries(), 0U ); - } + std::vector wkt; + wkt.push_back("POINT(1 2)"); + wkt.push_back("LINESTRING(0 0,1 1)"); + + for (std::vector::const_iterator it = wkt.begin(), + itE = wkt.end(); + it != itE; ++it) { + std::unique_ptr g(io::readWkt(*it)); + std::unique_ptr result(algorithm::straightSkeleton(*g)); + BOOST_CHECK_EQUAL(result->numGeometries(), 0U); + } } // See https://github.com/Oslandia/SFCGAL/issues/75 -BOOST_AUTO_TEST_CASE( testPostgisIssue3107 ) +BOOST_AUTO_TEST_CASE(testPostgisIssue3107) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON((1347259.25 7184745.94,1347273.17 7184758.16,1347280.39 7184749.95,1347278.04 7184747.88,1347281.66 7184743.76,1347284.01 7184745.83,1347293.5 7184735.05,1347279.61 7184722.85,1347269.29 7184734.6,1347270.35 7184735.51,1347267.31 7184738.96,1347266.22 7184738.01,1347259.25 7184745.94),(1347267.31 7184738.96,1347269.31 7184736.7,1347272.57 7184739.55,1347270.56 7184741.83,1347267.31 7184738.96))" ) ); - BOOST_CHECK_THROW( algorithm::straightSkeleton( *g ), NotImplementedException ); + std::unique_ptr g(io::readWkt( + "POLYGON((1347259.25 7184745.94,1347273.17 7184758.16,1347280.39 " + "7184749.95,1347278.04 7184747.88,1347281.66 7184743.76,1347284.01 " + "7184745.83,1347293.5 7184735.05,1347279.61 7184722.85,1347269.29 " + "7184734.6,1347270.35 7184735.51,1347267.31 7184738.96,1347266.22 " + "7184738.01,1347259.25 7184745.94),(1347267.31 7184738.96,1347269.31 " + "7184736.7,1347272.57 7184739.55,1347270.56 7184741.83,1347267.31 " + "7184738.96))")); + BOOST_CHECK_THROW(algorithm::straightSkeleton(*g), NotImplementedException); } // See https://github.com/Oslandia/SFCGAL/issues/91 -BOOST_AUTO_TEST_CASE( testMultiPolygonWithTouchingHoles ) +BOOST_AUTO_TEST_CASE(testMultiPolygonWithTouchingHoles) { - std::unique_ptr< Geometry > g( io::readWkt( "MULTIPOLYGON(((1347259.25 7184745.94,1347273.17 7184758.16,1347280.39 7184749.95,1347278.04 7184747.88,1347281.66 7184743.76,1347284.01 7184745.83,1347293.5 7184735.05,1347279.61 7184722.85,1347269.29 7184734.6,1347270.35 7184735.51,1347267.31 7184738.96,1347266.22 7184738.01,1347259.25 7184745.94),(1347267.31 7184738.96,1347269.31 7184736.7,1347272.57 7184739.55,1347270.56 7184741.83,1347267.31 7184738.96)))" ) ); - BOOST_CHECK_THROW( algorithm::straightSkeleton( *g ), NotImplementedException ); + std::unique_ptr g(io::readWkt( + "MULTIPOLYGON(((1347259.25 7184745.94,1347273.17 7184758.16,1347280.39 " + "7184749.95,1347278.04 7184747.88,1347281.66 7184743.76,1347284.01 " + "7184745.83,1347293.5 7184735.05,1347279.61 7184722.85,1347269.29 " + "7184734.6,1347270.35 7184735.51,1347267.31 7184738.96,1347266.22 " + "7184738.01,1347259.25 7184745.94),(1347267.31 7184738.96,1347269.31 " + "7184736.7,1347272.57 7184739.55,1347270.56 7184741.83,1347267.31 " + "7184738.96)))")); + BOOST_CHECK_THROW(algorithm::straightSkeleton(*g), NotImplementedException); } -BOOST_AUTO_TEST_CASE( testDistanceInM ) +BOOST_AUTO_TEST_CASE(testDistanceInM) { - std::unique_ptr< Geometry > g( io::readWkt( "POLYGON((0 0,1 0,1 1,0 1,0 0))" ) ); - std::unique_ptr out( algorithm::straightSkeleton( *g, /* autoOrientation */ true, /* innerOnly */ false, /* outputDistanceInM */ true ) ); - std::string expectedWKT( "MULTILINESTRING M((-0.0 -0.0 0.0,0.5 0.5 0.5),(1.0 -0.0 0.0,0.5 0.5 0.5),(1.0 1.0 0.0,0.5 0.5 0.5),(-0.0 1.0 0.0,0.5 0.5 0.5))" ); - BOOST_CHECK_EQUAL( out->asText( 1 ), expectedWKT ); + std::unique_ptr g(io::readWkt("POLYGON((0 0,1 0,1 1,0 1,0 0))")); + std::unique_ptr out(algorithm::straightSkeleton( + *g, /* autoOrientation */ true, /* innerOnly */ false, + /* outputDistanceInM */ true)); + std::string expectedWKT( + "MULTILINESTRING M((-0.0 -0.0 0.0,0.5 0.5 0.5),(1.0 -0.0 0.0,0.5 0.5 " + "0.5),(1.0 1.0 0.0,0.5 0.5 0.5),(-0.0 1.0 0.0,0.5 0.5 0.5))"); + BOOST_CHECK_EQUAL(out->asText(1), expectedWKT); } -BOOST_AUTO_TEST_CASE( testMultiEmptyEmpty ) +BOOST_AUTO_TEST_CASE(testMultiEmptyEmpty) { - std::unique_ptr< Geometry > g( io::readWkt( "MULTIPOLYGON(EMPTY,EMPTY)" ) ); - std::unique_ptr out( algorithm::straightSkeleton( *g ) ); - std::string expectedWKT( "MULTILINESTRING EMPTY" ); - BOOST_CHECK_EQUAL( out->asText( 1 ), expectedWKT ); + std::unique_ptr g(io::readWkt("MULTIPOLYGON(EMPTY,EMPTY)")); + std::unique_ptr out(algorithm::straightSkeleton(*g)); + std::string expectedWKT("MULTILINESTRING EMPTY"); + BOOST_CHECK_EQUAL(out->asText(1), expectedWKT); } // See https://gitlab.com/Oslandia/SFCGAL/-/issues/194 -BOOST_AUTO_TEST_CASE( testDegenerateMultiLineString ) +BOOST_AUTO_TEST_CASE(testDegenerateMultiLineString) +{ + std::unique_ptr g(io::readWkt( + "Polygon ((1294585.78643762995488942 200985.78643762698629871, 1294000 " + "202400, 1294000 212400, 1294585.78643762995488942 " + "213814.21356237301370129, 1296000 214400, 1297000 214400, " + "1298414.21356237004511058 213814.21356237301370129, 1299000 212400, " + "1299000 202400, 1298414.21356237004511058 200985.78643762698629871, " + "1297000 200400, 1296000 200400, 1294585.78643762995488942 " + "200985.78643762698629871),(1297000 202400, 1297000 212400, 1296000 " + "212400, 1296000 202400, 1297000 202400))")); + const double tolerance = 1e-8; + std::unique_ptr out(algorithm::straightSkeleton(*g, tolerance)); + for (size_t i = 0; i < out->numGeometries(); i++) { + BOOST_CHECK(algorithm::length(out->geometryN(i)) > tolerance); + } +} + +BOOST_AUTO_TEST_CASE(testExtrudeStraightSkeleton) +{ + + std::unique_ptr g( + io::readWkt("POLYGON (( 0 0, 5 0, 5 5, 4 5, 4 4, 0 4, 0 0 ))")); + std::unique_ptr out( + algorithm::extrudeStraightSkeleton(*g, 2.0)); + std::string expectedWKT( + "POLYHEDRALSURFACE Z(((4.00 5.00 0.00,5.00 5.00 0.00,4.00 4.00 0.00,4.00 " + "5.00 0.00)),((0.00 4.00 0.00,4.00 4.00 0.00,0.00 0.00 0.00,0.00 4.00 " + "0.00)),((4.00 4.00 0.00,5.00 0.00 0.00,0.00 0.00 0.00,4.00 4.00 " + "0.00)),((5.00 5.00 0.00,5.00 0.00 0.00,4.00 4.00 0.00,5.00 5.00 " + "0.00)),((0.00 4.00 0.00,0.00 0.00 0.00,2.00 2.00 2.00,0.00 4.00 " + "0.00)),((0.00 0.00 0.00,5.00 0.00 0.00,3.00 2.00 2.00,0.00 0.00 " + "0.00)),((2.00 2.00 2.00,0.00 0.00 0.00,3.00 2.00 2.00,2.00 2.00 " + "2.00)),((4.50 3.50 0.50,5.00 5.00 0.00,4.50 4.50 0.50,4.50 3.50 " + "0.50)),((3.00 2.00 2.00,5.00 0.00 0.00,4.50 3.50 0.50,3.00 2.00 " + "2.00)),((4.50 3.50 0.50,5.00 0.00 0.00,5.00 5.00 0.00,4.50 3.50 " + "0.50)),((5.00 5.00 0.00,4.00 5.00 0.00,4.50 4.50 0.50,5.00 5.00 " + "0.00)),((4.50 4.50 0.50,4.00 4.00 0.00,4.50 3.50 0.50,4.50 4.50 " + "0.50)),((4.50 4.50 0.50,4.00 5.00 0.00,4.00 4.00 0.00,4.50 4.50 " + "0.50)),((4.00 4.00 0.00,0.00 4.00 0.00,2.00 2.00 2.00,4.00 4.00 " + "0.00)),((4.50 3.50 0.50,4.00 4.00 0.00,3.00 2.00 2.00,4.50 3.50 " + "0.50)),((3.00 2.00 2.00,4.00 4.00 0.00,2.00 2.00 2.00,3.00 2.00 " + "2.00)))"); + BOOST_CHECK_EQUAL(out->asText(2), expectedWKT); +} + +BOOST_AUTO_TEST_CASE(testExtrudeStraightSkeletonPolygonWithHole) { - std::unique_ptr< Geometry > g( io::readWkt( "Polygon ((1294585.78643762995488942 200985.78643762698629871, 1294000 202400, 1294000 212400, 1294585.78643762995488942 213814.21356237301370129, 1296000 214400, 1297000 214400, 1298414.21356237004511058 213814.21356237301370129, 1299000 212400, 1299000 202400, 1298414.21356237004511058 200985.78643762698629871, 1297000 200400, 1296000 200400, 1294585.78643762995488942 200985.78643762698629871),(1297000 202400, 1297000 212400, 1296000 212400, 1296000 202400, 1297000 202400))")); - const double tolerance = 1e-8; - std::unique_ptr out( algorithm::straightSkeleton( *g, tolerance ) ); - for ( size_t i = 0; i < out->numGeometries(); i++ ) - { - BOOST_CHECK( algorithm::length( out->geometryN( i ) ) > tolerance ); - } + + std::unique_ptr g( + io::readWkt("POLYGON (( 0 0, 5 0, 5 5, 4 5, 4 4, 0 4, 0 0 ), (1 1, 1 2, " + "2 2, 2 1, 1 1))")); + std::unique_ptr out( + algorithm::extrudeStraightSkeleton(*g, 2.0)); + std::string expectedWKT( + "POLYHEDRALSURFACE Z(((4.00 5.00 0.00,5.00 5.00 0.00,4.00 4.00 0.00,4.00 " + "5.00 0.00)),((2.00 1.00 0.00,5.00 0.00 0.00,0.00 0.00 0.00,2.00 1.00 " + "0.00)),((5.00 5.00 0.00,5.00 0.00 0.00,4.00 4.00 0.00,5.00 5.00 " + "0.00)),((2.00 1.00 0.00,0.00 0.00 0.00,1.00 1.00 0.00,2.00 1.00 " + "0.00)),((1.00 2.00 0.00,1.00 1.00 0.00,0.00 0.00 0.00,1.00 2.00 " + "0.00)),((0.00 4.00 0.00,2.00 2.00 0.00,1.00 2.00 0.00,0.00 4.00 " + "0.00)),((0.00 4.00 0.00,1.00 2.00 0.00,0.00 0.00 0.00,0.00 4.00 " + "0.00)),((4.00 4.00 0.00,5.00 0.00 0.00,2.00 2.00 0.00,4.00 4.00 " + "0.00)),((4.00 4.00 0.00,2.00 2.00 0.00,0.00 4.00 0.00,4.00 4.00 " + "0.00)),((2.00 2.00 0.00,5.00 0.00 0.00,2.00 1.00 0.00,2.00 2.00 " + "0.00)),((0.50 2.50 0.50,0.00 0.00 0.00,0.50 0.50 0.50,0.50 2.50 " + "0.50)),((1.00 3.00 1.00,0.00 4.00 0.00,0.50 2.50 0.50,1.00 3.00 " + "1.00)),((0.50 2.50 0.50,0.00 4.00 0.00,0.00 0.00 0.00,0.50 2.50 " + "0.50)),((2.50 0.50 0.50,5.00 0.00 0.00,3.50 1.50 1.50,2.50 0.50 " + "0.50)),((0.00 0.00 0.00,5.00 0.00 0.00,2.50 0.50 0.50,0.00 0.00 " + "0.00)),((0.50 0.50 0.50,0.00 0.00 0.00,2.50 0.50 0.50,0.50 0.50 " + "0.50)),((4.50 3.50 0.50,5.00 5.00 0.00,4.50 4.50 0.50,4.50 3.50 " + "0.50)),((3.50 2.50 1.50,3.50 1.50 1.50,4.50 3.50 0.50,3.50 2.50 " + "1.50)),((4.50 3.50 0.50,5.00 0.00 0.00,5.00 5.00 0.00,4.50 3.50 " + "0.50)),((3.50 1.50 1.50,5.00 0.00 0.00,4.50 3.50 0.50,3.50 1.50 " + "1.50)),((5.00 5.00 0.00,4.00 5.00 0.00,4.50 4.50 0.50,5.00 5.00 " + "0.00)),((4.50 4.50 0.50,4.00 4.00 0.00,4.50 3.50 0.50,4.50 4.50 " + "0.50)),((4.50 4.50 0.50,4.00 5.00 0.00,4.00 4.00 0.00,4.50 4.50 " + "0.50)),((3.00 3.00 1.00,0.00 4.00 0.00,1.00 3.00 1.00,3.00 3.00 " + "1.00)),((3.50 2.50 1.50,4.50 3.50 0.50,3.00 3.00 1.00,3.50 2.50 " + "1.50)),((3.00 3.00 1.00,4.00 4.00 0.00,0.00 4.00 0.00,3.00 3.00 " + "1.00)),((4.50 3.50 0.50,4.00 4.00 0.00,3.00 3.00 1.00,4.50 3.50 " + "0.50)),((2.00 1.00 0.00,1.00 1.00 0.00,0.50 0.50 0.50,2.00 1.00 " + "0.00)),((2.50 0.50 0.50,2.00 1.00 0.00,0.50 0.50 0.50,2.50 0.50 " + "0.50)),((1.00 1.00 0.00,1.00 2.00 0.00,0.50 2.50 0.50,1.00 1.00 " + "0.00)),((0.50 0.50 0.50,1.00 1.00 0.00,0.50 2.50 0.50,0.50 0.50 " + "0.50)),((1.00 3.00 1.00,2.00 2.00 0.00,3.00 3.00 1.00,1.00 3.00 " + "1.00)),((0.50 2.50 0.50,1.00 2.00 0.00,1.00 3.00 1.00,0.50 2.50 " + "0.50)),((1.00 3.00 1.00,1.00 2.00 0.00,2.00 2.00 0.00,1.00 3.00 " + "1.00)),((2.00 2.00 0.00,2.00 1.00 0.00,2.50 0.50 0.50,2.00 2.00 " + "0.00)),((3.50 2.50 1.50,3.00 3.00 1.00,3.50 1.50 1.50,3.50 2.50 " + "1.50)),((3.50 1.50 1.50,2.00 2.00 0.00,2.50 0.50 0.50,3.50 1.50 " + "1.50)),((3.00 3.00 1.00,2.00 2.00 0.00,3.50 1.50 1.50,3.00 3.00 " + "1.00)))"); + BOOST_CHECK_EQUAL(out->asText(2), expectedWKT); +} + +BOOST_AUTO_TEST_CASE(testExtrudeStraightSkeletonGenerateBuilding) +{ + + std::unique_ptr g( + io::readWkt("POLYGON (( 0 0, 5 0, 5 5, 4 5, 4 4, 0 4, 0 0 ), (1 1, 1 2, " + "2 2, 2 1, 1 1))")); + std::unique_ptr out( + algorithm::extrudeStraightSkeleton(*g, 9.0, 2.0)); + std::string expectedWKT( + "POLYHEDRALSURFACE Z(((0.00 0.00 0.00,0.00 4.00 0.00,4.00 4.00 0.00,4.00 " + "5.00 0.00,5.00 5.00 0.00,5.00 0.00 0.00,0.00 0.00 0.00),(1.00 1.00 " + "0.00,2.00 1.00 0.00,2.00 2.00 0.00,1.00 2.00 0.00,1.00 1.00 " + "0.00)),((0.00 0.00 0.00,0.00 0.00 9.00,0.00 4.00 9.00,0.00 4.00 " + "0.00,0.00 0.00 0.00)),((0.00 4.00 0.00,0.00 4.00 9.00,4.00 4.00 " + "9.00,4.00 4.00 0.00,0.00 4.00 0.00)),((4.00 4.00 0.00,4.00 4.00 " + "9.00,4.00 5.00 9.00,4.00 5.00 0.00,4.00 4.00 0.00)),((4.00 5.00 " + "0.00,4.00 5.00 9.00,5.00 5.00 9.00,5.00 5.00 0.00,4.00 5.00 " + "0.00)),((5.00 5.00 0.00,5.00 5.00 9.00,5.00 0.00 9.00,5.00 0.00 " + "0.00,5.00 5.00 0.00)),((5.00 0.00 0.00,5.00 0.00 9.00,0.00 0.00 " + "9.00,0.00 0.00 0.00,5.00 0.00 0.00)),((1.00 1.00 0.00,1.00 1.00 " + "9.00,2.00 1.00 9.00,2.00 1.00 0.00,1.00 1.00 0.00)),((2.00 1.00 " + "0.00,2.00 1.00 9.00,2.00 2.00 9.00,2.00 2.00 0.00,2.00 1.00 " + "0.00)),((2.00 2.00 0.00,2.00 2.00 9.00,1.00 2.00 9.00,1.00 2.00 " + "0.00,2.00 2.00 0.00)),((1.00 2.00 0.00,1.00 2.00 9.00,1.00 1.00 " + "9.00,1.00 1.00 0.00,1.00 2.00 0.00)),((4.00 5.00 9.00,5.00 5.00 " + "9.00,4.00 4.00 9.00,4.00 5.00 9.00)),((2.00 1.00 9.00,5.00 0.00 " + "9.00,0.00 0.00 9.00,2.00 1.00 9.00)),((5.00 5.00 9.00,5.00 0.00 " + "9.00,4.00 4.00 9.00,5.00 5.00 9.00)),((2.00 1.00 9.00,0.00 0.00 " + "9.00,1.00 1.00 9.00,2.00 1.00 9.00)),((1.00 2.00 9.00,1.00 1.00 " + "9.00,0.00 0.00 9.00,1.00 2.00 9.00)),((0.00 4.00 9.00,2.00 2.00 " + "9.00,1.00 2.00 9.00,0.00 4.00 9.00)),((0.00 4.00 9.00,1.00 2.00 " + "9.00,0.00 0.00 9.00,0.00 4.00 9.00)),((4.00 4.00 9.00,5.00 0.00 " + "9.00,2.00 2.00 9.00,4.00 4.00 9.00)),((4.00 4.00 9.00,2.00 2.00 " + "9.00,0.00 4.00 9.00,4.00 4.00 9.00)),((2.00 2.00 9.00,5.00 0.00 " + "9.00,2.00 1.00 9.00,2.00 2.00 9.00)),((0.50 2.50 9.50,0.00 0.00 " + "9.00,0.50 0.50 9.50,0.50 2.50 9.50)),((1.00 3.00 10.00,0.00 4.00 " + "9.00,0.50 2.50 9.50,1.00 3.00 10.00)),((0.50 2.50 9.50,0.00 4.00 " + "9.00,0.00 0.00 9.00,0.50 2.50 9.50)),((2.50 0.50 9.50,5.00 0.00 " + "9.00,3.50 1.50 10.50,2.50 0.50 9.50)),((0.00 0.00 9.00,5.00 0.00 " + "9.00,2.50 0.50 9.50,0.00 0.00 9.00)),((0.50 0.50 9.50,0.00 0.00 " + "9.00,2.50 0.50 9.50,0.50 0.50 9.50)),((4.50 3.50 9.50,5.00 5.00 " + "9.00,4.50 4.50 9.50,4.50 3.50 9.50)),((3.50 2.50 10.50,3.50 1.50 " + "10.50,4.50 3.50 9.50,3.50 2.50 10.50)),((4.50 3.50 9.50,5.00 0.00 " + "9.00,5.00 5.00 9.00,4.50 3.50 9.50)),((3.50 1.50 10.50,5.00 0.00 " + "9.00,4.50 3.50 9.50,3.50 1.50 10.50)),((5.00 5.00 9.00,4.00 5.00 " + "9.00,4.50 4.50 9.50,5.00 5.00 9.00)),((4.50 4.50 9.50,4.00 4.00 " + "9.00,4.50 3.50 9.50,4.50 4.50 9.50)),((4.50 4.50 9.50,4.00 5.00 " + "9.00,4.00 4.00 9.00,4.50 4.50 9.50)),((3.00 3.00 10.00,0.00 4.00 " + "9.00,1.00 3.00 10.00,3.00 3.00 10.00)),((3.50 2.50 10.50,4.50 3.50 " + "9.50,3.00 3.00 10.00,3.50 2.50 10.50)),((3.00 3.00 10.00,4.00 4.00 " + "9.00,0.00 4.00 9.00,3.00 3.00 10.00)),((4.50 3.50 9.50,4.00 4.00 " + "9.00,3.00 3.00 10.00,4.50 3.50 9.50)),((2.00 1.00 9.00,1.00 1.00 " + "9.00,0.50 0.50 9.50,2.00 1.00 9.00)),((2.50 0.50 9.50,2.00 1.00 " + "9.00,0.50 0.50 9.50,2.50 0.50 9.50)),((1.00 1.00 9.00,1.00 2.00 " + "9.00,0.50 2.50 9.50,1.00 1.00 9.00)),((0.50 0.50 9.50,1.00 1.00 " + "9.00,0.50 2.50 9.50,0.50 0.50 9.50)),((1.00 3.00 10.00,2.00 2.00 " + "9.00,3.00 3.00 10.00,1.00 3.00 10.00)),((0.50 2.50 9.50,1.00 2.00 " + "9.00,1.00 3.00 10.00,0.50 2.50 9.50)),((1.00 3.00 10.00,1.00 2.00 " + "9.00,2.00 2.00 9.00,1.00 3.00 10.00)),((2.00 2.00 9.00,2.00 1.00 " + "9.00,2.50 0.50 9.50,2.00 2.00 9.00)),((3.50 2.50 10.50,3.00 3.00 " + "10.00,3.50 1.50 10.50,3.50 2.50 10.50)),((3.50 1.50 10.50,2.00 2.00 " + "9.00,2.50 0.50 9.50,3.50 1.50 10.50)),((3.00 3.00 10.00,2.00 2.00 " + "9.00,3.50 1.50 10.50,3.00 3.00 10.00)))"); + BOOST_CHECK_EQUAL(out->asText(2), expectedWKT); } BOOST_AUTO_TEST_SUITE_END()