intrinsic_delaunay_triangulation.cpp 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139
  1. #include <test_common.h>
  2. #include <igl/intrinsic_delaunay_triangulation.h>
  3. #include <igl/edge_lengths.h>
  4. #include <igl/triangulated_grid.h>
  5. #include <igl/is_delaunay.h>
  6. #include <igl/is_intrinsic_delaunay.h>
  7. #include <igl/is_edge_manifold.h>
  8. #include <igl/unique_simplices.h>
  9. #include <igl/get_seconds.h>
  10. #include <igl/matlab_format.h>
  11. TEST(intrinsic_delaunay_triangulation, two_triangles)
  12. {
  13. const Eigen::MatrixXd V =
  14. (Eigen::MatrixXd(4,2)<<
  15. 0,12,
  16. 1,0,
  17. 1,20,
  18. 2,9).finished();
  19. const Eigen::MatrixXi FN =
  20. (Eigen::MatrixXi(2,3)<<
  21. 0,1,2,
  22. 2,1,3).finished();
  23. Eigen::MatrixXd lN;
  24. igl::edge_lengths(V,FN,lN);
  25. Eigen::MatrixXd l;
  26. Eigen::MatrixXi F;
  27. igl::intrinsic_delaunay_triangulation( lN, FN, l, F);
  28. Eigen::MatrixXd lext;
  29. igl::edge_lengths(V,F,lext);
  30. test_common::assert_near(l,lext,1e-10);
  31. }
  32. TEST(intrinsic_delaunay_triangulation, skewed_grid)
  33. {
  34. Eigen::MatrixXd V;
  35. Eigen::MatrixXi F_in;
  36. igl::triangulated_grid(4,4,V,F_in);
  37. // Skew against diagonal direction
  38. V.col(0) -= 1.5*V.col(1);
  39. Eigen::MatrixXd l_in;
  40. igl::edge_lengths(V,F_in,l_in);
  41. Eigen::MatrixXd l;
  42. Eigen::MatrixXi F;
  43. igl::intrinsic_delaunay_triangulation( l_in, F_in, l, F);
  44. Eigen::MatrixXd lext;
  45. igl::edge_lengths(V,F,lext);
  46. test_common::assert_near(l,lext,1e-10);
  47. Eigen::Matrix<bool,Eigen::Dynamic,3> D;
  48. igl::is_delaunay(V,F,D);
  49. const Eigen::Matrix<bool,Eigen::Dynamic,3> Dtrue =
  50. Eigen::Matrix<bool,Eigen::Dynamic,3>::Constant(F.rows(),F.cols(),true);
  51. test_common::assert_eq(D,Dtrue);
  52. }
  53. TEST(intrinsic_delaunay_triangulation, peaks)
  54. {
  55. Eigen::MatrixXd V2;
  56. Eigen::MatrixXi F_in;
  57. igl::triangulated_grid(6,6,V2,F_in);
  58. Eigen::MatrixXd V(V2.rows(),3);
  59. for(int v=0;v<V.rows();v++)
  60. {
  61. const auto x = (V2(v,0)-0.5)*6.0;
  62. const auto y = (V2(v,1)-0.5)*6.0;
  63. // peaks.m
  64. const auto z = 3.*(1.-x)*(1.-x)*std::exp(-(x*x) - (y+1.)*(y+1.)) +
  65. - 10.*(x/5. - x*x*x - y*y*y*y*y)*std::exp(-x*x-y*y) +
  66. - 1./3.*std::exp(-(x+1.)*(x+1.) - y*y);
  67. V(v,0) = x;
  68. V(v,1) = y;
  69. V(v,2) = z;
  70. }
  71. Eigen::MatrixXd l_in;
  72. igl::edge_lengths(V,F_in,l_in);
  73. Eigen::MatrixXd l;
  74. Eigen::MatrixXi F;
  75. Eigen::MatrixXi E,uE;
  76. Eigen::VectorXi EMAP;
  77. std::vector<std::vector<int> > uE2E;
  78. igl::intrinsic_delaunay_triangulation(
  79. l_in, F_in, l, F, E, uE, EMAP, uE2E);
  80. //Eigen::MatrixXd lext;
  81. //igl::edge_lengths(V,F,lext);
  82. //std::cout<<igl::matlab_format(V,"V")<<std::endl;
  83. //std::cout<<igl::matlab_format(F_in.array()+1,"F_in")<<std::endl;
  84. //std::cout<<igl::matlab_format(F.array()+1,"F")<<std::endl;
  85. //std::cout<<igl::matlab_format(l,"l")<<std::endl;
  86. Eigen::Matrix<bool,Eigen::Dynamic,3> D;
  87. const Eigen::Matrix<bool,Eigen::Dynamic,3> D_gt =
  88. Eigen::Matrix<bool,Eigen::Dynamic,3>::Constant(F.rows(),F.cols(),true);
  89. igl::is_intrinsic_delaunay(l,F,uE2E,D);
  90. test_common::assert_eq(D,D_gt);
  91. //{
  92. // Eigen::MatrixXi Fu;
  93. // Eigen::VectorXi IA,IC;
  94. // igl::unique_simplices(F,Fu,IA,IC);
  95. // ASSERT_EQ(F.rows(),Fu.rows());
  96. //}
  97. // Input better have started manifold, otherwise this test doesn't make sense
  98. //assert(igl::is_edge_manifold(F_in));
  99. //ASSERT_TRUE(igl::is_edge_manifold(F));
  100. }
  101. //// Not sure if this is a good test... Even though the edge will exist twice
  102. //the intrinsic triangles on either edge are different...
  103. //TEST(intrinsic_delaunay_triangulation,unflippable_tet)
  104. //{
  105. // const Eigen::MatrixXd V = (Eigen::MatrixXd(4,3)<<
  106. // 10, 4,7,
  107. // 5, 9,0,
  108. // 8, 8,8,
  109. // 1,10,9).finished();
  110. // const Eigen::MatrixXi F_in = (Eigen::MatrixXi(4,3)<<
  111. // 0,1,2,
  112. // 0,2,3,
  113. // 0,3,1,
  114. // 1,3,2).finished();
  115. // const Eigen::Matrix<bool,Eigen::Dynamic,3> Dgt =
  116. // (Eigen::Matrix<bool,Eigen::Dynamic,3>(4,3)<<
  117. // 1,1,1,
  118. // 1,0,1,
  119. // 1,1,0,
  120. // 1,1,1).finished();
  121. // Eigen::Matrix<bool,Eigen::Dynamic,3> D;
  122. // Eigen::MatrixXd l_in;
  123. // igl::edge_lengths(V,F_in,l_in);
  124. // igl::is_intrinsic_delaunay(l_in,F_in,D);
  125. // test_common::assert_eq(D,Dgt);
  126. // Eigen::MatrixXd l;
  127. // Eigen::MatrixXi F;
  128. // igl::intrinsic_delaunay_triangulation( l_in, F_in, l, F);
  129. // // Nothing should have changed: no edges are flippable.
  130. // test_common::assert_eq(F,F_in);
  131. // test_common::assert_eq(l,l_in);
  132. //}