flip_avoiding_line_search.cpp 11 KB

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