arap.cpp 8.3 KB

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