arap.cpp 7.1 KB

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