cotmatrix_entries.cpp 7.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <[email protected]>
  4. // Copyright (C) 2018 Alec Jacobson <[email protected]>
  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 "cotmatrix_entries.h"
  10. #include "doublearea.h"
  11. #include "squared_edge_lengths.h"
  12. #include "edge_lengths.h"
  13. #include "face_areas.h"
  14. #include "volume.h"
  15. #include "dihedral_angles.h"
  16. #include "dihedral_angles_intrinsic.h"
  17. #include "verbose.h"
  18. template <typename DerivedV, typename DerivedF, typename DerivedC>
  19. IGL_INLINE void igl::cotmatrix_entries(
  20. const Eigen::MatrixBase<DerivedV>& V,
  21. const Eigen::MatrixBase<DerivedF>& F,
  22. Eigen::PlainObjectBase<DerivedC>& C)
  23. {
  24. // simplex size (3: triangles, 4: tetrahedra)
  25. int simplex_size = F.cols();
  26. // Number of elements
  27. int m = F.rows();
  28. // Law of cosines + law of sines
  29. switch(simplex_size)
  30. {
  31. case 3:
  32. {
  33. // Triangles
  34. //Compute Squared Edge lengths
  35. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,3> l2;
  36. igl::squared_edge_lengths(V,F,l2);
  37. //Compute Edge lengths
  38. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,3> l;
  39. l = l2.array().sqrt();
  40. // double area
  41. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,1> dblA;
  42. doublearea(l,0.,dblA);
  43. // cotangents and diagonal entries for element matrices
  44. // correctly divided by 4 (alec 2010)
  45. C.resize(m,3);
  46. for(int i = 0;i<m;i++)
  47. {
  48. // Alec: I'm doubtful that using l2 here is actually improving numerics.
  49. C(i,0) = (l2(i,1) + l2(i,2) - l2(i,0))/dblA(i)/4.0;
  50. C(i,1) = (l2(i,2) + l2(i,0) - l2(i,1))/dblA(i)/4.0;
  51. C(i,2) = (l2(i,0) + l2(i,1) - l2(i,2))/dblA(i)/4.0;
  52. }
  53. break;
  54. }
  55. case 4:
  56. {
  57. // edge lengths numbered same as opposite vertices
  58. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,6> l;
  59. edge_lengths(V,F,l);
  60. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,4> s;
  61. face_areas(l,s);
  62. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,6> cos_theta,theta;
  63. dihedral_angles_intrinsic(l,s,theta,cos_theta);
  64. // volume
  65. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,1> vol;
  66. volume(l,vol);
  67. // Law of sines
  68. // http://mathworld.wolfram.com/Tetrahedron.html
  69. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,6> sin_theta(m,6);
  70. sin_theta.col(0) = vol.array() / ((2./(3.*l.col(0).array())).array() * s.col(1).array() * s.col(2).array());
  71. sin_theta.col(1) = vol.array() / ((2./(3.*l.col(1).array())).array() * s.col(2).array() * s.col(0).array());
  72. sin_theta.col(2) = vol.array() / ((2./(3.*l.col(2).array())).array() * s.col(0).array() * s.col(1).array());
  73. sin_theta.col(3) = vol.array() / ((2./(3.*l.col(3).array())).array() * s.col(3).array() * s.col(0).array());
  74. sin_theta.col(4) = vol.array() / ((2./(3.*l.col(4).array())).array() * s.col(3).array() * s.col(1).array());
  75. sin_theta.col(5) = vol.array() / ((2./(3.*l.col(5).array())).array() * s.col(3).array() * s.col(2).array());
  76. // http://arxiv.org/pdf/1208.0354.pdf Page 18
  77. C = (1./6.) * l.array() * cos_theta.array() / sin_theta.array();
  78. break;
  79. }
  80. default:
  81. {
  82. fprintf(stderr,
  83. "cotmatrix_entries.h: Error: Simplex size (%d) not supported\n", simplex_size);
  84. assert(false);
  85. }
  86. }
  87. }
  88. template <typename Derivedl, typename DerivedC>
  89. IGL_INLINE void igl::cotmatrix_entries(
  90. const Eigen::MatrixBase<Derivedl>& l,
  91. Eigen::PlainObjectBase<DerivedC>& C)
  92. {
  93. const int m = l.rows();
  94. assert(l.cols() == 3 && "Only triangles accepted");
  95. //Compute squared Edge lengths
  96. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,3> l2;
  97. l2 = l.array().square();
  98. // Alec: It's a little annoying that there's duplicate code here. The
  99. // "extrinic" version above is first computing squared edge lengths, taking
  100. // the square root and calling this. We can't have a cotmatrix_entries(l,l2,C)
  101. // overload because it will confuse Eigen with the cotmatrix_entries(V,F,C)
  102. // overload. In the end, I'd like to be convinced that using l2 directly above
  103. // is actually better numerically (or significantly faster) than just calling
  104. // edge_lengths and this cotmatrix_entries(l,C);
  105. //
  106. // double area
  107. Eigen::Matrix<typename DerivedC::Scalar ,Eigen::Dynamic,1> dblA;
  108. doublearea(l,0.,dblA);
  109. // cotangents and diagonal entries for element matrices
  110. // correctly divided by 4 (alec 2010)
  111. C.resize(m,3);
  112. for(int i = 0;i<m;i++)
  113. {
  114. // Alec: I'm doubtful that using l2 here is actually improving numerics.
  115. C(i,0) = (l2(i,1) + l2(i,2) - l2(i,0))/dblA(i)/4.0;
  116. C(i,1) = (l2(i,2) + l2(i,0) - l2(i,1))/dblA(i)/4.0;
  117. C(i,2) = (l2(i,0) + l2(i,1) - l2(i,2))/dblA(i)/4.0;
  118. }
  119. }
  120. #ifdef IGL_STATIC_LIBRARY
  121. // Explicit template instantiation
  122. // generated by autoexplicit.sh
  123. template void igl::cotmatrix_entries<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::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> >&);
  124. // generated by autoexplicit.sh
  125. template void igl::cotmatrix_entries<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  126. // generated by autoexplicit.sh
  127. template void igl::cotmatrix_entries<Eigen::Matrix<double, -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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  128. template void igl::cotmatrix_entries<Eigen::Matrix<double, -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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  129. template void igl::cotmatrix_entries<Eigen::Matrix<double, -1, -1, 1, -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, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  130. template void igl::cotmatrix_entries<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> >&);
  131. template void igl::cotmatrix_entries<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 4, 0, -1, 4>, 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, 4, 0, -1, 4> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  132. template void igl::cotmatrix_entries<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> >&);
  133. #endif