slim.cpp 35 KB

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