intrinsic_delaunay_triangulation.cpp 3.5 KB

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