Browse Source

direct delta mush review fix and unit test

Dxyk 5 years ago
parent
commit
e4ab356031

+ 6 - 28
include/igl/direct_delta_mush.cpp

@@ -11,15 +11,11 @@
 template <
   typename DerivedV,
   typename DerivedF,
-  typename DerivedC,
-  typename DerivedE,
   typename DerivedOmega,
   typename DerivedU>
 IGL_INLINE void igl::direct_delta_mush(
   const Eigen::MatrixBase<DerivedV> & V,
   const Eigen::MatrixBase<DerivedF> & F,
-  const Eigen::MatrixBase<DerivedC> & C,
-  const Eigen::MatrixBase<DerivedE> & E,
   const std::vector<Eigen::Affine3d, Eigen::aligned_allocator<Eigen::Affine3d> > & T,
   const Eigen::MatrixBase<DerivedOmega> & Omega,
   Eigen::PlainObjectBase<DerivedU> & U)
@@ -29,9 +25,6 @@ IGL_INLINE void igl::direct_delta_mush(
   // Shape checks
   assert(V.cols() == 3 && "V should contain 3D positions.");
   assert(F.cols() == 3 && "F should contain triangles.");
-  assert(C.cols() == 3 && "C should contain 3D bone endpoint positions.");
-  assert(E.cols() == 2 && "E should contain 2 endpoint indices forming bone edges.");
-  assert(E.rows() == T.size() && "E.rows() should equal to T.size()");
   assert(Omega.rows() == V.rows() && "Omega contain the same number of rows as V.");
   assert(Omega.cols() == T.size() * 10 && "Omega should have #T*10 columns.");
 
@@ -93,15 +86,11 @@ IGL_INLINE void igl::direct_delta_mush(
 template <
   typename DerivedV,
   typename DerivedF,
-  typename DerivedC,
-  typename DerivedE,
   typename DerivedW,
   typename DerivedOmega>
 IGL_INLINE void igl::direct_delta_mush_precomputation(
   const Eigen::MatrixBase<DerivedV> & V,
   const Eigen::MatrixBase<DerivedF> & F,
-  const Eigen::MatrixBase<DerivedC> & C,
-  const Eigen::MatrixBase<DerivedE> & E,
   const Eigen::SparseMatrix<DerivedW> & W,
   const int p,
   const typename DerivedV::Scalar lambda,
@@ -114,10 +103,7 @@ IGL_INLINE void igl::direct_delta_mush_precomputation(
   // Shape checks
   assert(V.cols() == 3 && "V should contain 3D positions.");
   assert(F.cols() == 3 && "F should contain triangles.");
-  assert(C.cols() == 3 && "C should contain 3D bone endpoint positions.");
-  assert(E.cols() == 2 && "E should contain 2 endpoint indices forming bone edges.");
   assert(W.rows() == V.rows() && "W.rows() should be equal to V.rows().");
-  assert(W.cols() == E.rows() && "W.cols() should be equal to E.rows().");
 
   // Parameter checks
   assert(p > 0 && "Laplacian iteration p should be positive.");
@@ -151,7 +137,7 @@ IGL_INLINE void igl::direct_delta_mush_precomputation(
   };
 
   const int n = V.rows();
-  const int m = E.rows();
+  const int m = W.cols();
 
   // V_homogeneous: #V by 4, homogeneous version of V
   // Note:
@@ -188,11 +174,10 @@ IGL_INLINE void igl::direct_delta_mush_precomputation(
   // C positive semi-definite => ldlt solver
   SimplicialLDLT<SparseMatrix<DerivedW>> ldlt_W_prime;
   SparseMatrix<Scalar> c(I + kappa * L_bar);
-  SparseMatrix<DerivedW> W_prime(W);
+  Matrix<DerivedW, Dynamic, Dynamic> W_prime(W);
   ldlt_W_prime.compute(c.transpose());
   for (int iter = 0; iter < p; ++iter)
   {
-    W_prime.makeCompressed();
     W_prime = ldlt_W_prime.solve(W_prime);
   }
 
@@ -229,17 +214,14 @@ IGL_INLINE void igl::direct_delta_mush_precomputation(
   //   Psi = b.ldlt().solve(Psi);  // hangs here
   // }
   // Convert to sparse matrices and compute
-  Matrix<Scalar, Dynamic, Dynamic> Psi;
-  SparseMatrix<Scalar> Psi_sparse = U_prime.sparseView();
+  Matrix<Scalar, Dynamic, Dynamic> Psi = U_prime.sparseView();
   SparseMatrix<Scalar> b = (I + lambda * L_bar).transpose();
   SimplicialLDLT<SparseMatrix<Scalar>> ldlt_Psi;
   ldlt_Psi.compute(b);
   for (int iter = 0; iter < p; ++iter)
   {
-    Psi_sparse.makeCompressed();
-    Psi_sparse = ldlt_Psi.solve(Psi_sparse);
+    Psi = ldlt_Psi.solve(Psi);
   }
-  Psi = Psi_sparse.toDense();
 
   // P: #V by 10 precomputed upper triangle of
   //    p_i p_i^T , p_i
@@ -285,9 +267,7 @@ IGL_INLINE void igl::direct_delta_mush_precomputation(
 
 // Explicit template instantiation
 template void
-igl::direct_delta_mush<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>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(
-  Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &,
-  Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &,
+igl::direct_delta_mush<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<double, -1, -1, 0, -1, -1> >(
   Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &,
   Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &,
   std::vector<Eigen::Transform<double, 3, 2, 0>, Eigen::aligned_allocator<Eigen::Transform<double, 3, 2, 0> > > const &,
@@ -295,9 +275,7 @@ igl::direct_delta_mush<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<i
   Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &);
 
 template void
-igl::direct_delta_mush_precomputation<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>, double, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(
-  Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &,
-  Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &,
+igl::direct_delta_mush_precomputation<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(
   Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &,
   Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &,
   Eigen::SparseMatrix<double, 0, int> const &, int,

+ 1 - 13
include/igl/direct_delta_mush.h

@@ -22,8 +22,6 @@ namespace igl {
   // Inputs:
   //   V  #V by 3 list of rest pose vertex positions
   //   F  #F by 3 list of triangle indices into rows of V
-  //   C  #C by 3 list of rest pose bone endpoint positions
-  //   E  #T by 2 list of bone edge indices into rows of C
   //   T  #T list of bone pose transformations
   //   Omega #V by #T*10 list of precomputated matrix values
   // Outputs:
@@ -31,15 +29,11 @@ namespace igl {
   template <
     typename DerivedV,
     typename DerivedF,
-    typename DerivedC,
-    typename DerivedE,
     typename DerivedOmega,
     typename DerivedU>
   IGL_INLINE void direct_delta_mush(
     const Eigen::MatrixBase<DerivedV> & V,
     const Eigen::MatrixBase<DerivedF> & F,
-    const Eigen::MatrixBase<DerivedC> & C,
-    const Eigen::MatrixBase<DerivedE> & E,
     const std::vector<
       Eigen::Affine3d, Eigen::aligned_allocator<Eigen::Affine3d>
     > & T, /* should eventually be templated more generally than double */
@@ -51,9 +45,7 @@ namespace igl {
   // Inputs:
   //   V  #V by 3 list of rest pose vertex positions
   //   F  #F by 3 list of triangle indices into rows of V
-  //   C  #C by 3 list of rest pose bone endpoint positions
-  //   E  #E by 2 list of bone edge indices into rows of C
-  //   W  #V by #E list of weights
+  //   W  #V by #Edges list of weights
   //   p  number of smoothing iterations
   //   lambda  rotation smoothing step size
   //   kappa   translation smoothness step size
@@ -63,15 +55,11 @@ namespace igl {
   template <
     typename DerivedV,
     typename DerivedF,
-    typename DerivedC,
-    typename DerivedE,
     typename DerivedW,
     typename DerivedOmega>
   IGL_INLINE void direct_delta_mush_precomputation(
     const Eigen::MatrixBase<DerivedV> & V,
     const Eigen::MatrixBase<DerivedF> & F,
-    const Eigen::MatrixBase<DerivedC> & C,
-    const Eigen::MatrixBase<DerivedE> & E,
     const Eigen::SparseMatrix<DerivedW> & W,
     const int p,
     const typename DerivedV::Scalar lambda,

+ 25 - 0
tests/include/igl/direct_delta_mush.cpp

@@ -0,0 +1,25 @@
+#include <test_common.h>
+#include <igl/direct_delta_mush.cpp>
+#include <igl/adjacency_list.h>
+#include <iostream>
+
+TEST_CASE("direct_delta_mush: cube", "[igl]")
+{
+  Eigen::MatrixXd V, W, U, Omega;
+  Eigen::MatrixXi F;
+
+  // Test that direct delta mush with identity transform reproduces the rest state of the mesh.
+  // For simplicity, testing on a cube of dimensions 1.0 x 1.0 x 1.0,
+  // with all vertices bound strictly to one and only one imaginary bone (weight is a column of 1s)
+  igl::read_triangle_mesh(test_common::data_path("cube.off"), V, F);
+  Eigen::SparseMatrix<double> W_sparse = Eigen::MatrixXd::Ones(V.rows(), 1).sparseView();
+
+  // Parameters such as p, lambda, kappa, alpha do not matter given identity transform
+  igl::direct_delta_mush_precomputation(V, F, W_sparse, 1, 1., 0.5, 0.5, Omega);
+
+  std::vector<Eigen::Affine3d, Eigen::aligned_allocator<Eigen::Affine3d>> T_list;
+  T_list.push_back(Eigen::Affine3d::Identity());
+  igl::direct_delta_mush(V, F, T_list, Omega, U);
+
+  test_common::assert_near(U, V, 1e-4);
+}

+ 15 - 2
tutorial/408_DirectDeltaMush/main.cpp

@@ -15,6 +15,10 @@
 
 #include "tutorial_shared_path.h"
 
+// TODO:
+// This tutorial example shares a lot of code with 404_DualQuaternionSkinning.
+// We should either combine them or pull out the shared code.
+
 typedef std::vector<Eigen::Quaterniond, Eigen::aligned_allocator<Eigen::Quaterniond>>
   RotationList;
 
@@ -76,7 +80,7 @@ bool pre_draw(igl::opengl::glfw::Viewer & viewer)
     }
     if (use_ddm)
     {
-      igl::direct_delta_mush(V, F, C, BE, T_list, Omega, U);
+      igl::direct_delta_mush(V, F, T_list, Omega, U);
     }
     else
     {
@@ -139,6 +143,15 @@ int main(int argc, char *argv[])
   poses[3][2] = rest_pose[2] * bend * rest_pose[2].conjugate();
 
   igl::readDMAT(TUTORIAL_SHARED_PATH "/arm-weights.dmat", W);
+  // convert weight to piecewise-rigid weights
+  for (int i = 0; i < W.rows(); ++i)
+  {
+    const double Wi_max = W.row(i).maxCoeff();
+    for (int j = 0; j < W.cols(); j++)
+    {
+      W(i, j) = double(W(i, j) == Wi_max);
+    }
+  }
   SparseMatrix<double> W_sparse = W.sparseView();
   igl::lbs_matrix(V, W, M);
 
@@ -148,7 +161,7 @@ int main(int argc, char *argv[])
        << "- lambda: " << lambda << endl
        << "- kappa: " << kappa << endl
        << "- alpha: " << alpha << endl;
-  igl::direct_delta_mush_precomputation(V, F, C, BE, W_sparse, p, lambda, kappa, alpha, Omega);
+  igl::direct_delta_mush_precomputation(V, F, W_sparse, p, lambda, kappa, alpha, Omega);
 
   // Plot the mesh with pseudocolors
   igl::opengl::glfw::Viewer viewer;