mirror of https://github.com/CGAL/cgal
Stream Support: Return false if nothing was read (#9038)
## Summary of Changes The `read_..._WKT()` functions have `bool` as return type, but currently only check if the stream is bad, which is not the case if no input was read. This PR was triggered by the fact that CGALlab could not WTK files with polygons, but only polylines. It now can read polygons, but does only draw the contours and not the interior. ## Release Management * Affected package(s): Stream_support * License and copyright ownership: unchanged
This commit is contained in:
commit
39c7001a25
|
|
@ -2,6 +2,9 @@
|
|||
|
||||
#include <CGAL/Three/CGAL_Lab_io_plugin_interface.h>
|
||||
#include <CGAL/Three/Three.h>
|
||||
#include <CGAL/Projection_traits_xy_3.h>
|
||||
#include <CGAL/Polygon_with_holes_2.h>
|
||||
#include <CGAL/Kernel_traits.h>
|
||||
#include <QInputDialog>
|
||||
#include <QApplication>
|
||||
#include <fstream>
|
||||
|
|
@ -45,37 +48,51 @@ canLoad(QFileInfo) const {
|
|||
QList<Scene_item*>
|
||||
CGAL_Lab_wkt_plugin::
|
||||
load(QFileInfo fileinfo, bool& ok, bool add_to_scene) {
|
||||
typedef Scene_polylines_item::Point_3 Point_3;
|
||||
typedef CGAL::Kernel_traits<Point_3>::Kernel K;
|
||||
typedef std::vector<Point_3> Polyline;
|
||||
typedef CGAL::Projection_traits_xy_3<K> Kernel;
|
||||
typedef CGAL::Polygon_with_holes_2<Kernel> Polygon;
|
||||
|
||||
std::ifstream in(fileinfo.filePath().toUtf8());
|
||||
|
||||
if(!in)
|
||||
std::cerr << "Error!\n";
|
||||
|
||||
if(!in || fileinfo.size() == 0)
|
||||
{
|
||||
CGAL::Three::Three::warning( tr("The file you are trying to load does not exist or is empty."));
|
||||
ok = false;
|
||||
return QList<Scene_item*>();
|
||||
}
|
||||
|
||||
QApplication::setOverrideCursor(Qt::WaitCursor);
|
||||
|
||||
if(fileinfo.size() == 0)
|
||||
std::vector<Point_3> points;
|
||||
std::list<Polyline> polylines;
|
||||
std::vector<Polygon> polygons;
|
||||
bool success = CGAL::IO::read_WKT(in, points, polylines, polygons);
|
||||
for(const Polygon& p : polygons)
|
||||
{
|
||||
CGAL::Three::Three::warning( tr("The file you are trying to load is empty."));
|
||||
ok = false;
|
||||
QApplication::restoreOverrideCursor();
|
||||
return QList<Scene_item*>();
|
||||
Polyline polyline(p.outer_boundary().vertices_begin(), p.outer_boundary().vertices_end());
|
||||
polyline.push_back(polyline.front());
|
||||
polylines.push_back(polyline);
|
||||
for(auto hit = p.holes_begin(); hit != p.holes_end(); ++hit)
|
||||
{
|
||||
Polyline hole(hit->vertices_begin(), hit->vertices_end());
|
||||
hole.push_back(hole.front());
|
||||
polylines.push_back(hole);
|
||||
}
|
||||
}
|
||||
if(! polygons.empty()){
|
||||
CGAL::Three::Three::warning( tr("The polygons will be drawn as polylines"));
|
||||
}
|
||||
|
||||
std::list<std::vector<Scene_polylines_item::Point_3> > polylines;
|
||||
bool success = CGAL::IO::read_multi_linestring_WKT (in, polylines);
|
||||
if(! success){
|
||||
in.close();
|
||||
in.open(fileinfo.filePath().toUtf8());
|
||||
std::vector<Scene_polylines_item::Point_3> polyline;
|
||||
std::cout << " read" << std::endl;
|
||||
success = CGAL::IO::read_linestring_WKT (in, polyline);
|
||||
std::cout << " done " << std::boolalpha << success << " " << polyline.size() << std::endl;
|
||||
if(! success){
|
||||
if(!success || polylines.empty())
|
||||
{
|
||||
CGAL::Three::Three::warning( tr("The file you are trying to load is not a valid WKT file or is empty"));
|
||||
ok = false;
|
||||
QApplication::restoreOverrideCursor();
|
||||
return QList<Scene_item*>();
|
||||
}
|
||||
polylines.push_back(polyline);
|
||||
}
|
||||
|
||||
Scene_polylines_item* item = new Scene_polylines_item;
|
||||
item->polylines = polylines;
|
||||
|
|
|
|||
|
|
@ -1,15 +1,17 @@
|
|||
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
|
||||
#include <CGAL/Projection_traits_xy_3.h>
|
||||
#include <CGAL/IO/WKT.h>
|
||||
|
||||
#include <iostream>
|
||||
#include <fstream>
|
||||
#include <vector>
|
||||
|
||||
typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
|
||||
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
|
||||
typedef CGAL::Projection_traits_xy_3<K> Kernel;
|
||||
|
||||
int main(int argc, char* argv[])
|
||||
{
|
||||
typedef CGAL::Point_2<Kernel> Point;
|
||||
typedef Kernel::Point_2 Point;
|
||||
typedef std::vector<Point> MultiPoint;
|
||||
|
||||
typedef std::vector<Point> LineString;
|
||||
|
|
@ -23,6 +25,7 @@ int main(int argc, char* argv[])
|
|||
MultiPoint points;
|
||||
MultiLineString polylines;
|
||||
MultiPolygon polygons;
|
||||
|
||||
CGAL::IO::read_WKT(is, points,polylines,polygons);
|
||||
|
||||
for(Point p : points)
|
||||
|
|
@ -30,8 +33,9 @@ int main(int argc, char* argv[])
|
|||
for(LineString ls : polylines)
|
||||
for(Point p : ls)
|
||||
std::cout<<p<<std::endl;
|
||||
for(Polygon p : polygons)
|
||||
for(Polygon p : polygons){
|
||||
std::cout<<p<<std::endl;
|
||||
}
|
||||
|
||||
}
|
||||
return 0;
|
||||
|
|
|
|||
|
|
@ -133,16 +133,19 @@ bool read_multi_point_WKT(std::istream& in,
|
|||
MultiPoint& mp)
|
||||
{
|
||||
std::string line;
|
||||
bool found = false;
|
||||
while(internal::get_a_new_line(in, line))
|
||||
{
|
||||
if(line.substr(0, 10).compare("MULTIPOINT") == 0)
|
||||
{
|
||||
CGAL::internal::Geometry_container<MultiPoint, boost::geometry::multi_point_tag> gc(mp);
|
||||
internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
found = internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if(! found){
|
||||
return false;
|
||||
}
|
||||
return !in.fail();
|
||||
}
|
||||
|
||||
|
|
@ -168,16 +171,20 @@ bool read_linestring_WKT(std::istream& in,
|
|||
LineString& polyline)
|
||||
{
|
||||
std::string line;
|
||||
bool found = false;
|
||||
while(internal::get_a_new_line(in, line))
|
||||
{
|
||||
if(line.substr(0, 10).compare("LINESTRING") == 0)
|
||||
{
|
||||
CGAL::internal::Geometry_container<LineString, boost::geometry::linestring_tag> gc(polyline);
|
||||
internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
found = internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if(! found){
|
||||
return false;
|
||||
}
|
||||
return !in.fail();
|
||||
}
|
||||
|
||||
|
|
@ -199,6 +206,7 @@ bool read_multi_linestring_WKT(std::istream& in,
|
|||
MultiLineString& mls)
|
||||
{
|
||||
std::string line;
|
||||
bool found = false;
|
||||
while(internal::get_a_new_line(in, line))
|
||||
{
|
||||
if(line.substr(0, 15).compare("MULTILINESTRING") == 0)
|
||||
|
|
@ -209,7 +217,7 @@ bool read_multi_linestring_WKT(std::istream& in,
|
|||
std::vector<LineString> pr_range;
|
||||
CGAL::internal::Geometry_container<std::vector<LineString>, boost::geometry::multi_linestring_tag> gc(pr_range);
|
||||
|
||||
internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
found = internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
for(LineString& ls : gc) {
|
||||
mls.push_back(*ls.range);
|
||||
}
|
||||
|
|
@ -217,7 +225,9 @@ bool read_multi_linestring_WKT(std::istream& in,
|
|||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if(! found){
|
||||
return false;
|
||||
}
|
||||
return !in.fail();
|
||||
}
|
||||
|
||||
|
|
@ -237,15 +247,19 @@ bool read_polygon_WKT(std::istream& in,
|
|||
Polygon& polygon)
|
||||
{
|
||||
std::string line;
|
||||
bool found = false;
|
||||
while(internal::get_a_new_line(in, line))
|
||||
{
|
||||
if(line.substr(0, 7).compare("POLYGON") == 0)
|
||||
{
|
||||
internal::read_wkt_or_fail_stream(in, line, polygon);
|
||||
found = internal::read_wkt_or_fail_stream(in, line, polygon);
|
||||
internal::pop_back_if_equal_to_front(polygon);
|
||||
break;
|
||||
}
|
||||
}
|
||||
if(! found){
|
||||
return false;
|
||||
}
|
||||
return !in.fail();
|
||||
}
|
||||
|
||||
|
|
@ -268,12 +282,13 @@ bool read_multi_polygon_WKT(std::istream& in,
|
|||
MultiPolygon& polygons)
|
||||
{
|
||||
std::string line;
|
||||
bool found = false;
|
||||
while(internal::get_a_new_line(in, line))
|
||||
{
|
||||
if(line.substr(0, 12).compare("MULTIPOLYGON") == 0)
|
||||
{
|
||||
CGAL::internal::Geometry_container<MultiPolygon, boost::geometry::multi_polygon_tag> gc(polygons);
|
||||
internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
found = internal::read_wkt_or_fail_stream(in, line, gc);
|
||||
|
||||
for(auto& p : gc)
|
||||
internal::pop_back_if_equal_to_front(p);
|
||||
|
|
@ -281,7 +296,9 @@ bool read_multi_polygon_WKT(std::istream& in,
|
|||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if(! found){
|
||||
return false;
|
||||
}
|
||||
return !in.fail();
|
||||
}
|
||||
|
||||
|
|
@ -461,6 +478,7 @@ bool read_WKT(std::istream& is,
|
|||
{
|
||||
auto fail = [&is]() { is.clear(is.rdstate() | std::ios::failbit); return false; };
|
||||
|
||||
bool found = false;
|
||||
std::string line;
|
||||
while(is >> std::ws && is.good() && std::getline(is, line))
|
||||
{
|
||||
|
|
@ -490,44 +508,57 @@ bool read_WKT(std::istream& is,
|
|||
{
|
||||
Point p;
|
||||
if(!IO::read_point_WKT(iss, p) ) return fail();
|
||||
found = true;
|
||||
points.push_back(p);
|
||||
}
|
||||
else if(type == "LINESTRING")
|
||||
{
|
||||
LineString l;
|
||||
if(!IO::read_linestring_WKT(iss, l)) return fail();
|
||||
found = true;
|
||||
polylines.push_back(std::move(l));
|
||||
}
|
||||
else if(type == "POLYGON")
|
||||
{
|
||||
Polygon p;
|
||||
if(!IO::read_polygon_WKT(iss, p)) return fail();
|
||||
if(!p.outer_boundary().is_empty())
|
||||
if(!p.outer_boundary().is_empty()){
|
||||
found = true;
|
||||
polygons.push_back(std::move(p));
|
||||
}
|
||||
}
|
||||
else if(type == "MULTIPOINT")
|
||||
{
|
||||
MultiPoint mp;
|
||||
if(!IO::read_multi_point_WKT(iss, mp)) return fail();
|
||||
for(const Point& point : mp)
|
||||
for(const Point& point : mp){
|
||||
points.push_back(point);
|
||||
found = true;
|
||||
}
|
||||
}
|
||||
else if(type == "MULTILINESTRING")
|
||||
{
|
||||
MultiLineString mls;
|
||||
if(!IO::read_multi_linestring_WKT(iss, mls)) return fail();
|
||||
for(LineString& ls : mls)
|
||||
for(LineString& ls : mls){
|
||||
polylines.push_back(std::move(ls));
|
||||
found = true;
|
||||
}
|
||||
}
|
||||
else if(type == "MULTIPOLYGON")
|
||||
{
|
||||
MultiPolygon mp;
|
||||
if(!IO::read_multi_polygon_WKT(iss, mp)) return fail();
|
||||
for(Polygon& poly : mp)
|
||||
for(Polygon& poly : mp){
|
||||
polygons.push_back(std::move(poly));
|
||||
found = true;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
if(!found){
|
||||
return false;
|
||||
}
|
||||
|
||||
return !is.fail();
|
||||
}
|
||||
|
|
|
|||
|
|
@ -24,6 +24,42 @@ typedef std::vector<Point3> MultiPoint
|
|||
typedef std::vector<Linestring3> MultiLinestring3;
|
||||
////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
|
||||
bool test_mismatch()
|
||||
{
|
||||
Point3 p(1,2,3);
|
||||
MultiPoint3 mq;
|
||||
std::stringstream ss;
|
||||
CGAL::IO::write_point_WKT(ss, p);
|
||||
bool b = CGAL::IO::read_multi_point_WKT(ss, mq);
|
||||
assert(!b);
|
||||
return !b;
|
||||
}
|
||||
|
||||
bool test_read2Dinto3D()
|
||||
{
|
||||
Point p(10,11);
|
||||
Point3 p3(1,2,3);
|
||||
std::stringstream ss;
|
||||
CGAL::IO::write_point_WKT(ss, p);
|
||||
bool b = CGAL::IO::read_point_WKT(ss, p3);
|
||||
assert(p.x() == p3.x());
|
||||
assert(p.y() == p3.y());
|
||||
assert(p3.z() == 0);
|
||||
|
||||
Point q(12,13);
|
||||
Linestring ls;
|
||||
Linestring3 ls3;
|
||||
ls.push_back(p);
|
||||
ls.push_back(q);
|
||||
CGAL::IO::write_linestring_WKT(ss, ls);
|
||||
b = CGAL::IO::read_linestring_WKT(ss, ls3);
|
||||
assert(ls3[0]==p3);
|
||||
assert(ls3[1] == Point3(12,13,0));
|
||||
assert(b);
|
||||
return b;
|
||||
}
|
||||
|
||||
|
||||
bool test_WKT_3D()
|
||||
{
|
||||
{
|
||||
|
|
@ -337,6 +373,10 @@ int main()
|
|||
assert(ok);
|
||||
ok = test_WKT_3D();
|
||||
assert(ok);
|
||||
ok = test_mismatch();
|
||||
assert(ok);
|
||||
ok = test_read2Dinto3D();
|
||||
assert(ok);
|
||||
|
||||
return EXIT_SUCCESS;
|
||||
}
|
||||
|
|
|
|||
Loading…
Reference in New Issue