// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2018 Alec Jacobson // // This Source Code Form is subject to the terms of the Mozilla Public License // v. 2.0. If a copy of the MPL was not distributed with this file, You can // obtain one at http://mozilla.org/MPL/2.0/. #include "is_delaunay.h" #include "unique_edge_map.h" #include template < typename DerivedV, typename DerivedF, typename DerivedD> IGL_INLINE void igl::is_delaunay( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, Eigen::PlainObjectBase & D) { typedef typename DerivedV::Scalar Scalar; // Should use Shewchuk's predicates instead. const auto float_incircle = []( const Scalar pa[2], const Scalar pb[2], const Scalar pc[2], const Scalar pd[2])->short { // I acknowledge that I am cating to double const Eigen::Matrix3d A = (Eigen::Matrix3d(3,3)<< pa[0]-pd[0], pa[1]-pd[1],(pa[0]-pd[0])*(pa[0]-pd[0])+(pa[1]-pd[1])*(pa[1]-pd[1]), pb[0]-pd[0], pb[1]-pd[1],(pb[0]-pd[0])*(pb[0]-pd[0])+(pb[1]-pd[1])*(pb[1]-pd[1]), pc[0]-pd[0], pc[1]-pd[1],(pc[0]-pd[0])*(pc[0]-pd[0])+(pc[1]-pd[1])*(pc[1]-pd[1]) ).finished(); const Scalar detA = A.determinant(); return (Scalar(0) < detA) - (detA < Scalar(0)); }; typedef Eigen::Matrix MatrixX2I; typedef Eigen::Matrix VectorXI; MatrixX2I E,uE; VectorXI EMAP; std::vector > uE2E; igl::unique_edge_map(F, E, uE, EMAP, uE2E); const int num_faces = F.rows(); D.setConstant(F.rows(),F.cols(),false); // loop over all unique edges for(int ue = 0;ue < uE2E.size(); ue++) { const bool ue_is_d = is_delaunay(V,F,uE2E,float_incircle,ue); // Set for all instances for(int e = 0;e IGL_INLINE bool igl::is_delaunay( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, const std::vector > & uE2E, const InCircle incircle, const ueiType uei) { if(uE2E[uei].size() == 1) return true; if(uE2E[uei].size() > 2) return false; const int num_faces = F.rows(); typedef typename DerivedV::Scalar Scalar; const auto& half_edges = uE2E[uei]; assert((half_edges.size() == 2) && "uE2E[uei].size() should be 2"); const size_t f1 = half_edges[0] % num_faces; const size_t f2 = half_edges[1] % num_faces; const size_t c1 = half_edges[0] / num_faces; const size_t c2 = half_edges[1] / num_faces; assert(c1 < 3); assert(c2 < 3); assert(f1 != f2); const size_t v1 = F(f1, (c1+1)%3); const size_t v2 = F(f1, (c1+2)%3); const size_t v4 = F(f1, c1); const size_t v3 = F(f2, c2); const Scalar p1[] = {V(v1, 0), V(v1, 1)}; const Scalar p2[] = {V(v2, 0), V(v2, 1)}; const Scalar p3[] = {V(v3, 0), V(v3, 1)}; const Scalar p4[] = {V(v4, 0), V(v4, 1)}; auto orientation = incircle(p1, p2, p4, p3); return orientation <= 0; } #ifdef IGL_STATIC_LIBRARY // Explicit template instantiation // generated by autoexplicit.sh template void igl::is_delaunay, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template void igl::is_delaunay, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template bool igl::is_delaunay, Eigen::Matrix, int, short (*)(double const*, double const*, double const*, double const*), unsigned long>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, std::vector >, std::allocator > > > const&, short (*)(double const*, double const*, double const*, double const*), unsigned long); #ifdef WIN32 template bool igl::is_delaunay, class Eigen::Matrix, int, short(*)(double const *, double const *, double const *, double const *), unsigned __int64>(class Eigen::MatrixBase> const &, class Eigen::MatrixBase> const &, class std::vector>, class std::allocator>>> const &, short(*const)(double const *, double const *, double const *, double const *), unsigned __int64); #endif #endif