فهرست منبع

Merge pull request #101 from qnzhou/master

Exterior vertex/edge/facet

Former-commit-id: efb79bf38088c6dcb7ab5f00b59dd5f76fe422fa
Alec Jacobson 9 سال پیش
والد
کامیت
742e6f25cb
6فایلهای تغییر یافته به همراه354 افزوده شده و 176 حذف شده
  1. 8 8
      include/igl/cgal/outer_hull.cpp
  2. 1 0
      include/igl/facet_components.h
  3. 235 0
      include/igl/outer_element.cpp
  4. 110 0
      include/igl/outer_element.h
  5. 0 122
      include/igl/outer_facet.cpp
  6. 0 46
      include/igl/outer_facet.h

+ 8 - 8
include/igl/cgal/outer_hull.cpp

@@ -7,7 +7,7 @@
 // obtain one at http://mozilla.org/MPL/2.0/.
 #include "outer_hull.h"
 #include "order_facets_around_edges.h"
-#include "../outer_facet.h"
+#include "../outer_element.h"
 #include "../sortrows.h"
 #include "../facet_components.h"
 #include "../winding_number.h"
@@ -163,9 +163,9 @@ IGL_INLINE void igl::cgal::outer_hull(
     outer_facet(V,F,N,IM,f,f_flip);
 #ifdef IGL_OUTER_HULL_DEBUG
   cout<<"outer facet: "<<f<<endl;
-  cout << V.row(F(f, 0)) << std::endl;
-  cout << V.row(F(f, 1)) << std::endl;
-  cout << V.row(F(f, 2)) << std::endl;
+  //cout << V.row(F(f, 0)) << std::endl;
+  //cout << V.row(F(f, 1)) << std::endl;
+  //cout << V.row(F(f, 2)) << std::endl;
 #endif
     int FHcount = 1;
     FH[f] = true;
@@ -209,9 +209,9 @@ IGL_INLINE void igl::cgal::outer_hull(
       // edge valence
       const size_t val = uE2E[EMAP(e)].size();
 #ifdef IGL_OUTER_HULL_DEBUG
-      std::cout << "vd: " << V.row(fd) << std::endl;
-      std::cout << "vs: " << V.row(fs) << std::endl;
-      std::cout << "edge: " << V.row(fd) - V.row(fs) << std::endl;
+      //std::cout << "vd: " << V.row(fd) << std::endl;
+      //std::cout << "vs: " << V.row(fs) << std::endl;
+      //std::cout << "edge: " << V.row(fd) - V.row(fs) << std::endl;
       for (size_t i=0; i<val; i++) {
           if (i == diIM(e)) {
               std::cout << "* ";
@@ -536,7 +536,7 @@ IGL_INLINE void igl::cgal::outer_hull(
 // Explicit template specialization
 #undef IGL_STATIC_LIBRARY
 #include <igl/barycenter.cpp>
-#include <igl/outer_facet.cpp>
+#include <igl/outer_element.cpp>
 #include <igl/cgal/order_facets_around_edges.cpp>
 #define IGL_STATIC_LIBRARY
 template void igl::cgal::outer_hull<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<bool, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&);

+ 1 - 0
include/igl/facet_components.h

@@ -27,6 +27,7 @@ namespace igl
   //   triangle_triangle_adjacency.h)
   // Ouputs:
   //   C  #F list of connected component ids
+  //   counts #C list of number of facets in each components
   template <
     typename TTIndex, 
     typename DerivedC,

+ 235 - 0
include/igl/outer_element.cpp

@@ -0,0 +1,235 @@
+// This file is part of libigl, a simple c++ geometry processing library.
+// 
+// Copyright (C) 2015 Qingnan Zhou <qnzhou@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 "outer_element.h"
+#include <iostream>
+
+template <
+     typename DerivedV,
+     typename DerivedF,
+     typename DerivedI,
+     typename IndexType,
+     typename DerivedA
+     >
+IGL_INLINE void igl::outer_vertex(
+        const Eigen::PlainObjectBase<DerivedV> & V,
+        const Eigen::PlainObjectBase<DerivedF> & F,
+        const Eigen::PlainObjectBase<DerivedI> & I,
+        IndexType & v_index,
+        Eigen::PlainObjectBase<DerivedA> & A) {
+    // Algorithm: 
+    //    Find an outer vertex (i.e. vertex reachable from infinity)
+    //    Return the vertex with the largest X value.
+    //    If there is a tie, pick the one with largest Y value.
+    //    If there is still a tie, pick the one with the largest Z value.
+    //    If there is still a tie, then there are duplicated vertices within the
+    //    mesh, which violates the precondition.
+    const size_t INVALID = std::numeric_limits<size_t>::max();
+    const size_t num_selected_faces = I.rows();
+    std::vector<size_t> candidate_faces;
+    size_t outer_vid = INVALID;
+    typename DerivedV::Scalar outer_val = 0;
+    for (size_t i=0; i<num_selected_faces; i++) {
+        size_t f = I[i];
+        for (size_t j=0; j<3; j++) {
+            auto v = F(f, j);
+            auto vx = V(v, 0);
+            if (outer_vid == INVALID || vx > outer_val) {
+                outer_val = vx;
+                outer_vid = v;
+                candidate_faces = {f};
+            } else if (v == outer_vid) {
+                candidate_faces.push_back(f);
+            } else if (vx == outer_val) {
+                // Break tie.
+                auto vy = V(v,1);
+                auto vz = V(v, 2);
+                auto outer_y = V(outer_vid, 1);
+                auto outer_z = V(outer_vid, 2);
+                assert(!(vy == outer_y && vz == outer_z));
+                bool replace = (vy > outer_y) ||
+                    ((vy == outer_y) && (vz > outer_z));
+                if (replace) {
+                    outer_val = vx;
+                    outer_vid = v;
+                    candidate_faces = {f};
+                }
+            }
+        }
+    }
+
+    assert(outer_vid != INVALID);
+    assert(candidate_faces.size() > 0);
+    v_index = outer_vid;
+    A.resize(candidate_faces.size());
+    std::copy(candidate_faces.begin(), candidate_faces.end(), A.data());
+}
+
+template<
+    typename DerivedV,
+    typename DerivedF,
+    typename DerivedI,
+    typename IndexType,
+    typename DerivedA
+    >
+IGL_INLINE void igl::outer_edge(
+        const Eigen::PlainObjectBase<DerivedV> & V,
+        const Eigen::PlainObjectBase<DerivedF> & F,
+        const Eigen::PlainObjectBase<DerivedI> & I,
+        IndexType & v1,
+        IndexType & v2,
+        Eigen::PlainObjectBase<DerivedA> & A) {
+    // Algorithm:
+    //    Find an outer vertex first.
+    //    Find the incident edge with largest slope when projected onto XY plane.
+    //    If there is still a tie, break it using the projected slope onto ZX plane.
+    //    If there is still a tie, then there are multiple overlapping edges,
+    //    which violates the precondition.
+    typedef typename DerivedV::Scalar Scalar;
+    typedef typename DerivedV::Index Index;
+    typedef typename Eigen::Matrix<Scalar, 3, 1> ScalarArray3;
+    typedef typename Eigen::Matrix<typename DerivedF::Scalar, 3, 1> IndexArray3;
+    const size_t INVALID = std::numeric_limits<size_t>::max();
+
+    Index outer_vid;
+    DerivedI candidate_faces;
+    outer_vertex(V, F, I, outer_vid, candidate_faces);
+    const ScalarArray3& outer_v = V.row(outer_vid);
+    assert(candidate_faces.size() > 0);
+
+    auto get_vertex_index = [&](const IndexArray3& f,
+            Index vid) {
+        if (f[0] == vid) return 0;
+        if (f[1] == vid) return 1;
+        if (f[2] == vid) return 2;
+        assert(false);
+    };
+
+    Scalar outer_slope_YX = 0;
+    Scalar outer_slope_ZX = 0;
+    size_t outer_opp_vid = INVALID;
+    bool infinite_slope_detected = false;
+    std::vector<Index> incident_faces;
+    auto check_and_update_outer_edge = [&](Index opp_vid, Index fid) {
+        if (opp_vid == outer_opp_vid) {
+            incident_faces.push_back(fid);
+            return;
+        }
+
+        const ScalarArray3 opp_v = V.row(opp_vid);
+        if (!infinite_slope_detected && outer_v[0] != opp_v[0]) {
+            // Finite slope
+            const ScalarArray3 diff = opp_v - outer_v;
+            const Scalar slope_YX = diff[1] / diff[0];
+            const Scalar slope_ZX = diff[2] / diff[0];
+            if (outer_opp_vid == INVALID ||
+                    slope_YX > outer_slope_YX ||
+                    (slope_YX == outer_slope_YX &&
+                     slope_ZX > outer_slope_ZX)) {
+                outer_opp_vid = opp_vid;
+                outer_slope_YX = slope_YX;
+                outer_slope_ZX = slope_ZX;
+                incident_faces = {fid};
+            }
+        } else if (!infinite_slope_detected) {
+            // Infinite slope
+            outer_opp_vid = opp_vid;
+            infinite_slope_detected = true;
+            incident_faces = {fid};
+        }
+    };
+
+    const auto num_candidate_faces = candidate_faces.size();
+    for (size_t i=0; i<num_candidate_faces; i++) {
+        const Index fid = candidate_faces[i];
+        const IndexArray3& f = F.row(fid);
+        size_t id = get_vertex_index(f, outer_vid);
+        Index next_vid = f[(id+1)%3];
+        Index prev_vid = f[(id+2)%3];
+        check_and_update_outer_edge(next_vid, fid);
+        check_and_update_outer_edge(prev_vid, fid);
+    }
+
+    v1 = outer_vid;
+    v2 = outer_opp_vid;
+    A.resize(incident_faces.size());
+    std::copy(incident_faces.begin(), incident_faces.end(), A.data());
+}
+
+template<
+    typename DerivedV,
+    typename DerivedF,
+    typename DerivedN,
+    typename DerivedI,
+    typename IndexType
+    >
+IGL_INLINE void igl::outer_facet(
+        const Eigen::PlainObjectBase<DerivedV> & V,
+        const Eigen::PlainObjectBase<DerivedF> & F,
+        const Eigen::PlainObjectBase<DerivedN> & N,
+        const Eigen::PlainObjectBase<DerivedI> & I,
+        IndexType & f,
+        bool & flipped) {
+    // Algorithm:
+    //    Find an outer edge.
+    //    Find the incident facet with the largest absolute X normal component.
+    //    If there is a tie, keep the one with positive X component.
+    //    If there is still a tie, pick the face with the larger signed index
+    //    (flipped face has negative index).
+    typedef typename DerivedV::Scalar Scalar;
+    typedef typename DerivedV::Index Index;
+    const size_t INVALID = std::numeric_limits<size_t>::max();
+
+    Index v1,v2;
+    DerivedI incident_faces;
+    outer_edge(V, F, I, v1, v2, incident_faces);
+    assert(incident_faces.size() > 0);
+
+    auto generic_fabs = [&](const Scalar& val) -> const Scalar {
+        if (val >= 0) return val;
+        else return -val;
+    };
+
+    Scalar max_nx = 0;
+    size_t outer_fid = INVALID;
+    const size_t num_incident_faces = incident_faces.size();
+    for (size_t i=0; i<num_incident_faces; i++) {
+        const auto& fid = incident_faces[i];
+        const Scalar nx = N(fid, 0);
+        if (outer_fid == INVALID) {
+            max_nx = nx;
+            outer_fid = fid;
+        } else {
+            if (generic_fabs(nx) > generic_fabs(max_nx)) {
+                max_nx = nx;
+                outer_fid = fid;
+            } else if (nx == -max_nx && nx > 0) {
+                max_nx = nx;
+                outer_fid = fid;
+            } else if (nx == max_nx) {
+                if ((max_nx >= 0 && outer_fid < fid) ||
+                    (max_nx <  0 && outer_fid > fid)) {
+                    max_nx = nx;
+                    outer_fid = fid;
+                }
+            }
+        }
+    }
+
+    assert(outer_fid != INVALID);
+    f = outer_fid;
+    flipped = max_nx < 0;
+}
+
+
+#ifdef IGL_STATIC_LIBRARY
+// Explicit template specialization
+template void igl::outer_facet<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> > const&, int&, bool&);
+template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 1, -1, -1>, unsigned long>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 1, -1, -1> > const&, unsigned long&, bool&);
+template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int&, bool&);
+template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int&, bool&);
+#endif

+ 110 - 0
include/igl/outer_element.h

@@ -0,0 +1,110 @@
+// This file is part of libigl, a simple c++ geometry processing library.
+// 
+// Copyright (C) 2015 Qingan Zhou <qnzhou@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_OUTER_ELEMENT_H
+#define IGL_OUTER_ELEMENT_H
+#include "igl_inline.h"
+#include <Eigen/Core>
+namespace igl
+{
+  // Find a vertex that is reachable from infinite without crossing any faces.
+  // Such vertex is called "outer vertex."
+  //
+  // Precondition: The input mesh must have all self-intersection resolved and
+  // no duplicated vertices.  See cgal::remesh_self_intersections.h for how to
+  // obtain such input.
+  //
+  // Inputs:
+  //   V  #V by 3 list of vertex positions
+  //   F  #F by 3 list of triangle indices into V
+  //   I  #I list of facets to consider
+  // Outputs:
+  //   v_index  index of outer vertex
+  //   A  #A list of facets incident to the outer vertex
+  template <
+      typename DerivedV,
+      typename DerivedF,
+      typename DerivedI,
+      typename IndexType,
+      typename DerivedA
+      >
+  IGL_INLINE void outer_vertex(
+          const Eigen::PlainObjectBase<DerivedV> & V,
+          const Eigen::PlainObjectBase<DerivedF> & F,
+          const Eigen::PlainObjectBase<DerivedI> & I,
+          IndexType & v_index,
+          Eigen::PlainObjectBase<DerivedA> & A);
+
+
+  // Find an edge that is reachable from infinity without crossing any faces.
+  // Such edge is called "outer edge."
+  //
+  // Precondition: The input mesh must have all self-intersection resolved and
+  // no duplicated vertices.  The correctness of the output depends on the fact
+  // that there is no edge overlap.  See cgal::remesh_self_intersections.h for
+  // how to obtain such input.
+  //
+  // Inputs:
+  //   V  #V by 3 list of vertex positions
+  //   F  #F by 3 list of triangle indices into V
+  //   I  #I list of facets to consider
+  // Outputs:
+  //   v1 index of the first end point of outer edge
+  //   v2 index of the second end point of outer edge
+  //   A  #A list of facets incident to the outer edge
+  template<
+      typename DerivedV,
+      typename DerivedF,
+      typename DerivedI,
+      typename IndexType,
+      typename DerivedA
+      >
+  IGL_INLINE void outer_edge(
+          const Eigen::PlainObjectBase<DerivedV> & V,
+          const Eigen::PlainObjectBase<DerivedF> & F,
+          const Eigen::PlainObjectBase<DerivedI> & I,
+          IndexType & v1,
+          IndexType & v2,
+          Eigen::PlainObjectBase<DerivedA> & A);
+
+
+  // Find a facet that is reachable from infinity without crossing any faces.
+  // Such facet is called "outer facet."
+  //
+  // Precondition: The input mesh must have all self-intersection resolved.  I.e
+  // there is no duplicated vertices, no overlapping edge and no intersecting
+  // faces (the only exception is there could be topologically duplicated faces).
+  // See cgal::remesh_self_intersections.h for how to obtain such input.
+  //
+  // Inputs:
+  //   V  #V by 3 list of vertex positions
+  //   F  #F by 3 list of triangle indices into V
+  //   N  #N by 3 list of face normals
+  //   I  #I list of facets to consider
+  // Outputs:
+  //   f  Index of the outer facet.
+  //   flipped  true iff the normal of f points inwards.
+  template<
+      typename DerivedV,
+      typename DerivedF,
+      typename DerivedN,
+      typename DerivedI,
+      typename IndexType
+      >
+  IGL_INLINE void outer_facet(
+          const Eigen::PlainObjectBase<DerivedV> & V,
+          const Eigen::PlainObjectBase<DerivedF> & F,
+          const Eigen::PlainObjectBase<DerivedN> & N,
+          const Eigen::PlainObjectBase<DerivedI> & I,
+          IndexType & f,
+          bool & flipped);
+}
+
+#ifndef IGL_STATIC_LIBRARY
+#  include "outer_element.cpp"
+#endif
+#endif

+ 0 - 122
include/igl/outer_facet.cpp

@@ -1,122 +0,0 @@
-// 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 "outer_facet.h"
-#include "sort.h"
-#include "vertex_triangle_adjacency.h"
-#include <iostream>
-#define IGL_OUTER_FACET_DEBUG
-
-template <
-  typename DerivedV,
-  typename DerivedF,
-  typename DerivedN,
-  typename DerivedI,
-  typename f_type>
-IGL_INLINE void igl::outer_facet(
-  const Eigen::PlainObjectBase<DerivedV> & V,
-  const Eigen::PlainObjectBase<DerivedF> & F,
-  const Eigen::PlainObjectBase<DerivedN> & N,
-  const Eigen::PlainObjectBase<DerivedI> & I,
-  f_type & max_f,
-  bool & flip)
-{
-  using namespace std;
-  typedef typename DerivedV::Scalar Scalar;
-  typedef typename DerivedV::Index Index;
-  typedef 
-    typename Eigen::Matrix<Index, DerivedV::RowsAtCompileTime,1> VectorXI;
-  typedef 
-    typename Eigen::Matrix<Scalar, DerivedV::RowsAtCompileTime,1> VectorXS;
-  // "Direct repair of self-intersecting meshes" [Attene 14]
-  const Index mi = I.size();
-  assert(V.cols() == 3);
-  assert(N.cols() == 3);
-  Index max_v = -1;
-  auto generic_fabs = [&](Scalar val) { 
-      if (val >= 0) return val;
-      else return -val;
-  };
-  for(size_t d = 0;d<(size_t)V.cols();d++)
-  {
-    Scalar max_d = -1e26;
-    Scalar max_nd = -1e26;
-    for(Index i = 0;i<mi;i++)
-    {
-      const Index f = I(i);
-      const Scalar nd = N(f,d);
-      if(generic_fabs(nd)>0)
-      {
-        for(Index c = 0;c<3;c++)
-        {
-          const Index v = F(f,c);
-          if(v == max_v)
-          {
-            if(generic_fabs(nd) > max_nd)
-            {
-              // Just update max face and normal
-              max_f = f;
-              max_nd = generic_fabs(nd);
-              flip = nd<0;
-            } else if (generic_fabs(nd) == max_nd) {
-                if (nd == max_nd) {
-                    if (flip) {
-                        max_f = f;
-                        max_nd = nd;
-                        flip = false;
-                    } else if (f > (Index)max_f){
-                        max_f = f;
-                        max_nd = nd;
-                        flip = false;
-                    }
-                } else {
-                    if (flip && f < (Index)max_f) {
-                        max_f = f;
-                        max_nd = generic_fabs(nd);
-                        flip = true;
-                    }
-                }
-            }
-          }else
-          {
-            const Scalar vd = V(v,d);
-            if(vd>max_d)
-            {
-              // update max vertex, face and normal
-              max_v = v;
-              max_d = vd;
-              max_f = f;
-              max_nd = generic_fabs(nd);
-              flip = nd<0;
-            }
-          }
-        }
-      }
-    }
-    if(max_v >= 0)
-    {
-      break;
-    }
-    // if we get here and max_v is still -1 then there were no faces with
-    // |N(d)| > 0
-  }
-#ifdef IGL_OUTER_FACET_DEBUG
-  if(max_v <0)
-  {
-    cerr<<"Very degenerate case, no suitable face found."<<endl;
-  }
-#endif
-  assert(max_v >=0 && "Very degenerate case, no suitable face found.");
-}
-
-#ifdef IGL_STATIC_LIBRARY
-// Explicit template specialization
-template void igl::outer_facet<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> > const&, int&, bool&);
-template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 1, -1, -1>, unsigned long>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 1, -1, -1> > const&, unsigned long&, bool&);
-template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int&, bool&);
-template void igl::outer_facet<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, int>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int&, bool&);
-#endif

+ 0 - 46
include/igl/outer_facet.h

@@ -1,46 +0,0 @@
-// 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_OUTER_FACET_H
-#define IGL_OUTER_FACET_H
-#include "igl_inline.h"
-#include <Eigen/Core>
-namespace igl
-{
-  // Compute a single facet which is guaranteed to be part of the "outer hull of
-  // a mesh (V,F). This implementation follows Section 3.6 of "Direct repair of
-  // self-intersecting meshes" [Attene 2014].
-  //
-  // Inputs:
-  //   V  #V by 3 list of vertex positions
-  //   F  #F by 3 list of triangle indices into V
-  //   N  #F by 3 list of face normals
-  //   I  #I list of facets to actually consider
-  // Outputs:
-  //   f  index of facet into V
-  //   flip  whether facet's orientation should be flipped so that
-  //     counter-clockwise normal points outward.
-  //
-  // See also: cgal/outer_hull.h
-  template <
-    typename DerivedV,
-    typename DerivedF,
-    typename DerivedN,
-    typename DerivedI,
-    typename f_type>
-  IGL_INLINE void outer_facet(
-    const Eigen::PlainObjectBase<DerivedV> & V,
-    const Eigen::PlainObjectBase<DerivedF> & F,
-    const Eigen::PlainObjectBase<DerivedN> & N,
-    const Eigen::PlainObjectBase<DerivedI> & I,
-    f_type & f,
-    bool & flip);
-}
-#ifndef IGL_STATIC_LIBRARY
-#  include "outer_facet.cpp"
-#endif
-#endif