slim.cpp 32 KB

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