delaunay_triangulation.cpp 2.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Qingnan Zhou <qnzhou@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "delaunay_triangulation.h"
  9. #include "flip_edge.h"
  10. #include "lexicographic_triangulation.h"
  11. #include "unique_edge_map.h"
  12. #include "is_delaunay.h"
  13. #include <vector>
  14. #include <sstream>
  15. template<
  16. typename DerivedV,
  17. typename Orient2D,
  18. typename InCircle,
  19. typename DerivedF>
  20. IGL_INLINE void igl::delaunay_triangulation(
  21. const Eigen::MatrixBase<DerivedV>& V,
  22. Orient2D orient2D,
  23. InCircle incircle,
  24. Eigen::PlainObjectBase<DerivedF>& F)
  25. {
  26. assert(V.cols() == 2);
  27. typedef typename DerivedF::Scalar Index;
  28. typedef typename DerivedV::Scalar Scalar;
  29. igl::lexicographic_triangulation(V, orient2D, F);
  30. const size_t num_faces = F.rows();
  31. if (num_faces == 0) {
  32. // Input points are degenerate. No faces will be generated.
  33. return;
  34. }
  35. assert(F.cols() == 3);
  36. typedef Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,2> MatrixX2I;
  37. MatrixX2I E,uE;
  38. Eigen::VectorXi EMAP;
  39. std::vector<std::vector<Index> > uE2E;
  40. igl::unique_edge_map(F, E, uE, EMAP, uE2E);
  41. bool all_delaunay = false;
  42. while(!all_delaunay) {
  43. all_delaunay = true;
  44. for (size_t i=0; i<uE2E.size(); i++) {
  45. if (uE2E[i].size() == 2) {
  46. if (!is_delaunay(V,F,uE2E,incircle,i)) {
  47. all_delaunay = false;
  48. flip_edge(F, E, uE, EMAP, uE2E, i);
  49. }
  50. }
  51. }
  52. }
  53. }
  54. #ifdef IGL_STATIC_LIBRARY
  55. // Explicit template instantiation
  56. template void igl::delaunay_triangulation<Eigen::Matrix<double, -1, -1, 0, -1, -1>, short (*)(double const*, double const*, double const*), short (*)(double const*, 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*), short (*)(double const*, double const*, double const*, double const*), Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  57. #ifdef WIN32
  58. template void igl::delaunay_triangulation<class Eigen::Matrix<double, -1, -1, 0, -1, -1>, short(*)(double const *, double const *, double const *), short(*)(double const *, double const *, double const *, double const *), class Eigen::Matrix<int, -1, -1, 0, -1, -1>>(class Eigen::MatrixBase<class Eigen::Matrix<double, -1, -1, 0, -1, -1>> const &, short(*const)(double const *, double const *, double const *), short(*const)(double const *, double const *, double const *, double const *), class Eigen::PlainObjectBase<class Eigen::Matrix<int, -1, -1, 0, -1, -1>> &);
  59. #endif
  60. #endif