cotmatrix.cpp 5.4 KB

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