internal_angles.cpp 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
  4. // Copyright (C) 2015 Daniele Panozzo <daniele.panozzo@gmail.com>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla Public License
  7. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  8. // obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "internal_angles.h"
  10. #include "squared_edge_lengths.h"
  11. #include "parallel_for.h"
  12. #include "get_seconds.h"
  13. template <typename DerivedV, typename DerivedF, typename DerivedK>
  14. IGL_INLINE void igl::internal_angles(
  15. const Eigen::MatrixBase<DerivedV>& V,
  16. const Eigen::MatrixBase<DerivedF>& F,
  17. Eigen::PlainObjectBase<DerivedK> & K)
  18. {
  19. using namespace Eigen;
  20. using namespace std;
  21. typedef typename DerivedV::Scalar Scalar;
  22. if(F.cols() == 3)
  23. {
  24. // Edge lengths
  25. Matrix<
  26. Scalar,
  27. DerivedF::RowsAtCompileTime,
  28. DerivedF::ColsAtCompileTime> L_sq;
  29. igl::squared_edge_lengths(V,F,L_sq);
  30. assert(F.cols() == 3 && "F should contain triangles");
  31. igl::internal_angles_using_squared_edge_lengths(L_sq,K);
  32. }else
  33. {
  34. assert(V.cols() == 3 && "If F contains non-triangle facets, V must be 3D");
  35. K.resizeLike(F);
  36. auto corner = [](
  37. const typename DerivedV::ConstRowXpr & x,
  38. const typename DerivedV::ConstRowXpr & y,
  39. const typename DerivedV::ConstRowXpr & z)
  40. {
  41. typedef Eigen::Matrix<Scalar,1,3> RowVector3S;
  42. RowVector3S v1 = (x-y).normalized();
  43. RowVector3S v2 = (z-y).normalized();
  44. // http://stackoverflow.com/questions/10133957/signed-angle-between-two-vectors-without-a-reference-plane
  45. Scalar s = v1.cross(v2).norm();
  46. Scalar c = v1.dot(v2);
  47. return atan2(s, c);
  48. };
  49. for(unsigned i=0; i<F.rows(); ++i)
  50. {
  51. for(unsigned j=0; j<F.cols(); ++j)
  52. {
  53. K(i,j) = corner(
  54. V.row(F(i,int(j-1+F.cols())%F.cols())),
  55. V.row(F(i,j)),
  56. V.row(F(i,(j+1+F.cols())%F.cols()))
  57. );
  58. }
  59. }
  60. }
  61. }
  62. template <typename DerivedL, typename DerivedK>
  63. IGL_INLINE void igl::internal_angles_using_squared_edge_lengths(
  64. const Eigen::MatrixBase<DerivedL>& L_sq,
  65. Eigen::PlainObjectBase<DerivedK> & K)
  66. {
  67. typedef typename DerivedL::Index Index;
  68. assert(L_sq.cols() == 3 && "Edge-lengths should come from triangles");
  69. const Index m = L_sq.rows();
  70. K.resize(m,3);
  71. parallel_for(
  72. m,
  73. [&L_sq,&K](const Index f)
  74. {
  75. for(size_t d = 0;d<3;d++)
  76. {
  77. const auto & s1 = L_sq(f,d);
  78. const auto & s2 = L_sq(f,(d+1)%3);
  79. const auto & s3 = L_sq(f,(d+2)%3);
  80. K(f,d) = acos((s3 + s2 - s1)/(2.*sqrt(s3*s2)));
  81. }
  82. },
  83. 1000l);
  84. }
  85. template <typename DerivedL, typename DerivedK>
  86. IGL_INLINE void igl::internal_angles_using_edge_lengths(
  87. const Eigen::MatrixBase<DerivedL>& L,
  88. Eigen::PlainObjectBase<DerivedK> & K)
  89. {
  90. // Note:
  91. // Usage of internal_angles_using_squared_edge_lengths() is preferred to internal_angles_using_squared_edge_lengths()
  92. // This function is deprecated and probably will be removed in future versions
  93. typedef typename DerivedL::Index Index;
  94. assert(L.cols() == 3 && "Edge-lengths should come from triangles");
  95. const Index m = L.rows();
  96. K.resize(m,3);
  97. parallel_for(
  98. m,
  99. [&L,&K](const Index f)
  100. {
  101. for(size_t d = 0;d<3;d++)
  102. {
  103. const auto & s1 = L(f,d);
  104. const auto & s2 = L(f,(d+1)%3);
  105. const auto & s3 = L(f,(d+2)%3);
  106. K(f,d) = acos((s3*s3 + s2*s2 - s1*s1)/(2.*s3*s2));
  107. }
  108. },
  109. 1000l);
  110. }
  111. #ifdef IGL_STATIC_LIBRARY
  112. // Explicit template instantiation
  113. // generated by autoexplicit.sh
  114. template void igl::internal_angles<Eigen::Matrix<float, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<float, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> >&);
  115. // generated by autoexplicit.sh
  116. template void igl::internal_angles<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, Eigen::Matrix<float, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> >&);
  117. // generated by autoexplicit.sh
  118. template void igl::internal_angles<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, 3, 1, -1, 3>, Eigen::Matrix<float, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> >&);
  119. template void igl::internal_angles<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::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> >&);
  120. template void igl::internal_angles<Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  121. template void igl::internal_angles<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  122. template void igl::internal_angles<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -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<double, -1, 3, 0, -1, 3> >&);
  123. template void igl::internal_angles<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&);
  124. template void igl::internal_angles<Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&);
  125. #endif