arap_linear_block.cpp 6.5 KB

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