#include "is_inside.h" #include #include #include #include #include #include #include #include #include "order_facets_around_edge.h" #include "assign_scalar.h" #include "intersect_other.h" #include "RemeshSelfIntersectionsParam.h" namespace igl { namespace cgal { namespace is_inside_helper { typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel; typedef Kernel::Ray_3 Ray_3; typedef Kernel::Point_3 Point_3; typedef Kernel::Vector_3 Vector_3; typedef Kernel::Triangle_3 Triangle; typedef Kernel::Plane_3 Plane_3; typedef std::vector::iterator Iterator; typedef CGAL::AABB_triangle_primitive Primitive; typedef CGAL::AABB_traits AABB_triangle_traits; typedef CGAL::AABB_tree Tree; template bool intersect_each_other( const Eigen::PlainObjectBase& V1, const Eigen::PlainObjectBase& F1, const Eigen::PlainObjectBase& I1, const Eigen::PlainObjectBase& V2, const Eigen::PlainObjectBase& F2, const Eigen::PlainObjectBase& I2) { const size_t num_faces_1 = I1.rows(); DerivedF F1_selected(num_faces_1, F1.cols()); for (size_t i=0; i ElementType determine_element_type( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const size_t fid, const Point_3& p, size_t& element_index) { const Eigen::Vector3i f = F.row(I(fid, 0)); const Point_3 p0(V(f[0], 0), V(f[0], 1), V(f[0], 2)); const Point_3 p1(V(f[1], 0), V(f[1], 1), V(f[1], 2)); const Point_3 p2(V(f[2], 0), V(f[2], 1), V(f[2], 2)); if (p == p0) { element_index = 0; return VERTEX; } if (p == p1) { element_index = 1; return VERTEX; } if (p == p2) { element_index = 2; return VERTEX; } if (CGAL::collinear(p0, p1, p)) { element_index = 2; return EDGE; } if (CGAL::collinear(p1, p2, p)) { element_index = 0; return EDGE; } if (CGAL::collinear(p2, p0, p)) { element_index = 1; return EDGE; } element_index = 0; return FACE; } template void extract_adj_faces( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const size_t s, const size_t d, std::vector& adj_faces) { const size_t num_faces = I.rows(); for (size_t i=0; i void extract_adj_vertices( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const size_t v, std::vector& adj_vertices) { std::set unique_adj_vertices; const size_t num_faces = I.rows(); for (size_t i=0; i bool determine_point_edge_orientation( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const Point_3& query, size_t s, size_t d) { // Algorithm: // // If the query point is projected onto an edge, all adjacent // faces of that edge must be on or belong to a single half // space (i.e. there exists a plane passing through the edge and // all adjacent faces are either on the plane or on the same // side of that plane). // // If these adjacent faces are not coplanar, query is inside iff // the edge is concave. // // If two or more faces are coplanar, the query point is // definitely outside of the std::vector adj_faces; extract_adj_faces(V, F, I, s, d, adj_faces); const size_t num_adj_faces = adj_faces.size(); assert(num_adj_faces > 0); //std::cout << "adj faces: "; //for (size_t i=0; i size_t{ // Eigen::Vector3i f = F.row(abs(fid)-1); // if (f[0] != s && f[0] != d) return f[0]; // if (f[1] != s && f[1] != d) return f[1]; // if (f[2] != s && f[2] != d) return f[2]; // return -1; // }; // Point_3 p_s(V(s,0), V(s,1), V(s,2)); // Point_3 p_d(V(d,0), V(d,1), V(d,2)); // //std::cout << "s: " // // << CGAL::to_double(V(s,0)) << " " // // << CGAL::to_double(V(s,1)) << " " // // << CGAL::to_double(V(s,2)) << std::endl; // //std::cout << "d: " // // << CGAL::to_double(V(d,0)) << " " // // << CGAL::to_double(V(d,1)) << " " // // << CGAL::to_double(V(d,2)) << std::endl; // for (size_t i=0; i 0 && adj_faces[order[num_adj_faces-1] < 0]) { return true; } else if (adj_faces[order[0]] < 0 && adj_faces[order[num_adj_faces-1] > 0]) { return false; } else { assert(false); } assert(false); return false; } template bool determine_point_vertex_orientation( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const Point_3& query, size_t s) { std::vector adj_vertices; extract_adj_vertices(V, F, I, s, adj_vertices); const size_t num_adj_vertices = adj_vertices.size(); //std::cout << "Q: " // << CGAL::to_double(query.x()) << " " // << CGAL::to_double(query.y()) << " " // << CGAL::to_double(query.z()) << " " // << std::endl; std::vector adj_points; for (size_t i=0; i::max(); //std::cout << "P: " // << CGAL::to_double(V(s,0)) << " " // << CGAL::to_double(V(s,1)) << " " // << CGAL::to_double(V(s,2)) << " " // << std::endl; Point_3 p(V(s,0), V(s,1), V(s,2)); for (size_t i=0; i V.rows()) { // All adj faces are coplanar, use the first edge. d = adj_vertices[0]; //std::cout << "all adj faces are coplanar" << std::endl; //return false; } //std::cout << "s: " << s << " d: " << d << std::endl; return determine_point_edge_orientation(V, F, I, query, s, d); } template bool determine_point_face_orientation( const Eigen::PlainObjectBase& V, const Eigen::PlainObjectBase& F, const Eigen::PlainObjectBase& I, const Point_3& query, size_t fid) { // Algorithm: A point is on the inside of a face if the // tetrahedron formed by them is negatively oriented. Eigen::Vector3i f = F.row(I(fid, 0)); const Point_3 v0(V(f[0], 0), V(f[0], 1), V(f[0], 2)); const Point_3 v1(V(f[1], 0), V(f[1], 1), V(f[1], 2)); const Point_3 v2(V(f[2], 0), V(f[2], 1), V(f[2], 2)); auto result = CGAL::orientation(v0, v1, v2, query); assert(result != CGAL::COPLANAR); return result == CGAL::NEGATIVE; } } } } template IGL_INLINE bool igl::cgal::is_inside( const Eigen::PlainObjectBase& V1, const Eigen::PlainObjectBase& F1, const Eigen::PlainObjectBase& I1, const Eigen::PlainObjectBase& V2, const Eigen::PlainObjectBase& F2, const Eigen::PlainObjectBase& I2) { using namespace igl::cgal::is_inside_helper; assert(F1.rows() > 0); assert(I1.rows() > 0); assert(F2.rows() > 0); assert(I2.rows() > 0); //assert(!intersect_each_other(V1, F1, I1, V2, F2, I2)); const size_t num_faces = I2.rows(); std::vector triangles; for (size_t i=0; i IGL_INLINE bool igl::cgal::is_inside( const Eigen::PlainObjectBase& V1, const Eigen::PlainObjectBase& F1, const Eigen::PlainObjectBase& V2, const Eigen::PlainObjectBase& F2) { Eigen::VectorXi I1(F1.rows()), I2(F2.rows()); I1.setLinSpaced(F1.rows(), 0, F1.rows()-1); I2.setLinSpaced(F2.rows(), 0, F2.rows()-1); return igl::cgal::is_inside(V1, F1, I1, V2, F2, I2); }