arap.cpp 5.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200
  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/group_sum_matrix.h>
  12. #include <igl/covariance_scatter_matrix.h>
  13. #include <igl/speye.h>
  14. #include <igl/mode.h>
  15. #include <igl/slice.h>
  16. #include <igl/arap_rhs.h>
  17. #include <igl/repdiag.h>
  18. #include <igl/columnize.h>
  19. #include "fit_rotations.h"
  20. #include <cassert>
  21. #include <iostream>
  22. template <
  23. typename DerivedV,
  24. typename DerivedF,
  25. typename Derivedb>
  26. IGL_INLINE bool igl::arap_precomputation(
  27. const Eigen::PlainObjectBase<DerivedV> & V,
  28. const Eigen::PlainObjectBase<DerivedF> & F,
  29. const Eigen::PlainObjectBase<Derivedb> & b,
  30. ARAPData & data)
  31. {
  32. using namespace igl;
  33. using namespace Eigen;
  34. typedef typename DerivedV::Scalar Scalar;
  35. // number of vertices
  36. const int n = V.rows();
  37. data.n = n;
  38. assert((b.size() == 0 || b.maxCoeff() < n) && "b out of bounds");
  39. assert((b.size() == 0 || b.minCoeff() >=0) && "b out of bounds");
  40. // remember b
  41. data.b = b;
  42. assert(F.cols() == 3 && "For now only triangles");
  43. // dimension
  44. const int dim = V.cols();
  45. assert(dim == 3 && "Only 3d supported");
  46. // Defaults
  47. data.f_ext = Eigen::MatrixXd::Zero(n,dim);
  48. SparseMatrix<Scalar> L;
  49. cotmatrix(V,F,L);
  50. ARAPEnergyType eff_energy = data.energy;
  51. if(eff_energy == ARAP_ENERGY_TYPE_DEFAULT)
  52. {
  53. switch(F.cols())
  54. {
  55. case 3:
  56. if(dim == 3)
  57. {
  58. eff_energy = ARAP_ENERGY_TYPE_SPOKES_AND_RIMS;
  59. }else
  60. {
  61. eff_energy = ARAP_ENERGY_TYPE_ELEMENTS;
  62. }
  63. break;
  64. case 4:
  65. eff_energy = ARAP_ENERGY_TYPE_ELEMENTS;
  66. default:
  67. assert(false);
  68. }
  69. }
  70. // Get covariance scatter matrix, when applied collects the covariance matrices
  71. // used to fit rotations to during optimization
  72. covariance_scatter_matrix(V,F,eff_energy,data.CSM);
  73. // Get group sum scatter matrix, when applied sums all entries of the same
  74. // group according to G
  75. SparseMatrix<double> G_sum;
  76. if(data.G.size() == 0)
  77. {
  78. speye(n,G_sum);
  79. }else
  80. {
  81. // groups are defined per vertex, convert to per face using mode
  82. Eigen::Matrix<int,Eigen::Dynamic,1> GG;
  83. if(eff_energy == ARAP_ENERGY_TYPE_ELEMENTS)
  84. {
  85. MatrixXi GF(F.rows(),F.cols());
  86. for(int j = 0;j<F.cols();j++)
  87. {
  88. Matrix<int,Eigen::Dynamic,1> GFj;
  89. slice(data.G,F.col(j),GFj);
  90. GF.col(j) = GFj;
  91. }
  92. mode<int>(GF,2,GG);
  93. }else
  94. {
  95. GG=data.G;
  96. }
  97. //printf("group_sum_matrix()\n");
  98. group_sum_matrix(GG,G_sum);
  99. }
  100. SparseMatrix<double> G_sum_dim;
  101. repdiag(G_sum,dim,G_sum_dim);
  102. data.CSM = (G_sum_dim * data.CSM).eval();
  103. arap_rhs(V,F,eff_energy,data.K);
  104. SparseMatrix<double> Q = (-0.5*L).eval();
  105. return min_quad_with_fixed_precompute(
  106. Q,b,SparseMatrix<double>(),true,data.solver_data);
  107. }
  108. template <
  109. typename Derivedbc,
  110. typename DerivedU>
  111. IGL_INLINE bool igl::arap_solve(
  112. const Eigen::PlainObjectBase<Derivedbc> & bc,
  113. ARAPData & data,
  114. Eigen::PlainObjectBase<DerivedU> & U)
  115. {
  116. using namespace igl;
  117. using namespace Eigen;
  118. using namespace std;
  119. assert(data.b.size() == bc.rows());
  120. const int dim = bc.cols();
  121. const int n = data.n;
  122. int iter = 0;
  123. if(U.size() == 0)
  124. {
  125. // terrible initial guess.. should at least copy input mesh
  126. U = MatrixXd::Zero(data.n,dim);
  127. }
  128. MatrixXd U_prev = U;
  129. while(iter < data.max_iter)
  130. {
  131. U_prev = U;
  132. // enforce boundary conditions exactly
  133. for(int bi = 0;bi<bc.rows();bi++)
  134. {
  135. U.row(data.b(bi)) = bc.row(bi);
  136. }
  137. MatrixXd S = data.CSM * U.replicate(dim,1);
  138. MatrixXd R(dim,data.CSM.rows());
  139. #ifdef __SSE__ // fit_rotations_SSE will convert to float if necessary
  140. fit_rotations_SSE(S,R);
  141. #else
  142. fit_rotations(S,R);
  143. #endif
  144. //for(int k = 0;k<(data.CSM.rows()/dim);k++)
  145. //{
  146. // R.block(0,dim*k,dim,dim) = MatrixXd::Identity(dim,dim);
  147. //}
  148. // distribute group rotations to vertices in each group
  149. MatrixXd eff_R;
  150. if(data.G.size() == 0)
  151. {
  152. // copy...
  153. eff_R = R;
  154. }else
  155. {
  156. eff_R.resize(dim,dim*n);
  157. for(int v = 0;v<n;v++)
  158. {
  159. eff_R.block(0,dim*v,dim,dim) =
  160. R.block(0,dim*data.G(v),dim,dim);
  161. }
  162. }
  163. VectorXd Rcol;
  164. columnize(eff_R,n,2,Rcol);
  165. VectorXd Bcol = -data.K * Rcol;
  166. for(int c = 0;c<dim;c++)
  167. {
  168. VectorXd Uc,Bc,bcc,Beq;
  169. Bc = Bcol.block(c*n,0,n,1);
  170. bcc = bc.col(c);
  171. min_quad_with_fixed_solve(
  172. data.solver_data,
  173. Bc,bcc,Beq,
  174. Uc);
  175. U.col(c) = Uc;
  176. }
  177. iter++;
  178. }
  179. return true;
  180. }
  181. #ifndef IGL_HEADER_ONLY
  182. 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&);
  183. 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> >&);
  184. #endif