#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/write_MEDIT.h>
#include <CGAL/Surface_mesh/Surface_mesh.h>
#include <CGAL/Polygon_mesh_processing/self_intersections.h>
#include <CGAL/Polygon_mesh_processing/autorefinement.h>
#include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
#include <CGAL/boost/graph/copy_face_graph.h>
#include <CGAL/make_conforming_constrained_Delaunay_triangulation_3.h>
#include <CGAL/draw_constrained_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/mpi_and_sphere.off");
std::cerr << "Error: cannot read file " << filename << std::endl;
return EXIT_FAILURE;
}
std::cout << "Number of facets in " << filename << ": "
if(!triangle_mesh)
{
std::cout << "Mesh is not a triangle mesh, triangulate faces"
<< " to check self-intersections...\n";
PMP::triangulate_faces(trimesh);
if(PMP::does_self_intersect(trimesh))
{
std::cout << "Mesh self-intersects, let's keep the triangulated version"
<< " for future autorefinement\n";
mesh = std::move(trimesh);
triangle_mesh = true;
}
}
if(triangle_mesh && PMP::does_self_intersect(mesh))
{
std::cout << "Mesh is a self-intersecting triangle mesh, perform autorefinement...\n";
std::vector<K::Point_3> 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 facets after preprocessing: "
<< polygons.size() << "\n";
if(PMP::does_polygon_soup_self_intersect(points, polygons))
{
std::cerr << "Error: mesh still self-intersects after autorefine\n";
return EXIT_FAILURE;
}
}
else
{
std::cout << "Number of facets after preprocessing: "
}
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 is_triangle_mesh(const FaceGraph &g)
void copy_face_graph(const SourceMesh &sm, TargetMesh &tm, const NamedParameters1 &np1=parameters::default_values(), const NamedParameters2 &np2=parameters::default_values())
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())