Przeglądaj źródła

arap parameterization working, really

Former-commit-id: 1d0ab487b9cdab0a60e2be847e8c1d4b18189039
Alec Jacobson 11 lat temu
rodzic
commit
49dadb340c

+ 26 - 8
include/igl/arap_rhs.cpp

@@ -20,8 +20,8 @@ IGL_INLINE void igl::arap_rhs(
 {
   using namespace igl;
   using namespace Eigen;
-  //// Number of dimensions
-  //int dim = V.cols();
+  // Number of dimensions
+  int Vdim = V.cols();
   //// Number of mesh vertices
   //int n = V.rows();
   //// Number of mesh elements
@@ -42,7 +42,7 @@ IGL_INLINE void igl::arap_rhs(
     default:
       fprintf(
         stderr,
-        "covariance_scatter_matrix.h: Error: Unsupported arap energy %d\n",
+        "arap_rhs.h: Error: Unsupported arap energy %d\n",
         energy);
       return;
   }
@@ -50,19 +50,37 @@ IGL_INLINE void igl::arap_rhs(
   SparseMatrix<double> KX,KY,KZ;
   arap_linear_block(V,F,0,energy,KX);
   arap_linear_block(V,F,1,energy,KY);
-  if(dim == 2)
+  if(Vdim == 2)
   {
     K = cat(2,repdiag(KX,dim),repdiag(KY,dim));
-  }else if(dim == 3)
+  }else if(Vdim == 3)
   {
     arap_linear_block(V,F,2,energy,KZ);
-    K = cat(2,cat(2,repdiag(KX,dim),repdiag(KY,dim)),repdiag(KZ,dim));
+    if(dim == 3)
+    {
+      K = cat(2,cat(2,repdiag(KX,dim),repdiag(KY,dim)),repdiag(KZ,dim));
+    }else if(dim ==2)
+    {
+      SparseMatrix<double> ZZ(KX.rows()*2,KX.cols());
+      K = cat(2,cat(2,
+            cat(2,repdiag(KX,dim),ZZ),
+            cat(2,repdiag(KY,dim),ZZ)),
+            cat(2,repdiag(KZ,dim),ZZ));
+    }else
+    {
+      assert(false);
+      fprintf(
+      stderr,
+      "arap_rhs.h: Error: Unsupported dimension %d\n",
+      dim);
+    }
   }else
   {
+    assert(false);
     fprintf(
      stderr,
-     "covariance_scatter_matrix.h: Error: Unsupported dimension %d\n",
-     dim);
+     "arap_rhs.h: Error: Unsupported dimension %d\n",
+     Vdim);
     return;
   }
   

+ 2 - 3
include/igl/covariance_scatter_matrix.cpp

@@ -17,7 +17,6 @@
 IGL_INLINE void igl::covariance_scatter_matrix(
   const Eigen::MatrixXd & V, 
   const Eigen::MatrixXi & F,
-  const int dim,
   const ARAPEnergyType energy,
   Eigen::SparseMatrix<double>& CSM)
 {
@@ -26,8 +25,8 @@ IGL_INLINE void igl::covariance_scatter_matrix(
   // number of mesh vertices
   int n = V.rows();
   assert(n > F.maxCoeff());
-  //// dimension of mesh
-  //int dim = V.cols();
+  // dimension of mesh
+  int dim = V.cols();
   // Number of mesh elements
   int m = F.rows();
 

+ 2 - 5
include/igl/covariance_scatter_matrix.h

@@ -17,20 +17,17 @@ namespace igl
 {
   // Construct the covariance scatter matrix for a given arap energy
   // Inputs:
-  //   V  #V by dim list of initial domain positions
+  //   V  #V by Vdim list of initial domain positions
   //   F  #F by 3 list of triangle indices into V
-  //   dim  dimension being used at solve time. For deformation usually dim =
-  //     V.cols(), for surface parameterization V.cols() = 3 and dim = 2
   //   energy  ARAPEnergyType enum value defining which energy is being used.
   //     See ARAPEnergyType.h for valid options and explanations.
   // Outputs:
-  //   CSM dim*#V by dim*#V sparse matrix containing special laplacians along
+  //   CSM dim*#V/#F by dim*#V sparse matrix containing special laplacians along
   //     the diagonal so that when multiplied by V gives covariance matrix
   //     elements, can be used to speed up covariance matrix computation
   IGL_INLINE void covariance_scatter_matrix(
     const Eigen::MatrixXd & V, 
     const Eigen::MatrixXi & F,
-    const int dim,
     const ARAPEnergyType energy,
     Eigen::SparseMatrix<double>& CSM);
 }

+ 1 - 0
include/igl/matlab_format.cpp

@@ -118,4 +118,5 @@ template Eigen::WithFormat<Eigen::Matrix<double, 3, 1, 0, 3, 1> > const igl::mat
 template Eigen::WithFormat<Eigen::Matrix<float, 2, 1, 0, 2, 1> > const igl::matlab_format<Eigen::Matrix<float, 2, 1, 0, 2, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, 2, 1, 0, 2, 1> > const&, std::basic_string<char, std::char_traits<char>, std::allocator<char> >);
 template Eigen::WithFormat<Eigen::Matrix<float, 2, 2, 1, 2, 2> > const igl::matlab_format<Eigen::Matrix<float, 2, 2, 1, 2, 2> >(Eigen::PlainObjectBase<Eigen::Matrix<float, 2, 2, 1, 2, 2> > const&, std::basic_string<char, std::char_traits<char>, std::allocator<char> >);
 template Eigen::WithFormat<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const igl::matlab_format<Eigen::Matrix<float, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, std::basic_string<char, std::char_traits<char>, std::allocator<char> >);
+template Eigen::WithFormat<Eigen::Matrix<double, 3, 3, 1, 3, 3> > const igl::matlab_format<Eigen::Matrix<double, 3, 3, 1, 3, 3> >(Eigen::PlainObjectBase<Eigen::Matrix<double, 3, 3, 1, 3, 3> > const&, std::basic_string<char, std::char_traits<char>, std::allocator<char> >);
 #endif

+ 9 - 0
include/igl/svd3x3/Makefile

@@ -27,6 +27,15 @@ SINGULAR_VALUE_DECOMPOSITION_INC=\
 EIGEN3_INC=-I$(DEFAULT_PREFIX)/include/eigen3 -I$(DEFAULT_PREFIX)/include/eigen3/unsupported
 INC+=$(EIGEN3_INC) $(SINGULAR_VALUE_DECOMPOSITION_INC)
 
+# Matlab dependency
+ifndef MATLAB
+  MATLAB=/Applications/MATLAB_R2013b.app/
+  $(warning MATLAB undefined. Setting MATLAB=${MATLAB})
+endif
+MATLAB_INC=-I$(MATLAB)/extern/include/
+MATLAB_LIB=-L$(MATLAB)/bin/maci64 -lmx -leng
+INC+=$(MATLAB_INC)
+
 obj: 
 	mkdir -p obj
 

+ 25 - 14
include/igl/svd3x3/arap.cpp

@@ -21,7 +21,6 @@
 #include <cassert>
 #include <iostream>
 
-
 template <
   typename DerivedV,
   typename DerivedF,
@@ -49,6 +48,7 @@ IGL_INLINE bool igl::arap_precomputation(
   //const int dim = V.cols();
   assert((dim == 3 || dim ==2) && "dim should be 2 or 3");
   data.dim = dim;
+  data.Vdim = V.cols();
   //assert(dim == 3 && "Only 3d supported");
   // Defaults
   data.f_ext = Eigen::MatrixXd::Zero(n,data.dim);
@@ -81,7 +81,7 @@ IGL_INLINE bool igl::arap_precomputation(
 
   // Get covariance scatter matrix, when applied collects the covariance
   // matrices used to fit rotations to during optimization
-  covariance_scatter_matrix(V,F,data.dim,eff_energy,data.CSM);
+  covariance_scatter_matrix(V,F,eff_energy,data.CSM);
 
   // Get group sum scatter matrix, when applied sums all entries of the same
   // group according to G
@@ -115,12 +115,13 @@ IGL_INLINE bool igl::arap_precomputation(
     group_sum_matrix(data.G,G_sum);
   }
   SparseMatrix<double> G_sum_dim;
-  repdiag(G_sum,data.dim,G_sum_dim);
+  repdiag(G_sum,data.Vdim,G_sum_dim);
   assert(G_sum_dim.cols() == data.CSM.rows());
   data.CSM = (G_sum_dim * data.CSM).eval();
 
 
   arap_rhs(V,F,data.dim,eff_energy,data.K);
+  assert(data.K.rows() == data.n*data.dim);
 
   SparseMatrix<double> Q = (-0.5*L).eval();
 
@@ -187,19 +188,25 @@ IGL_INLINE bool igl::arap_solve(
       U.row(data.b(bi)) = bc.row(bi);
     }
 
-    const auto & Udim = U.replicate(data.dim,1);
-    MatrixXd S = data.CSM * Udim;
+    const auto & Udim = U.replicate(data.Vdim,1);
+    assert(U.cols() == data.dim);
+    assert(data.Vdim >= data.dim);
+    // As if U.col(2) was 0
+    MatrixXd S = MatrixXd::Zero(data.CSM.rows(),data.Vdim);
+    S.leftCols(data.dim) = data.CSM * Udim;
 
-    MatrixXd R(data.dim,data.CSM.rows());
-    if(data.dim == 2)
+    const int Rdim = data.Vdim;
+    MatrixXd R(Rdim,data.CSM.rows());
+    if(R.rows() == 2)
     {
       fit_rotations_planar(S,R);
+      cerr<<"PLANAR"<<endl;
     }else
     {
 #ifdef __SSE__ // fit_rotations_SSE will convert to float if necessary
-      fit_rotations_SSE(S,R);
+    fit_rotations_SSE(S,R);
 #else
-      fit_rotations(S,R);
+    fit_rotations(S,true,R);
 #endif
     }
     //for(int k = 0;k<(data.CSM.rows()/dim);k++)
@@ -209,7 +216,7 @@ IGL_INLINE bool igl::arap_solve(
 
 
     // Number of rotations: #vertices or #elements
-    int num_rots = data.K.cols()/data.dim/data.dim;
+    int num_rots = data.K.cols()/Rdim/Rdim;
     // distribute group rotations to vertices in each group
     MatrixXd eff_R;
     if(data.G.size() == 0)
@@ -218,11 +225,11 @@ IGL_INLINE bool igl::arap_solve(
       eff_R = R;
     }else
     {
-      eff_R.resize(data.dim,num_rots*data.dim);
+      eff_R.resize(Rdim,num_rots*Rdim);
       for(int r = 0;r<num_rots;r++)
       {
-        eff_R.block(0,data.dim*r,data.dim,data.dim) = 
-          R.block(0,data.dim*data.G(r),data.dim,data.dim);
+        eff_R.block(0,Rdim*r,Rdim,Rdim) = 
+          R.block(0,Rdim*data.G(r),Rdim,Rdim);
       }
     }
 
@@ -244,6 +251,7 @@ IGL_INLINE bool igl::arap_solve(
     VectorXd Rcol;
     columnize(eff_R,num_rots,2,Rcol);
     VectorXd Bcol = -data.K * Rcol;
+    assert(Bcol.size() == data.n*data.dim);
     for(int c = 0;c<data.dim;c++)
     {
       VectorXd Uc,Bc,bcc,Beq;
@@ -252,7 +260,10 @@ IGL_INLINE bool igl::arap_solve(
       {
         Bc += Dl.col(c);
       }
-      bcc = bc.col(c);
+      if(bc.size()>0)
+      {
+        bcc = bc.col(c);
+      }
       min_quad_with_fixed_solve(
         data.solver_data,
         Bc,bcc,Beq,

+ 4 - 1
include/igl/svd3x3/arap.h

@@ -32,6 +32,7 @@ namespace igl
     // solver_data  quadratic solver data
     // b  list of boundary indices into V
     // dim  dimension being used for solving
+    // Vdim  dimension of original V (determines dimension of rotation fitting)
     int n;
     Eigen::VectorXi G;
     ARAPEnergyType energy;
@@ -44,6 +45,7 @@ namespace igl
     min_quad_with_fixed_data<double> solver_data;
     Eigen::VectorXi b;
     int dim;
+    int Vdim;
       ARAPData():
         n(0),
         G(),
@@ -56,7 +58,8 @@ namespace igl
         CSM(),
         solver_data(),
         b(),
-        dim(-1) // force this to be set by _precomputation
+        dim(-1), // force this to be set by _precomputation
+        Vdim(-1) // force this to be set by _precomputation
     {
     };
   };

+ 1 - 1
include/igl/svd3x3/arap_dof.cpp

@@ -123,7 +123,7 @@ IGL_INLINE bool igl::arap_dof_precomputation(
   // used to fit rotations to during optimization
   SparseMatrix<double> CSM;
   //printf("covariance_scatter_matrix()\n");
-  covariance_scatter_matrix(V,F,V.cols(),data.energy,CSM);
+  covariance_scatter_matrix(V,F,data.energy,CSM);
 #ifdef EXTREME_VERBOSE
   cout<<"CSMIJV=["<<endl;print_ijv(CSM,1);cout<<endl<<"];"<<
     endl<<"CSM=sparse(CSMIJV(:,1),CSMIJV(:,2),CSMIJV(:,3),"<<

+ 24 - 5
include/igl/svd3x3/fit_rotations.cpp

@@ -12,17 +12,20 @@
 #include <igl/polar_dec.h>
 #include <igl/polar_svd.h>
 #include <igl/matlab_format.h>
+#include <igl/C_STR.h>
 #include <iostream>
 
 template <typename DerivedS, typename DerivedD>
 IGL_INLINE void igl::fit_rotations(
   const Eigen::PlainObjectBase<DerivedS> & S,
+  const bool single_precision,
         Eigen::PlainObjectBase<DerivedD> & R)
 {
   using namespace std;
   const int dim = S.cols();
   const int nr = S.rows()/dim;
   assert(nr * dim == S.rows());
+  assert(dim == 3);
 
   // resize output
   R.resize(dim,dim*nr); // hopefully no op (should be already allocated)
@@ -41,12 +44,29 @@ IGL_INLINE void igl::fit_rotations(
         si(i,j) = S(i*nr+r,j);
       }
     }
-    Eigen::Matrix<typename DerivedD::Scalar,3,3> ri;
-    Eigen::Matrix<typename DerivedD::Scalar,3,3> ti;
-    polar_svd3x3(si, ri);
+    typedef Eigen::Matrix<typename DerivedD::Scalar,3,3> Mat3;
+    typedef Eigen::Matrix<typename DerivedD::Scalar,3,1> Vec3;
+    Mat3 ri;
+    if(single_precision)
+    {
+      polar_svd3x3(si, ri);
+    }else
+    {
+      Mat3 ti,ui,vi;
+      Vec3 _;
+      igl::polar_svd(si,ri,ti,ui,_,vi);
+      // Check for reflection
+      if(ri.determinant() < 0)
+      {
+        vi.col(1) *= -1.;
+        ri = ui * vi.transpose();
+      }
+    }
     assert(ri.determinant() >= 0);
     // Not sure why polar_dec computes transpose...
     R.block(0,r*dim,dim,dim) = ri.block(0,0,dim,dim).transpose();
+    //cout<<matlab_format(si,C_STR("si_"<<r))<<endl;
+    //cout<<matlab_format(ri.transpose().eval(),C_STR("ri_"<<r))<<endl;
   }
 }
 
@@ -81,7 +101,6 @@ IGL_INLINE void igl::fit_rotations_planar(
     Mat2 ri,ti,ui,vi;
     Vec2 _;
     igl::polar_svd(si,ri,ti,ui,_,vi);
-
 #ifndef FIT_ROTATIONS_ALLOW_FLIPS
     // Check for reflection
     if(ri.determinant() < 0)
@@ -206,7 +225,7 @@ IGL_INLINE void igl::fit_rotations_AVX(
 
 #ifndef IGL_HEADER_ONLY
 // Explicit template instanciation
-template void igl::fit_rotations<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
+template void igl::fit_rotations<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
 template void igl::fit_rotations_planar<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
 template void igl::fit_rotations_planar<Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&);
 #endif

+ 2 - 0
include/igl/svd3x3/fit_rotations.h

@@ -17,12 +17,14 @@ namespace igl
   // 
   // Inputs:
   //   S  nr*dim by dim stack of covariance matrices
+  //   single_precision  whether to use single precision (faster)
   // Outputs:
   //   R  dim by dim * nr list of rotations
   //
   template <typename DerivedS, typename DerivedD>
   IGL_INLINE void fit_rotations(
     const Eigen::PlainObjectBase<DerivedS> & S,
+    const bool single_precision,
           Eigen::PlainObjectBase<DerivedD> & R);
   
   // FIT_ROTATIONS Given an input mesh and new positions find 2D rotations for

+ 1 - 0
include/igl/writeOBJ.cpp

@@ -120,4 +120,5 @@ IGL_INLINE bool igl::writeOBJ(
 template bool igl::writeOBJ<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::basic_string<char, std::char_traits<char>, std::allocator<char> >, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&);
 template bool igl::writeOBJ<Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, 2, 1, -1, 2> >(std::basic_string<char, std::char_traits<char>, std::allocator<char> >, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 2, 1, -1, 2> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&);
 template bool igl::writeOBJ<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<float, -1, 2, 1, -1, 2> >(std::basic_string<char, std::char_traits<char>, std::allocator<char> >, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 2, 1, -1, 2> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&);
+template bool igl::writeOBJ<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(std::basic_string<char, std::char_traits<char>, std::allocator<char> >, 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&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&);
 #endif

+ 5 - 4
tutorial/503_ARAPParam/main.cpp

@@ -34,6 +34,7 @@ bool key_down(igl::Viewer& viewer, unsigned char key, int modifier)
 
 int main(int argc, char *argv[])
 {
+  using namespace std;
   // Load a mesh in OFF format
   igl::readOFF("../shared/camelhead.off", V, F);
 
@@ -52,9 +53,9 @@ int main(int argc, char *argv[])
   Eigen::MatrixXd bc = Eigen::MatrixXd::Zero(0,0);
   
   // Initialize ARAP
-  arap_data.max_iter = 1000;
-  Eigen::MatrixXd V_2d = V.block(0,0,V.rows(),2);
-  arap_precomputation(V_2d,F,b,arap_data);
+  arap_data.max_iter = 100;
+  // 2 means that we're going to *solve* in 2d
+  arap_precomputation(V,F,2,b,arap_data);
 
 
   // Solve arap using the harmonic map as initial guess
@@ -64,7 +65,7 @@ int main(int argc, char *argv[])
 
 
   // Scale UV to make the texture more clear
-  // V_uv *= 5;
+  V_uv *= 20;
 
   // Plot the mesh
   igl::Viewer viewer;

+ 1 - 0
tutorial/CMakeLists.shared

@@ -65,6 +65,7 @@ endif(WIN32)
 
 
 set(SHARED_LIBRARIES
+		${LIBIGL_LIBRARIES}
 		${OPENGL_LIBRARIES}
 		${ANT_TWEAK_BAR_LIBRARY}
 		${EIGEN_LIBRARIES}

+ 44 - 1
tutorial/cmake/FindLIBIGL.cmake

@@ -20,12 +20,55 @@ FIND_PATH(LIBIGL_INCLUDE_DIR igl/readOBJ.h
 if(LIBIGL_INCLUDE_DIR)
    set(LIBIGL_FOUND TRUE)
    set(LIBIGL_INCLUDE_DIR ${LIBIGL_INCLUDE_DIR}  ${LIBIGL_INCLUDE_DIR}/../external/Singular_Value_Decomposition)
-   add_definitions(-DIGL_HEADER_ONLY)
+   if(NOT LIBIGL_USE_STATIC_LIBRARY)
+      add_definitions(-DIGL_HEADER_ONLY)
+   endif(NOT LIBIGL_USE_STATIC_LIBRARY)
    #set(LIBIGL_SOURCES
    #   ${LIBIGL_INCLUDE_DIR}/igl/viewer/Viewer.cpp
    #)
 endif(LIBIGL_INCLUDE_DIR)
 
+if(LIBIGL_USE_STATIC_LIBRARY)
+  set(LIBIGL_LIB_DIRS 
+   /usr/lib
+   /usr/local/lib
+   $ENV{LIBIGLROOT}/lib
+   $ENV{LIBIGL_ROOT}/lib
+   $ENV{LIBIGL_DIR}/lib
+   ${PROJECT_SOURCE_DIR}/../libigl/lib
+   ${PROJECT_SOURCE_DIR}/../../libigl/lib
+   ${PROJECT_SOURCE_DIR}/../../lib)
+  FIND_LIBRARY( LIBIGL_LIBRARY NAMES igl PATHS ${LIBIGL_LIB_DIRS})
+  if(NOT LIBIGL_LIBRARY)
+    set(LIBIGL_FOUND FALSE)
+  endif(NOT LIBIGL_LIBRARY)
+  set(LIBIGL_LIBRARIES ${LIBIGL_LIBRARIES}  ${LIBIGL_LIBRARY})
+  FIND_LIBRARY( LIBIGLMATLAB_LIBRARY NAMES iglmatlab PATHS ${LIBIGL_LIB_DIRS})
+  if(NOT LIBIGLMATLAB_LIBRARY)
+    set(LIBIGL_FOUND FALSE)
+  endif(NOT LIBIGLMATLAB_LIBRARY)
+  set(LIBIGL_LIBRARIES ${LIBIGL_LIBRARIES}  ${LIBIGLMATLAB_LIBRARY})
+  FIND_LIBRARY( LIBIGLSVD3X3_LIBRARY NAMES iglsvd3x3 PATHS ${LIBIGL_LIB_DIRS})
+  if(NOT LIBIGLSVD3X3_LIBRARY)
+    set(LIBIGL_FOUND FALSE)
+  endif(NOT LIBIGLSVD3X3_LIBRARY)
+  set(LIBIGL_LIBRARIES ${LIBIGL_LIBRARIES}  ${LIBIGLSVD3X3_LIBRARY})
+  FIND_LIBRARY( LIBIGLVIEWER_LIBRARY NAMES iglviewer PATHS ${LIBIGL_LIB_DIRS})
+  if(NOT LIBIGLVIEWER_LIBRARY)
+    set(LIBIGL_FOUND FALSE)
+  endif(NOT LIBIGLVIEWER_LIBRARY)
+  set(LIBIGL_LIBRARIES ${LIBIGL_LIBRARIES}  ${LIBIGLVIEWER_LIBRARY})
+  find_package(OpenMP)
+  if (OPENMP_FOUND)
+    set (CMAKE_C_FLAGS "${CMAKE_C_FLAGS} ${OpenMP_C_FLAGS}")
+    set (CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} ${OpenMP_CXX_FLAGS}")
+  else(OPENMP_FOUND)
+    set(LIBIGL_FOUND false)
+  endif(OPENMP_FOUND)
+endif(LIBIGL_USE_STATIC_LIBRARY)
+
+
+
 if(LIBIGL_FOUND)
    if(NOT LIBIGL_FIND_QUIETLY)
       message(STATUS "Found LIBIGL: ${LIBIGL_INCLUDE_DIR}")

+ 1 - 1
tutorial/compile_example.sh

@@ -1,5 +1,5 @@
 rm -fr build
 mkdir build
 cd build
-cmake ../
+cmake ../ $@
 make

+ 6 - 0
tutorial/readme.md

@@ -18,6 +18,8 @@ All examples depends on glfw, glew and anttweakbar. A copy
 of the sourcecode of each library is provided together with libigl
 and they can be precompiled using:
 
+** Alec: Is this just compiling the dependencies? Then perhaps rename `compile_dependencies_*`**
+
     sh compile_macosx.sh (MACOSX)
     sh compile_linux.sh (LINUX)
     compile_windows.bat (Visual Studio 2012)
@@ -27,6 +29,10 @@ On Linux and MacOSX, you can use the provided bash script:
 
     sh ../compile_example.sh
 
+## Optional compilation with libigl as static library
+
+    ../compile_example.sh -DLIBIGL_USE_STATIC_LIBRARY=ON
+
 # Chapter 2: Discrete Geometric Quantities and Operators
 This chapter illustrates a few discrete quantities that libigl can compute on a
 mesh. This also provides an introduction to basic drawing and coloring routines