arap.cpp 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@gmail.com>
  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 "arap.h"
  9. #include <igl/colon.h>
  10. #include <igl/cotmatrix.h>
  11. #include <igl/massmatrix.h>
  12. #include <igl/group_sum_matrix.h>
  13. #include <igl/covariance_scatter_matrix.h>
  14. #include <igl/speye.h>
  15. #include <igl/mode.h>
  16. #include <igl/slice.h>
  17. #include <igl/arap_rhs.h>
  18. #include <igl/repdiag.h>
  19. #include <igl/columnize.h>
  20. #include "fit_rotations.h"
  21. #include <cassert>
  22. #include <iostream>
  23. template <
  24. typename DerivedV,
  25. typename DerivedF,
  26. typename Derivedb>
  27. IGL_INLINE bool igl::arap_precomputation(
  28. const Eigen::PlainObjectBase<DerivedV> & V,
  29. const Eigen::PlainObjectBase<DerivedF> & F,
  30. const Eigen::PlainObjectBase<Derivedb> & b,
  31. ARAPData & data)
  32. {
  33. using namespace igl;
  34. using namespace Eigen;
  35. typedef typename DerivedV::Scalar Scalar;
  36. // number of vertices
  37. const int n = V.rows();
  38. data.n = n;
  39. assert((b.size() == 0 || b.maxCoeff() < n) && "b out of bounds");
  40. assert((b.size() == 0 || b.minCoeff() >=0) && "b out of bounds");
  41. // remember b
  42. data.b = b;
  43. //assert(F.cols() == 3 && "For now only triangles");
  44. // dimension
  45. const int dim = V.cols();
  46. assert(dim == 3 && "Only 3d supported");
  47. // Defaults
  48. data.f_ext = Eigen::MatrixXd::Zero(n,dim);
  49. SparseMatrix<Scalar> L;
  50. cotmatrix(V,F,L);
  51. ARAPEnergyType eff_energy = data.energy;
  52. if(eff_energy == ARAP_ENERGY_TYPE_DEFAULT)
  53. {
  54. switch(F.cols())
  55. {
  56. case 3:
  57. if(dim == 3)
  58. {
  59. eff_energy = ARAP_ENERGY_TYPE_SPOKES_AND_RIMS;
  60. }else
  61. {
  62. eff_energy = ARAP_ENERGY_TYPE_ELEMENTS;
  63. }
  64. break;
  65. case 4:
  66. eff_energy = ARAP_ENERGY_TYPE_ELEMENTS;
  67. break;
  68. default:
  69. assert(false);
  70. }
  71. }
  72. // Get covariance scatter matrix, when applied collects the covariance
  73. // matrices used to fit rotations to during optimization
  74. covariance_scatter_matrix(V,F,eff_energy,data.CSM);
  75. // Get group sum scatter matrix, when applied sums all entries of the same
  76. // group according to G
  77. SparseMatrix<double> G_sum;
  78. if(data.G.size() == 0)
  79. {
  80. if(eff_energy == ARAP_ENERGY_TYPE_ELEMENTS)
  81. {
  82. speye(F.rows(),G_sum);
  83. }else
  84. {
  85. speye(n,G_sum);
  86. }
  87. }else
  88. {
  89. // groups are defined per vertex, convert to per face using mode
  90. if(eff_energy == ARAP_ENERGY_TYPE_ELEMENTS)
  91. {
  92. Eigen::Matrix<int,Eigen::Dynamic,1> GG;
  93. MatrixXi GF(F.rows(),F.cols());
  94. for(int j = 0;j<F.cols();j++)
  95. {
  96. Matrix<int,Eigen::Dynamic,1> GFj;
  97. slice(data.G,F.col(j),GFj);
  98. GF.col(j) = GFj;
  99. }
  100. mode<int>(GF,2,GG);
  101. data.G=GG;
  102. }
  103. //printf("group_sum_matrix()\n");
  104. group_sum_matrix(data.G,G_sum);
  105. }
  106. SparseMatrix<double> G_sum_dim;
  107. repdiag(G_sum,dim,G_sum_dim);
  108. assert(G_sum_dim.cols() == data.CSM.rows());
  109. data.CSM = (G_sum_dim * data.CSM).eval();
  110. arap_rhs(V,F,eff_energy,data.K);
  111. SparseMatrix<double> Q = (-0.5*L).eval();
  112. if(data.with_dynamics)
  113. {
  114. const double h = data.h;
  115. assert(h != 0);
  116. SparseMatrix<double> M;
  117. massmatrix(V,F,MASSMATRIX_DEFAULT,data.M);
  118. SparseMatrix<double> DQ = 0.5/(h*h)*data.M;
  119. Q += DQ;
  120. // Dummy external forces
  121. data.f_ext = MatrixXd::Zero(n,dim);
  122. data.vel = MatrixXd::Zero(n,dim);
  123. }
  124. return min_quad_with_fixed_precompute(
  125. Q,b,SparseMatrix<double>(),true,data.solver_data);
  126. }
  127. template <
  128. typename Derivedbc,
  129. typename DerivedU>
  130. IGL_INLINE bool igl::arap_solve(
  131. const Eigen::PlainObjectBase<Derivedbc> & bc,
  132. ARAPData & data,
  133. Eigen::PlainObjectBase<DerivedU> & U)
  134. {
  135. using namespace igl;
  136. using namespace Eigen;
  137. using namespace std;
  138. assert(data.b.size() == bc.rows());
  139. const int dim = bc.cols();
  140. const int n = data.n;
  141. int iter = 0;
  142. if(U.size() == 0)
  143. {
  144. // terrible initial guess.. should at least copy input mesh
  145. U = MatrixXd::Zero(data.n,dim);
  146. }
  147. // changes each arap iteration
  148. MatrixXd U_prev = U;
  149. // doesn't change for fixed with_dynamics timestep
  150. MatrixXd U0;
  151. if(data.with_dynamics)
  152. {
  153. U0 = U_prev;
  154. }
  155. while(iter < data.max_iter)
  156. {
  157. U_prev = U;
  158. // enforce boundary conditions exactly
  159. for(int bi = 0;bi<bc.rows();bi++)
  160. {
  161. U.row(data.b(bi)) = bc.row(bi);
  162. }
  163. MatrixXd S = data.CSM * U.replicate(dim,1);
  164. MatrixXd R(dim,data.CSM.rows());
  165. #ifdef __SSE__ // fit_rotations_SSE will convert to float if necessary
  166. fit_rotations_SSE(S,R);
  167. #else
  168. fit_rotations(S,R);
  169. #endif
  170. //for(int k = 0;k<(data.CSM.rows()/dim);k++)
  171. //{
  172. // R.block(0,dim*k,dim,dim) = MatrixXd::Identity(dim,dim);
  173. //}
  174. // Number of rotations: #vertices or #elements
  175. int num_rots = data.K.cols()/dim/dim;
  176. // distribute group rotations to vertices in each group
  177. MatrixXd eff_R;
  178. if(data.G.size() == 0)
  179. {
  180. // copy...
  181. eff_R = R;
  182. }else
  183. {
  184. eff_R.resize(dim,num_rots*dim);
  185. for(int r = 0;r<num_rots;r++)
  186. {
  187. eff_R.block(0,dim*r,dim,dim) =
  188. R.block(0,dim*data.G(r),dim,dim);
  189. }
  190. }
  191. MatrixXd Dl;
  192. if(data.with_dynamics)
  193. {
  194. assert(M.rows() == n &&
  195. "No mass matrix. Call arap_precomputation if changing with_dynamics");
  196. const double h = data.h;
  197. assert(h != 0);
  198. //Dl = 1./(h*h*h)*M*(-2.*V0 + Vm1) - fext;
  199. // data.vel = (V0-Vm1)/h
  200. // h*data.vel = (V0-Vm1)
  201. // -h*data.vel = -V0+Vm1)
  202. // -V0-h*data.vel = -2V0+Vm1
  203. Dl = 1./(h*h)*data.M*(-U0 - h*data.vel) - data.f_ext;
  204. }
  205. VectorXd Rcol;
  206. columnize(eff_R,num_rots,2,Rcol);
  207. VectorXd Bcol = -data.K * Rcol;
  208. for(int c = 0;c<dim;c++)
  209. {
  210. VectorXd Uc,Bc,bcc,Beq;
  211. Bc = Bcol.block(c*n,0,n,1);
  212. if(data.with_dynamics)
  213. {
  214. Bc += Dl.col(c);
  215. }
  216. bcc = bc.col(c);
  217. min_quad_with_fixed_solve(
  218. data.solver_data,
  219. Bc,bcc,Beq,
  220. Uc);
  221. U.col(c) = Uc;
  222. }
  223. iter++;
  224. }
  225. if(data.with_dynamics)
  226. {
  227. // Keep track of velocity for next time
  228. data.vel = (U-U0)/data.h;
  229. }
  230. return true;
  231. }
  232. #ifndef IGL_HEADER_ONLY
  233. template bool igl::arap_precomputation<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(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<int, -1, 1, 0, -1, 1> > const&, igl::ARAPData&);
  234. template bool igl::arap_solve<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&, igl::ARAPData&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  235. #endif