cotmatrix.cpp 5.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176
  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. for(int f = 0;f<L1.rows();f++)
  46. {
  47. #ifdef IGL_EDGE_LENGTHS_SQUARED_H
  48. //Hard assert if we have edge_lenght_squared
  49. REQUIRE (L1.coeff(f,f) == -3.0);
  50. REQUIRE (L1.row(f).sum() == 0.0);
  51. REQUIRE (L1.col(f).sum() == 0.0);
  52. #else
  53. //Soft assert if we have not edge_lenght_squared
  54. REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  55. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  56. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  57. #endif
  58. }
  59. //Same for huge cube.
  60. igl::cotmatrix(V_huge,F,L1);
  61. REQUIRE (L1.rows() == V.rows());
  62. REQUIRE (L1.cols() == V.rows());
  63. for(int f = 0;f<L1.rows();f++)
  64. {
  65. REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  66. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  67. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  68. }
  69. //Same for tiny cube. we need to use a tolerance this time...
  70. igl::cotmatrix(V_tiny,F,L1);
  71. REQUIRE (L1.rows() == V.rows());
  72. REQUIRE (L1.cols() == V.rows());
  73. for(int f = 0;f<L1.rows();f++)
  74. {
  75. REQUIRE (L1.coeff(f,f) == Approx (-3.0).margin( epsilon));
  76. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  77. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  78. }
  79. }
  80. TEST_CASE("cotmatrix: tetrahedron", "[igl]")
  81. {
  82. //The allowed error for this test
  83. const double epsilon = 1e-15;
  84. Eigen::MatrixXd V;
  85. Eigen::MatrixXi F;
  86. //This is a cube of dimensions 1.0x1.0x1.0
  87. test_common::load_mesh("cube.obj", V, F);
  88. //Prepare another mesh with triangles along side diagonals of the cube
  89. //These triangles are form a regular tetrahedron of side sqrt(2)
  90. Eigen::MatrixXi F_equi(4,3);
  91. F_equi << 4,6,1,
  92. 6,4,3,
  93. 4,1,3,
  94. 1,6,3;
  95. //Scale the cube to have huge sides
  96. Eigen::MatrixXd V_huge = V * 1.0e8;
  97. //Scale the cube to have tiny sides
  98. Eigen::MatrixXd V_tiny = V * 1.0e-8;
  99. //Check cotmatrix (Laplacian)
  100. //The laplacian for the cube is quite singular.
  101. //Each edge in a diagonal has two opposite angles of 90, with cotangent 0.0 each
  102. //Each edge in a side has two opposite angle of 45, with (half)cotangen 0.5 each
  103. //So the cotangent matrix always are (0+0) or (0.5+0.5)
  104. Eigen::SparseMatrix<double> L1;
  105. //Check the regular tetrahedron of side sqrt(2)
  106. igl::cotmatrix(V,F_equi,L1);
  107. REQUIRE (L1.rows() == V.rows());
  108. REQUIRE (L1.cols() == V.rows());
  109. for(int f = 0;f<L1.rows();f++)
  110. {
  111. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  112. if (L1.coeff(f,f) < -0.1)
  113. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(M_PI / 3.0)).margin( epsilon));
  114. else
  115. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  116. #ifdef IGL_EDGE_LENGTHS_SQUARED_H
  117. //Hard assert if we have edge_lenght_squared
  118. REQUIRE (L1.row(f).sum() == 0.0);
  119. REQUIRE (L1.col(f).sum() == 0.0);
  120. #else
  121. //Soft assert if we have not edge_lenght_squared
  122. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  123. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  124. #endif
  125. }
  126. //Check the huge regular tetrahedron
  127. igl::cotmatrix(V_huge,F_equi,L1);
  128. REQUIRE (L1.rows() == V.rows());
  129. REQUIRE (L1.cols() == V.rows());
  130. for(int f = 0;f<L1.rows();f++)
  131. {
  132. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  133. if (L1.coeff(f,f) < -0.1)
  134. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(M_PI / 3.0)).margin( epsilon));
  135. else
  136. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  137. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  138. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  139. }
  140. //Check the tiny regular tetrahedron
  141. igl::cotmatrix(V_tiny,F_equi,L1);
  142. REQUIRE (L1.rows() == V.rows());
  143. REQUIRE (L1.cols() == V.rows());
  144. for(int f = 0;f<L1.rows();f++)
  145. {
  146. //Check the diagonal. Only can value 0.0 for unused vertex or -3 / tan(60)
  147. if (L1.coeff(f,f) < -0.1)
  148. REQUIRE (L1.coeff(f,f) == Approx (-3 / tan(M_PI / 3.0)).margin( epsilon));
  149. else
  150. REQUIRE (L1.coeff(f,f) == Approx (0.0).margin( epsilon));
  151. REQUIRE (L1.row(f).sum() == Approx (0.0).margin( epsilon));
  152. REQUIRE (L1.col(f).sum() == Approx (0.0).margin( epsilon));
  153. }
  154. }