CGAL 5.1 - 3D Generalized Barycentric Coordinates
Barycentric_coordinates_3/wachspress_coordinates.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/point_generators_3.h>
#include <CGAL/convex_hull_3.h>
#include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
#include <CGAL/Barycentric_coordinates_3/Wachspress_coordinates_3.h>
using FT = Kernel::FT;
using Point_3 = Kernel::Point_3;
using Surface_mesh = CGAL::Surface_mesh<Point_3>;
namespace PMP = CGAL::Polygon_mesh_processing;
int main()
{
std::vector<Point_3> points;
const std::size_t number_of_points = 250;
std::copy_n(gen, number_of_points, std::back_inserter(points));
Surface_mesh sm;
CGAL::convex_hull_3(points.begin(), points.end(), sm);
const std::size_t number_of_vertices = num_vertices(sm);
WP wp(sm, CP3::FAST_WITH_EDGE_CASES);
std::cout << "Computed Wachspress coordinates: " << std::endl << std::endl;
for(std::size_t i = 0; i < number_of_points; i++){
std::vector<FT> coordinates;
coordinates.reserve(number_of_vertices);
wp(points[i], std::back_inserter(coordinates));
std::cout << "Point " << i + 1 << ": " << std::endl;
for(std::size_t j = 0; j < number_of_vertices; j++)
std::cout << "Coordinate " << j + 1 << " = " << coordinates[j] << "; " << std::endl;
std::cout << std::endl;
}
return EXIT_SUCCESS;
}
CGAL::Exact_predicates_inexact_constructions_kernel
CGAL::Random_points_in_sphere_3
CGAL::Barycentric_coordinates::Computation_policy_3
Computation_policy_3
Computation_policy_3 provides a way to choose an asymptotic time complexity of the algorithm and its ...
Definition: barycentric_enum_3.h:36
copy_n
OutputIterator copy_n(InputIterator first, Size n, OutputIterator result)
Kernel::Point_3
unspecified_type Point_3
Kernel
CGAL::Barycentric_coordinates::Wachspress_coordinates_3
3D Wachspress coordinates.
Definition: Wachspress_coordinates_3.h:55
Kernel::FT
unspecified_type FT