mirror of https://github.com/CGAL/cgal
add test using Mesh_3
This commit is contained in:
parent
266ab1e4f7
commit
38970cae1a
|
|
@ -25,6 +25,9 @@ if(TARGET CGAL::Eigen3_support)
|
|||
create_single_source_cgal_program("poisson_reconstruction_test_surface_mesher.cpp")
|
||||
target_link_libraries(poisson_reconstruction_test_surface_mesher PUBLIC CGAL::Eigen3_support)
|
||||
|
||||
create_single_source_cgal_program("poisson_reconstruction_test_mesh_3.cpp")
|
||||
target_link_libraries(poisson_reconstruction_test_mesh_3 PUBLIC CGAL::Eigen3_support)
|
||||
|
||||
find_package(TBB QUIET)
|
||||
include(CGAL_TBB_support)
|
||||
if (TBB_FOUND)
|
||||
|
|
|
|||
|
|
@ -0,0 +1 @@
|
|||
${CGAL_DATA_DIR}/meshes/ChineseDragon-10kv.off ${CGAL_DATA_DIR}/points_3/oni.pwn data/robocat_deci.off ${CGAL_DATA_DIR}/points_3/sphere_20k.xyz
|
||||
|
|
@ -0,0 +1,300 @@
|
|||
// poisson_reconstruction_test_mesh_3.cpp
|
||||
|
||||
//----------------------------------------------------------
|
||||
// Test the Poisson Delaunay Reconstruction method:
|
||||
// For each input point set or mesh's set of vertices, reconstruct a surface.
|
||||
// No output.
|
||||
//----------------------------------------------------------
|
||||
// poisson_reconstruction_test_mesh_3 mesh1.off point_set2.xyz...
|
||||
|
||||
// CGAL
|
||||
|
||||
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
|
||||
#include <CGAL/Timer.h>
|
||||
#include <CGAL/Memory_sizer.h>
|
||||
#include <CGAL/Polyhedron_3.h>
|
||||
|
||||
#include <CGAL/Mesh_triangulation_3.h>
|
||||
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
|
||||
#include <CGAL/Mesh_criteria_3.h>
|
||||
#include <CGAL/Labeled_mesh_domain_3.h>
|
||||
#include <CGAL/make_mesh_3.h>
|
||||
#include <CGAL/facets_in_complex_3_to_triangle_mesh.h>
|
||||
|
||||
#include <CGAL/Poisson_reconstruction_function.h>
|
||||
#include <CGAL/Point_with_normal_3.h>
|
||||
#include <CGAL/property_map.h>
|
||||
#include <CGAL/IO/read_points.h>
|
||||
#include <CGAL/compute_average_spacing.h>
|
||||
#include <CGAL/Polygon_mesh_processing/compute_normal.h>
|
||||
|
||||
#include <deque>
|
||||
#include <cstdlib>
|
||||
#include <fstream>
|
||||
#include <math.h>
|
||||
|
||||
|
||||
#include <CGAL/disable_warnings.h>
|
||||
|
||||
// ----------------------------------------------------------------------------
|
||||
// Types
|
||||
// ----------------------------------------------------------------------------
|
||||
|
||||
// kernel
|
||||
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
|
||||
|
||||
// Simple geometric types
|
||||
typedef Kernel::FT FT;
|
||||
typedef Kernel::Point_3 Point;
|
||||
typedef Kernel::Vector_3 Vector;
|
||||
typedef CGAL::Point_with_normal_3<Kernel> Point_with_normal;
|
||||
typedef Kernel::Sphere_3 Sphere;
|
||||
typedef std::deque<Point_with_normal> PointList;
|
||||
|
||||
// polyhedron
|
||||
typedef CGAL::Polyhedron_3<Kernel> Polyhedron;
|
||||
|
||||
// Poisson implicit function
|
||||
typedef CGAL::Poisson_reconstruction_function<Kernel> Poisson_reconstruction_function;
|
||||
|
||||
// Mesh_3
|
||||
typedef CGAL::Labeled_mesh_domain_3<Kernel> Mesh_domain;
|
||||
typedef CGAL::Mesh_triangulation_3<Mesh_domain>::type Tr;
|
||||
typedef CGAL::Mesh_complex_3_in_triangulation_3<Tr> C3t3;
|
||||
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
|
||||
|
||||
namespace params = CGAL::parameters;
|
||||
|
||||
// ----------------------------------------------------------------------------
|
||||
// main()
|
||||
// ----------------------------------------------------------------------------
|
||||
|
||||
int main(int argc, char * argv[])
|
||||
{
|
||||
std::cerr << "Test the Poisson Delaunay Reconstruction method" << std::endl;
|
||||
|
||||
//***************************************
|
||||
// decode parameters
|
||||
//***************************************
|
||||
|
||||
// usage
|
||||
if (argc-1 == 0)
|
||||
{
|
||||
std::cerr << "For each input point set or mesh's set of vertices, reconstruct a surface.\n";
|
||||
std::cerr << "\n";
|
||||
std::cerr << "Usage: " << argv[0] << " mesh1.off point_set2.xyz..." << std::endl;
|
||||
std::cerr << "Input file formats are .off (mesh) and .xyz or .pwn (point set).\n";
|
||||
std::cerr << "No output" << std::endl;
|
||||
return EXIT_FAILURE;
|
||||
}
|
||||
|
||||
// Poisson options
|
||||
FT sm_angle = 20.0; // Min triangle angle (degrees).
|
||||
FT sm_radius = 100; // Max triangle size w.r.t. point set average spacing.
|
||||
FT sm_distance = 0.5; // Approximation error w.r.t. point set average spacing.
|
||||
|
||||
// Accumulated errors
|
||||
int accumulated_fatal_err = EXIT_SUCCESS;
|
||||
|
||||
// Process each input file
|
||||
for (int i = 1; i <= argc-1; i++)
|
||||
{
|
||||
CGAL::Timer task_timer; task_timer.start();
|
||||
|
||||
std::cerr << std::endl;
|
||||
|
||||
//***************************************
|
||||
// Loads mesh/point set
|
||||
//***************************************
|
||||
|
||||
// File name is:
|
||||
std::string input_filename = argv[i];
|
||||
|
||||
PointList points;
|
||||
|
||||
// If OFF file format
|
||||
std::cerr << "Open " << input_filename << " for reading..." << std::endl;
|
||||
std::string extension = input_filename.substr(input_filename.find_last_of('.'));
|
||||
if (extension == ".off" || extension == ".OFF")
|
||||
{
|
||||
// Reads the mesh file in a polyhedron
|
||||
std::ifstream stream(input_filename.c_str());
|
||||
Polyhedron input_mesh;
|
||||
CGAL::scan_OFF(stream, input_mesh, true /* verbose */);
|
||||
if(!stream || !input_mesh.is_valid() || input_mesh.empty())
|
||||
{
|
||||
std::cerr << "Error: cannot read file " << input_filename << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
// Converts Polyhedron vertices to point set.
|
||||
// Computes vertices normal from connectivity.
|
||||
for(boost::graph_traits<Polyhedron>::vertex_descriptor v :
|
||||
vertices(input_mesh)){
|
||||
const Point& p = v->point();
|
||||
Vector n = CGAL::Polygon_mesh_processing::compute_vertex_normal(v,input_mesh);
|
||||
points.push_back(Point_with_normal(p,n));
|
||||
}
|
||||
}
|
||||
// If XYZ file format
|
||||
else if (extension == ".xyz" || extension == ".XYZ" ||
|
||||
extension == ".pwn" || extension == ".PWN")
|
||||
{
|
||||
// Reads the point set file in points[].
|
||||
// Note: read_points() requires an iterator over points
|
||||
// + property maps to access each point's position and normal.
|
||||
// The position property map can be omitted here as we use iterators over Point_3 elements.
|
||||
if (!CGAL::IO::read_points(
|
||||
input_filename.c_str(),
|
||||
std::back_inserter(points),
|
||||
CGAL::parameters::normal_map
|
||||
(CGAL::make_normal_of_point_with_normal_map(PointList::value_type()))
|
||||
))
|
||||
{
|
||||
std::cerr << "Error: cannot read file " << input_filename << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
std::cerr << "Error: cannot read file " << input_filename << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
// Prints status
|
||||
std::size_t memory = CGAL::Memory_sizer().virtual_size();
|
||||
std::size_t nb_points = points.size();
|
||||
std::cerr << "Reads file " << input_filename << ": " << nb_points << " points, "
|
||||
<< task_timer.time() << " seconds, "
|
||||
<< (memory>>20) << " Mb allocated"
|
||||
<< std::endl;
|
||||
task_timer.reset();
|
||||
|
||||
//***************************************
|
||||
// Checks requirements
|
||||
//***************************************
|
||||
|
||||
if (nb_points == 0)
|
||||
{
|
||||
std::cerr << "Error: empty point set" << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
bool points_have_normals = (points.begin()->normal() != CGAL::NULL_VECTOR);
|
||||
if ( ! points_have_normals )
|
||||
{
|
||||
std::cerr << "Input point set not supported: this reconstruction method requires oriented normals" << std::endl;
|
||||
// this is not a bug => do not set accumulated_fatal_err
|
||||
continue;
|
||||
}
|
||||
|
||||
CGAL::Timer reconstruction_timer; reconstruction_timer.start();
|
||||
|
||||
//***************************************
|
||||
// Computes implicit function
|
||||
//***************************************
|
||||
|
||||
std::cerr << "Computes Poisson implicit function...\n";
|
||||
|
||||
// Creates implicit function from the read points.
|
||||
// Note: this method requires an iterator over points
|
||||
// + property maps to access each point's position and normal.
|
||||
// The position property map can be omitted here as we use iterators over Point_3 elements.
|
||||
Poisson_reconstruction_function function(
|
||||
points.begin(), points.end(),
|
||||
CGAL::make_normal_of_point_with_normal_map(PointList::value_type())
|
||||
);
|
||||
|
||||
// Computes the Poisson indicator function f()
|
||||
// at each vertex of the triangulation.
|
||||
if ( ! function.compute_implicit_function() )
|
||||
{
|
||||
std::cerr << "Error: cannot compute implicit function" << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
// Prints status
|
||||
std::cerr << "Total implicit function (triangulation+refinement+solver): " << task_timer.time() << " seconds\n";
|
||||
task_timer.reset();
|
||||
|
||||
//***************************************
|
||||
// Surface mesh generation
|
||||
//***************************************
|
||||
|
||||
std::cerr << "Surface meshing...\n";
|
||||
|
||||
// Computes average spacing
|
||||
FT average_spacing = CGAL::compute_average_spacing<CGAL::Sequential_tag>(points, 6 /* knn = 1 ring */);
|
||||
|
||||
// Gets one point inside the implicit surface
|
||||
Point inner_point = function.get_inner_point();
|
||||
FT inner_point_value = function(inner_point);
|
||||
if(inner_point_value >= 0.0)
|
||||
{
|
||||
std::cerr << "Error: unable to seed (" << inner_point_value << " at inner_point)" << std::endl;
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
// Gets implicit function's radius
|
||||
Sphere bsphere = function.bounding_sphere();
|
||||
FT radius = std::sqrt(bsphere.squared_radius());
|
||||
|
||||
// Defines the implicit surface: requires defining a
|
||||
// conservative bounding sphere centered at inner point.
|
||||
FT sm_sphere_radius = 5.0 * radius;
|
||||
FT sm_dichotomy_error = sm_distance*average_spacing/1000.0; // Dichotomy error must be << sm_distance
|
||||
|
||||
// Defines surface mesh generation criteria
|
||||
Mesh_criteria criteria(params::facet_angle = sm_angle,
|
||||
params::facet_size = sm_radius*average_spacing,
|
||||
params::facet_distance = sm_distance*average_spacing);
|
||||
|
||||
std::cerr << " make_mesh_3 with sphere center=("<<inner_point << "),\n"
|
||||
<< " sphere radius="<<sm_sphere_radius<<",\n"
|
||||
<< " angle="<<sm_angle << " degrees,\n"
|
||||
<< " triangle size="<<sm_radius<<" * average spacing="<<sm_radius*average_spacing<<",\n"
|
||||
<< " distance="<<sm_distance<<" * average spacing="<<sm_distance*average_spacing<<",\n"
|
||||
<< " dichotomy = distance/"<<sm_distance*average_spacing/sm_dichotomy_error<<",\n"
|
||||
<< " manifold_with_boundary()\n";
|
||||
|
||||
// Generates surface mesh with manifold option
|
||||
Mesh_domain domain = Mesh_domain::create_implicit_mesh_domain(function, bsphere,
|
||||
params::relative_error_bound(sm_dichotomy_error / sm_sphere_radius));
|
||||
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria,
|
||||
params::no_exude().no_perturb().manifold_with_boundary());
|
||||
|
||||
// Prints status
|
||||
/*long*/ memory = CGAL::Memory_sizer().virtual_size();
|
||||
const Tr& tr = c3t3.triangulation();
|
||||
std::cerr << "Surface meshing: " << task_timer.time() << " seconds, "
|
||||
<< tr.number_of_vertices() << " output vertices, "
|
||||
<< (memory>>20) << " Mb allocated"
|
||||
<< std::endl;
|
||||
task_timer.reset();
|
||||
|
||||
if(tr.number_of_vertices() == 0) {
|
||||
accumulated_fatal_err = EXIT_FAILURE;
|
||||
continue;
|
||||
}
|
||||
|
||||
// Converts to polyhedron
|
||||
Polyhedron output_mesh;
|
||||
CGAL::facets_in_complex_3_to_triangle_mesh(c3t3, output_mesh);
|
||||
|
||||
// Prints total reconstruction duration
|
||||
std::cerr << "Total reconstruction (implicit function + meshing): " << reconstruction_timer.time() << " seconds\n";
|
||||
|
||||
} // for each input file
|
||||
|
||||
std::cerr << std::endl;
|
||||
|
||||
// Returns accumulated fatal error
|
||||
std::cerr << "Tool returned " << accumulated_fatal_err << std::endl;
|
||||
return accumulated_fatal_err;
|
||||
}
|
||||
Loading…
Reference in New Issue