1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86
|
#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/Polyhedral_mesh_domain_with_features_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/Mesh_3/experimental/Lipschitz_sizing_polyhedron.h>
//Kernel
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef K::FT FT;
// Domain
typedef CGAL::Mesh_polyhedron_3<K>::type Polyhedron;
typedef CGAL::Polyhedral_mesh_domain_with_features_3<K> Mesh_domain;
#ifdef CGAL_CONCURRENT_MESH_3
typedef CGAL::Parallel_tag Concurrency_tag;
#else
typedef CGAL::Sequential_tag Concurrency_tag;
#endif
// Triangulation
typedef CGAL::Mesh_triangulation_3<Mesh_domain,CGAL::Default,Concurrency_tag>::type Tr;
typedef CGAL::Mesh_complex_3_in_triangulation_3<
Tr,Mesh_domain::Corner_index,Mesh_domain::Curve_index> C3t3;
// Criteria
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
// Sizing field
typedef CGAL::Mesh_3::Lipschitz_sizing<K, Mesh_domain, Mesh_domain::AABB_tree> Lip_sizing;
namespace params = CGAL::parameters;
int main(int argc, char*argv[])
{
const std::string fname = (argc>1) ? argv[1] : CGAL::data_file_path("meshes/fandisk.off");
std::ifstream input(fname);
Polyhedron polyhedron;
input >> polyhedron;
if (input.fail()){
std::cerr << "Error: Cannot read file " << fname << std::endl;
return EXIT_FAILURE;
}
if (!CGAL::is_triangle_mesh(polyhedron)){
std::cerr << "Input geometry is not triangulated." << std::endl;
return EXIT_FAILURE;
}
// Create domain
Mesh_domain domain(polyhedron);
// Get sharp features
domain.detect_features();
// Create Lipschitz sizing field
Lip_sizing lip_sizing(domain, &domain.aabb_tree());
FT min_size = 0.02;
lip_sizing.add_parameters_for_subdomain(1, //subdomain id
0.3, //k
min_size,//min_size
0.5); //max_size
// Mesh criteria
Mesh_criteria criteria(params::edge_size(min_size).
facet_angle(25).
facet_size(min_size).
facet_distance(0.005).
cell_radius_edge_ratio(3).
cell_size(lip_sizing));
// Mesh generation
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria);
// Output
std::ofstream medit_file("out.mesh");
CGAL::IO::write_MEDIT(medit_file, c3t3);
medit_file.close();
return EXIT_SUCCESS;
}
|