intrinsic_delaunay_triangulation.cpp 8.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199
  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 "EPS.h"
  14. #include "matlab_format.h"
  15. #include <iostream>
  16. #include <queue>
  17. #include <map>
  18. template <
  19. typename Derivedl_in,
  20. typename DerivedF_in,
  21. typename Derivedl,
  22. typename DerivedF>
  23. IGL_INLINE void igl::intrinsic_delaunay_triangulation(
  24. const Eigen::MatrixBase<Derivedl_in> & l_in,
  25. const Eigen::MatrixBase<DerivedF_in> & F_in,
  26. Eigen::PlainObjectBase<Derivedl> & l,
  27. Eigen::PlainObjectBase<DerivedF> & F)
  28. {
  29. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,2> MatrixX2I;
  30. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,1> VectorXI;
  31. MatrixX2I E,uE;
  32. VectorXI EMAP;
  33. std::vector<std::vector<typename DerivedF::Scalar> > uE2E;
  34. return intrinsic_delaunay_triangulation(l_in,F_in,l,F,E,uE,EMAP,uE2E);
  35. }
  36. template <
  37. typename Derivedl_in,
  38. typename DerivedF_in,
  39. typename Derivedl,
  40. typename DerivedF,
  41. typename DerivedE,
  42. typename DeriveduE,
  43. typename DerivedEMAP,
  44. typename uE2EType>
  45. IGL_INLINE void igl::intrinsic_delaunay_triangulation(
  46. const Eigen::MatrixBase<Derivedl_in> & l_in,
  47. const Eigen::MatrixBase<DerivedF_in> & F_in,
  48. Eigen::PlainObjectBase<Derivedl> & l,
  49. Eigen::PlainObjectBase<DerivedF> & F,
  50. Eigen::PlainObjectBase<DerivedE> & E,
  51. Eigen::PlainObjectBase<DeriveduE> & uE,
  52. Eigen::PlainObjectBase<DerivedEMAP> & EMAP,
  53. std::vector<std::vector<uE2EType> > & uE2E)
  54. {
  55. igl::unique_edge_map(F_in, E, uE, EMAP, uE2E);
  56. // We're going to work in place
  57. l = l_in;
  58. F = F_in;
  59. typedef typename DerivedF::Scalar Index;
  60. typedef typename Derivedl::Scalar Scalar;
  61. const Index num_faces = F.rows();
  62. // Vector is faster than queue...
  63. std::vector<Index> Q;
  64. Q.reserve(uE2E.size());
  65. for (size_t uei=0; uei<uE2E.size(); uei++)
  66. {
  67. Q.push_back(uei);
  68. }
  69. // I tried using a "delaunay_since = iter" flag to avoid duplicates, but there
  70. // was no speed up.
  71. while(!Q.empty())
  72. {
  73. #ifdef IGL_INTRINSIC_DELAUNAY_TRIANGULATION_DEBUG
  74. // Expensive sanity check
  75. {
  76. Eigen::Matrix<bool,Eigen::Dynamic,1> inQ(uE2E.size(),1);
  77. inQ.setConstant(false);
  78. for(const auto uei : Q)
  79. {
  80. inQ(uei) = true;
  81. }
  82. for (Index uei=0; uei<uE2E.size(); uei++)
  83. {
  84. if(!inQ(uei) && !is_intrinsic_delaunay(l,uE2E,num_faces,uei))
  85. {
  86. std::cout<<" "<<uei<<" might never be fixed!"<<std::endl;
  87. }
  88. }
  89. }
  90. #endif
  91. const Index uei = Q.back();
  92. Q.pop_back();
  93. if (uE2E[uei].size() == 2)
  94. {
  95. if(!is_intrinsic_delaunay(l,uE2E,num_faces,uei))
  96. {
  97. // update l just before flipping edge
  98. // . //
  99. // /|\ //
  100. // a/ | \d //
  101. // / e \ //
  102. // / | \ //
  103. // .----|-f--. //
  104. // \ | / //
  105. // \ | / //
  106. // b\α|δ/c //
  107. // \|/ //
  108. // . //
  109. // Annotated from flip_edge:
  110. // Edge to flip [v1,v2] --> [v3,v4]
  111. // Before:
  112. // F(f1,:) = [v1,v2,v4] // in some cyclic order
  113. // F(f2,:) = [v1,v3,v2] // in some cyclic order
  114. // After:
  115. // F(f1,:) = [v1,v3,v4] // in *this* order
  116. // F(f2,:) = [v2,v4,v3] // in *this* order
  117. //
  118. // v1 v1
  119. // /|\ / \
  120. // c/ | \b c/f1 \b
  121. // v3 /f2|f1\ v4 => v3 /__f__\ v4
  122. // \ e / \ f2 /
  123. // d\ | /a d\ /a
  124. // \|/ \ /
  125. // v2 v2
  126. //
  127. // Compute intrinsic length of oppposite edge
  128. assert(uE2E[uei].size() == 2 && "edge should have 2 incident faces");
  129. const Index f1 = uE2E[uei][0]%num_faces;
  130. const Index f2 = uE2E[uei][1]%num_faces;
  131. const Index c1 = uE2E[uei][0]/num_faces;
  132. const Index c2 = uE2E[uei][1]/num_faces;
  133. assert(c1 < 3);
  134. assert(c2 < 3);
  135. assert(f1 != f2);
  136. const Index v1 = F(f1, (c1+1)%3);
  137. const Index v2 = F(f1, (c1+2)%3);
  138. const Index v4 = F(f1, c1);
  139. const Index v3 = F(f2, c2);
  140. assert(F(f2, (c2+2)%3) == v1);
  141. assert(F(f2, (c2+1)%3) == v2);
  142. assert( std::abs(l(f1,c1)-l(f2,c2)) < igl::EPS<Scalar>() );
  143. const Scalar e = l(f1,c1);
  144. const Scalar a = l(f1,(c1+1)%3);
  145. const Scalar b = l(f1,(c1+2)%3);
  146. const Scalar c = l(f2,(c2+1)%3);
  147. const Scalar d = l(f2,(c2+2)%3);
  148. // tan(α/2)
  149. const Scalar tan_a_2= tan_half_angle(a,b,e);
  150. // tan(δ/2)
  151. const Scalar tan_d_2 = tan_half_angle(d,e,c);
  152. // tan((α+δ)/2)
  153. const Scalar tan_a_d_2 = (tan_a_2 + tan_d_2)/(1.0-tan_a_2*tan_d_2);
  154. // cos(α+δ)
  155. const Scalar cos_a_d =
  156. (1.0 - tan_a_d_2*tan_a_d_2)/(1.0+tan_a_d_2*tan_a_d_2);
  157. const Scalar f = sqrt(b*b + c*c - 2.0*b*c*cos_a_d);
  158. l(f1,0) = f;
  159. l(f1,1) = b;
  160. l(f1,2) = c;
  161. l(f2,0) = f;
  162. l(f2,1) = d;
  163. l(f2,2) = a;
  164. // Important to grab these indices _before_ calling flip_edges (they
  165. // will be correct after)
  166. const size_t e_24 = f1 + ((c1 + 1) % 3) * num_faces;
  167. const size_t e_41 = f1 + ((c1 + 2) % 3) * num_faces;
  168. const size_t e_13 = f2 + ((c2 + 1) % 3) * num_faces;
  169. const size_t e_32 = f2 + ((c2 + 2) % 3) * num_faces;
  170. const size_t ue_24 = EMAP(e_24);
  171. const size_t ue_41 = EMAP(e_41);
  172. const size_t ue_13 = EMAP(e_13);
  173. const size_t ue_32 = EMAP(e_32);
  174. flip_edge(F, E, uE, EMAP, uE2E, uei);
  175. Q.push_back(ue_24);
  176. Q.push_back(ue_41);
  177. Q.push_back(ue_13);
  178. Q.push_back(ue_32);
  179. }
  180. }
  181. }
  182. }
  183. #ifdef IGL_STATIC_LIBRARY
  184. // Explicit template instantiation
  185. // generated by autoexplicit.sh
  186. template void igl::intrinsic_delaunay_triangulation<Eigen::Matrix<double, -1, 3, 0, -1, 3>, 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, 3, 0, -1, 3> > 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> >&);
  187. // generated by autoexplicit.sh
  188. template void igl::intrinsic_delaunay_triangulation<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  189. // generated by autoexplicit.sh
  190. 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::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, int>(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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > >&);
  191. // generated by autoexplicit.sh
  192. 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> >&);
  193. #endif