direct_delta_mush.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265
  1. // This file is part of libigl, a simple C++ geometry processing library.
  2. //
  3. // Copyright (C) 2020 Xiangyu Kong <[email protected]>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "direct_delta_mush.h"
  9. #include "cotmatrix.h"
  10. #include "PlainMatrix.h"
  11. template <
  12. typename DerivedV,
  13. typename DerivedOmega,
  14. typename DerivedU>
  15. IGL_INLINE void igl::direct_delta_mush(
  16. const Eigen::MatrixBase<DerivedV> & V,
  17. const std::vector<Eigen::Affine3d, Eigen::aligned_allocator<Eigen::Affine3d> > & T,
  18. const Eigen::MatrixBase<DerivedOmega> & Omega,
  19. Eigen::PlainObjectBase<DerivedU> & U)
  20. {
  21. // Shape checks
  22. assert(V.cols() == 3 && "V should contain 3D positions.");
  23. assert(Omega.rows() == V.rows() && "Omega contain the same number of rows as V.");
  24. assert(Omega.cols() == T.size() * 10 && "Omega should have #T*10 columns.");
  25. typedef typename DerivedV::Scalar Scalar;
  26. int n = V.rows();
  27. int m = T.size();
  28. // V_homogeneous: #V by 4, homogeneous version of V
  29. // Note:
  30. // In the paper, the rest pose vertices are represented in U \in R^{4 x #V}
  31. // Thus the formulae involving U would differ from the paper by a transpose.
  32. Eigen::Matrix<Scalar, Eigen::Dynamic, 4> V_homogeneous(n, 4);
  33. V_homogeneous << V, Eigen::Matrix<Scalar, Eigen::Dynamic, 1>::Ones(n, 1);
  34. U.resize(n, 3);
  35. for (int i = 0; i < n; ++i)
  36. {
  37. // Construct Q matrix using Omega and Transformations
  38. Eigen::Matrix<Scalar, 4, 4> Q_mat(4, 4);
  39. Q_mat = Eigen::Matrix<Scalar, 4, 4>::Zero(4, 4);
  40. for (int j = 0; j < m; ++j)
  41. {
  42. Eigen::Matrix<typename DerivedOmega::Scalar, 4, 4> Omega_curr(4, 4);
  43. Eigen::Matrix<typename DerivedOmega::Scalar, 10, 1> curr = Omega.block(i, j * 10, 1, 10).transpose();
  44. Omega_curr << curr(0), curr(1), curr(2), curr(3),
  45. curr(1), curr(4), curr(5), curr(6),
  46. curr(2), curr(5), curr(7), curr(8),
  47. curr(3), curr(6), curr(8), curr(9);
  48. Eigen::Affine3d M_curr = T[j];
  49. Q_mat += M_curr.matrix() * Omega_curr;
  50. }
  51. // Normalize so that the last element is 1
  52. Q_mat /= Q_mat(Q_mat.rows() - 1, Q_mat.cols() - 1);
  53. Eigen::Matrix<Scalar, 3, 3> Q_i = Q_mat.block(0, 0, 3, 3);
  54. Eigen::Matrix<Scalar, 3, 1> q_i = Q_mat.block(0, 3, 3, 1);
  55. Eigen::Matrix<Scalar, 3, 1> p_i = Q_mat.block(3, 0, 1, 3).transpose();
  56. // Get rotation and translation matrices using SVD
  57. Eigen::Matrix<Scalar, 3, 3> SVD_i = Q_i - q_i * p_i.transpose();
  58. Eigen::JacobiSVD<Eigen::Matrix<Scalar, 3, 3>> svd;
  59. svd.compute(SVD_i, Eigen::ComputeFullU | Eigen::ComputeFullV);
  60. Eigen::Matrix<Scalar, 3, 3> R_i = svd.matrixU() * svd.matrixV().transpose();
  61. Eigen::Matrix<Scalar, 3, 1> t_i = q_i - R_i * p_i;
  62. // Gamma final transformation matrix
  63. Eigen::Matrix<Scalar, 3, 4> Gamma_i(3, 4);
  64. Gamma_i.block(0, 0, 3, 3) = R_i;
  65. Gamma_i.block(0, 3, 3, 1) = t_i;
  66. // Final deformed position
  67. Eigen::Matrix<Scalar, 4, 1> v_i = V_homogeneous.row(i);
  68. U.row(i) = Gamma_i * v_i;
  69. }
  70. }
  71. template <
  72. typename DerivedV,
  73. typename DerivedF,
  74. typename DerivedW,
  75. typename DerivedOmega>
  76. IGL_INLINE void igl::direct_delta_mush_precomputation(
  77. const Eigen::MatrixBase<DerivedV> & V,
  78. const Eigen::MatrixBase<DerivedF> & F,
  79. const Eigen::MatrixBase<DerivedW> & W,
  80. const int p,
  81. const typename DerivedV::Scalar lambda,
  82. const typename DerivedV::Scalar kappa,
  83. const typename DerivedV::Scalar alpha,
  84. Eigen::PlainObjectBase<DerivedOmega> & Omega)
  85. {
  86. // Shape checks
  87. assert(V.cols() == 3 && "V should contain 3D positions.");
  88. assert(F.cols() == 3 && "F should contain triangles.");
  89. assert(W.rows() == V.rows() && "W.rows() should be equal to V.rows().");
  90. // Parameter checks
  91. assert(p > 0 && "Laplacian iteration p should be positive.");
  92. assert(lambda > 0 && "lambda should be positive.");
  93. assert(kappa > 0 && kappa < lambda && "kappa should be positive and less than lambda.");
  94. assert(alpha >= 0 && alpha < 1 && "alpha should be non-negative and less than 1.");
  95. typedef typename DerivedV::Scalar Scalar;
  96. // lambda helper
  97. // Given a square matrix, extract the upper triangle (including diagonal) to an array.
  98. // E.g. 1 2 3 4
  99. // 5 6 7 8 -> [1, 2, 3, 4, 6, 7, 8, 11, 12, 16]
  100. // 9 10 11 12 0 1 2 3 4 5 6 7 8 9
  101. // 13 14 15 16
  102. auto extract_upper_triangle = [](
  103. const Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> & full) -> Eigen::Matrix<Scalar, Eigen::Dynamic, 1>
  104. {
  105. int dims = full.rows();
  106. Eigen::Matrix<Scalar, Eigen::Dynamic, 1> upper_triangle((dims * (dims + 1)) / 2);
  107. int vector_idx = 0;
  108. for (int i = 0; i < dims; ++i)
  109. {
  110. for (int j = i; j < dims; ++j)
  111. {
  112. upper_triangle(vector_idx) = full(i, j);
  113. vector_idx++;
  114. }
  115. }
  116. return upper_triangle;
  117. };
  118. const int n = V.rows();
  119. const int m = W.cols();
  120. // V_homogeneous: #V by 4, homogeneous version of V
  121. // Note:
  122. // in the paper, the rest pose vertices are represented in U \in R^{4 \times #V}
  123. // Thus the formulae involving U would differ from the paper by a transpose.
  124. Eigen::Matrix<Scalar, Eigen::Dynamic, 4> V_homogeneous(n, 4);
  125. V_homogeneous << V, Eigen::Matrix<Scalar, Eigen::Dynamic, 1>::Ones(n);
  126. // Identity matrix of #V by #V
  127. Eigen::SparseMatrix<Scalar> I(n, n);
  128. I.setIdentity();
  129. // Laplacian matrix of #V by #V
  130. // L_bar = L \times D_L^{-1}
  131. Eigen::SparseMatrix<Scalar> L;
  132. igl::cotmatrix(V, F, L);
  133. L = -L;
  134. // Inverse of diagonal matrix = reciprocal elements in diagonal
  135. Eigen::Matrix<Scalar, Eigen::Dynamic, 1> D_L = L.diagonal();
  136. // D_L = D_L.array().pow(-1); // Not using this since not sure if diagonal contains 0
  137. for (int i = 0; i < D_L.size(); ++i)
  138. {
  139. if (D_L(i) != 0)
  140. {
  141. D_L(i) = 1 / D_L(i);
  142. }
  143. }
  144. Eigen::SparseMatrix<Scalar> D_L_inv = D_L.asDiagonal().toDenseMatrix().sparseView();
  145. Eigen::SparseMatrix<Scalar> L_bar = L * D_L_inv;
  146. // Implicitly and iteratively solve for W'
  147. // w'_{ij} = \sum_{k=1}^{n}{C_{ki} w_{kj}} where C = (I + kappa L_bar)^{-p}:
  148. // W' = C^T \times W => c^T W_k = W_{k-1} where c = (I + kappa L_bar)
  149. // C positive semi-definite => ldlt solver
  150. Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>> ldlt_W_prime;
  151. Eigen::SparseMatrix<Scalar> c(I + kappa * L_bar);
  152. // working copy
  153. PlainMatrix<DerivedW ,Eigen::Dynamic ,Eigen::Dynamic> W_prime(W);
  154. ldlt_W_prime.compute(c.transpose());
  155. for (int iter = 0; iter < p; ++iter)
  156. {
  157. W_prime = ldlt_W_prime.solve(W_prime);
  158. }
  159. // U_precomputed: #V by 10
  160. // Cache u_i^T \dot u_i \in R^{4 x 4} to reduce computation time.
  161. Eigen::Matrix<Scalar, Eigen::Dynamic, 10> U_precomputed(n, 10);
  162. for (int k = 0; k < n; ++k)
  163. {
  164. Eigen::Matrix<Scalar, 4, 4> u_full = V_homogeneous.row(k).transpose() * V_homogeneous.row(k);
  165. U_precomputed.row(k) = extract_upper_triangle(u_full);
  166. }
  167. // U_prime: #V by #T*10 of u_{jx}
  168. // Each column of U_prime (u_{jx}) is the element-wise product of
  169. // W_j and U_precomputed_x where j \in {1...m}, x \in {1...10}
  170. Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> U_prime(n, m * 10);
  171. for (int j = 0; j < m; ++j)
  172. {
  173. Eigen::Matrix<Scalar, Eigen::Dynamic, 1> w_j = W.col(j);
  174. for (int x = 0; x < 10; ++x)
  175. {
  176. Eigen::Matrix<Scalar, Eigen::Dynamic, 1> u_x = U_precomputed.col(x);
  177. U_prime.col(10 * j + x) = w_j.array() * u_x.array();
  178. }
  179. }
  180. // Implicitly and iteratively solve for Psi: #V by #T*10 of \Psi_{ij}s.
  181. // Note: Using dense matrices to solve for Psi will cause the program to hang.
  182. // The following won't work
  183. // Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Psi(U_prime);
  184. // Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> b((I + lambda * L_bar).transpose());
  185. // for (int iter = 0; iter < p; ++iter)
  186. // {
  187. // Psi = b.ldlt().solve(Psi); // hangs here
  188. // }
  189. // Convert to sparse matrices and compute
  190. Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Psi = U_prime.sparseView();
  191. Eigen::SparseMatrix<Scalar> b = (I + lambda * L_bar).transpose();
  192. Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>> ldlt_Psi;
  193. ldlt_Psi.compute(b);
  194. for (int iter = 0; iter < p; ++iter)
  195. {
  196. Psi = ldlt_Psi.solve(Psi);
  197. }
  198. // P: #V by 10 precomputed upper triangle of
  199. // p_i p_i^T , p_i
  200. // p_i^T , 1
  201. // where p_i = (\sum_{j=1}^{n} Psi_{ij})'s top right 3 by 1 column
  202. Eigen::Matrix<Scalar, Eigen::Dynamic, 10> P(n, 10);
  203. for (int i = 0; i < n; ++i)
  204. {
  205. Eigen::Matrix<Scalar, 3, 1> p_i = Eigen::Matrix<Scalar, 3, 1>::Zero(3);
  206. Scalar last = 0;
  207. for (int j = 0; j < m; ++j)
  208. {
  209. Eigen::Matrix<Scalar, 3, 1> p_i_curr(3);
  210. p_i_curr << Psi(i, j * 10 + 3), Psi(i, j * 10 + 6), Psi(i, j * 10 + 8);
  211. p_i += p_i_curr;
  212. last += Psi(i, j * 10 + 9);
  213. }
  214. p_i /= last; // normalize
  215. Eigen::Matrix<Scalar, 4, 4> p_matrix(4, 4);
  216. p_matrix.block(0, 0, 3, 3) = p_i * p_i.transpose();
  217. p_matrix.block(0, 3, 3, 1) = p_i;
  218. p_matrix.block(3, 0, 1, 3) = p_i.transpose();
  219. p_matrix(3, 3) = 1;
  220. P.row(i) = extract_upper_triangle(p_matrix);
  221. }
  222. // Omega
  223. Omega.resize(n, m * 10);
  224. for (int i = 0; i < n; ++i)
  225. {
  226. Eigen::Matrix<Scalar, 10, 1> p_vector = P.row(i);
  227. for (int j = 0; j < m; ++j)
  228. {
  229. Eigen::Matrix<Scalar, 10, 1> Omega_curr(10);
  230. Eigen::Matrix<Scalar, 10, 1> Psi_curr = Psi.block(i, j * 10, 1, 10).transpose();
  231. Omega_curr = (1. - alpha) * Psi_curr + alpha * W_prime(i, j) * p_vector;
  232. Omega.block(i, j * 10, 1, 10) = Omega_curr.transpose();
  233. }
  234. }
  235. }
  236. #ifdef IGL_STATIC_LIBRARY
  237. // Explicit template instantiation
  238. template void igl::direct_delta_mush<Eigen::Matrix<double, -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&, std::vector<Eigen::Transform<double, 3, 2, 0>, Eigen::aligned_allocator<Eigen::Transform<double, 3, 2, 0> > > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, 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<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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  239. #endif