is_intrinsic_delaunay.cpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146
  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 "is_intrinsic_delaunay.h"
  9. #include "unique_edge_map.h"
  10. #include "tan_half_angle.h"
  11. #include "EPS.h"
  12. #include "cotmatrix_entries.h"
  13. #include <iostream>
  14. #include <cassert>
  15. template <
  16. typename Derivedl,
  17. typename DerivedF,
  18. typename DerivedD>
  19. IGL_INLINE void igl::is_intrinsic_delaunay(
  20. const Eigen::MatrixBase<Derivedl> & l,
  21. const Eigen::MatrixBase<DerivedF> & F,
  22. Eigen::PlainObjectBase<DerivedD> & D)
  23. {
  24. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,2> MatrixX2I;
  25. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,1> VectorXI;
  26. MatrixX2I E,uE;
  27. VectorXI EMAP;
  28. std::vector<std::vector<typename DerivedF::Scalar> > uE2E;
  29. igl::unique_edge_map(F, E, uE, EMAP, uE2E);
  30. return is_intrinsic_delaunay(l,F,uE2E,D);
  31. }
  32. template <
  33. typename Derivedl,
  34. typename DerivedF,
  35. typename uE2EType,
  36. typename DerivedD>
  37. IGL_INLINE void igl::is_intrinsic_delaunay(
  38. const Eigen::MatrixBase<Derivedl> & l,
  39. const Eigen::MatrixBase<DerivedF> & F,
  40. const std::vector<std::vector<uE2EType> > & uE2E,
  41. Eigen::PlainObjectBase<DerivedD> & D)
  42. {
  43. const int num_faces = F.rows();
  44. D.setConstant(F.rows(),F.cols(),false);
  45. // loop over all unique edges
  46. for(int ue = 0;ue < uE2E.size(); ue++)
  47. {
  48. const bool ue_is_d = is_intrinsic_delaunay(l,uE2E,num_faces,ue);
  49. // Set for all instances
  50. for(int e = 0;e<uE2E[ue].size();e++)
  51. {
  52. D( uE2E[ue][e]%num_faces, uE2E[ue][e]/num_faces) = ue_is_d;
  53. }
  54. }
  55. };
  56. template <
  57. typename Derivedl,
  58. typename uE2EType,
  59. typename Index>
  60. IGL_INLINE bool igl::is_intrinsic_delaunay(
  61. const Eigen::MatrixBase<Derivedl> & l,
  62. const std::vector<std::vector<uE2EType> > & uE2E,
  63. const Index num_faces,
  64. const Index uei)
  65. {
  66. if(uE2E[uei].size() == 1) return true;
  67. if(uE2E[uei].size() > 2) return false;
  68. typedef typename Derivedl::Scalar Scalar;
  69. const auto cot_alpha = [](
  70. const Scalar & a,
  71. const Scalar & b,
  72. const Scalar & c)->Scalar
  73. {
  74. // Fisher 2007
  75. const Scalar t = tan_half_angle(a,b,c);
  76. return (1.0-t*t)/(2*t);
  77. };
  78. // 2 //
  79. // /|\ //
  80. // a/ | \d //
  81. // / e \ //
  82. // / | \ //
  83. //0.α---|-f-β.3 //
  84. // \ | / //
  85. // \ | / //
  86. // b\ | /c //
  87. // \|/ //
  88. // . //
  89. // 1
  90. // Fisher 2007
  91. assert(uE2E[uei].size() == 2 && "edge should have 2 incident faces");
  92. const Index he1 = uE2E[uei][0];
  93. const Index he2 = uE2E[uei][1];
  94. const Index f1 = he1%num_faces;
  95. const Index c1 = he1/num_faces;
  96. const Index f2 = he2%num_faces;
  97. const Index c2 = he2/num_faces;
  98. assert( std::abs(l(f1,c1)-l(f2,c2)) < igl::EPS<Scalar>());
  99. const Scalar e = l(f1,c1);
  100. const Scalar a = l(f1,(c1+1)%3);
  101. const Scalar b = l(f1,(c1+2)%3);
  102. const Scalar c = l(f2,(c2+1)%3);
  103. const Scalar d = l(f2,(c2+2)%3);
  104. const Scalar w = cot_alpha(e,a,b) + cot_alpha(e,c,d);
  105. //// Test
  106. //{
  107. // Eigen::MatrixXd l(2,3);
  108. // l<<e,a,b,
  109. // d,e,c;
  110. // //Eigen::MatrixXi F(2,3);
  111. // //F<<0,1,2,
  112. // // 1,3,2;
  113. // Eigen::MatrixXd C;
  114. // cotmatrix_entries(l,C);
  115. // if(std::abs(w-(2.0*(C(0,0)+C(1,1))))>1e-10)
  116. // {
  117. // std::cout<<matlab_format(l,"l")<<std::endl;
  118. // std::cout<<w<<std::endl;
  119. // std::cout<<(2.0*(C(0,0)+C(1,1)))<<std::endl;
  120. // std::cout<<w-(2.0*(C(0,0)+C(1,1)))<<std::endl;
  121. // std::cout<<std::endl;
  122. // }
  123. //}
  124. return w >= 0;
  125. }
  126. #ifdef IGL_STATIC_LIBRARY
  127. // Explicit template instantiation
  128. // generated by autoexplicit.sh
  129. template bool igl::is_intrinsic_delaunay<Eigen::Matrix<double, -1, 3, 0, -1, 3>, int, int>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > const&, int, int);
  130. // generated by autoexplicit.sh
  131. template bool igl::is_intrinsic_delaunay<Eigen::Matrix<double, -1, -1, 0, -1, -1>, int, int>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > const&, int, int);
  132. // generated by autoexplicit.sh
  133. template void igl::is_intrinsic_delaunay<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, int, Eigen::Matrix<bool, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 3, 0, -1, 3> >&);
  134. // generated by autoexplicit.sh
  135. template void igl::is_intrinsic_delaunay<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<bool, -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<bool, -1, -1, 0, -1, -1> >&);
  136. // generated by autoexplicit.sh
  137. template void igl::is_intrinsic_delaunay<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<bool, -1, 3, 0, -1, 3> >(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<bool, -1, 3, 0, -1, 3> >&);
  138. #endif