lexicographic_triangulation.cpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Alec Jacobson <[email protected]>
  4. // Qingan Zhou <[email protected]>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla Public License
  7. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  8. // obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "lexicographic_triangulation.h"
  10. #include "sortrows.h"
  11. #include "PlainMatrix.h"
  12. #include <vector>
  13. #include <list>
  14. template<
  15. typename DerivedP,
  16. typename Orient2D,
  17. typename DerivedF
  18. >
  19. IGL_INLINE void igl::lexicographic_triangulation(
  20. const Eigen::MatrixBase<DerivedP>& P,
  21. Orient2D orient2D,
  22. Eigen::PlainObjectBase<DerivedF>& F)
  23. {
  24. typedef typename DerivedP::Scalar Scalar;
  25. const size_t num_pts = P.rows();
  26. if (num_pts < 3) {
  27. throw "At least 3 points are required for triangulation!";
  28. }
  29. // Sort points in lexicographic order.
  30. PlainMatrix<DerivedP> ordered_P;
  31. Eigen::VectorXi order;
  32. igl::sortrows(P, true, ordered_P, order);
  33. std::vector<Eigen::Vector3i> faces;
  34. std::list<int> boundary;
  35. const Scalar p0[] = {ordered_P(0, 0), ordered_P(0, 1)};
  36. const Scalar p1[] = {ordered_P(1, 0), ordered_P(1, 1)};
  37. for (size_t i=2; i<num_pts; i++) {
  38. const Scalar curr_p[] = {ordered_P(i, 0), ordered_P(i, 1)};
  39. if (faces.size() == 0) {
  40. // All points processed so far are collinear.
  41. // Check if the current point is collinear with every points before it.
  42. auto orientation = orient2D(p0, p1, curr_p);
  43. if (orientation != 0) {
  44. // Add a fan of triangles eminating from curr_p.
  45. if (orientation > 0) {
  46. for (size_t j=0; j<=i-2; j++) {
  47. faces.push_back({order[j], order[j+1], order[i]});
  48. }
  49. } else if (orientation < 0) {
  50. for (size_t j=0; j<=i-2; j++) {
  51. faces.push_back({order[j+1], order[j], order[i]});
  52. }
  53. }
  54. // Initialize current boundary.
  55. boundary.insert(boundary.end(), order.data(), order.data()+i+1);
  56. if (orientation < 0) {
  57. boundary.reverse();
  58. }
  59. }
  60. } else {
  61. const size_t bd_size = boundary.size();
  62. assert(bd_size >= 3);
  63. std::vector<short> orientations;
  64. for (auto itr=boundary.begin(); itr!=boundary.end(); itr++) {
  65. auto next_itr = std::next(itr, 1);
  66. if (next_itr == boundary.end()) {
  67. next_itr = boundary.begin();
  68. }
  69. const Scalar bd_p0[] = {P(*itr, 0), P(*itr, 1)};
  70. const Scalar bd_p1[] = {P(*next_itr, 0), P(*next_itr, 1)};
  71. auto orientation = orient2D(bd_p0, bd_p1, curr_p);
  72. if (orientation < 0) {
  73. faces.push_back({*next_itr, *itr, order[i]});
  74. }
  75. orientations.push_back(orientation);
  76. }
  77. auto left_itr = boundary.begin();
  78. auto right_itr = boundary.begin();
  79. auto curr_itr = boundary.begin();
  80. for (size_t j=0; j<bd_size; j++, curr_itr++) {
  81. size_t prev = (j+bd_size-1) % bd_size;
  82. if (orientations[j] >= 0 && orientations[prev] < 0) {
  83. right_itr = curr_itr;
  84. } else if (orientations[j] < 0 && orientations[prev] >= 0) {
  85. left_itr = curr_itr;
  86. }
  87. }
  88. assert(left_itr != right_itr);
  89. for (auto itr=left_itr; itr!=right_itr; itr++) {
  90. if (itr == boundary.end()) itr = boundary.begin();
  91. if (itr == right_itr) break;
  92. if (itr == left_itr) continue;
  93. itr = boundary.erase(itr);
  94. if (itr == boundary.begin()) {
  95. itr = boundary.end();
  96. }
  97. itr--;
  98. }
  99. if (right_itr == boundary.begin()) {
  100. assert(std::next(left_itr, 1) == boundary.end());
  101. boundary.insert(boundary.end(), order[i]);
  102. } else {
  103. assert(std::next(left_itr, 1) == right_itr);
  104. boundary.insert(right_itr, order[i]);
  105. }
  106. }
  107. }
  108. const size_t num_faces = faces.size();
  109. if (num_faces == 0) {
  110. // All input points are collinear.
  111. // Do nothing here.
  112. } else {
  113. F.resize(num_faces, 3);
  114. for (size_t i=0; i<num_faces; i++) {
  115. F.row(i) = faces[i];
  116. }
  117. }
  118. }
  119. #ifdef IGL_STATIC_LIBRARY
  120. template void igl::lexicographic_triangulation<Eigen::Matrix<double, -1, -1, 0, -1, -1>, short (*)(double const*, double const*, double const*), Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, short (*)(double const*, double const*, double const*), Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  121. #endif