intrinsic_delaunay_triangulation.cpp 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2018 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 "intrinsic_delaunay_triangulation.h"
  9. #include "is_intrinsic_delaunay.h"
  10. #include "tan_half_angle.h"
  11. #include "unique_edge_map.h"
  12. #include "flip_edge.h"
  13. #include <iostream>
  14. template <
  15. typename Derivedl_in,
  16. typename DerivedF_in,
  17. typename Derivedl,
  18. typename DerivedF>
  19. IGL_INLINE void igl::intrinsic_delaunay_triangulation(
  20. const Eigen::MatrixBase<Derivedl_in> & l_in,
  21. const Eigen::MatrixBase<DerivedF_in> & F_in,
  22. Eigen::PlainObjectBase<Derivedl> & l,
  23. Eigen::PlainObjectBase<DerivedF> & F)
  24. {
  25. // We're going to work in place
  26. l = l_in;
  27. F = F_in;
  28. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,2> MatrixX2I;
  29. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,1> VectorXI;
  30. MatrixX2I E,uE;
  31. VectorXI EMAP;
  32. std::vector<std::vector<typename DerivedF::Scalar> > uE2E;
  33. igl::unique_edge_map(F, E, uE, EMAP, uE2E);
  34. typedef typename DerivedF::Scalar Index;
  35. typedef typename Derivedl::Scalar Scalar;
  36. const Index num_faces = F.rows();
  37. // Copied from delaunay_triangulation
  38. bool all_delaunay = false;
  39. while(!all_delaunay)
  40. {
  41. all_delaunay = true;
  42. for (size_t uei=0; uei<uE2E.size(); uei++)
  43. {
  44. if (uE2E[uei].size() == 2)
  45. {
  46. if(!is_intrinsic_delaunay(l,F,uE2E,uei))
  47. {
  48. all_delaunay = false;
  49. // update l just before flipping edge
  50. // . //
  51. // /|\ //
  52. // a/ | \d //
  53. // / e \ //
  54. // / | \ //
  55. // .----|-f--. //
  56. // \ | / //
  57. // \ | / //
  58. // b\α|δ/c //
  59. // \|/ //
  60. // . //
  61. // Compute intrinsic length of oppposite edge
  62. assert(uE2E[uei].size() == 2 && "edge should have 2 incident faces");
  63. const Index he1 = uE2E[uei][0];
  64. const Index he2 = uE2E[uei][1];
  65. const Index f1 = he1%num_faces;
  66. const Index c1 = he1/num_faces;
  67. const Index f2 = he2%num_faces;
  68. const Index c2 = he2/num_faces;
  69. assert( std::abs(l(f1,c1)-l(f2,c2) < igl::EPS<Scalar>()) );
  70. const Scalar e = l(f1,c1);
  71. const Scalar a = l(f1,(c1+1)%3);
  72. const Scalar b = l(f1,(c1+2)%3);
  73. const Scalar c = l(f2,(c2+1)%3);
  74. const Scalar d = l(f2,(c2+2)%3);
  75. // tan(α/2)
  76. const Scalar tan_a_2= tan_half_angle(a,b,e);
  77. // tan(δ/2)
  78. const Scalar tan_d_2 = tan_half_angle(d,e,c);
  79. // tan((α+δ)/2)
  80. const Scalar tan_a_d_2 = (tan_a_2 + tan_d_2)/(1.0-tan_a_2*tan_d_2);
  81. // cos(α+δ)
  82. const Scalar cos_a_d =
  83. (1.0 - tan_a_d_2*tan_a_d_2)/(1.0+tan_a_d_2*tan_a_d_2);
  84. const Scalar f = sqrt(b*b + c*c - 2.0*b*c*cos_a_d);
  85. // Annotated from flip_edge:
  86. // Edge to flip [v1,v2] --> [v3,v4]
  87. // Before:
  88. // F(f1,:) = [v1,v2,v4] // in some cyclic order
  89. // F(f2,:) = [v1,v3,v2] // in some cyclic order
  90. // After:
  91. // F(f1,:) = [v1,v3,v4] // in *this* order
  92. // F(f2,:) = [v2,v4,v3] // in *this* order
  93. //
  94. // v1 v1
  95. // /|\ / \
  96. // c/ | \b c/f1 \b
  97. // v3 /f2|f1\ v4 => v3 /__f__\ v4
  98. // \ e / \ f2 /
  99. // d\ | /a d\ /a
  100. // \|/ \ /
  101. // v2 v2
  102. //
  103. l(f1,0) = f;
  104. l(f1,1) = b;
  105. l(f1,2) = c;
  106. l(f2,0) = f;
  107. l(f2,1) = d;
  108. l(f2,2) = a;
  109. flip_edge(F, E, uE, EMAP, uE2E, uei);
  110. }
  111. }
  112. }
  113. }
  114. }
  115. #ifdef IGL_STATIC_LIBRARY
  116. // Explicit template instantiation
  117. // generated by autoexplicit.sh
  118. template void igl::intrinsic_delaunay_triangulation<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  119. #endif