cotmatrix_entries.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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 "cotmatrix_entries.h"
  9. #include "doublearea.h"
  10. #include "squared_edge_lengths.h"
  11. #include "edge_lengths.h"
  12. #include "face_areas.h"
  13. #include "volume.h"
  14. #include "dihedral_angles.h"
  15. #include "verbose.h"
  16. template <typename DerivedV, typename DerivedF, typename DerivedC>
  17. IGL_INLINE void igl::cotmatrix_entries(
  18. const Eigen::MatrixBase<DerivedV>& V,
  19. const Eigen::MatrixBase<DerivedF>& F,
  20. Eigen::PlainObjectBase<DerivedC>& C)
  21. {
  22. using namespace std;
  23. using namespace Eigen;
  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 lenghts
  35. Matrix<typename DerivedC::Scalar,Dynamic,3> l2;
  36. igl::squared_edge_lengths(V,F,l2);
  37. //Compute Edge lenghts
  38. Matrix<typename DerivedC::Scalar,Dynamic,3> l;
  39. l = l2.array().sqrt();
  40. // double area
  41. Matrix<typename DerivedC::Scalar,Dynamic,1> dblA;
  42. doublearea(l,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. C(i,0) = (l2(i,1) + l2(i,2) - l2(i,0))/dblA(i)/4.0;
  49. C(i,1) = (l2(i,2) + l2(i,0) - l2(i,1))/dblA(i)/4.0;
  50. C(i,2) = (l2(i,0) + l2(i,1) - l2(i,2))/dblA(i)/4.0;
  51. }
  52. break;
  53. }
  54. case 4:
  55. {
  56. // edge lengths numbered same as opposite vertices
  57. Matrix<typename DerivedC::Scalar,Dynamic,6> l;
  58. edge_lengths(V,F,l);
  59. Matrix<typename DerivedC::Scalar,Dynamic,4> s;
  60. face_areas(l,s);
  61. Matrix<typename DerivedC::Scalar,Dynamic,6> cos_theta,theta;
  62. dihedral_angles_intrinsic(l,s,theta,cos_theta);
  63. // volume
  64. Matrix<typename DerivedC::Scalar,Dynamic,1> vol;
  65. volume(l,vol);
  66. // Law of sines
  67. // http://mathworld.wolfram.com/Tetrahedron.html
  68. Matrix<typename DerivedC::Scalar,Dynamic,6> sin_theta(m,6);
  69. sin_theta.col(0) = vol.array() / ((2./(3.*l.col(0).array())).array() * s.col(1).array() * s.col(2).array());
  70. sin_theta.col(1) = vol.array() / ((2./(3.*l.col(1).array())).array() * s.col(2).array() * s.col(0).array());
  71. sin_theta.col(2) = vol.array() / ((2./(3.*l.col(2).array())).array() * s.col(0).array() * s.col(1).array());
  72. sin_theta.col(3) = vol.array() / ((2./(3.*l.col(3).array())).array() * s.col(3).array() * s.col(0).array());
  73. sin_theta.col(4) = vol.array() / ((2./(3.*l.col(4).array())).array() * s.col(3).array() * s.col(1).array());
  74. sin_theta.col(5) = vol.array() / ((2./(3.*l.col(5).array())).array() * s.col(3).array() * s.col(2).array());
  75. // http://arxiv.org/pdf/1208.0354.pdf Page 18
  76. C = (1./6.) * l.array() * cos_theta.array() / sin_theta.array();
  77. break;
  78. }
  79. default:
  80. {
  81. fprintf(stderr,
  82. "cotmatrix_entries.h: Error: Simplex size (%d) not supported\n", simplex_size);
  83. assert(false);
  84. }
  85. }
  86. }
  87. #ifdef IGL_STATIC_LIBRARY
  88. // Explicit template specialization
  89. // generated by autoexplicit.sh
  90. 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> >&);
  91. // generated by autoexplicit.sh
  92. 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> >&);
  93. // generated by autoexplicit.sh
  94. 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> >&);
  95. // generated by autoexplicit.sh
  96. 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> >&);
  97. #endif