flip_avoiding_line_search.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317
  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 "flip_avoiding_line_search.h"
  9. #include "line_search.h"
  10. #include <Eigen/Dense>
  11. #include <vector>
  12. namespace igl
  13. {
  14. namespace flip_avoiding
  15. {
  16. //---------------------------------------------------------------------------
  17. // x - array of size 3
  18. // In case 3 real roots: => x[0], x[1], x[2], return 3
  19. // 2 real roots: x[0], x[1], return 2
  20. // 1 real root : x[0], x[1] ± i*x[2], return 1
  21. // http://math.ivanovo.ac.ru/dalgebra/Khashin/poly/index.html
  22. IGL_INLINE int SolveP3(std::vector<double>& x,double a,double b,double c)
  23. { // solve cubic equation x^3 + a*x^2 + b*x + c
  24. using namespace std;
  25. double a2 = a*a;
  26. double q = (a2 - 3*b)/9;
  27. double r = (a*(2*a2-9*b) + 27*c)/54;
  28. double r2 = r*r;
  29. double q3 = q*q*q;
  30. double A,B;
  31. if(r2<q3)
  32. {
  33. double t=r/sqrt(q3);
  34. if( t<-1) t=-1;
  35. if( t> 1) t= 1;
  36. t=acos(t);
  37. a/=3; q=-2*sqrt(q);
  38. x[0]=q*cos(t/3)-a;
  39. x[1]=q*cos((t+(2*M_PI))/3)-a;
  40. x[2]=q*cos((t-(2*M_PI))/3)-a;
  41. return(3);
  42. }
  43. else
  44. {
  45. A =-pow(fabs(r)+sqrt(r2-q3),1./3);
  46. if( r<0 ) A=-A;
  47. B = A==0? 0 : B=q/A;
  48. a/=3;
  49. x[0] =(A+B)-a;
  50. x[1] =-0.5*(A+B)-a;
  51. x[2] = 0.5*sqrt(3.)*(A-B);
  52. if(fabs(x[2])<1e-14)
  53. {
  54. x[2]=x[1]; return(2);
  55. }
  56. return(1);
  57. }
  58. }
  59. IGL_INLINE double get_smallest_pos_quad_zero(double a,double b, double c)
  60. {
  61. using namespace std;
  62. double t1,t2;
  63. if (a != 0)
  64. {
  65. double delta_in = pow(b,2) - 4*a*c;
  66. if (delta_in < 0)
  67. {
  68. return INFINITY;
  69. }
  70. double delta = sqrt(delta_in);
  71. t1 = (-b + delta)/ (2*a);
  72. t2 = (-b - delta)/ (2*a);
  73. }
  74. else
  75. {
  76. t1 = t2 = -b/c;
  77. }
  78. assert (std::isfinite(t1));
  79. assert (std::isfinite(t2));
  80. double tmp_n = min(t1,t2);
  81. t1 = max(t1,t2); t2 = tmp_n;
  82. if (t1 == t2)
  83. {
  84. return INFINITY; // means the orientation flips twice = doesn't flip
  85. }
  86. // return the smallest negative root if it exists, otherwise return infinity
  87. if (t1 > 0)
  88. {
  89. if (t2 > 0)
  90. {
  91. return t2;
  92. }
  93. else
  94. {
  95. return t1;
  96. }
  97. }
  98. else
  99. {
  100. return INFINITY;
  101. }
  102. }
  103. IGL_INLINE double get_min_pos_root_2D(const Eigen::MatrixXd& uv,
  104. const Eigen::MatrixXi& F,
  105. Eigen::MatrixXd& d,
  106. int f)
  107. {
  108. using namespace std;
  109. /*
  110. Finding the smallest timestep t s.t a triangle get degenerated (<=> det = 0)
  111. The following code can be derived by a symbolic expression in matlab:
  112. Symbolic matlab:
  113. U11 = sym('U11');
  114. U12 = sym('U12');
  115. U21 = sym('U21');
  116. U22 = sym('U22');
  117. U31 = sym('U31');
  118. U32 = sym('U32');
  119. V11 = sym('V11');
  120. V12 = sym('V12');
  121. V21 = sym('V21');
  122. V22 = sym('V22');
  123. V31 = sym('V31');
  124. V32 = sym('V32');
  125. t = sym('t');
  126. U1 = [U11,U12];
  127. U2 = [U21,U22];
  128. U3 = [U31,U32];
  129. V1 = [V11,V12];
  130. V2 = [V21,V22];
  131. V3 = [V31,V32];
  132. A = [(U2+V2*t) - (U1+ V1*t)];
  133. B = [(U3+V3*t) - (U1+ V1*t)];
  134. C = [A;B];
  135. solve(det(C), t);
  136. cf = coeffs(det(C),t); % Now cf(1),cf(2),cf(3) holds the coefficients for the polynom. at order c,b,a
  137. */
  138. int v1 = F(f,0); int v2 = F(f,1); int v3 = F(f,2);
  139. // get quadratic coefficients (ax^2 + b^x + c)
  140. const double& U11 = uv(v1,0);
  141. const double& U12 = uv(v1,1);
  142. const double& U21 = uv(v2,0);
  143. const double& U22 = uv(v2,1);
  144. const double& U31 = uv(v3,0);
  145. const double& U32 = uv(v3,1);
  146. const double& V11 = d(v1,0);
  147. const double& V12 = d(v1,1);
  148. const double& V21 = d(v2,0);
  149. const double& V22 = d(v2,1);
  150. const double& V31 = d(v3,0);
  151. const double& V32 = d(v3,1);
  152. double a = V11*V22 - V12*V21 - V11*V32 + V12*V31 + V21*V32 - V22*V31;
  153. double b = U11*V22 - U12*V21 - U21*V12 + U22*V11 - U11*V32 + U12*V31 + U31*V12 - U32*V11 + U21*V32 - U22*V31 - U31*V22 + U32*V21;
  154. double c = U11*U22 - U12*U21 - U11*U32 + U12*U31 + U21*U32 - U22*U31;
  155. return get_smallest_pos_quad_zero(a,b,c);
  156. }
  157. IGL_INLINE double get_min_pos_root_3D(const Eigen::MatrixXd& uv,
  158. const Eigen::MatrixXi& F,
  159. Eigen::MatrixXd& direc,
  160. int f)
  161. {
  162. using namespace std;
  163. /*
  164. Searching for the roots of:
  165. +-1/6 * |ax ay az 1|
  166. |bx by bz 1|
  167. |cx cy cz 1|
  168. |dx dy dz 1|
  169. Every point ax,ay,az has a search direction a_dx,a_dy,a_dz, and so we add those to the matrix, and solve the cubic to find the step size t for a 0 volume
  170. Symbolic matlab:
  171. syms a_x a_y a_z a_dx a_dy a_dz % tetrahedera point and search direction
  172. syms b_x b_y b_z b_dx b_dy b_dz
  173. syms c_x c_y c_z c_dx c_dy c_dz
  174. syms d_x d_y d_z d_dx d_dy d_dz
  175. syms t % Timestep var, this is what we're looking for
  176. a_plus_t = [a_x,a_y,a_z] + t*[a_dx,a_dy,a_dz];
  177. b_plus_t = [b_x,b_y,b_z] + t*[b_dx,b_dy,b_dz];
  178. c_plus_t = [c_x,c_y,c_z] + t*[c_dx,c_dy,c_dz];
  179. d_plus_t = [d_x,d_y,d_z] + t*[d_dx,d_dy,d_dz];
  180. vol_mat = [a_plus_t,1;b_plus_t,1;c_plus_t,1;d_plus_t,1]
  181. //cf = coeffs(det(vol_det),t); % Now cf(1),cf(2),cf(3),cf(4) holds the coefficients for the polynom
  182. [coefficients,terms] = coeffs(det(vol_det),t); % terms = [ t^3, t^2, t, 1], Coefficients hold the coeff we seek
  183. */
  184. int v1 = F(f,0); int v2 = F(f,1); int v3 = F(f,2); int v4 = F(f,3);
  185. const double& a_x = uv(v1,0);
  186. const double& a_y = uv(v1,1);
  187. const double& a_z = uv(v1,2);
  188. const double& b_x = uv(v2,0);
  189. const double& b_y = uv(v2,1);
  190. const double& b_z = uv(v2,2);
  191. const double& c_x = uv(v3,0);
  192. const double& c_y = uv(v3,1);
  193. const double& c_z = uv(v3,2);
  194. const double& d_x = uv(v4,0);
  195. const double& d_y = uv(v4,1);
  196. const double& d_z = uv(v4,2);
  197. const double& a_dx = direc(v1,0);
  198. const double& a_dy = direc(v1,1);
  199. const double& a_dz = direc(v1,2);
  200. const double& b_dx = direc(v2,0);
  201. const double& b_dy = direc(v2,1);
  202. const double& b_dz = direc(v2,2);
  203. const double& c_dx = direc(v3,0);
  204. const double& c_dy = direc(v3,1);
  205. const double& c_dz = direc(v3,2);
  206. const double& d_dx = direc(v4,0);
  207. const double& d_dy = direc(v4,1);
  208. const double& d_dz = direc(v4,2);
  209. // Find solution for: a*t^3 + b*t^2 + c*d +d = 0
  210. double a = a_dx*b_dy*c_dz - a_dx*b_dz*c_dy - a_dy*b_dx*c_dz + a_dy*b_dz*c_dx + a_dz*b_dx*c_dy - a_dz*b_dy*c_dx - a_dx*b_dy*d_dz + a_dx*b_dz*d_dy + a_dy*b_dx*d_dz - a_dy*b_dz*d_dx - a_dz*b_dx*d_dy + a_dz*b_dy*d_dx + a_dx*c_dy*d_dz - a_dx*c_dz*d_dy - a_dy*c_dx*d_dz + a_dy*c_dz*d_dx + a_dz*c_dx*d_dy - a_dz*c_dy*d_dx - b_dx*c_dy*d_dz + b_dx*c_dz*d_dy + b_dy*c_dx*d_dz - b_dy*c_dz*d_dx - b_dz*c_dx*d_dy + b_dz*c_dy*d_dx;
  211. double b = a_dy*b_dz*c_x - a_dy*b_x*c_dz - a_dz*b_dy*c_x + a_dz*b_x*c_dy + a_x*b_dy*c_dz - a_x*b_dz*c_dy - a_dx*b_dz*c_y + a_dx*b_y*c_dz + a_dz*b_dx*c_y - a_dz*b_y*c_dx - a_y*b_dx*c_dz + a_y*b_dz*c_dx + a_dx*b_dy*c_z - a_dx*b_z*c_dy - a_dy*b_dx*c_z + a_dy*b_z*c_dx + a_z*b_dx*c_dy - a_z*b_dy*c_dx - a_dy*b_dz*d_x + a_dy*b_x*d_dz + a_dz*b_dy*d_x - a_dz*b_x*d_dy - a_x*b_dy*d_dz + a_x*b_dz*d_dy + a_dx*b_dz*d_y - a_dx*b_y*d_dz - a_dz*b_dx*d_y + a_dz*b_y*d_dx + a_y*b_dx*d_dz - a_y*b_dz*d_dx - a_dx*b_dy*d_z + a_dx*b_z*d_dy + a_dy*b_dx*d_z - a_dy*b_z*d_dx - a_z*b_dx*d_dy + a_z*b_dy*d_dx + a_dy*c_dz*d_x - a_dy*c_x*d_dz - a_dz*c_dy*d_x + a_dz*c_x*d_dy + a_x*c_dy*d_dz - a_x*c_dz*d_dy - a_dx*c_dz*d_y + a_dx*c_y*d_dz + a_dz*c_dx*d_y - a_dz*c_y*d_dx - a_y*c_dx*d_dz + a_y*c_dz*d_dx + a_dx*c_dy*d_z - a_dx*c_z*d_dy - a_dy*c_dx*d_z + a_dy*c_z*d_dx + a_z*c_dx*d_dy - a_z*c_dy*d_dx - b_dy*c_dz*d_x + b_dy*c_x*d_dz + b_dz*c_dy*d_x - b_dz*c_x*d_dy - b_x*c_dy*d_dz + b_x*c_dz*d_dy + b_dx*c_dz*d_y - b_dx*c_y*d_dz - b_dz*c_dx*d_y + b_dz*c_y*d_dx + b_y*c_dx*d_dz - b_y*c_dz*d_dx - b_dx*c_dy*d_z + b_dx*c_z*d_dy + b_dy*c_dx*d_z - b_dy*c_z*d_dx - b_z*c_dx*d_dy + b_z*c_dy*d_dx;
  212. double c = a_dz*b_x*c_y - a_dz*b_y*c_x - a_x*b_dz*c_y + a_x*b_y*c_dz + a_y*b_dz*c_x - a_y*b_x*c_dz - a_dy*b_x*c_z + a_dy*b_z*c_x + a_x*b_dy*c_z - a_x*b_z*c_dy - a_z*b_dy*c_x + a_z*b_x*c_dy + a_dx*b_y*c_z - a_dx*b_z*c_y - a_y*b_dx*c_z + a_y*b_z*c_dx + a_z*b_dx*c_y - a_z*b_y*c_dx - a_dz*b_x*d_y + a_dz*b_y*d_x + a_x*b_dz*d_y - a_x*b_y*d_dz - a_y*b_dz*d_x + a_y*b_x*d_dz + a_dy*b_x*d_z - a_dy*b_z*d_x - a_x*b_dy*d_z + a_x*b_z*d_dy + a_z*b_dy*d_x - a_z*b_x*d_dy - a_dx*b_y*d_z + a_dx*b_z*d_y + a_y*b_dx*d_z - a_y*b_z*d_dx - a_z*b_dx*d_y + a_z*b_y*d_dx + a_dz*c_x*d_y - a_dz*c_y*d_x - a_x*c_dz*d_y + a_x*c_y*d_dz + a_y*c_dz*d_x - a_y*c_x*d_dz - a_dy*c_x*d_z + a_dy*c_z*d_x + a_x*c_dy*d_z - a_x*c_z*d_dy - a_z*c_dy*d_x + a_z*c_x*d_dy + a_dx*c_y*d_z - a_dx*c_z*d_y - a_y*c_dx*d_z + a_y*c_z*d_dx + a_z*c_dx*d_y - a_z*c_y*d_dx - b_dz*c_x*d_y + b_dz*c_y*d_x + b_x*c_dz*d_y - b_x*c_y*d_dz - b_y*c_dz*d_x + b_y*c_x*d_dz + b_dy*c_x*d_z - b_dy*c_z*d_x - b_x*c_dy*d_z + b_x*c_z*d_dy + b_z*c_dy*d_x - b_z*c_x*d_dy - b_dx*c_y*d_z + b_dx*c_z*d_y + b_y*c_dx*d_z - b_y*c_z*d_dx - b_z*c_dx*d_y + b_z*c_y*d_dx;
  213. double d = a_x*b_y*c_z - a_x*b_z*c_y - a_y*b_x*c_z + a_y*b_z*c_x + a_z*b_x*c_y - a_z*b_y*c_x - a_x*b_y*d_z + a_x*b_z*d_y + a_y*b_x*d_z - a_y*b_z*d_x - a_z*b_x*d_y + a_z*b_y*d_x + a_x*c_y*d_z - a_x*c_z*d_y - a_y*c_x*d_z + a_y*c_z*d_x + a_z*c_x*d_y - a_z*c_y*d_x - b_x*c_y*d_z + b_x*c_z*d_y + b_y*c_x*d_z - b_y*c_z*d_x - b_z*c_x*d_y + b_z*c_y*d_x;
  214. if (a==0)
  215. {
  216. return get_smallest_pos_quad_zero(b,c,d);
  217. }
  218. b/=a; c/=a; d/=a; // normalize it all
  219. std::vector<double> res(3);
  220. int real_roots_num = SolveP3(res,b,c,d);
  221. switch (real_roots_num)
  222. {
  223. case 1:
  224. return (res[0] >= 0) ? res[0]:INFINITY;
  225. case 2:
  226. {
  227. double max_root = max(res[0],res[1]); double min_root = min(res[0],res[1]);
  228. if (min_root > 0) return min_root;
  229. if (max_root > 0) return max_root;
  230. return INFINITY;
  231. }
  232. case 3:
  233. default:
  234. {
  235. std::sort(res.begin(),res.end());
  236. if (res[0] > 0) return res[0];
  237. if (res[1] > 0) return res[1];
  238. if (res[2] > 0) return res[2];
  239. return INFINITY;
  240. }
  241. }
  242. }
  243. IGL_INLINE double compute_max_step_from_singularities(const Eigen::MatrixXd& uv,
  244. const Eigen::MatrixXi& F,
  245. Eigen::MatrixXd& d)
  246. {
  247. using namespace std;
  248. double max_step = INFINITY;
  249. // The if statement is outside the for loops to avoid branching/ease parallelizing
  250. if (uv.cols() == 2)
  251. {
  252. for (int f = 0; f < F.rows(); f++)
  253. {
  254. double min_positive_root = get_min_pos_root_2D(uv,F,d,f);
  255. max_step = min(max_step, min_positive_root);
  256. }
  257. }
  258. else
  259. { // volumetric deformation
  260. for (int f = 0; f < F.rows(); f++)
  261. {
  262. double min_positive_root = get_min_pos_root_3D(uv,F,d,f);
  263. max_step = min(max_step, min_positive_root);
  264. }
  265. }
  266. return max_step;
  267. }
  268. }
  269. }
  270. IGL_INLINE double igl::flip_avoiding_line_search(
  271. const Eigen::MatrixXi F,
  272. Eigen::MatrixXd& cur_v,
  273. Eigen::MatrixXd& dst_v,
  274. std::function<double(Eigen::MatrixXd&)> energy,
  275. double cur_energy)
  276. {
  277. using namespace std;
  278. Eigen::MatrixXd d = dst_v - cur_v;
  279. double min_step_to_singularity = igl::flip_avoiding::compute_max_step_from_singularities(cur_v,F,d);
  280. double max_step_size = min(1., min_step_to_singularity*0.8);
  281. return igl::line_search(cur_v,d,max_step_size, energy, cur_energy);
  282. }
  283. #ifdef IGL_STATIC_LIBRARY
  284. #endif