biharmonic_coordinates.cpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151
  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 "massmatrix.h"
  11. #include "min_quad_with_fixed.h"
  12. #include "normal_derivative.h"
  13. #include "on_boundary.h"
  14. #include <Eigen/Sparse>
  15. template <
  16. typename DerivedV,
  17. typename DerivedT,
  18. typename SType,
  19. typename DerivedW>
  20. IGL_INLINE bool igl::biharmonic_coordinates(
  21. const Eigen::PlainObjectBase<DerivedV> & V,
  22. const Eigen::PlainObjectBase<DerivedT> & T,
  23. const std::vector<std::vector<SType> > & S,
  24. Eigen::PlainObjectBase<DerivedW> & W)
  25. {
  26. return biharmonic_coordinates(V,T,S,2,W);
  27. }
  28. template <
  29. typename DerivedV,
  30. typename DerivedT,
  31. typename SType,
  32. typename DerivedW>
  33. IGL_INLINE bool igl::biharmonic_coordinates(
  34. const Eigen::PlainObjectBase<DerivedV> & V,
  35. const Eigen::PlainObjectBase<DerivedT> & T,
  36. const std::vector<std::vector<SType> > & S,
  37. const int k,
  38. Eigen::PlainObjectBase<DerivedW> & W)
  39. {
  40. using namespace Eigen;
  41. using namespace std;
  42. // This is not the most efficient way to build A, but follows "Linear
  43. // Subspace Design for Real-Time Shape Deformation" [Wang et al. 2015].
  44. SparseMatrix<double> A;
  45. {
  46. SparseMatrix<double> N,Z,L,K,M;
  47. normal_derivative(V,T,N);
  48. Array<bool,Dynamic,1> I;
  49. Array<bool,Dynamic,Dynamic> C;
  50. on_boundary(T,I,C);
  51. {
  52. std::vector<Triplet<double> >ZIJV;
  53. for(int t =0;t<T.rows();t++)
  54. {
  55. for(int f =0;f<T.cols();f++)
  56. {
  57. if(C(t,f))
  58. {
  59. const int i = t+f*T.rows();
  60. for(int c = 1;c<T.cols();c++)
  61. {
  62. ZIJV.emplace_back(T(t,(f+c)%T.cols()),i,1);
  63. }
  64. }
  65. }
  66. }
  67. Z.resize(V.rows(),N.rows());
  68. Z.setFromTriplets(ZIJV.begin(),ZIJV.end());
  69. N = (Z*N).eval();
  70. }
  71. cotmatrix(V,T,L);
  72. K = N+L;
  73. massmatrix(V,T,MASSMATRIX_TYPE_DEFAULT,M);
  74. // normalize
  75. M /= ((VectorXd)M.diagonal()).array().abs().maxCoeff();
  76. DiagonalMatrix<double,Dynamic> Minv =
  77. ((VectorXd)M.diagonal().array().inverse()).asDiagonal();
  78. switch(k)
  79. {
  80. default:
  81. assert(false && "unsupported");
  82. case 2:
  83. // For C1 smoothness in 2D, one should use bi-harmonic
  84. A = K.transpose() * (Minv * K);
  85. break;
  86. case 3:
  87. // For C1 smoothness in 3D, one should use tri-harmonic
  88. A = K.transpose() * (Minv * (-L * (Minv * K)));
  89. break;
  90. }
  91. }
  92. // Vertices in point handles
  93. const size_t mp =
  94. count_if(S.begin(),S.end(),[](const vector<int> & h){return h.size()==1;});
  95. // number of region handles
  96. const size_t r = S.size()-mp;
  97. // Vertices in region handles
  98. size_t mr = 0;
  99. for(const auto & h : S)
  100. {
  101. if(h.size() > 1)
  102. {
  103. mr += h.size();
  104. }
  105. }
  106. const size_t dim = T.cols()-1;
  107. // Might as well be dense... I think...
  108. MatrixXd J = MatrixXd::Zero(mp+mr,mp+r*(dim+1));
  109. VectorXi b(mp+mr);
  110. MatrixXd H(mp+r*(dim+1),dim);
  111. {
  112. int v = 0;
  113. int c = 0;
  114. for(int h = 0;h<S.size();h++)
  115. {
  116. if(S[h].size()==1)
  117. {
  118. H.row(c) = V.block(S[h][0],0,1,dim);
  119. J(v,c++) = 1;
  120. b(v) = S[h][0];
  121. v++;
  122. }else
  123. {
  124. assert(S[h].size() >= dim+1);
  125. for(int p = 0;p<S[h].size();p++)
  126. {
  127. for(int d = 0;d<dim;d++)
  128. {
  129. J(v,c+d) = V(S[h][p],d);
  130. }
  131. J(v,c+dim) = 1;
  132. b(v) = S[h][p];
  133. v++;
  134. }
  135. H.block(c,0,dim+1,dim).setIdentity();
  136. c+=dim+1;
  137. }
  138. }
  139. }
  140. // minimize ½ W' A W'
  141. // subject to W(b,:) = J
  142. return min_quad_with_fixed(
  143. A,VectorXd::Zero(A.rows()).eval(),b,J,SparseMatrix<double>(),VectorXd(),true,W);
  144. }
  145. #ifdef IGL_STATIC_LIBRARY
  146. // Explicit template specialization
  147. 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> >&);
  148. #endif