slim.cpp 28 KB

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