slim.cpp 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Michael Rabinovich
  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 "slim.h"
  9. #include "boundary_loop.h"
  10. #include "cotmatrix.h"
  11. #include "edge_lengths.h"
  12. #include "grad.h"
  13. #include "local_basis.h"
  14. #include "repdiag.h"
  15. #include "vector_area_matrix.h"
  16. #include "arap.h"
  17. #include "cat.h"
  18. #include "doublearea.h"
  19. #include "grad.h"
  20. #include "local_basis.h"
  21. #include "per_face_normals.h"
  22. #include "slice_into.h"
  23. #include "volume.h"
  24. #include "polar_svd.h"
  25. #include "flip_avoiding_line_search.h"
  26. #include <iostream>
  27. #include <map>
  28. #include <set>
  29. #include <vector>
  30. #include <Eigen/IterativeLinearSolvers>
  31. #include <Eigen/SparseCholesky>
  32. #include <Eigen/IterativeLinearSolvers>
  33. namespace igl
  34. {
  35. namespace slim
  36. {
  37. // Definitions of internal functions
  38. IGL_INLINE void compute_surface_gradient_matrix(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F,
  39. const Eigen::MatrixXd &F1, const Eigen::MatrixXd &F2,
  40. Eigen::SparseMatrix<double> &D1, Eigen::SparseMatrix<double> &D2);
  41. IGL_INLINE void buildA(igl::SLIMData& s, Eigen::SparseMatrix<double> &A);
  42. IGL_INLINE void buildRhs(igl::SLIMData& s, const Eigen::SparseMatrix<double> &At);
  43. IGL_INLINE void add_soft_constraints(igl::SLIMData& s, Eigen::SparseMatrix<double> &L);
  44. IGL_INLINE double compute_energy(igl::SLIMData& s, Eigen::MatrixXd &V_new);
  45. IGL_INLINE double compute_soft_const_energy(igl::SLIMData& s,
  46. const Eigen::MatrixXd &V,
  47. const Eigen::MatrixXi &F,
  48. Eigen::MatrixXd &V_o);
  49. IGL_INLINE double compute_energy_with_jacobians(igl::SLIMData& s,
  50. const Eigen::MatrixXd &V,
  51. const Eigen::MatrixXi &F, const Eigen::MatrixXd &Ji,
  52. Eigen::MatrixXd &uv, Eigen::VectorXd &areas);
  53. IGL_INLINE void solve_weighted_arap(igl::SLIMData& s,
  54. const Eigen::MatrixXd &V,
  55. const Eigen::MatrixXi &F,
  56. Eigen::MatrixXd &uv,
  57. Eigen::VectorXi &soft_b_p,
  58. Eigen::MatrixXd &soft_bc_p);
  59. IGL_INLINE void update_weights_and_closest_rotations( igl::SLIMData& s,
  60. const Eigen::MatrixXd &V,
  61. const Eigen::MatrixXi &F,
  62. Eigen::MatrixXd &uv);
  63. IGL_INLINE void compute_jacobians(igl::SLIMData& s, const Eigen::MatrixXd &uv);
  64. IGL_INLINE void build_linear_system(igl::SLIMData& s, Eigen::SparseMatrix<double> &L);
  65. IGL_INLINE void pre_calc(igl::SLIMData& s);
  66. // Implementation
  67. IGL_INLINE void compute_surface_gradient_matrix(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F,
  68. const Eigen::MatrixXd &F1, const Eigen::MatrixXd &F2,
  69. Eigen::SparseMatrix<double> &D1, Eigen::SparseMatrix<double> &D2)
  70. {
  71. Eigen::SparseMatrix<double> G;
  72. igl::grad(V, F, G);
  73. Eigen::SparseMatrix<double> Dx = G.block(0, 0, F.rows(), V.rows());
  74. Eigen::SparseMatrix<double> Dy = G.block(F.rows(), 0, F.rows(), V.rows());
  75. Eigen::SparseMatrix<double> Dz = G.block(2 * F.rows(), 0, F.rows(), V.rows());
  76. D1 = F1.col(0).asDiagonal() * Dx + F1.col(1).asDiagonal() * Dy + F1.col(2).asDiagonal() * Dz;
  77. D2 = F2.col(0).asDiagonal() * Dx + F2.col(1).asDiagonal() * Dy + F2.col(2).asDiagonal() * Dz;
  78. }
  79. IGL_INLINE void compute_jacobians(igl::SLIMData& s, const Eigen::MatrixXd &uv)
  80. {
  81. if (s.F.cols() == 3)
  82. {
  83. // Ji=[D1*u,D2*u,D1*v,D2*v];
  84. s.Ji.col(0) = s.Dx * uv.col(0);
  85. s.Ji.col(1) = s.Dy * uv.col(0);
  86. s.Ji.col(2) = s.Dx * uv.col(1);
  87. s.Ji.col(3) = s.Dy * uv.col(1);
  88. }
  89. else /*tet mesh*/{
  90. // Ji=[D1*u,D2*u,D3*u, D1*v,D2*v, D3*v, D1*w,D2*w,D3*w];
  91. s.Ji.col(0) = s.Dx * uv.col(0);
  92. s.Ji.col(1) = s.Dy * uv.col(0);
  93. s.Ji.col(2) = s.Dz * uv.col(0);
  94. s.Ji.col(3) = s.Dx * uv.col(1);
  95. s.Ji.col(4) = s.Dy * uv.col(1);
  96. s.Ji.col(5) = s.Dz * uv.col(1);
  97. s.Ji.col(6) = s.Dx * uv.col(2);
  98. s.Ji.col(7) = s.Dy * uv.col(2);
  99. s.Ji.col(8) = s.Dz * uv.col(2);
  100. }
  101. }
  102. IGL_INLINE void update_weights_and_closest_rotations(igl::SLIMData& s,
  103. const Eigen::MatrixXd &V,
  104. const Eigen::MatrixXi &F,
  105. Eigen::MatrixXd &uv)
  106. {
  107. compute_jacobians(s, uv);
  108. const double eps = 1e-8;
  109. double exp_f = s.exp_factor;
  110. if (s.dim == 2)
  111. {
  112. for (int i = 0; i < s.Ji.rows(); ++i)
  113. {
  114. typedef Eigen::Matrix<double, 2, 2> Mat2;
  115. typedef Eigen::Matrix<double, 2, 1> Vec2;
  116. Mat2 ji, ri, ti, ui, vi;
  117. Vec2 sing;
  118. Vec2 closest_sing_vec;
  119. Mat2 mat_W;
  120. Vec2 m_sing_new;
  121. double s1, s2;
  122. ji(0, 0) = s.Ji(i, 0);
  123. ji(0, 1) = s.Ji(i, 1);
  124. ji(1, 0) = s.Ji(i, 2);
  125. ji(1, 1) = s.Ji(i, 3);
  126. igl::polar_svd(ji, ri, ti, ui, sing, vi);
  127. s1 = sing(0);
  128. s2 = sing(1);
  129. // Update Weights according to energy
  130. switch (s.slim_energy)
  131. {
  132. case igl::SLIMData::ARAP:
  133. {
  134. m_sing_new << 1, 1;
  135. break;
  136. }
  137. case igl::SLIMData::SYMMETRIC_DIRICHLET:
  138. {
  139. double s1_g = 2 * (s1 - pow(s1, -3));
  140. double s2_g = 2 * (s2 - pow(s2, -3));
  141. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1)));
  142. break;
  143. }
  144. case igl::SLIMData::LOG_ARAP:
  145. {
  146. double s1_g = 2 * (log(s1) / s1);
  147. double s2_g = 2 * (log(s2) / s2);
  148. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1)));
  149. break;
  150. }
  151. case igl::SLIMData::CONFORMAL:
  152. {
  153. double s1_g = 1 / (2 * s2) - s2 / (2 * pow(s1, 2));
  154. double s2_g = 1 / (2 * s1) - s1 / (2 * pow(s2, 2));
  155. double geo_avg = sqrt(s1 * s2);
  156. double s1_min = geo_avg;
  157. double s2_min = geo_avg;
  158. m_sing_new << sqrt(s1_g / (2 * (s1 - s1_min))), sqrt(s2_g / (2 * (s2 - s2_min)));
  159. // change local step
  160. closest_sing_vec << s1_min, s2_min;
  161. ri = ui * closest_sing_vec.asDiagonal() * vi.transpose();
  162. break;
  163. }
  164. case igl::SLIMData::EXP_CONFORMAL:
  165. {
  166. double s1_g = 2 * (s1 - pow(s1, -3));
  167. double s2_g = 2 * (s2 - pow(s2, -3));
  168. double geo_avg = sqrt(s1 * s2);
  169. double s1_min = geo_avg;
  170. double s2_min = geo_avg;
  171. double in_exp = exp_f * ((pow(s1, 2) + pow(s2, 2)) / (2 * s1 * s2));
  172. double exp_thing = exp(in_exp);
  173. s1_g *= exp_thing * exp_f;
  174. s2_g *= exp_thing * exp_f;
  175. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1)));
  176. break;
  177. }
  178. case igl::SLIMData::EXP_SYMMETRIC_DIRICHLET:
  179. {
  180. double s1_g = 2 * (s1 - pow(s1, -3));
  181. double s2_g = 2 * (s2 - pow(s2, -3));
  182. double in_exp = exp_f * (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2));
  183. double exp_thing = exp(in_exp);
  184. s1_g *= exp_thing * exp_f;
  185. s2_g *= exp_thing * exp_f;
  186. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1)));
  187. break;
  188. }
  189. }
  190. if (std::abs(s1 - 1) < eps) m_sing_new(0) = 1;
  191. if (std::abs(s2 - 1) < eps) m_sing_new(1) = 1;
  192. mat_W = ui * m_sing_new.asDiagonal() * ui.transpose();
  193. s.W_11(i) = mat_W(0, 0);
  194. s.W_12(i) = mat_W(0, 1);
  195. s.W_21(i) = mat_W(1, 0);
  196. s.W_22(i) = mat_W(1, 1);
  197. // 2) Update local step (doesn't have to be a rotation, for instance in case of conformal energy)
  198. s.Ri(i, 0) = ri(0, 0);
  199. s.Ri(i, 1) = ri(1, 0);
  200. s.Ri(i, 2) = ri(0, 1);
  201. s.Ri(i, 3) = ri(1, 1);
  202. }
  203. }
  204. else
  205. {
  206. typedef Eigen::Matrix<double, 3, 1> Vec3;
  207. typedef Eigen::Matrix<double, 3, 3> Mat3;
  208. Mat3 ji;
  209. Vec3 m_sing_new;
  210. Vec3 closest_sing_vec;
  211. const double sqrt_2 = sqrt(2);
  212. for (int i = 0; i < s.Ji.rows(); ++i)
  213. {
  214. ji(0, 0) = s.Ji(i, 0);
  215. ji(0, 1) = s.Ji(i, 1);
  216. ji(0, 2) = s.Ji(i, 2);
  217. ji(1, 0) = s.Ji(i, 3);
  218. ji(1, 1) = s.Ji(i, 4);
  219. ji(1, 2) = s.Ji(i, 5);
  220. ji(2, 0) = s.Ji(i, 6);
  221. ji(2, 1) = s.Ji(i, 7);
  222. ji(2, 2) = s.Ji(i, 8);
  223. Mat3 ri, ti, ui, vi;
  224. Vec3 sing;
  225. igl::polar_svd(ji, ri, ti, ui, sing, vi);
  226. double s1 = sing(0);
  227. double s2 = sing(1);
  228. double s3 = sing(2);
  229. // 1) Update Weights
  230. switch (s.slim_energy)
  231. {
  232. case igl::SLIMData::ARAP:
  233. {
  234. m_sing_new << 1, 1, 1;
  235. break;
  236. }
  237. case igl::SLIMData::LOG_ARAP:
  238. {
  239. double s1_g = 2 * (log(s1) / s1);
  240. double s2_g = 2 * (log(s2) / s2);
  241. double s3_g = 2 * (log(s3) / s3);
  242. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1))), sqrt(s3_g / (2 * (s3 - 1)));
  243. break;
  244. }
  245. case igl::SLIMData::SYMMETRIC_DIRICHLET:
  246. {
  247. double s1_g = 2 * (s1 - pow(s1, -3));
  248. double s2_g = 2 * (s2 - pow(s2, -3));
  249. double s3_g = 2 * (s3 - pow(s3, -3));
  250. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1))), sqrt(s3_g / (2 * (s3 - 1)));
  251. break;
  252. }
  253. case igl::SLIMData::EXP_SYMMETRIC_DIRICHLET:
  254. {
  255. double s1_g = 2 * (s1 - pow(s1, -3));
  256. double s2_g = 2 * (s2 - pow(s2, -3));
  257. double s3_g = 2 * (s3 - pow(s3, -3));
  258. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1))), sqrt(s3_g / (2 * (s3 - 1)));
  259. double in_exp = exp_f * (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2) + pow(s3, 2) + pow(s3, -2));
  260. double exp_thing = exp(in_exp);
  261. s1_g *= exp_thing * exp_f;
  262. s2_g *= exp_thing * exp_f;
  263. s3_g *= exp_thing * exp_f;
  264. m_sing_new << sqrt(s1_g / (2 * (s1 - 1))), sqrt(s2_g / (2 * (s2 - 1))), sqrt(s3_g / (2 * (s3 - 1)));
  265. break;
  266. }
  267. case igl::SLIMData::CONFORMAL:
  268. {
  269. double common_div = 9 * (pow(s1 * s2 * s3, 5. / 3.));
  270. double s1_g = (-2 * s2 * s3 * (pow(s2, 2) + pow(s3, 2) - 2 * pow(s1, 2))) / common_div;
  271. double s2_g = (-2 * s1 * s3 * (pow(s1, 2) + pow(s3, 2) - 2 * pow(s2, 2))) / common_div;
  272. double s3_g = (-2 * s1 * s2 * (pow(s1, 2) + pow(s2, 2) - 2 * pow(s3, 2))) / common_div;
  273. double closest_s = sqrt(pow(s1, 2) + pow(s3, 2)) / sqrt_2;
  274. double s1_min = closest_s;
  275. double s2_min = closest_s;
  276. double s3_min = closest_s;
  277. m_sing_new << sqrt(s1_g / (2 * (s1 - s1_min))), sqrt(s2_g / (2 * (s2 - s2_min))), sqrt(
  278. s3_g / (2 * (s3 - s3_min)));
  279. // change local step
  280. closest_sing_vec << s1_min, s2_min, s3_min;
  281. ri = ui * closest_sing_vec.asDiagonal() * vi.transpose();
  282. break;
  283. }
  284. case igl::SLIMData::EXP_CONFORMAL:
  285. {
  286. // E_conf = (s1^2 + s2^2 + s3^2)/(3*(s1*s2*s3)^(2/3) )
  287. // dE_conf/ds1 = (-2*(s2*s3)*(s2^2+s3^2 -2*s1^2) ) / (9*(s1*s2*s3)^(5/3))
  288. // Argmin E_conf(s1): s1 = sqrt(s1^2+s2^2)/sqrt(2)
  289. double common_div = 9 * (pow(s1 * s2 * s3, 5. / 3.));
  290. double s1_g = (-2 * s2 * s3 * (pow(s2, 2) + pow(s3, 2) - 2 * pow(s1, 2))) / common_div;
  291. double s2_g = (-2 * s1 * s3 * (pow(s1, 2) + pow(s3, 2) - 2 * pow(s2, 2))) / common_div;
  292. double s3_g = (-2 * s1 * s2 * (pow(s1, 2) + pow(s2, 2) - 2 * pow(s3, 2))) / common_div;
  293. double in_exp = exp_f * ((pow(s1, 2) + pow(s2, 2) + pow(s3, 2)) / (3 * pow((s1 * s2 * s3), 2. / 3)));;
  294. double exp_thing = exp(in_exp);
  295. double closest_s = sqrt(pow(s1, 2) + pow(s3, 2)) / sqrt_2;
  296. double s1_min = closest_s;
  297. double s2_min = closest_s;
  298. double s3_min = closest_s;
  299. s1_g *= exp_thing * exp_f;
  300. s2_g *= exp_thing * exp_f;
  301. s3_g *= exp_thing * exp_f;
  302. m_sing_new << sqrt(s1_g / (2 * (s1 - s1_min))), sqrt(s2_g / (2 * (s2 - s2_min))), sqrt(
  303. s3_g / (2 * (s3 - s3_min)));
  304. // change local step
  305. closest_sing_vec << s1_min, s2_min, s3_min;
  306. ri = ui * closest_sing_vec.asDiagonal() * vi.transpose();
  307. }
  308. }
  309. if (std::abs(s1 - 1) < eps) m_sing_new(0) = 1;
  310. if (std::abs(s2 - 1) < eps) m_sing_new(1) = 1;
  311. if (std::abs(s3 - 1) < eps) m_sing_new(2) = 1;
  312. Mat3 mat_W;
  313. mat_W = ui * m_sing_new.asDiagonal() * ui.transpose();
  314. s.W_11(i) = mat_W(0, 0);
  315. s.W_12(i) = mat_W(0, 1);
  316. s.W_13(i) = mat_W(0, 2);
  317. s.W_21(i) = mat_W(1, 0);
  318. s.W_22(i) = mat_W(1, 1);
  319. s.W_23(i) = mat_W(1, 2);
  320. s.W_31(i) = mat_W(2, 0);
  321. s.W_32(i) = mat_W(2, 1);
  322. s.W_33(i) = mat_W(2, 2);
  323. // 2) Update closest rotations (not rotations in case of conformal energy)
  324. s.Ri(i, 0) = ri(0, 0);
  325. s.Ri(i, 1) = ri(1, 0);
  326. s.Ri(i, 2) = ri(2, 0);
  327. s.Ri(i, 3) = ri(0, 1);
  328. s.Ri(i, 4) = ri(1, 1);
  329. s.Ri(i, 5) = ri(2, 1);
  330. s.Ri(i, 6) = ri(0, 2);
  331. s.Ri(i, 7) = ri(1, 2);
  332. s.Ri(i, 8) = ri(2, 2);
  333. } // for loop end
  334. } // if dim end
  335. }
  336. IGL_INLINE void solve_weighted_arap(igl::SLIMData& s,
  337. const Eigen::MatrixXd &V,
  338. const Eigen::MatrixXi &F,
  339. Eigen::MatrixXd &uv,
  340. Eigen::VectorXi &soft_b_p,
  341. Eigen::MatrixXd &soft_bc_p)
  342. {
  343. using namespace Eigen;
  344. Eigen::SparseMatrix<double> L;
  345. build_linear_system(s,L);
  346. // solve
  347. Eigen::VectorXd Uc;
  348. if (s.dim == 2)
  349. {
  350. SimplicialLDLT<Eigen::SparseMatrix<double> > solver;
  351. Uc = solver.compute(L).solve(s.rhs);
  352. }
  353. else
  354. { // seems like CG performs much worse for 2D and way better for 3D
  355. Eigen::VectorXd guess(uv.rows() * s.dim);
  356. for (int i = 0; i < s.dim; i++) for (int j = 0; j < s.dim; j++) guess(uv.rows() * i + j) = uv(i, j); // flatten vector
  357. ConjugateGradient<Eigen::SparseMatrix<double>, Eigen::Lower | Upper> solver;
  358. solver.setTolerance(1e-8);
  359. Uc = solver.compute(L).solveWithGuess(s.rhs, guess);
  360. }
  361. for (int i = 0; i < s.dim; i++)
  362. uv.col(i) = Uc.block(i * s.v_n, 0, s.v_n, 1);
  363. }
  364. IGL_INLINE void pre_calc(igl::SLIMData& s)
  365. {
  366. if (!s.has_pre_calc)
  367. {
  368. s.v_n = s.v_num;
  369. s.f_n = s.f_num;
  370. if (s.F.cols() == 3)
  371. {
  372. s.dim = 2;
  373. Eigen::MatrixXd F1, F2, F3;
  374. igl::local_basis(s.V, s.F, F1, F2, F3);
  375. compute_surface_gradient_matrix(s.V, s.F, F1, F2, s.Dx, s.Dy);
  376. s.W_11.resize(s.f_n);
  377. s.W_12.resize(s.f_n);
  378. s.W_21.resize(s.f_n);
  379. s.W_22.resize(s.f_n);
  380. }
  381. else
  382. {
  383. s.dim = 3;
  384. Eigen::SparseMatrix<double> G;
  385. igl::grad(s.V, s.F, G,
  386. s.mesh_improvement_3d /*use normal gradient, or one from a "regular" tet*/);
  387. s.Dx = G.block(0, 0, s.F.rows(), s.V.rows());
  388. s.Dy = G.block(s.F.rows(), 0, s.F.rows(), s.V.rows());
  389. s.Dz = G.block(2 * s.F.rows(), 0, s.F.rows(), s.V.rows());
  390. s.W_11.resize(s.f_n);
  391. s.W_12.resize(s.f_n);
  392. s.W_13.resize(s.f_n);
  393. s.W_21.resize(s.f_n);
  394. s.W_22.resize(s.f_n);
  395. s.W_23.resize(s.f_n);
  396. s.W_31.resize(s.f_n);
  397. s.W_32.resize(s.f_n);
  398. s.W_33.resize(s.f_n);
  399. }
  400. s.Dx.makeCompressed();
  401. s.Dy.makeCompressed();
  402. s.Dz.makeCompressed();
  403. s.Ri.resize(s.f_n, s.dim * s.dim);
  404. s.Ji.resize(s.f_n, s.dim * s.dim);
  405. s.rhs.resize(s.dim * s.v_num);
  406. // flattened weight matrix
  407. s.WGL_M.resize(s.dim * s.dim * s.f_n);
  408. for (int i = 0; i < s.dim * s.dim; i++)
  409. for (int j = 0; j < s.f_n; j++)
  410. s.WGL_M(i * s.f_n + j) = s.M(j);
  411. s.first_solve = true;
  412. s.has_pre_calc = true;
  413. }
  414. }
  415. IGL_INLINE void build_linear_system(igl::SLIMData& s, Eigen::SparseMatrix<double> &L)
  416. {
  417. // formula (35) in paper
  418. Eigen::SparseMatrix<double> A(s.dim * s.dim * s.f_n, s.dim * s.v_n);
  419. buildA(s,A);
  420. Eigen::SparseMatrix<double> At = A.transpose();
  421. At.makeCompressed();
  422. Eigen::SparseMatrix<double> id_m(At.rows(), At.rows());
  423. id_m.setIdentity();
  424. // add proximal penalty
  425. L = At * s.WGL_M.asDiagonal() * A + s.proximal_p * id_m; //add also a proximal term
  426. L.makeCompressed();
  427. buildRhs(s, At);
  428. Eigen::SparseMatrix<double> OldL = L;
  429. add_soft_constraints(s,L);
  430. L.makeCompressed();
  431. }
  432. IGL_INLINE void add_soft_constraints(igl::SLIMData& s, Eigen::SparseMatrix<double> &L)
  433. {
  434. int v_n = s.v_num;
  435. for (int d = 0; d < s.dim; d++)
  436. {
  437. for (int i = 0; i < s.b.rows(); i++)
  438. {
  439. int v_idx = s.b(i);
  440. s.rhs(d * v_n + v_idx) += s.soft_const_p * s.bc(i, d); // rhs
  441. L.coeffRef(d * v_n + v_idx, d * v_n + v_idx) += s.soft_const_p; // diagonal of matrix
  442. }
  443. }
  444. }
  445. IGL_INLINE double compute_energy(igl::SLIMData& s, Eigen::MatrixXd &V_new)
  446. {
  447. compute_jacobians(s,V_new);
  448. return compute_energy_with_jacobians(s, s.V, s.F, s.Ji, V_new, s.M) +
  449. compute_soft_const_energy(s, s.V, s.F, V_new);
  450. }
  451. IGL_INLINE double compute_soft_const_energy(igl::SLIMData& s,
  452. const Eigen::MatrixXd &V,
  453. const Eigen::MatrixXi &F,
  454. Eigen::MatrixXd &V_o)
  455. {
  456. double e = 0;
  457. for (int i = 0; i < s.b.rows(); i++)
  458. {
  459. e += s.soft_const_p * (s.bc.row(i) - V_o.row(s.b(i))).squaredNorm();
  460. }
  461. return e;
  462. }
  463. IGL_INLINE double compute_energy_with_jacobians(igl::SLIMData& s,
  464. const Eigen::MatrixXd &V,
  465. const Eigen::MatrixXi &F, const Eigen::MatrixXd &Ji,
  466. Eigen::MatrixXd &uv, Eigen::VectorXd &areas)
  467. {
  468. double energy = 0;
  469. if (s.dim == 2)
  470. {
  471. Eigen::Matrix<double, 2, 2> ji;
  472. for (int i = 0; i < s.f_n; i++)
  473. {
  474. ji(0, 0) = Ji(i, 0);
  475. ji(0, 1) = Ji(i, 1);
  476. ji(1, 0) = Ji(i, 2);
  477. ji(1, 1) = Ji(i, 3);
  478. typedef Eigen::Matrix<double, 2, 2> Mat2;
  479. typedef Eigen::Matrix<double, 2, 1> Vec2;
  480. Mat2 ri, ti, ui, vi;
  481. Vec2 sing;
  482. igl::polar_svd(ji, ri, ti, ui, sing, vi);
  483. double s1 = sing(0);
  484. double s2 = sing(1);
  485. switch (s.slim_energy)
  486. {
  487. case igl::SLIMData::ARAP:
  488. {
  489. energy += areas(i) * (pow(s1 - 1, 2) + pow(s2 - 1, 2));
  490. break;
  491. }
  492. case igl::SLIMData::SYMMETRIC_DIRICHLET:
  493. {
  494. energy += areas(i) * (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2));
  495. break;
  496. }
  497. case igl::SLIMData::EXP_SYMMETRIC_DIRICHLET:
  498. {
  499. energy += areas(i) * exp(s.exp_factor * (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2)));
  500. break;
  501. }
  502. case igl::SLIMData::LOG_ARAP:
  503. {
  504. energy += areas(i) * (pow(log(s1), 2) + pow(log(s2), 2));
  505. break;
  506. }
  507. case igl::SLIMData::CONFORMAL:
  508. {
  509. energy += areas(i) * ((pow(s1, 2) + pow(s2, 2)) / (2 * s1 * s2));
  510. break;
  511. }
  512. case igl::SLIMData::EXP_CONFORMAL:
  513. {
  514. energy += areas(i) * exp(s.exp_factor * ((pow(s1, 2) + pow(s2, 2)) / (2 * s1 * s2)));
  515. break;
  516. }
  517. }
  518. }
  519. }
  520. else
  521. {
  522. Eigen::Matrix<double, 3, 3> ji;
  523. for (int i = 0; i < s.f_n; i++)
  524. {
  525. ji(0, 0) = Ji(i, 0);
  526. ji(0, 1) = Ji(i, 1);
  527. ji(0, 2) = Ji(i, 2);
  528. ji(1, 0) = Ji(i, 3);
  529. ji(1, 1) = Ji(i, 4);
  530. ji(1, 2) = Ji(i, 5);
  531. ji(2, 0) = Ji(i, 6);
  532. ji(2, 1) = Ji(i, 7);
  533. ji(2, 2) = Ji(i, 8);
  534. typedef Eigen::Matrix<double, 3, 3> Mat3;
  535. typedef Eigen::Matrix<double, 3, 1> Vec3;
  536. Mat3 ri, ti, ui, vi;
  537. Vec3 sing;
  538. igl::polar_svd(ji, ri, ti, ui, sing, vi);
  539. double s1 = sing(0);
  540. double s2 = sing(1);
  541. double s3 = sing(2);
  542. switch (s.slim_energy)
  543. {
  544. case igl::SLIMData::ARAP:
  545. {
  546. energy += areas(i) * (pow(s1 - 1, 2) + pow(s2 - 1, 2) + pow(s3 - 1, 2));
  547. break;
  548. }
  549. case igl::SLIMData::SYMMETRIC_DIRICHLET:
  550. {
  551. energy += areas(i) * (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2) + pow(s3, 2) + pow(s3, -2));
  552. break;
  553. }
  554. case igl::SLIMData::EXP_SYMMETRIC_DIRICHLET:
  555. {
  556. energy += areas(i) * exp(s.exp_factor *
  557. (pow(s1, 2) + pow(s1, -2) + pow(s2, 2) + pow(s2, -2) + pow(s3, 2) + pow(s3, -2)));
  558. break;
  559. }
  560. case igl::SLIMData::LOG_ARAP:
  561. {
  562. energy += areas(i) * (pow(log(s1), 2) + pow(log(std::abs(s2)), 2) + pow(log(std::abs(s3)), 2));
  563. break;
  564. }
  565. case igl::SLIMData::CONFORMAL:
  566. {
  567. energy += areas(i) * ((pow(s1, 2) + pow(s2, 2) + pow(s3, 2)) / (3 * pow(s1 * s2 * s3, 2. / 3.)));
  568. break;
  569. }
  570. case igl::SLIMData::EXP_CONFORMAL:
  571. {
  572. energy += areas(i) * exp((pow(s1, 2) + pow(s2, 2) + pow(s3, 2)) / (3 * pow(s1 * s2 * s3, 2. / 3.)));
  573. break;
  574. }
  575. }
  576. }
  577. }
  578. return energy;
  579. }
  580. IGL_INLINE void buildA(igl::SLIMData& s, Eigen::SparseMatrix<double> &A)
  581. {
  582. // formula (35) in paper
  583. std::vector<Eigen::Triplet<double> > IJV;
  584. if (s.dim == 2)
  585. {
  586. IJV.reserve(4 * (s.Dx.outerSize() + s.Dy.outerSize()));
  587. /*A = [W11*Dx, W12*Dx;
  588. W11*Dy, W12*Dy;
  589. W21*Dx, W22*Dx;
  590. W21*Dy, W22*Dy];*/
  591. for (int k = 0; k < s.Dx.outerSize(); ++k)
  592. {
  593. for (Eigen::SparseMatrix<double>::InnerIterator it(s.Dx, k); it; ++it)
  594. {
  595. int dx_r = it.row();
  596. int dx_c = it.col();
  597. double val = it.value();
  598. IJV.push_back(Eigen::Triplet<double>(dx_r, dx_c, val * s.W_11(dx_r)));
  599. IJV.push_back(Eigen::Triplet<double>(dx_r, s.v_n + dx_c, val * s.W_12(dx_r)));
  600. IJV.push_back(Eigen::Triplet<double>(2 * s.f_n + dx_r, dx_c, val * s.W_21(dx_r)));
  601. IJV.push_back(Eigen::Triplet<double>(2 * s.f_n + dx_r, s.v_n + dx_c, val * s.W_22(dx_r)));
  602. }
  603. }
  604. for (int k = 0; k < s.Dy.outerSize(); ++k)
  605. {
  606. for (Eigen::SparseMatrix<double>::InnerIterator it(s.Dy, k); it; ++it)
  607. {
  608. int dy_r = it.row();
  609. int dy_c = it.col();
  610. double val = it.value();
  611. IJV.push_back(Eigen::Triplet<double>(s.f_n + dy_r, dy_c, val * s.W_11(dy_r)));
  612. IJV.push_back(Eigen::Triplet<double>(s.f_n + dy_r, s.v_n + dy_c, val * s.W_12(dy_r)));
  613. IJV.push_back(Eigen::Triplet<double>(3 * s.f_n + dy_r, dy_c, val * s.W_21(dy_r)));
  614. IJV.push_back(Eigen::Triplet<double>(3 * s.f_n + dy_r, s.v_n + dy_c, val * s.W_22(dy_r)));
  615. }
  616. }
  617. }
  618. else
  619. {
  620. /*A = [W11*Dx, W12*Dx, W13*Dx;
  621. W11*Dy, W12*Dy, W13*Dy;
  622. W11*Dz, W12*Dz, W13*Dz;
  623. W21*Dx, W22*Dx, W23*Dx;
  624. W21*Dy, W22*Dy, W23*Dy;
  625. W21*Dz, W22*Dz, W23*Dz;
  626. W31*Dx, W32*Dx, W33*Dx;
  627. W31*Dy, W32*Dy, W33*Dy;
  628. W31*Dz, W32*Dz, W33*Dz;];*/
  629. IJV.reserve(9 * (s.Dx.outerSize() + s.Dy.outerSize() + s.Dz.outerSize()));
  630. for (int k = 0; k < s.Dx.outerSize(); k++)
  631. {
  632. for (Eigen::SparseMatrix<double>::InnerIterator it(s.Dx, k); it; ++it)
  633. {
  634. int dx_r = it.row();
  635. int dx_c = it.col();
  636. double val = it.value();
  637. IJV.push_back(Eigen::Triplet<double>(dx_r, dx_c, val * s.W_11(dx_r)));
  638. IJV.push_back(Eigen::Triplet<double>(dx_r, s.v_n + dx_c, val * s.W_12(dx_r)));
  639. IJV.push_back(Eigen::Triplet<double>(dx_r, 2 * s.v_n + dx_c, val * s.W_13(dx_r)));
  640. IJV.push_back(Eigen::Triplet<double>(3 * s.f_n + dx_r, dx_c, val * s.W_21(dx_r)));
  641. IJV.push_back(Eigen::Triplet<double>(3 * s.f_n + dx_r, s.v_n + dx_c, val * s.W_22(dx_r)));
  642. IJV.push_back(Eigen::Triplet<double>(3 * s.f_n + dx_r, 2 * s.v_n + dx_c, val * s.W_23(dx_r)));
  643. IJV.push_back(Eigen::Triplet<double>(6 * s.f_n + dx_r, dx_c, val * s.W_31(dx_r)));
  644. IJV.push_back(Eigen::Triplet<double>(6 * s.f_n + dx_r, s.v_n + dx_c, val * s.W_32(dx_r)));
  645. IJV.push_back(Eigen::Triplet<double>(6 * s.f_n + dx_r, 2 * s.v_n + dx_c, val * s.W_33(dx_r)));
  646. }
  647. }
  648. for (int k = 0; k < s.Dy.outerSize(); k++)
  649. {
  650. for (Eigen::SparseMatrix<double>::InnerIterator it(s.Dy, k); it; ++it)
  651. {
  652. int dy_r = it.row();
  653. int dy_c = it.col();
  654. double val = it.value();
  655. IJV.push_back(Eigen::Triplet<double>(s.f_n + dy_r, dy_c, val * s.W_11(dy_r)));
  656. IJV.push_back(Eigen::Triplet<double>(s.f_n + dy_r, s.v_n + dy_c, val * s.W_12(dy_r)));
  657. IJV.push_back(Eigen::Triplet<double>(s.f_n + dy_r, 2 * s.v_n + dy_c, val * s.W_13(dy_r)));
  658. IJV.push_back(Eigen::Triplet<double>(4 * s.f_n + dy_r, dy_c, val * s.W_21(dy_r)));
  659. IJV.push_back(Eigen::Triplet<double>(4 * s.f_n + dy_r, s.v_n + dy_c, val * s.W_22(dy_r)));
  660. IJV.push_back(Eigen::Triplet<double>(4 * s.f_n + dy_r, 2 * s.v_n + dy_c, val * s.W_23(dy_r)));
  661. IJV.push_back(Eigen::Triplet<double>(7 * s.f_n + dy_r, dy_c, val * s.W_31(dy_r)));
  662. IJV.push_back(Eigen::Triplet<double>(7 * s.f_n + dy_r, s.v_n + dy_c, val * s.W_32(dy_r)));
  663. IJV.push_back(Eigen::Triplet<double>(7 * s.f_n + dy_r, 2 * s.v_n + dy_c, val * s.W_33(dy_r)));
  664. }
  665. }
  666. for (int k = 0; k < s.Dz.outerSize(); k++)
  667. {
  668. for (Eigen::SparseMatrix<double>::InnerIterator it(s.Dz, k); it; ++it)
  669. {
  670. int dz_r = it.row();
  671. int dz_c = it.col();
  672. double val = it.value();
  673. IJV.push_back(Eigen::Triplet<double>(2 * s.f_n + dz_r, dz_c, val * s.W_11(dz_r)));
  674. IJV.push_back(Eigen::Triplet<double>(2 * s.f_n + dz_r, s.v_n + dz_c, val * s.W_12(dz_r)));
  675. IJV.push_back(Eigen::Triplet<double>(2 * s.f_n + dz_r, 2 * s.v_n + dz_c, val * s.W_13(dz_r)));
  676. IJV.push_back(Eigen::Triplet<double>(5 * s.f_n + dz_r, dz_c, val * s.W_21(dz_r)));
  677. IJV.push_back(Eigen::Triplet<double>(5 * s.f_n + dz_r, s.v_n + dz_c, val * s.W_22(dz_r)));
  678. IJV.push_back(Eigen::Triplet<double>(5 * s.f_n + dz_r, 2 * s.v_n + dz_c, val * s.W_23(dz_r)));
  679. IJV.push_back(Eigen::Triplet<double>(8 * s.f_n + dz_r, dz_c, val * s.W_31(dz_r)));
  680. IJV.push_back(Eigen::Triplet<double>(8 * s.f_n + dz_r, s.v_n + dz_c, val * s.W_32(dz_r)));
  681. IJV.push_back(Eigen::Triplet<double>(8 * s.f_n + dz_r, 2 * s.v_n + dz_c, val * s.W_33(dz_r)));
  682. }
  683. }
  684. }
  685. A.setFromTriplets(IJV.begin(), IJV.end());
  686. }
  687. IGL_INLINE void buildRhs(igl::SLIMData& s, const Eigen::SparseMatrix<double> &At)
  688. {
  689. Eigen::VectorXd f_rhs(s.dim * s.dim * s.f_n);
  690. f_rhs.setZero();
  691. if (s.dim == 2)
  692. {
  693. /*b = [W11*R11 + W12*R21; (formula (36))
  694. W11*R12 + W12*R22;
  695. W21*R11 + W22*R21;
  696. W21*R12 + W22*R22];*/
  697. for (int i = 0; i < s.f_n; i++)
  698. {
  699. f_rhs(i + 0 * s.f_n) = s.W_11(i) * s.Ri(i, 0) + s.W_12(i) * s.Ri(i, 1);
  700. f_rhs(i + 1 * s.f_n) = s.W_11(i) * s.Ri(i, 2) + s.W_12(i) * s.Ri(i, 3);
  701. f_rhs(i + 2 * s.f_n) = s.W_21(i) * s.Ri(i, 0) + s.W_22(i) * s.Ri(i, 1);
  702. f_rhs(i + 3 * s.f_n) = s.W_21(i) * s.Ri(i, 2) + s.W_22(i) * s.Ri(i, 3);
  703. }
  704. }
  705. else
  706. {
  707. /*b = [W11*R11 + W12*R21 + W13*R31;
  708. W11*R12 + W12*R22 + W13*R32;
  709. W11*R13 + W12*R23 + W13*R33;
  710. W21*R11 + W22*R21 + W23*R31;
  711. W21*R12 + W22*R22 + W23*R32;
  712. W21*R13 + W22*R23 + W23*R33;
  713. W31*R11 + W32*R21 + W33*R31;
  714. W31*R12 + W32*R22 + W33*R32;
  715. W31*R13 + W32*R23 + W33*R33;];*/
  716. for (int i = 0; i < s.f_n; i++)
  717. {
  718. f_rhs(i + 0 * s.f_n) = s.W_11(i) * s.Ri(i, 0) + s.W_12(i) * s.Ri(i, 1) + s.W_13(i) * s.Ri(i, 2);
  719. f_rhs(i + 1 * s.f_n) = s.W_11(i) * s.Ri(i, 3) + s.W_12(i) * s.Ri(i, 4) + s.W_13(i) * s.Ri(i, 5);
  720. f_rhs(i + 2 * s.f_n) = s.W_11(i) * s.Ri(i, 6) + s.W_12(i) * s.Ri(i, 7) + s.W_13(i) * s.Ri(i, 8);
  721. f_rhs(i + 3 * s.f_n) = s.W_21(i) * s.Ri(i, 0) + s.W_22(i) * s.Ri(i, 1) + s.W_23(i) * s.Ri(i, 2);
  722. f_rhs(i + 4 * s.f_n) = s.W_21(i) * s.Ri(i, 3) + s.W_22(i) * s.Ri(i, 4) + s.W_23(i) * s.Ri(i, 5);
  723. f_rhs(i + 5 * s.f_n) = s.W_21(i) * s.Ri(i, 6) + s.W_22(i) * s.Ri(i, 7) + s.W_23(i) * s.Ri(i, 8);
  724. f_rhs(i + 6 * s.f_n) = s.W_31(i) * s.Ri(i, 0) + s.W_32(i) * s.Ri(i, 1) + s.W_33(i) * s.Ri(i, 2);
  725. f_rhs(i + 7 * s.f_n) = s.W_31(i) * s.Ri(i, 3) + s.W_32(i) * s.Ri(i, 4) + s.W_33(i) * s.Ri(i, 5);
  726. f_rhs(i + 8 * s.f_n) = s.W_31(i) * s.Ri(i, 6) + s.W_32(i) * s.Ri(i, 7) + s.W_33(i) * s.Ri(i, 8);
  727. }
  728. }
  729. Eigen::VectorXd uv_flat(s.dim *s.v_n);
  730. for (int i = 0; i < s.dim; i++)
  731. for (int j = 0; j < s.v_n; j++)
  732. uv_flat(s.v_n * i + j) = s.V_o(j, i);
  733. s.rhs = (At * s.WGL_M.asDiagonal() * f_rhs + s.proximal_p * uv_flat);
  734. }
  735. }
  736. }
  737. /// Slim Implementation
  738. IGL_INLINE void igl::slim_precompute(Eigen::MatrixXd &V, Eigen::MatrixXi &F, Eigen::MatrixXd &V_init, SLIMData &data,
  739. SLIMData::SLIM_ENERGY slim_energy, Eigen::VectorXi &b, Eigen::MatrixXd &bc,
  740. double soft_p)
  741. {
  742. data.V = V;
  743. data.F = F;
  744. data.V_o = V_init;
  745. data.v_num = V.rows();
  746. data.f_num = F.rows();
  747. data.slim_energy = slim_energy;
  748. data.b = b;
  749. data.bc = bc;
  750. data.soft_const_p = soft_p;
  751. data.proximal_p = 0.0001;
  752. igl::doublearea(V, F, data.M);
  753. data.M /= 2.;
  754. data.mesh_area = data.M.sum();
  755. data.mesh_improvement_3d = false; // whether to use a jacobian derived from a real mesh or an abstract regular mesh (used for mesh improvement)
  756. data.exp_factor = 1.0; // param used only for exponential energies (e.g exponential symmetric dirichlet)
  757. assert (F.cols() == 3 || F.cols() == 4);
  758. igl::slim::pre_calc(data);
  759. data.energy = igl::slim::compute_energy(data,data.V_o) / data.mesh_area;
  760. }
  761. IGL_INLINE Eigen::MatrixXd igl::slim_solve(SLIMData &data, int iter_num)
  762. {
  763. for (int i = 0; i < iter_num; i++)
  764. {
  765. Eigen::MatrixXd dest_res;
  766. dest_res = data.V_o;
  767. // Solve Weighted Proxy
  768. igl::slim::update_weights_and_closest_rotations(data,data.V, data.F, dest_res);
  769. igl::slim::solve_weighted_arap(data,data.V, data.F, dest_res, data.b, data.bc);
  770. double old_energy = data.energy;
  771. std::function<double(Eigen::MatrixXd &)> compute_energy = [&](
  772. Eigen::MatrixXd &aaa) { return igl::slim::compute_energy(data,aaa); };
  773. data.energy = igl::flip_avoiding_line_search(data.F, data.V_o, dest_res, compute_energy,
  774. data.energy * data.mesh_area) / data.mesh_area;
  775. }
  776. return data.V_o;
  777. }