CGAL 6.1 - 3D Constrained Triangulations
Loading...
Searching...
No Matches
Constrained_triangulation_3/conforming_constrained_Delaunay_triangulation_3.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 3D constrained Delaunay triangulation from a polygon mesh.

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/polygon_mesh_io.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/make_conforming_constrained_Delaunay_triangulation_3.h>
#include <CGAL/IO/write_MEDIT.h>
#include <cassert>
int main(int argc, char* argv[])
{
auto filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/mpi.off");
if(!CGAL::IO::read_polygon_mesh(filename, mesh)) {
std::cerr << "Error: cannot read file " << filename << std::endl;
return EXIT_FAILURE;
}
std::cout << "Read " << mesh.number_of_vertices() << " vertices and "
<< mesh.number_of_faces() << " faces" << std::endl;
std::cout << "Number of vertices in the CDT: "
<< ccdt.triangulation().number_of_vertices() << '\n';
std::cout << "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);
auto tr = std::move(ccdt).triangulation();
// Now `tr` is a valid `CGAL::Triangulation_3` object that can be used for further processing.
// and the triangulation of `ccdt` is empty.
std::cout << "Number of vertices in the triangulation `tr`: "
<< tr.number_of_vertices() << '\n';
std::cout << "Number of vertices in `ccdt`: "
<< ccdt.triangulation().number_of_vertices() << '\n';
assert(ccdt.triangulation().number_of_vertices() == 0);
return EXIT_SUCCESS;
}
size_type number_of_vertices() const
size_type number_of_faces() const
bool read_polygon_mesh(const std::string &fname, Graph &g, 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 write_MEDIT(std::ostream &os, const T3 &t3, const NamedParameters &np=parameters::default_values())