mirror of https://github.com/CGAL/cgal
reorganize example and use non-manifold input
This commit is contained in:
parent
b6bb23d81f
commit
9b99d0e754
|
|
@ -10,6 +10,7 @@
|
||||||
|
|
||||||
#include <CGAL/draw_constrained_triangulation_3.h>
|
#include <CGAL/draw_constrained_triangulation_3.h>
|
||||||
#include <CGAL/IO/polygon_mesh_io.h>
|
#include <CGAL/IO/polygon_mesh_io.h>
|
||||||
|
#include <CGAL/IO/polygon_soup_io.h>
|
||||||
#include <CGAL/IO/write_MEDIT.h>
|
#include <CGAL/IO/write_MEDIT.h>
|
||||||
|
|
||||||
#include <vector>
|
#include <vector>
|
||||||
|
|
@ -45,38 +46,47 @@ bool verify_preconditions_soup(const PointRange& points, const PolygonRange& pol
|
||||||
int main(int argc, char* argv[])
|
int main(int argc, char* argv[])
|
||||||
{
|
{
|
||||||
const auto filename = (argc > 1) ? argv[1]
|
const auto filename = (argc > 1) ? argv[1]
|
||||||
: CGAL::data_file_path("meshes/mpi_and_sphere.off");
|
: CGAL::data_file_path("meshes/cubes.off");
|
||||||
|
|
||||||
CGAL::Surface_mesh<K::Point_3> mesh;
|
// Read polygon soup
|
||||||
if(!CGAL::IO::read_polygon_mesh(filename, mesh))
|
using Points = std::vector<Point>;
|
||||||
|
using PLC_face = std::vector<std::size_t>;
|
||||||
|
using PLC_faces = std::vector<PLC_face>;
|
||||||
|
Points points;
|
||||||
|
PLC_faces faces;
|
||||||
|
if(!CGAL::IO::read_polygon_soup(filename, points, faces))
|
||||||
{
|
{
|
||||||
std::cerr << "Error: cannot read file " << filename << std::endl;
|
std::cerr << "Error: cannot read file " << filename << std::endl;
|
||||||
return EXIT_FAILURE;
|
return EXIT_FAILURE;
|
||||||
}
|
}
|
||||||
|
|
||||||
std::cout << "Number of facets in " << filename << ": "
|
// When possible, convert polygon soup to polygon mesh
|
||||||
<< mesh.number_of_faces() << "\n";
|
CGAL::Surface_mesh<K::Point_3> mesh;
|
||||||
|
if(CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh(faces))
|
||||||
|
{
|
||||||
|
CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh(points, faces, mesh);
|
||||||
|
std::cout << "Number of facets in " << filename << ": " << mesh.number_of_faces() << "\n";
|
||||||
|
}
|
||||||
|
|
||||||
// Verify preconditions for the mesh input
|
// Verify preconditions for the mesh input
|
||||||
if(!verify_preconditions_mesh(mesh))
|
const bool is_polygon_mesh = !mesh.is_empty();
|
||||||
|
if(is_polygon_mesh && !verify_preconditions_mesh(mesh))
|
||||||
{
|
{
|
||||||
// If the mesh is not a valid triangle mesh or has self-intersections,
|
|
||||||
// convert it to a polygon soup and verify the preconditions for the soup
|
|
||||||
std::vector<K::Point_3> points;
|
|
||||||
std::vector<std::vector<std::size_t>> polygons;
|
|
||||||
CGAL::Polygon_mesh_processing::polygon_mesh_to_polygon_soup(mesh, points, polygons);
|
|
||||||
if(!verify_preconditions_soup(points, polygons))
|
|
||||||
{
|
|
||||||
std::cerr << "Error: input polygon soup is not a valid input for CCDT_3\n";
|
|
||||||
return EXIT_FAILURE;
|
|
||||||
}
|
|
||||||
|
|
||||||
std::cerr << "Error: input mesh is not a valid input for CCDT_3\n";
|
std::cerr << "Error: input mesh is not a valid input for CCDT_3\n";
|
||||||
return EXIT_FAILURE;
|
return EXIT_FAILURE;
|
||||||
}
|
}
|
||||||
|
else if(!verify_preconditions_soup(points, faces))
|
||||||
|
{
|
||||||
|
std::cerr << "Error: input polygon soup is not a valid input for CCDT_3\n";
|
||||||
|
return EXIT_FAILURE;
|
||||||
|
}
|
||||||
|
|
||||||
auto ccdt = CGAL::make_conforming_constrained_Delaunay_triangulation_3(mesh);
|
// Build conforming constrained Delaunay triangulation
|
||||||
|
auto ccdt = is_polygon_mesh
|
||||||
|
? CGAL::make_conforming_constrained_Delaunay_triangulation_3(mesh)
|
||||||
|
: CGAL::make_conforming_constrained_Delaunay_triangulation_3(points, faces);
|
||||||
|
|
||||||
|
// Print mesh details
|
||||||
if(ccdt.number_of_constrained_facets() == 0)
|
if(ccdt.number_of_constrained_facets() == 0)
|
||||||
{
|
{
|
||||||
std::cerr << "Error: no constrained facets in the CDT.\n";
|
std::cerr << "Error: no constrained facets in the CDT.\n";
|
||||||
|
|
@ -87,6 +97,7 @@ int main(int argc, char* argv[])
|
||||||
std::cout << "Number of constrained facets in the CDT: "
|
std::cout << "Number of constrained facets in the CDT: "
|
||||||
<< ccdt.number_of_constrained_facets() << '\n';
|
<< ccdt.number_of_constrained_facets() << '\n';
|
||||||
|
|
||||||
|
// Output
|
||||||
std::ofstream ofs(argc > 2 ? argv[2] : "out.mesh");
|
std::ofstream ofs(argc > 2 ? argv[2] : "out.mesh");
|
||||||
ofs.precision(17);
|
ofs.precision(17);
|
||||||
CGAL::IO::write_MEDIT(ofs, ccdt);
|
CGAL::IO::write_MEDIT(ofs, ccdt);
|
||||||
|
|
|
||||||
Loading…
Reference in New Issue