CGAL 6.1 - 3D Constrained Triangulations
Loading...
Searching...
No Matches
Constrained_triangulation_3/ccdt_3_from_soup.cpp

Simple example demonstrating the usage of the Constrained_triangulation_3 package.

Simple example demonstrating the usage of the Constrained_triangulation_3 package.

It constructs a constrained Delaunay triangulation from a polygon soup.

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/polygon_soup_io.h>
#include <CGAL/IO/write_MEDIT.h>
#include <CGAL/draw_constrained_triangulation_3.h>
#include <CGAL/make_conforming_constrained_Delaunay_triangulation_3.h>
#include <vector>
int main(int argc, char* argv[])
{
auto filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/cubes.off");
std::vector<K::Point_3> points;
std::vector<std::vector<std::size_t>> polygons;
if(!CGAL::IO::read_polygon_soup(filename, points, polygons)) {
std::cerr << "Error: cannot read file " << filename << std::endl;
return EXIT_FAILURE;
}
std::cout << "Read " << points.size() << " vertices and "
<< polygons.size() << " polygons" << std::endl;
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';
std::ofstream ofs(argc > 2 ? argv[2] : "out.mesh");
ofs.precision(17);
CGAL::draw(ccdt);
return EXIT_SUCCESS;
}
bool read_polygon_soup(const std::string &fname, PointRange &points, PolygonRange &polygons, 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
void draw(const T3 &at3, const GSOptions &gso)
void write_MEDIT(std::ostream &os, const T3 &t3, const NamedParameters &np=parameters::default_values())