refine.cpp 2.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118
  1. #include <test_common.h>
  2. #include <igl/triangle/refine.h>
  3. #include <igl/triangle/triangulate.h>
  4. #include <igl/matlab_format.h>
  5. #include <igl/find.h>
  6. #include <igl/barycenter.h>
  7. #include <igl/winding_number.h>
  8. #include <igl/doublearea.h>
  9. #include <iostream>
  10. TEST_CASE("refine", "[igl][triangle]") {
  11. Eigen::MatrixXd V(3,2);
  12. V <<
  13. 0,0,
  14. 1,0,
  15. 0,1;
  16. Eigen::MatrixXi F(1,3);
  17. F<<0,1,2;
  18. std::string flags = "a0.0625Q";
  19. Eigen::MatrixXd V2;
  20. Eigen::MatrixXi F2;
  21. Eigen::MatrixXi E;
  22. igl::triangle::refine(V,E,F,flags,V2,F2);
  23. REQUIRE(F2.rows()>1);
  24. Eigen::VectorXd A2;
  25. igl::doublearea(V2,F2,A2);
  26. REQUIRE((A2.array()<2*0.0625).all());
  27. }
  28. TEST_CASE("refine-inner-pinch", "[igl][triangle]") {
  29. // Mesh an outline that looks like:
  30. //
  31. // o-----------------o
  32. // | |
  33. // | o o |
  34. // | |\ /| |
  35. // | | \ / | |
  36. // | | \ / | |
  37. // | | \ / | |
  38. // | | o | |
  39. // | o---------o |
  40. // | |
  41. // o-----------------o
  42. //
  43. // If we mesh the whole convex hull and then remove the inner cavity we'll get
  44. // over-refinement near the sharp corner.
  45. //
  46. // We could tell triangle a point inside this cavity but that's not often easy
  47. // to do robustly/for broken input. Instead we can use a trivial mesh of the
  48. // whole domain, remove the cavity with the winding number and then refine
  49. // what's left.
  50. //
  51. Eigen::MatrixXd V(9,2);
  52. V <<
  53. 0,0,
  54. 3,0,
  55. 3,3,
  56. 0,3,
  57. 1,1,
  58. 2,1,
  59. 2,2,
  60. 1.5,1.00000001,
  61. 1,2;
  62. Eigen::MatrixXi E(9,2);
  63. E <<
  64. 0,1,
  65. 1,2,
  66. 2,3,
  67. 3,0,
  68. 8,7,
  69. 7,6,
  70. 6,5,
  71. 5,4,
  72. 4,8;
  73. Eigen::MatrixXd H(0,2);
  74. const auto filter = [&V,&E](
  75. const Eigen::MatrixXd & Vc,
  76. Eigen::MatrixXi & Fc)
  77. {
  78. Eigen::MatrixXd BC;
  79. igl::barycenter(Vc,Fc,BC);
  80. Eigen::VectorXd W;
  81. // This is sadly O(|E| ⋅ |BC|) for 2D inputs.
  82. igl::winding_number(V,E,BC,W);
  83. auto keep = igl::find((W.array().abs() > 0.5).eval());
  84. Fc = Fc(keep,Eigen::all).eval();
  85. };
  86. Eigen::MatrixXd Vq;
  87. Eigen::MatrixXi Fq;
  88. igl::triangle::triangulate(V,E,H,"cq33a0.0625Q",Vq,Fq);
  89. filter(Vq,Fq);
  90. //std::cout<<igl::matlab_format(Vq,"Vq")<<std::endl;
  91. //std::cout<<igl::matlab_format_index(Fq,"Fq")<<std::endl;
  92. Eigen::MatrixXd Vc;
  93. Eigen::MatrixXi Fc;
  94. igl::triangle::triangulate(V,E,H,"cQ",Vc,Fc);
  95. REQUIRE(V.rows() == Vc.rows());
  96. REQUIRE(Fc.rows() == 12);
  97. filter(Vc,Fc);
  98. REQUIRE(Fc.rows() == 9);
  99. //std::cout<<igl::matlab_format(Vc,"Vc")<<std::endl;
  100. //std::cout<<igl::matlab_format_index(Fc,"Fc")<<std::endl;
  101. Eigen::MatrixXd Vr;
  102. Eigen::MatrixXi Fr;
  103. Eigen::MatrixXi _;
  104. igl::triangle::refine(Vc,_,Fc,"q33a0.0625Q",Vr,Fr);
  105. //std::cout<<igl::matlab_format(Vr,"Vr")<<std::endl;
  106. //std::cout<<igl::matlab_format_index(Fr,"Fr")<<std::endl;
  107. REQUIRE(Fr.rows() > Fc.rows());
  108. REQUIRE(Fr.rows() < 0.15 * Fq.rows());
  109. }