cotmatrix.h 3.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129
  1. #ifndef IGL_COTMATRIX_H
  2. #define IGL_COTMATRIX_H
  3. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  4. #include <Eigen/Dense>
  5. #include <Eigen/Sparse>
  6. // History:
  7. // Used const references rather than copying the entire mesh
  8. // Alec 9 October 2011
  9. // removed cotan (uniform weights) optional parameter it was building a buggy
  10. // half of the uniform laplacian, please see adjacency_matrix istead
  11. // Alec 9 October 2011
  12. namespace igl
  13. {
  14. // Constructs the cotangent stiffness matrix (discrete laplacian) for a given
  15. // mesh (V,F).
  16. // Inputs:
  17. // V #V by dim list of mesh vertex positions
  18. // F #F by 3 list of mesh faces (must be triangles)
  19. // Outputs:
  20. // L #V by #V cotangent matrix, each row i corresponding to V(i,:)
  21. //
  22. // See also: adjacency_matrix
  23. void cotmatrix(
  24. const Eigen::MatrixXd & V,
  25. const Eigen::MatrixXi & F,
  26. Eigen::SparseMatrix<double>& L);
  27. // Helper function that computes the cotangent weights for each corner of a
  28. // given triangle
  29. // Inputs:
  30. // v1 position of corner #1 of triangle
  31. // v2 position of corner #2 of triangle
  32. // v3 position of corner #3 of triangle
  33. // Outputs:
  34. // cot1 cotangent of angle at corner #1
  35. // cot2 cotangent of angle at corner #2
  36. // cot3 cotangent of angle at corner #3
  37. //
  38. void computeCotWeights(
  39. const Eigen::Vector3d& v1,
  40. const Eigen::Vector3d& v2,
  41. const Eigen::Vector3d& v3,
  42. double& cot1,
  43. double& cot2,
  44. double& cot3);
  45. }
  46. // Implementation
  47. // For error printing
  48. #include <cstdio>
  49. void igl::computeCotWeights(
  50. const Eigen::Vector3d& v1,
  51. const Eigen::Vector3d& v2,
  52. const Eigen::Vector3d& v3,
  53. double& cot1,
  54. double& cot2,
  55. double& cot3)
  56. {
  57. Eigen::Vector3d v12 = v2-v1;
  58. Eigen::Vector3d v13 = v3-v1;
  59. Eigen::Vector3d v23 = v3-v2;
  60. double halfArea = (v12.cross(v13)).norm();//squaredNorm();
  61. //improve numerical stability
  62. const double cotTolerance = 1e-10;
  63. if(halfArea < cotTolerance)
  64. {
  65. fprintf(stderr,"Cot weights are close to singular!\n");
  66. halfArea = cotTolerance;
  67. }
  68. cot1 = (v12.dot(v13)) / halfArea /2;
  69. cot2 = (v23.dot(-v12)) / halfArea /2;
  70. cot3 = (-v23.dot(-v13)) / halfArea /2;
  71. }
  72. void igl::cotmatrix(
  73. const Eigen::MatrixXd & V,
  74. const Eigen::MatrixXi & F,
  75. Eigen::SparseMatrix<double>& L)
  76. {
  77. // Assumes vertices are 3D
  78. assert(V.cols() == 3);
  79. // Assumes faces are triangles
  80. assert(F.cols() == 3);
  81. Eigen::DynamicSparseMatrix<double, Eigen::RowMajor> dyn_L (V.rows(), V.rows());
  82. // Loop over triangles
  83. for (unsigned i = 0; i < F.rows(); i++)
  84. {
  85. // Corner indices of this triangle
  86. int vi1 = F(i,0);
  87. int vi2 = F(i,1);
  88. int vi3 = F(i,2);
  89. // Grab corner positions of this triangle
  90. Eigen::Vector3d v1(V(vi1,0), V(vi1,1), V(vi1,2));
  91. Eigen::Vector3d v2(V(vi2,0), V(vi2,1), V(vi2,2));
  92. Eigen::Vector3d v3(V(vi3,0), V(vi3,1), V(vi3,2));
  93. // Compute cotangent of angles at each corner
  94. double cot1, cot2, cot3;
  95. computeCotWeights(v1, v2, v3, cot1, cot2, cot3);
  96. // Throw each corner's cotangent at opposite edge (in both directions)
  97. dyn_L.coeffRef(vi1, vi2) += cot3;
  98. dyn_L.coeffRef(vi2, vi1) += cot3;
  99. dyn_L.coeffRef(vi2, vi3) += cot1;
  100. dyn_L.coeffRef(vi3, vi2) += cot1;
  101. dyn_L.coeffRef(vi3, vi1) += cot2;
  102. dyn_L.coeffRef(vi1, vi3) += cot2;
  103. }
  104. for (int k=0; k < dyn_L.outerSize(); ++k)
  105. {
  106. double tmp = 0.0f;
  107. for(Eigen::DynamicSparseMatrix<double, Eigen::RowMajor>::InnerIterator it (dyn_L, k); it; ++it)
  108. {
  109. tmp += it.value();
  110. }
  111. dyn_L.coeffRef(k,k) = -tmp;
  112. }
  113. L = Eigen::SparseMatrix<double>(dyn_L);
  114. }
  115. #endif