projected_delaunay.cpp 4.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "projected_delaunay.h"
  9. #include "../REDRUM.h"
  10. #include <iostream>
  11. #include <cassert>
  12. template <typename Kernel>
  13. IGL_INLINE void igl::cgal::projected_delaunay(
  14. const CGAL::Triangle_3<Kernel> & A,
  15. const std::vector<CGAL::Object> & A_objects_3,
  16. CGAL::Constrained_triangulation_plus_2<
  17. CGAL::Constrained_Delaunay_triangulation_2<
  18. Kernel,
  19. CGAL::Triangulation_data_structure_2<
  20. CGAL::Triangulation_vertex_base_2<Kernel>,
  21. CGAL::Constrained_triangulation_face_base_2<Kernel> >,
  22. CGAL::Exact_intersections_tag> > & cdt)
  23. {
  24. using namespace std;
  25. // 3D Primitives
  26. typedef CGAL::Point_3<Kernel> Point_3;
  27. typedef CGAL::Segment_3<Kernel> Segment_3;
  28. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  29. typedef CGAL::Plane_3<Kernel> Plane_3;
  30. typedef CGAL::Tetrahedron_3<Kernel> Tetrahedron_3;
  31. typedef CGAL::Point_2<Kernel> Point_2;
  32. typedef CGAL::Segment_2<Kernel> Segment_2;
  33. typedef CGAL::Triangle_2<Kernel> Triangle_2;
  34. typedef CGAL::Triangulation_vertex_base_2<Kernel> TVB_2;
  35. typedef CGAL::Constrained_triangulation_face_base_2<Kernel> CTFB_2;
  36. typedef CGAL::Triangulation_data_structure_2<TVB_2,CTFB_2> TDS_2;
  37. typedef CGAL::Exact_intersections_tag Itag;
  38. typedef CGAL::Constrained_Delaunay_triangulation_2<Kernel,TDS_2,Itag>
  39. CDT_2;
  40. typedef CGAL::Constrained_triangulation_plus_2<CDT_2> CDT_plus_2;
  41. // http://www.cgal.org/Manual/3.2/doc_html/cgal_manual/Triangulation_2/Chapter_main.html#Section_2D_Triangulations_Constrained_Plus
  42. // Plane of triangle A
  43. Plane_3 P(A.vertex(0),A.vertex(1),A.vertex(2));
  44. // Insert triangle into vertices
  45. typename CDT_plus_2::Vertex_handle corners[3];
  46. typedef size_t Index;
  47. for(Index i = 0;i<3;i++)
  48. {
  49. const Point_3 & p3 = A.vertex(i);
  50. const Point_2 & p2 = P.to_2d(p3);
  51. typename CDT_plus_2::Vertex_handle corner = cdt.insert(p2);
  52. corners[i] = corner;
  53. }
  54. // Insert triangle edges as constraints
  55. for(Index i = 0;i<3;i++)
  56. {
  57. cdt.insert_constraint( corners[(i+1)%3], corners[(i+2)%3]);
  58. }
  59. // Insert constraints for intersection objects
  60. for( const auto & obj : A_objects_3)
  61. {
  62. if(const Segment_3 *iseg = CGAL::object_cast<Segment_3 >(&obj))
  63. {
  64. // Add segment constraint
  65. cdt.insert_constraint(P.to_2d(iseg->vertex(0)),P.to_2d(iseg->vertex(1)));
  66. }else if(const Point_3 *ipoint = CGAL::object_cast<Point_3 >(&obj))
  67. {
  68. // Add point
  69. cdt.insert(P.to_2d(*ipoint));
  70. } else if(const Triangle_3 *itri = CGAL::object_cast<Triangle_3 >(&obj))
  71. {
  72. // Add 3 segment constraints
  73. cdt.insert_constraint(P.to_2d(itri->vertex(0)),P.to_2d(itri->vertex(1)));
  74. cdt.insert_constraint(P.to_2d(itri->vertex(1)),P.to_2d(itri->vertex(2)));
  75. cdt.insert_constraint(P.to_2d(itri->vertex(2)),P.to_2d(itri->vertex(0)));
  76. } else if(const std::vector<Point_3 > *polyp =
  77. CGAL::object_cast< std::vector<Point_3 > >(&obj))
  78. {
  79. //cerr<<REDRUM("Poly...")<<endl;
  80. const std::vector<Point_3 > & poly = *polyp;
  81. const Index m = poly.size();
  82. assert(m>=2);
  83. for(Index p = 0;p<m;p++)
  84. {
  85. const Index np = (p+1)%m;
  86. cdt.insert_constraint(P.to_2d(poly[p]),P.to_2d(poly[np]));
  87. }
  88. }else
  89. {
  90. cerr<<REDRUM("What is this object?!")<<endl;
  91. assert(false);
  92. }
  93. }
  94. }
  95. #ifdef IGL_STATIC_LIBRARY
  96. // Explicit template specialization
  97. template void igl::cgal::projected_delaunay<CGAL::Epeck>(CGAL::Triangle_3<CGAL::Epeck> const&, std::__1::vector<CGAL::Object, std::__1::allocator<CGAL::Object> > const&, CGAL::Constrained_triangulation_plus_2<CGAL::Constrained_Delaunay_triangulation_2<CGAL::Epeck, CGAL::Triangulation_data_structure_2<CGAL::Triangulation_vertex_base_2<CGAL::Epeck, CGAL::Triangulation_ds_vertex_base_2<void> >, CGAL::Constrained_triangulation_face_base_2<CGAL::Epeck, CGAL::Triangulation_face_base_2<CGAL::Epeck, CGAL::Triangulation_ds_face_base_2<void> > > >, CGAL::Exact_intersections_tag> >&);
  98. template void igl::cgal::projected_delaunay<CGAL::Epick>(CGAL::Triangle_3<CGAL::Epick> const&, std::__1::vector<CGAL::Object, std::__1::allocator<CGAL::Object> > const&, CGAL::Constrained_triangulation_plus_2<CGAL::Constrained_Delaunay_triangulation_2<CGAL::Epick, CGAL::Triangulation_data_structure_2<CGAL::Triangulation_vertex_base_2<CGAL::Epick, CGAL::Triangulation_ds_vertex_base_2<void> >, CGAL::Constrained_triangulation_face_base_2<CGAL::Epick, CGAL::Triangulation_face_base_2<CGAL::Epick, CGAL::Triangulation_ds_face_base_2<void> > > >, CGAL::Exact_intersections_tag> >&);
  99. #endif