Преглед изворни кода

mesh to polyhedron missing

Former-commit-id: 2fe2399930be185441999014dbb9160c90cc6399
Alec Jacobson пре 10 година
родитељ
комит
144127cde1

+ 3 - 4
include/igl/cgal/CGAL_includes.hpp

@@ -31,10 +31,9 @@
 #include <CGAL/AABB_traits.h>
 #include <CGAL/AABB_triangle_primitive.h>
 
-// ALEC: IS ANY IGL FUNCTION USING THESE?
-//// Boolean operations
-//#include <CGAL/Polyhedron_3.h>
-//#include <CGAL/Nef_polyhedron_3.h>
+// Boolean operations
+#include <CGAL/Polyhedron_3.h>
+#include <CGAL/Nef_polyhedron_3.h>
 
 // Delaunay Triangulation in 3D
 #include <CGAL/Delaunay_triangulation_3.h>

+ 1 - 1
include/igl/cgal/SelfIntersectMesh.h

@@ -911,7 +911,7 @@ inline bool igl::SelfIntersectMesh<
     try
     {
       CGAL::Object result = CGAL::intersection(A,B);
-      if(result)
+      if(!result.empty())
       {
         if(CGAL::object_cast<Segment_3 >(&result))
         {

+ 1 - 0
include/igl/cgal/complex_to_mesh.cpp

@@ -11,6 +11,7 @@
 #include <igl/remove_unreferenced.h>
 
 #include <CGAL/Surface_mesh_default_triangulation_3.h>
+#include <CGAL/Triangulation_cell_base_with_circumcenter_3.h>
 #include <set>
 #include <stack>
 

+ 45 - 0
include/igl/cgal/mesh_to_polyhedron.cpp

@@ -0,0 +1,45 @@
+// This file is part of libigl, a simple c++ geometry processing library.
+// 
+// Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
+// 
+// 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 "mesh_to_polyhedron.h"
+#include <CGAL/Polyhedron_3.h>
+#include <CGAL/Polyhedron_incremental_builder_3.h>
+template <typename Polyhedron>
+IGL_INLINE bool igl::mesh_to_polyhedron(
+  const Eigen::MatrixXd & V,
+  const Eigen::MatrixXi & F,
+  Polyhedron & poly)
+{
+  typedef typename Polyhedron::HalfedgeDS HalfedgeDS;
+  // Postcondition: hds is a valid polyhedral surface.
+  CGAL::Polyhedron_incremental_builder_3<HalfedgeDS> B(poly.hds());
+  B.begin_surface(V.rows(),F.rows());
+  typedef typename HalfedgeDS::Vertex   Vertex;
+  typedef typename Vertex::Point Point;
+  assert(V.cols() == 3 && "V must be #V by 3");
+  for(int v = 0;v<V.rows();v++)
+  {
+    B.add_vertex(Point(V(v,0),V(v,1),V(v,2)));
+  }
+  assert(F.cols() == 3 && "F must be #F by 3");
+  for(int f=0;f<F.rows();f++)
+  {
+    B.begin_facet();
+    for(int c = 0;c<3;c++)
+    {
+      B.add_vertex_to_facet(F(f,c));
+    }
+    B.end_facet();
+  }
+  if(B.error())
+  {
+    B.rollback();
+    return false;
+  }
+  B.end_surface();
+  return poly.is_valid();
+}

+ 36 - 0
include/igl/cgal/mesh_to_polyhedron.h

@@ -0,0 +1,36 @@
+// This file is part of libigl, a simple c++ geometry processing library.
+// 
+// Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
+// 
+// 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/.
+#ifndef IGL_MESH_TO_POLYHEDRON_H
+#define IGL_MESH_TO_POLYHEDRON_H
+#include <igl/igl_inline.h>
+#include <Eigen/Core>
+
+namespace igl
+{
+  // Convert a mesh (V,F) to a CGAL Polyhedron
+  //
+  // Templates:
+  //   Polyhedron  CGAL Polyhedron type (e.g. Polyhedron_3)
+  // Inputs:
+  //   V  #V by 3 list of vertex positions
+  //   F  #F by 3 list of triangle indices
+  // Outputs:
+  //   poly  cgal polyhedron
+  // Returns true only if (V,F) can be converted to a valid polyhedron (i.e. if
+  // (V,F) is vertex and edge manifold).
+  template <typename Polyhedron>
+  IGL_INLINE bool mesh_to_polyhedron(
+    const Eigen::MatrixXd & V,
+    const Eigen::MatrixXi & F,
+    Polyhedron & poly);
+}
+#ifndef IGL_STATIC_LIBRARY
+#  include "mesh_to_polyhedron.cpp"
+#endif
+
+#endif