arap_linear_block.cpp 6.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257
  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_linear_block.h"
  9. #include "verbose.h"
  10. #include "cotmatrix_entries.h"
  11. #include <Eigen/Dense>
  12. template <typename MatV, typename MatF, typename Scalar>
  13. IGL_INLINE void igl::arap_linear_block(
  14. const MatV & V,
  15. const MatF & F,
  16. const int d,
  17. const igl::ARAPEnergyType energy,
  18. Eigen::SparseMatrix<Scalar> & Kd)
  19. {
  20. using namespace igl;
  21. switch(energy)
  22. {
  23. case ARAP_ENERGY_TYPE_SPOKES:
  24. return igl::arap_linear_block_spokes(V,F,d,Kd);
  25. break;
  26. case ARAP_ENERGY_TYPE_SPOKES_AND_RIMS:
  27. return igl::arap_linear_block_spokes_and_rims(V,F,d,Kd);
  28. break;
  29. case ARAP_ENERGY_TYPE_ELEMENTS:
  30. return igl::arap_linear_block_elements(V,F,d,Kd);
  31. break;
  32. default:
  33. verbose("Unsupported energy type: %d\n",energy);
  34. assert(false);
  35. }
  36. }
  37. template <typename MatV, typename MatF, typename Scalar>
  38. IGL_INLINE void igl::arap_linear_block_spokes(
  39. const MatV & V,
  40. const MatF & F,
  41. const int d,
  42. Eigen::SparseMatrix<Scalar> & Kd)
  43. {
  44. using namespace igl;
  45. using namespace std;
  46. using namespace Eigen;
  47. // simplex size (3: triangles, 4: tetrahedra)
  48. int simplex_size = F.cols();
  49. // Number of elements
  50. int m = F.rows();
  51. // Temporary output
  52. Matrix<int,Dynamic,2> edges;
  53. Kd.resize(V.rows(), V.rows());
  54. vector<Triplet<Scalar> > Kd_IJV;
  55. if(simplex_size == 3)
  56. {
  57. // triangles
  58. Kd.reserve(7*V.rows());
  59. Kd_IJV.reserve(7*V.rows());
  60. edges.resize(3,2);
  61. edges <<
  62. 1,2,
  63. 2,0,
  64. 0,1;
  65. }else if(simplex_size == 4)
  66. {
  67. // tets
  68. Kd.reserve(17*V.rows());
  69. Kd_IJV.reserve(17*V.rows());
  70. edges.resize(6,2);
  71. edges <<
  72. 1,2,
  73. 2,0,
  74. 0,1,
  75. 3,0,
  76. 3,1,
  77. 3,2;
  78. }
  79. // gather cotangent weights
  80. Matrix<Scalar,Dynamic,Dynamic> C;
  81. cotmatrix_entries(V,F,C);
  82. // should have weights for each edge
  83. assert(C.cols() == edges.rows());
  84. // loop over elements
  85. for(int i = 0;i<m;i++)
  86. {
  87. // loop over edges of element
  88. for(int e = 0;e<edges.rows();e++)
  89. {
  90. int source = F(i,edges(e,0));
  91. int dest = F(i,edges(e,1));
  92. double v = 0.5*C(i,e)*(V(source,d)-V(dest,d));
  93. Kd_IJV.push_back(Triplet<Scalar>(source,dest,v));
  94. Kd_IJV.push_back(Triplet<Scalar>(dest,source,-v));
  95. Kd_IJV.push_back(Triplet<Scalar>(source,source,v));
  96. Kd_IJV.push_back(Triplet<Scalar>(dest,dest,-v));
  97. }
  98. }
  99. Kd.setFromTriplets(Kd_IJV.begin(),Kd_IJV.end());
  100. Kd.makeCompressed();
  101. }
  102. template <typename MatV, typename MatF, typename Scalar>
  103. IGL_INLINE void igl::arap_linear_block_spokes_and_rims(
  104. const MatV & V,
  105. const MatF & F,
  106. const int d,
  107. Eigen::SparseMatrix<Scalar> & Kd)
  108. {
  109. using namespace igl;
  110. using namespace std;
  111. using namespace Eigen;
  112. // simplex size (3: triangles, 4: tetrahedra)
  113. int simplex_size = F.cols();
  114. // Number of elements
  115. int m = F.rows();
  116. // Temporary output
  117. Kd.resize(V.rows(), V.rows());
  118. vector<Triplet<Scalar> > Kd_IJV;
  119. Matrix<int,Dynamic,2> edges;
  120. if(simplex_size == 3)
  121. {
  122. // triangles
  123. Kd.reserve(7*V.rows());
  124. Kd_IJV.reserve(7*V.rows());
  125. edges.resize(3,2);
  126. edges <<
  127. 1,2,
  128. 2,0,
  129. 0,1;
  130. }else if(simplex_size == 4)
  131. {
  132. // tets
  133. Kd.reserve(17*V.rows());
  134. Kd_IJV.reserve(17*V.rows());
  135. edges.resize(6,2);
  136. edges <<
  137. 1,2,
  138. 2,0,
  139. 0,1,
  140. 3,0,
  141. 3,1,
  142. 3,2;
  143. // Not implemented yet for tets
  144. assert(false);
  145. }
  146. // gather cotangent weights
  147. Matrix<Scalar,Dynamic,Dynamic> C;
  148. cotmatrix_entries(V,F,C);
  149. // should have weights for each edge
  150. assert(C.cols() == edges.rows());
  151. // loop over elements
  152. for(int i = 0;i<m;i++)
  153. {
  154. // loop over edges of element
  155. for(int e = 0;e<edges.rows();e++)
  156. {
  157. int source = F(i,edges(e,0));
  158. int dest = F(i,edges(e,1));
  159. double v = C(i,e)*(V(source,d)-V(dest,d))/3.0;
  160. // loop over edges again
  161. for(int f = 0;f<edges.rows();f++)
  162. {
  163. int Rs = F(i,edges(f,0));
  164. int Rd = F(i,edges(f,1));
  165. if(Rs == source && Rd == dest)
  166. {
  167. Kd_IJV.push_back(Triplet<Scalar>(Rs,Rd,v));
  168. Kd_IJV.push_back(Triplet<Scalar>(Rd,Rs,-v));
  169. }else if(Rd == source)
  170. {
  171. Kd_IJV.push_back(Triplet<Scalar>(Rd,Rs,v));
  172. }else if(Rs == dest)
  173. {
  174. Kd_IJV.push_back(Triplet<Scalar>(Rs,Rd,-v));
  175. }
  176. }
  177. Kd_IJV.push_back(Triplet<Scalar>(source,source,v));
  178. Kd_IJV.push_back(Triplet<Scalar>(dest,dest,-v));
  179. }
  180. }
  181. Kd.setFromTriplets(Kd_IJV.begin(),Kd_IJV.end());
  182. Kd.makeCompressed();
  183. }
  184. template <typename MatV, typename MatF, typename Scalar>
  185. IGL_INLINE void igl::arap_linear_block_elements(
  186. const MatV & V,
  187. const MatF & F,
  188. const int d,
  189. Eigen::SparseMatrix<Scalar> & Kd)
  190. {
  191. using namespace igl;
  192. using namespace std;
  193. using namespace Eigen;
  194. // simplex size (3: triangles, 4: tetrahedra)
  195. int simplex_size = F.cols();
  196. // Number of elements
  197. int m = F.rows();
  198. // Temporary output
  199. Kd.resize(V.rows(), F.rows());
  200. vector<Triplet<Scalar> > Kd_IJV;
  201. Matrix<int,Dynamic,2> edges;
  202. if(simplex_size == 3)
  203. {
  204. // triangles
  205. Kd.reserve(7*V.rows());
  206. Kd_IJV.reserve(7*V.rows());
  207. edges.resize(3,2);
  208. edges <<
  209. 1,2,
  210. 2,0,
  211. 0,1;
  212. }else if(simplex_size == 4)
  213. {
  214. // tets
  215. Kd.reserve(17*V.rows());
  216. Kd_IJV.reserve(17*V.rows());
  217. edges.resize(6,2);
  218. edges <<
  219. 1,2,
  220. 2,0,
  221. 0,1,
  222. 3,0,
  223. 3,1,
  224. 3,2;
  225. }
  226. // gather cotangent weights
  227. Matrix<Scalar,Dynamic,Dynamic> C;
  228. cotmatrix_entries(V,F,C);
  229. // should have weights for each edge
  230. assert(C.cols() == edges.rows());
  231. // loop over elements
  232. for(int i = 0;i<m;i++)
  233. {
  234. // loop over edges of element
  235. for(int e = 0;e<edges.rows();e++)
  236. {
  237. int source = F(i,edges(e,0));
  238. int dest = F(i,edges(e,1));
  239. double v = C(i,e)*(V(source,d)-V(dest,d));
  240. Kd_IJV.push_back(Triplet<Scalar>(source,i,v));
  241. Kd_IJV.push_back(Triplet<Scalar>(dest,i,-v));
  242. }
  243. }
  244. Kd.setFromTriplets(Kd_IJV.begin(),Kd_IJV.end());
  245. Kd.makeCompressed();
  246. }
  247. #ifdef IGL_STATIC_LIBRARY
  248. // Explicit template specialization
  249. template IGL_INLINE void igl::arap_linear_block<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, int, igl::ARAPEnergyType, Eigen::SparseMatrix<double, 0, int>&);
  250. #endif