\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.11.3 - 3D Mesh Generation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Mesh_3/mesh_polyhedral_domain.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Mesh_triangulation_3.h>
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
#include <CGAL/Mesh_criteria_3.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/boost/graph/helpers.h>
#include <CGAL/Polyhedral_mesh_domain_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/refine_mesh_3.h>
// Domain
typedef CGAL::Polyhedron_3<K> Polyhedron;
#ifdef CGAL_CONCURRENT_MESH_3
typedef CGAL::Parallel_tag Concurrency_tag;
#else
typedef CGAL::Sequential_tag Concurrency_tag;
#endif
// Triangulation
// Criteria
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
// To avoid verbose function and named parameters call
using namespace CGAL::parameters;
int main(int argc, char*argv[])
{
const char* fname = (argc>1)?argv[1]:"data/elephant.off";
// Create input polyhedron
Polyhedron polyhedron;
std::ifstream input(fname);
input >> polyhedron;
if(input.fail()){
std::cerr << "Error: Cannot read file " << fname << std::endl;
return EXIT_FAILURE;
}
input.close();
if (!CGAL::is_triangle_mesh(polyhedron)){
std::cerr << "Input geometry is not triangulated." << std::endl;
return EXIT_FAILURE;
}
// Create domain
Mesh_domain domain(polyhedron);
// Mesh criteria (no cell_size set)
Mesh_criteria criteria(facet_angle=25, facet_size=0.15, facet_distance=0.008,
cell_radius_edge_ratio=3);
// Mesh generation
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria, no_perturb(), no_exude());
// Output
std::ofstream medit_file("out_1.mesh");
c3t3.output_to_medit(medit_file);
medit_file.close();
// Set tetrahedron size (keep cell_radius_edge_ratio), ignore facets
Mesh_criteria new_criteria(cell_radius_edge_ratio=3, cell_size=0.03);
// Mesh refinement
CGAL::refine_mesh_3(c3t3, domain, new_criteria);
// Output
medit_file.open("out_2.mesh");
c3t3.output_to_medit(medit_file);
return EXIT_SUCCESS;
}