cotmatrix.cpp 6.5 KB

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