From a self-intersecting and non-triangulated surface in an OFF file, construct the constrained Delaunay triangulation after preprocessing by autorefinement.
From a self-intersecting and non-triangulated surface in an OFF file, construct the constrained Delaunay triangulation after preprocessing by autorefinement.
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/polygon_mesh_io.h>
#include <CGAL/IO/write_MEDIT.h>
#include <CGAL/Polygon_mesh_processing/autorefinement.h>
#include <CGAL/Polygon_mesh_processing/polygon_soup_self_intersections.h>
#include <CGAL/Polygon_mesh_processing/self_intersections.h>
#include <CGAL/Surface_mesh/Surface_mesh.h>
#include <CGAL/draw_constrained_triangulation_3.h>
#include <CGAL/make_conforming_constrained_Delaunay_triangulation_3.h>
using Point = K::Point_3;
int main(int argc, char* argv[])
{
const auto filename = (argc > 1) ? argv[1]
: CGAL::data_file_path("meshes/spheres_intersecting.off");
std::cerr << "Error: cannot read file " << filename << std::endl;
return EXIT_FAILURE;
}
std::cout << "Number of facets in " << filename << ": "
if(PMP::does_self_intersect(mesh))
{
std::cout << "Mesh self-intersects, performing autorefine...\n";
std::vector<Point> points;
std::vector<std::vector<std::size_t>> polygons;
PMP::polygon_mesh_to_polygon_soup(mesh, points, polygons);
PMP::autorefine_triangle_soup(points, polygons);
std::cout << "Number of input triangles after autorefine: "
<< polygons.size() << "\n";
if(PMP::does_polygon_soup_self_intersect(points, polygons))
{
std::cerr << "Error: mesh still self-intersects after autorefine\n";
std::cerr << "Run autorefinement again with an exact kernel ";
std::cerr << "such as CGAL::Exact_predicates_exact_constructions_kernel \n";
return EXIT_FAILURE;
}
}
else
{
}
std::cout << "Number of constrained facets in the CDT: "
std::ofstream ofs(argc > 2 ? argv[2] : "out.mesh");
ofs.precision(17);
return EXIT_SUCCESS;
}
size_type number_of_faces() const
bool read_polygon_mesh(const std::string &fname, Graph &g, const NamedParameters &np=parameters::default_values())
auto make_conforming_constrained_Delaunay_triangulation_3(const PolygonMesh &mesh, const NamedParameters &np=parameters::default_values())
creates a 3D constrained Delaunay triangulation conforming to the faces of a polygon mesh.
Definition: make_conforming_constrained_Delaunay_triangulation_3.h:146
void draw(const T3 &at3, const GSOptions &gso)
void write_MEDIT(std::ostream &os, const T3 &t3, const NamedParameters &np=parameters::default_values())