#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/write_MEDIT.h>
#include <CGAL/Polygon_mesh_processing/IO/polygon_mesh_io.h>
#include <CGAL/Polygon_mesh_processing/bbox.h>
#include <CGAL/Polygon_mesh_processing/region_growing.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/make_conforming_constrained_Delaunay_triangulation_3.h>
using face_descriptor = boost::graph_traits<Mesh>::face_descriptor;
int main(int argc, char* argv[])
{
std::string filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/cross_quad.off");
Mesh mesh;
if(!PMP::IO::read_polygon_mesh(filename, mesh)) {
std::cerr << "Invalid input." << std::endl;
return 1;
}
std::cout << "Read " << mesh.number_of_vertices() << " vertices and "
<< mesh.number_of_faces() << " facets\n";
auto [face_patch_map, _] =mesh.add_property_map<face_descriptor, std::size_t>("f:patch_id");
const double bbox_max_span = (std::max)({bbox.x_span(), bbox.y_span(), bbox.z_span()});
std::cout << "Merging facets into coplanar patches..." << std::endl;
mesh,
face_patch_map,
CGAL::parameters::maximum_distance(bbox_max_span * 1.e-6)
.maximum_angle(5.));
for(auto f: faces(mesh))
{
if(get(face_patch_map, f) == static_cast<std::size_t>(-1)) {
put(face_patch_map, f, number_of_patches++);
}
}
std::cout << "Number of patches: " << number_of_patches << std::endl;
filename = argc > 2 ? argv[2] : "mesh.ply";
CGAL::parameters::stream_precision(17)
.use_binary_mode(false)
.face_patch_map(face_patch_map));
std::cout << "-- Wrote segmented mesh to \"" << filename << "\"\n";
std::cout << "Creating a conforming constrained Delaunay triangulation...\n";
CGAL::parameters::plc_face_id(face_patch_map));
std::cout << "Number of vertices in the CDT: "
<< ccdt.triangulation().number_of_vertices() << '\n'
<< "Number of constrained facets in the CDT: "
<< ccdt.number_of_constrained_facets() << '\n';
filename = argc > 3 ? argv[3] : "out.mesh";
std::ofstream out(filename);
out.precision(17);
std::cout << "-- Wrote CDT to \"" << filename << "\"\n";
return EXIT_SUCCESS;
}
bool write_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
std::size_t region_growing_of_planes_on_faces(const PolygonMesh &mesh, RegionMap region_map, const NamedParameters &np=parameters::default_values())
CGAL::Bbox_3 bbox(const PolygonMesh &pmesh, const NamedParameters &np=parameters::default_values())
void write_MEDIT(std::ostream &os, const T3 &t3, const NamedParameters &np=parameters::default_values())