cotmatrix.cpp 7.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215
  1. #include <test_common.h>
  2. #include <igl/PI.h>
  3. #include <igl/cotmatrix.h>
  4. #include <igl/matrix_to_list.h>
  5. #include <igl/polygon_corners.h>
  6. TEST_CASE("cotmatrix: poly", "[igl]" )
  7. {
  8. const auto test_case = [](const std::string &param)
  9. {
  10. Eigen::MatrixXd V;
  11. Eigen::MatrixXi F;
  12. // Load example mesh: GetParam() will be name of mesh file
  13. igl::read_triangle_mesh(test_common::data_path(param), V, F);
  14. Eigen::SparseMatrix<double> tL,pL,pM,pP;
  15. igl::cotmatrix(V,F,tL);
  16. std::vector<std::vector<int> > vF;
  17. igl::matrix_to_list(F,vF);
  18. // trivial polygon mesh
  19. Eigen::VectorXi I,C;
  20. igl::polygon_corners(vF,I,C);
  21. igl::cotmatrix(V,I,C,pL,pM,pP);
  22. REQUIRE (tL.cols() == pL.cols());
  23. REQUIRE (tL.rows() == pL.rows());
  24. REQUIRE ( tL.isApprox(pL,1e-7) );
  25. };
  26. test_common::run_test_cases(test_common::all_meshes(), test_case);
  27. }
  28. TEST_CASE("cotmatrix: constant_in_null_space", "[igl]" "[slow]")
  29. {
  30. const auto test_case = [](const std::string &param)
  31. {
  32. Eigen::MatrixXd V;
  33. Eigen::MatrixXi F;
  34. Eigen::SparseMatrix<double> L;
  35. // Load example mesh: GetParam() will be name of mesh file
  36. igl::read_triangle_mesh(test_common::data_path(param), V, F);
  37. igl::cotmatrix(V,F,L);
  38. REQUIRE (L.rows() == V.rows());
  39. REQUIRE (L.cols() == L.rows());
  40. Eigen::VectorXd C = Eigen::VectorXd::Ones(L.rows());
  41. Eigen::VectorXd Z = Eigen::VectorXd::Zero(L.rows());
  42. // REQUIRE (b == a);
  43. // REQUIRE (a==b);
  44. // ASSERT_NEAR(a,b,1e-15)
  45. REQUIRE (1e-12 > ((L*C)-(Z)).norm());
  46. };
  47. test_common::run_test_cases(test_common::all_meshes(), test_case);
  48. }
  49. TEST_CASE("cotmatrix: cube", "[igl]")
  50. {
  51. //The allowed error for this test
  52. const double epsilon = 1e-15;
  53. Eigen::MatrixXd V;
  54. Eigen::MatrixXi F;
  55. //This is a cube of dimensions 1.0x1.0x1.0
  56. igl::read_triangle_mesh(test_common::data_path("cube.obj"), V, F);
  57. //Scale the cube to have huge sides
  58. Eigen::MatrixXd V_huge = V * 1.0e8;
  59. //Scale the cube to have tiny sides
  60. Eigen::MatrixXd V_tiny = V * 1.0e-8;
  61. //Check cotmatrix (Laplacian)
  62. //The laplacian for the cube is quite singular.
  63. //Each edge in a diagonal has two opposite angles of 90, with cotangent 0.0 each
  64. //Each edge in a side has two opposite angle of 45, with (half)cotangen 0.5 each
  65. //So the cotangent matrix always are (0+0) or (0.5+0.5)
  66. Eigen::SparseMatrix<double> L1;
  67. igl::cotmatrix(V,F,L1);
  68. REQUIRE (L1.rows() == V.rows());
  69. REQUIRE (L1.cols() == V.rows());
  70. //// This is hitting an Eigen bug. https://github.com/libigl/libigl/pull/1064
  71. // for(int f = 0;f<L1.rows();f++)
  72. // {
  73. //#ifdef IGL_EDGE_LENGTHS_SQUARED_H
  74. // //Hard assert if we have edge_lenght_squared
  75. // REQUIRE (L1.coeff(f,f) == -3.0);
  76. // REQUIRE (L1.row(f).sum() == 0.0);
  77. // REQUIRE (L1.col(f).sum() == 0.0);
  78. //#else
  79. // //Soft assert if we have not edge_lenght_squared
  80. // REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  81. // REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  82. // REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  83. //#endif
  84. // }
  85. Eigen::VectorXd row_sum = L1 * Eigen::VectorXd::Constant(L1.rows(),1,1);
  86. Eigen::RowVectorXd col_sum = Eigen::RowVectorXd::Constant(1,L1.rows(),1) * L1;
  87. Eigen::VectorXd diag = L1.diagonal();
  88. #ifdef IGL_EDGE_LENGTHS_SQUARED_H
  89. test_common::assert_eq( row_sum, Eigen::VectorXd::Zero(L1.rows()) );
  90. test_common::assert_eq( col_sum, Eigen::RowVectorXd::Zero(L1.rows()) );
  91. test_common::assert_eq( diag, Eigen::VectorXd::Constant(L1.rows(),1,-3) );
  92. #else
  93. test_common::assert_near( row_sum, Eigen::VectorXd::Zero(L1.rows()) , epsilon);
  94. test_common::assert_near( col_sum, Eigen::RowVectorXd::Zero(L1.rows()) , epsilon);
  95. test_common::assert_near( diag, Eigen::VectorXd::Constant(L1.rows(),1,-3) , epsilon);
  96. #endif
  97. //Same for huge cube.
  98. igl::cotmatrix(V_huge,F,L1);
  99. REQUIRE (L1.rows() == V.rows());
  100. REQUIRE (L1.cols() == V.rows());
  101. for(int f = 0;f<L1.rows();f++)
  102. {
  103. REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  104. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  105. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  106. }
  107. //Same for tiny cube. we need to use a tolerance this time...
  108. igl::cotmatrix(V_tiny,F,L1);
  109. REQUIRE (L1.rows() == V.rows());
  110. REQUIRE (L1.cols() == V.rows());
  111. for(int f = 0;f<L1.rows();f++)
  112. {
  113. REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  114. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  115. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  116. }
  117. }
  118. TEST_CASE("cotmatrix: tetrahedron", "[igl]")
  119. {
  120. //The allowed error for this test
  121. const double epsilon = 1e-15;
  122. Eigen::MatrixXd V;
  123. Eigen::MatrixXi F;
  124. //This is a cube of dimensions 1.0x1.0x1.0
  125. igl::read_triangle_mesh(test_common::data_path("cube.obj"), V, F);
  126. //Prepare another mesh with triangles along side diagonals of the cube
  127. //These triangles are form a regular tetrahedron of side sqrt(2)
  128. Eigen::MatrixXi F_equi(4,3);
  129. F_equi << 4,6,1,
  130. 6,4,3,
  131. 4,1,3,
  132. 1,6,3;
  133. //Scale the cube to have huge sides
  134. Eigen::MatrixXd V_huge = V * 1.0e8;
  135. //Scale the cube to have tiny sides
  136. Eigen::MatrixXd V_tiny = V * 1.0e-8;
  137. //Check cotmatrix (Laplacian)
  138. //The laplacian for the cube is quite singular.
  139. //Each edge in a diagonal has two opposite angles of 90, with cotangent 0.0 each
  140. //Each edge in a side has two opposite angle of 45, with (half)cotangen 0.5 each
  141. //So the cotangent matrix always are (0+0) or (0.5+0.5)
  142. Eigen::SparseMatrix<double> L1;
  143. //Check the regular tetrahedron of side sqrt(2)
  144. igl::cotmatrix(V,F_equi,L1);
  145. REQUIRE (L1.rows() == V.rows());
  146. REQUIRE (L1.cols() == V.rows());
  147. for(int f = 0;f<L1.rows();f++)
  148. {
  149. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  150. if (L1.coeff(f,f) < -0.1)
  151. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(igl::PI / 3.0)).margin( epsilon));
  152. else
  153. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  154. #ifdef IGL_EDGE_LENGTHS_SQUARED_H
  155. //Hard assert if we have edge_lenght_squared
  156. REQUIRE (L1.row(f).sum() == 0.0);
  157. REQUIRE (L1.col(f).sum() == 0.0);
  158. #else
  159. //Soft assert if we have not edge_lenght_squared
  160. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  161. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  162. #endif
  163. }
  164. //Check the huge regular tetrahedron
  165. igl::cotmatrix(V_huge,F_equi,L1);
  166. REQUIRE (L1.rows() == V.rows());
  167. REQUIRE (L1.cols() == V.rows());
  168. for(int f = 0;f<L1.rows();f++)
  169. {
  170. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  171. if (L1.coeff(f,f) < -0.1)
  172. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(igl::PI / 3.0)).margin( epsilon));
  173. else
  174. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  175. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  176. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  177. }
  178. //Check the tiny regular tetrahedron
  179. igl::cotmatrix(V_tiny,F_equi,L1);
  180. REQUIRE (L1.rows() == V.rows());
  181. REQUIRE (L1.cols() == V.rows());
  182. for(int f = 0;f<L1.rows();f++)
  183. {
  184. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  185. if (L1.coeff(f,f) < -0.1)
  186. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(igl::PI / 3.0)).margin( epsilon));
  187. else
  188. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  189. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  190. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  191. }
  192. }