表面メッシュからテトラメッシュの作成(CGALライブラリ) †
表面メッシュからテトラメッシュを作成します。ここで使用するCGALは、テトラメッシュ分割法としてデローニー分割法を使用しています。
入力となる表面メッシュは、2本の直管を含む画像データから、直管表面をマーチンキューブス法を利用して取り出したものです。
テトラメッシュの生成は、計算幾何学ライブラリCGALの例題3D Polyhedral Domainsを修正して実施しました。出力ファイル形式はParaViewで可視化できるように、VTKライブラリの非構造メッシュ(Unstructured grid)形式(*.vtu)にしました。結果は以下の通り、表面は入力の表面と同様、滑らかではないですが内部にテトラメッシュが生成されました。
mesh_polyhedral_domain.cpp †
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
|
-
!
-
!
-
!
-
!
-
!
-
-
!
|
|
|
-
!
|
-
!
|
|
-
!
|
-
!
|
-
!
|
-
!
|
|
|
|
|
|
|
|
|
|
!
| #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_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/refine_mesh_3.h>
#include <CGAL/IO/Complex_3_in_triangulation_3_to_vtk.h>
#include <vtkSmartPointer.h>
#include <vtkUnstructuredGrid.h>
#include <vtkXMLUnstructuredGridWriter.h>
#include <CGAL/IO/Polyhedron_iostream.h>
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Polyhedron_3<K> Polyhedron;
typedef CGAL::Polyhedral_mesh_domain_3<Polyhedron, K> Mesh_domain;
typedef CGAL::Mesh_triangulation_3<Mesh_domain>::type Tr;
typedef CGAL::Mesh_complex_3_in_triangulation_3<Tr> C3t3;
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
using namespace CGAL::parameters;
int main()
{
Polyhedron polyhedron;
std::ifstream input("data/surface.off");
input >> polyhedron;
Mesh_domain domain(polyhedron);
Mesh_criteria criteria(facet_angle=25, facet_size=0.1, facet_distance=0.008,
cell_radius_edge_ratio=3);
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria, no_perturb(), no_exude());
Mesh_criteria new_criteria(cell_radius_edge_ratio=3, cell_size=0.05);
CGAL::refine_mesh_3(c3t3, domain, new_criteria);
vtkSmartPointer<vtkUnstructuredGrid> output =
vtkSmartPointer<vtkUnstructuredGrid>::New();
CGAL::output_c3t3_to_vtk_unstructured_grid(c3t3, output);
vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer =
vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
writer->SetFileName("out.vtu");
writer->SetInput(output);
writer->Update();
return 0;
}
|
ダウンロードとビルド †
ソースコードとCMakeLists.txtファイル: