File: mesh_implicit_ellipsoid.cpp

package info (click to toggle)
cgal 4.13-1
  • links: PTS
  • area: main
  • in suites: buster
  • size: 101,504 kB
  • sloc: cpp: 703,154; ansic: 163,044; sh: 674; fortran: 616; python: 411; makefile: 115
file content (73 lines) | stat: -rw-r--r-- 2,060 bytes parent folder | download
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
#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/Labeled_mesh_domain_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/perturb_mesh_3.h>
#include <CGAL/exude_mesh_3.h>

// Domain
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef K::FT FT;
typedef K::Point_3 Point;
typedef FT (Function)(const Point&);
typedef CGAL::Labeled_mesh_domain_3<K> Mesh_domain;

// Triangulation
typedef CGAL::Mesh_triangulation_3<Mesh_domain>::type Tr;
typedef CGAL::Mesh_complex_3_in_triangulation_3<Tr> C3t3;

// Criteria
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
typedef Mesh_criteria::Facet_criteria Facet_criteria;
typedef Mesh_criteria::Cell_criteria Cell_criteria;

// Function
FT ellipsoid_function (const Point& p)
{
  const FT x2=p.x()*p.x();
  const FT y2=p.y()*p.y();
  const FT z2=p.z()*p.z();
  
  return x2+2*y2+4*z2-1;
}

// To avoid verbose function and named parameters call
using namespace CGAL::parameters;

int main()
{
  // Domain (Warning: Sphere_3 constructor uses square radius !)
  Mesh_domain domain =
    Mesh_domain::create_implicit_mesh_domain(ellipsoid_function,
                                             K::Sphere_3(CGAL::ORIGIN, 2.));
  
  // Criteria
  Facet_criteria facet_criteria(30, 0.08, 0.025); // angle, size, approximation
  Cell_criteria cell_criteria(2, 0.1); // radius-edge ratio, size
  Mesh_criteria criteria(facet_criteria, cell_criteria);
  
  // Mesh generation (without optimization)
  C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria, no_perturb(), no_exude());
  
  // Output
  std::ofstream medit_file("out_wo.mesh");
  c3t3.output_to_medit(medit_file);
  medit_file.close();
  
  // Perturbation (5s, 12degree)
  CGAL::perturb_mesh_3(c3t3, domain, time_limit=5, sliver_bound=12);
  
  // Exudation
  CGAL::exude_mesh_3(c3t3);
  
  // Output
  medit_file.open("out_optimized.mesh");
  c3t3.output_to_medit(medit_file);
  
  return 0;
}