File: AABB_indexed_triangle_2_example.cpp

package info (click to toggle)
cgal 6.0.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 141,840 kB
  • sloc: cpp: 797,081; ansic: 203,398; sh: 490; python: 411; makefile: 286; javascript: 174
file content (91 lines) | stat: -rw-r--r-- 3,250 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
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
#include <iostream>
#include <vector>
#include <array>

#include <CGAL/Simple_cartesian.h>
#include <CGAL/AABB_tree.h>
#include <CGAL/AABB_traits_2.h>
#include <CGAL/AABB_indexed_triangle_primitive_2.h>

typedef CGAL::Simple_cartesian<double> K;
typedef K::Point_3 Point_3;
typedef K::Point_2 Point_2;
typedef K::Ray_2 Ray;

template <class GeomTraits>
struct Projection_xy_point_map {

  typedef typename GeomTraits::Point_3 key_type;
  typedef typename GeomTraits::Point_2 value_type;
  typedef value_type reference;

  typedef boost::readable_property_map_tag category;
  typedef Projection_xy_point_map<GeomTraits> Self;

  Projection_xy_point_map() {}

  inline friend value_type get(Self, key_type p)
  {
    return value_type(p.x(), p.y());
  }
};

typedef std::vector<std::array<uint8_t, 3> >::const_iterator IndexIterator;
typedef std::vector<Point_3> PointRange;
typedef CGAL::AABB_indexed_triangle_primitive_2<K, IndexIterator, PointRange, CGAL::Tag_false, Projection_xy_point_map<K>> Primitive;
typedef CGAL::AABB_traits_2<K, Primitive> AABB_triangle_traits;
typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
typedef Tree::Point_and_primitive_id Point_and_primitive_id;
typedef std::optional<Tree::Intersection_and_primitive_id<Ray>::Type> Ray_intersection;

int main()
{
  Point_3 a(0.0, 0.0, 0.0);
  Point_3 b(0.0, 1.0, 0.0);
  Point_3 c(1.0, 0.0, 0.0);
  Point_3 d(1.0, 1.0, 0.0);
  Point_3 e(2.0, 0.0, 0.0);
  Point_3 f(2.0, 1.0, 0.0);

  std::vector<Point_3> points = { a, b, c, d, e, f };

  std::vector<std::array<uint8_t, 3> > triangles;
  triangles.push_back({ 0, 2, 1 });
  triangles.push_back({ 1, 2, 3 });
  triangles.push_back({ 3, 2, 4 });
  triangles.push_back({ 3, 4, 5 });

  // constructs AABB tree
  Tree tree(triangles.begin(), triangles.end(), points);

  // point sampling
  Point_and_primitive_id id;
  id = tree.closest_point_and_primitive(Point_2(0.5, 0.4));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;
  id = tree.closest_point_and_primitive(Point_2(0.5, 0.6));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;
  id = tree.closest_point_and_primitive(Point_2(1.5, 0.5));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;
  id = tree.closest_point_and_primitive(Point_2(1.5, 0.6));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;
  id = tree.closest_point_and_primitive(Point_2(1.0, 0.0));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;
  id = tree.closest_point_and_primitive(Point_2(3.0, 0.5));
  std::cout << std::distance(triangles.cbegin(), id.second) << ". triangle" << std::endl;

  Ray ray(Point_2(5.5, 0.5), Point_2(1.5, 0.5));
  Ray_intersection intersection = tree.first_intersection(ray);

  if (!intersection) {
    std::cout << "Ray does not intersect with triangles although it should!" << std::endl;
    return EXIT_FAILURE;
  }
  else {
    std::cout << std::distance(triangles.cbegin(), intersection->second) << ". triangle" << std::endl;
  }

  std::list<Ray_intersection> intersections;
  tree.all_intersections(ray, std::back_inserter(intersections));

  return EXIT_SUCCESS;
}