biharmonic_coordinates.cpp 5.9 KB

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