1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45
|
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/Polygon_mesh_processing/remesh.h>
#include <CGAL/Polygon_mesh_processing/IO/polygon_mesh_io.h>
#include <CGAL/Polygon_mesh_processing/Adaptive_sizing_field.h>
#include <fstream>
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Surface_mesh<K::Point_3> Mesh;
namespace PMP = CGAL::Polygon_mesh_processing;
int main(int argc, char* argv[])
{
const std::string filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/nefertiti.off");
Mesh mesh;
if (!PMP::IO::read_polygon_mesh(filename, mesh) || !CGAL::is_triangle_mesh(mesh)) {
std::cerr << "Not a valid input file." << std::endl;
return 1;
}
std::cout << "Start remeshing of " << filename
<< " (" << num_faces(mesh) << " faces)..." << std::endl;
const double tol = 0.001;
const std::pair edge_min_max{0.001, 0.5};
PMP::Adaptive_sizing_field<Mesh> sizing_field(tol, edge_min_max, faces(mesh), mesh);
unsigned int nb_iter = 5;
PMP::isotropic_remeshing(
faces(mesh),
sizing_field,
mesh,
CGAL::parameters::number_of_iterations(nb_iter)
.number_of_relaxation_steps(3)
);
CGAL::IO::write_polygon_mesh("out.off", mesh, CGAL::parameters::stream_precision(17));
std::cout << "Remeshing done." << std::endl;
return 0;
}
|