CGAL 6.2 - Convex Decomposition of Polyhedra
Loading...
Searching...
No Matches
Convex_decomposition_3/approximate_convex_decomposition.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/Polygon_mesh_processing/IO/polygon_mesh_io.h>
#include <CGAL/approximate_convex_decomposition.h>
#include <iostream>
#include <string>
#include <vector>
using Point = K::Point_3;
using Convex_hull = std::pair<std::vector<Point>, std::vector<std::array<unsigned int, 3> > >;
using Mesh = CGAL::Surface_mesh<Point>;
int main(int argc, char* argv[])
{
const std::string filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/knot2.off");
std::cout << filename << std::endl;
Mesh mesh;
if(!PMP::IO::read_polygon_mesh(filename, mesh)) {
std::cerr << "Invalid input." << std::endl;
return 1;
}
std::vector<Convex_hull> convex_volumes;
convex_volumes.reserve(9);
CGAL::approximate_convex_decomposition(mesh, std::back_inserter(convex_volumes),
CGAL::parameters::maximum_depth(10)
.volume_error(0.1)
.maximum_number_of_convex_volumes(9)
.split_at_concavity(true)
.refitting(true)
.maximum_number_of_voxels(1000000)
.concurrency_tag(CGAL::Parallel_if_available_tag()));
for (std::size_t i = 0;i<convex_volumes.size();i++) {
const Convex_hull& ch = convex_volumes[i];
CGAL::IO::write_polygon_soup(std::to_string(i) + ".off", ch.first, ch.second);
}
return 0;
}
bool write_polygon_soup(const std::string &fname, const PointRange &points, const PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
std::size_t approximate_convex_decomposition(const FaceGraph &tmesh, OutputIterator out_volumes, const NamedParameters &np=parameters::default_values())
approximates the closed input mesh by a number of convex volumes.
Definition: approximate_convex_decomposition.h:1940
std::string data_file_path(const std::string &filename)