biharmonic_coordinates.cpp 6.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 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 "biharmonic_coordinates.h"
  9. #include "cotmatrix.h"
  10. #include "sum.h"
  11. #include "massmatrix.h"
  12. #include "min_quad_with_fixed.h"
  13. #include "crouzeix_raviart_massmatrix.h"
  14. #include "crouzeix_raviart_cotmatrix.h"
  15. #include "normal_derivative.h"
  16. #include "on_boundary.h"
  17. #include <Eigen/Sparse>
  18. template <
  19. typename DerivedV,
  20. typename DerivedT,
  21. typename SType,
  22. typename DerivedW>
  23. IGL_INLINE bool igl::biharmonic_coordinates(
  24. const Eigen::MatrixBase<DerivedV> & V,
  25. const Eigen::MatrixBase<DerivedT> & T,
  26. const std::vector<std::vector<SType> > & S,
  27. Eigen::PlainObjectBase<DerivedW> & W)
  28. {
  29. return biharmonic_coordinates(V,T,S,2,W);
  30. }
  31. template <
  32. typename DerivedV,
  33. typename DerivedT,
  34. typename SType,
  35. typename DerivedW>
  36. IGL_INLINE bool igl::biharmonic_coordinates(
  37. const Eigen::MatrixBase<DerivedV> & V,
  38. const Eigen::MatrixBase<DerivedT> & T,
  39. const std::vector<std::vector<SType> > & S,
  40. const int k,
  41. Eigen::PlainObjectBase<DerivedW> & W)
  42. {
  43. using namespace Eigen;
  44. using namespace std;
  45. typedef typename DerivedV::Scalar Scalar;
  46. typedef typename DerivedT::Scalar Integer;
  47. // This is not the most efficient way to build A, but follows "Linear
  48. // Subspace Design for Real-Time Shape Deformation" [Wang et al. 2015].
  49. SparseMatrix<Scalar> A;
  50. {
  51. DiagonalMatrix<Scalar, Dynamic> Minv;
  52. SparseMatrix<Scalar> L, K;
  53. Array<bool,Dynamic,Dynamic> C;
  54. {
  55. Array<bool,Dynamic,1> I;
  56. on_boundary(T,I,C);
  57. }
  58. #ifdef false
  59. // Version described in paper is "wrong"
  60. // http://www.cs.toronto.edu/~jacobson/images/error-in-linear-subspace-design-for-real-time-shape-deformation-2017-wang-et-al.pdf
  61. SparseMatrix<Scalar> N, Z, M;
  62. normal_derivative(V,T,N);
  63. {
  64. std::vector<Triplet<Scalar>> ZIJV;
  65. for(int t =0;t<T.rows();t++)
  66. {
  67. for(int f =0;f<T.cols();f++)
  68. {
  69. if(C(t,f))
  70. {
  71. const int i = t+f*T.rows();
  72. for(int c = 1;c<T.cols();c++)
  73. {
  74. ZIJV.emplace_back(T(t,(f+c)%T.cols()),i,1);
  75. }
  76. }
  77. }
  78. }
  79. Z.resize(V.rows(),N.rows());
  80. Z.setFromTriplets(ZIJV.begin(),ZIJV.end());
  81. N = (Z*N).eval();
  82. }
  83. cotmatrix(V,T,L);
  84. K = N+L;
  85. massmatrix(V,T,MASSMATRIX_TYPE_DEFAULT,M);
  86. // normalize
  87. M /= ((Matrix<Scalar, Dynamic, 1>)M.diagonal()).array().abs().maxCoeff();
  88. Minv =
  89. ((Matrix<Scalar, Dynamic, 1>)M.diagonal().array().inverse()).asDiagonal();
  90. #else
  91. Eigen::SparseMatrix<Scalar> M;
  92. Eigen::Matrix<Integer, Dynamic, Dynamic> E;
  93. Eigen::Matrix<Integer, Dynamic, 1> EMAP;
  94. crouzeix_raviart_massmatrix(V,T,M,E,EMAP);
  95. crouzeix_raviart_cotmatrix(V,T,E,EMAP,L);
  96. // Ad #E by #V facet-vertex incidence matrix
  97. Eigen::SparseMatrix<Scalar> Ad(E.rows(),V.rows());
  98. {
  99. std::vector<Eigen::Triplet<Scalar>> AIJV(E.size());
  100. for(int e = 0;e<E.rows();e++)
  101. {
  102. for(int c = 0;c<E.cols();c++)
  103. {
  104. AIJV[e + c * E.rows()] = Eigen::Triplet<Scalar>(e, E(e, c), 1);
  105. }
  106. }
  107. Ad.setFromTriplets(AIJV.begin(),AIJV.end());
  108. }
  109. // Degrees
  110. Eigen::Matrix<Scalar, Dynamic, 1> De;
  111. sum(Ad,2,De);
  112. Eigen::DiagonalMatrix<Scalar,Eigen::Dynamic> De_diag =
  113. De.array().inverse().matrix().asDiagonal();
  114. K = L*(De_diag*Ad);
  115. // normalize
  116. M /= ((Matrix<Scalar, Dynamic, 1>)M.diagonal()).array().abs().maxCoeff();
  117. Minv = ((Matrix<Scalar, Dynamic, 1>)M.diagonal().array().inverse()).asDiagonal();
  118. // kill boundary edges
  119. for(int f = 0;f<T.rows();f++)
  120. {
  121. for(int c = 0;c<T.cols();c++)
  122. {
  123. if(C(f,c))
  124. {
  125. const int e = EMAP(f+T.rows()*c);
  126. Minv.diagonal()(e) = 0;
  127. }
  128. }
  129. }
  130. #endif
  131. switch(k)
  132. {
  133. default:
  134. assert(false && "unsupported");
  135. case 2:
  136. // For C1 smoothness in 2D, one should use bi-harmonic
  137. A = K.transpose() * (Minv * K);
  138. break;
  139. case 3:
  140. // For C1 smoothness in 3D, one should use tri-harmonic
  141. A = K.transpose() * (Minv * (-L * (Minv * K)));
  142. break;
  143. }
  144. }
  145. // Vertices in point handles
  146. const size_t mp =
  147. count_if(S.begin(),S.end(),[](const vector<int> & h){return h.size()==1;});
  148. // number of region handles
  149. const size_t r = S.size()-mp;
  150. // Vertices in region handles
  151. size_t mr = 0;
  152. for(const auto & h : S)
  153. {
  154. if(h.size() > 1)
  155. {
  156. mr += h.size();
  157. }
  158. }
  159. const size_t dim = T.cols()-1;
  160. // Might as well be dense... I think...
  161. Matrix<Scalar, Dynamic, Dynamic> J = Matrix<Scalar, Dynamic, Dynamic>::Zero(mp+mr,mp+r*(dim+1));
  162. Matrix<Integer, Dynamic, 1> b(mp+mr);
  163. Matrix<Scalar, Dynamic, Dynamic> H(mp+r*(dim+1),dim);
  164. {
  165. int v = 0;
  166. int c = 0;
  167. for(int h = 0;h<S.size();h++)
  168. {
  169. if(S[h].size()==1)
  170. {
  171. H.row(c) = V.block(S[h][0],0,1,dim);
  172. J(v,c++) = 1;
  173. b(v) = S[h][0];
  174. v++;
  175. }else
  176. {
  177. assert(S[h].size() >= dim+1);
  178. for(int p = 0;p<S[h].size();p++)
  179. {
  180. for(int d = 0;d<dim;d++)
  181. {
  182. J(v,c+d) = V(S[h][p],d);
  183. }
  184. J(v,c+dim) = 1;
  185. b(v) = S[h][p];
  186. v++;
  187. }
  188. H.block(c,0,dim+1,dim).setIdentity();
  189. c+=dim+1;
  190. }
  191. }
  192. }
  193. // minimize ½ W' A W'
  194. // subject to W(b,:) = J
  195. return min_quad_with_fixed(
  196. A,Matrix<Scalar, Dynamic, 1>::Zero(A.rows()).eval(),b,J,SparseMatrix<Scalar>(),Matrix<Scalar, Dynamic, 1>(),true,W);
  197. }
  198. #ifdef IGL_STATIC_LIBRARY
  199. // Explicit template instantiation
  200. template bool igl::biharmonic_coordinates<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, int, 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&, std::vector<std::vector<int, std::allocator<int> >, std::allocator<std::vector<int, std::allocator<int> > > > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  201. #endif