fit_cubic_bezier.cpp 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2020 Alec Jacobson <[email protected]>
  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 "fit_cubic_bezier.h"
  9. #include "bezier.h"
  10. #include "EPS.h"
  11. // Adapted from main.c accompanying
  12. // An Algorithm for Automatically Fitting Digitized Curves
  13. // by Philip J. Schneider
  14. // from "Graphics Gems", Academic Press, 1990
  15. IGL_INLINE void igl::fit_cubic_bezier(
  16. const Eigen::MatrixXd & d,
  17. const double error,
  18. std::vector<Eigen::MatrixXd> & cubics)
  19. {
  20. const int nPts = d.rows();
  21. // Don't attempt to fit curve to single point
  22. if(nPts==1) { return; }
  23. // Avoid using zero tangent
  24. const static auto tangent = [](
  25. const Eigen::MatrixXd & d,
  26. const int i,const int dir)->Eigen::RowVectorXd
  27. {
  28. int j = i;
  29. const int nPts = d.rows();
  30. Eigen::RowVectorXd t;
  31. while(true)
  32. {
  33. // look at next point
  34. j += dir;
  35. if(j < 0 || j>=nPts)
  36. {
  37. // All points are coincident?
  38. // give up and use zero tangent...
  39. return Eigen::RowVectorXd::Zero(1,d.cols());
  40. }
  41. t = d.row(j)-d.row(i);
  42. if(t.squaredNorm() > igl::DOUBLE_EPS)
  43. {
  44. break;
  45. }
  46. }
  47. return t.normalized();
  48. };
  49. Eigen::RowVectorXd tHat1 = tangent(d,0,+1);
  50. Eigen::RowVectorXd tHat2 = tangent(d,nPts-1,-1);
  51. // If first and last points are identically equal, then consider closed
  52. const bool closed = (d.row(0) - d.row(d.rows()-1)).squaredNorm() == 0;
  53. // If closed loop make tangents match
  54. if(closed)
  55. {
  56. tHat1 = (tHat1 - tHat2).eval().normalized();
  57. tHat2 = -tHat1;
  58. }
  59. cubics.clear();
  60. fit_cubic_bezier_substring(d,0,nPts-1,tHat1,tHat2,error,closed,cubics);
  61. };
  62. IGL_INLINE void igl::fit_cubic_bezier_substring(
  63. const Eigen::MatrixXd & d,
  64. const int first,
  65. const int last,
  66. const Eigen::RowVectorXd & tHat1,
  67. const Eigen::RowVectorXd & tHat2,
  68. const double error,
  69. const bool force_split,
  70. std::vector<Eigen::MatrixXd> & cubics)
  71. {
  72. // Helper functions
  73. // Evaluate a Bezier curve at a particular parameter value
  74. const static auto bezier_eval = [](const Eigen::MatrixXd & V, const double t)
  75. { Eigen::RowVectorXd P; bezier(V,t,P); return P; };
  76. //
  77. // Use Newton-Raphson iteration to find better root.
  78. const static auto NewtonRaphsonRootFind = [](
  79. const Eigen::MatrixXd & Q,
  80. const Eigen::RowVectorXd & P,
  81. const double u)->double
  82. {
  83. /* Compute Q(u) */
  84. Eigen::RowVectorXd Q_u = bezier_eval(Q, u);
  85. Eigen::MatrixXd Q1(3,Q.cols());
  86. Eigen::MatrixXd Q2(2,Q.cols());
  87. /* Generate control vertices for Q' */
  88. for (int i = 0; i <= 2; i++)
  89. {
  90. Q1.row(i) = (Q.row(i+1) - Q.row(i)) * 3.0;
  91. }
  92. /* Generate control vertices for Q'' */
  93. for (int i = 0; i <= 1; i++)
  94. {
  95. Q2.row(i) = (Q1.row(i+1) - Q1.row(i)) * 2.0;
  96. }
  97. /* Compute Q'(u) and Q''(u) */
  98. const Eigen::RowVectorXd Q1_u = bezier_eval(Q1, u);
  99. const Eigen::RowVectorXd Q2_u = bezier_eval(Q2, u);
  100. /* Compute f(u)/f'(u) */
  101. const double numerator = ((Q_u-P).array() * Q1_u.array()).array().sum();
  102. const double denominator =
  103. Q1_u.squaredNorm() + ((Q_u-P).array() * Q2_u.array()).array().sum();
  104. /* u = u - f(u)/f'(u) */
  105. return u - (numerator/denominator);
  106. };
  107. const static auto ComputeMaxError = [](
  108. const Eigen::MatrixXd & d,
  109. const int first,
  110. const int last,
  111. const Eigen::MatrixXd & bezCurve,
  112. const Eigen::VectorXd & u,
  113. int & splitPoint)->double
  114. {
  115. Eigen::VectorXd E(last - (first+1));
  116. splitPoint = (last-first + 1)/2;
  117. double maxDist = 0.0;
  118. for (int i = first + 1; i < last; i++)
  119. {
  120. Eigen::RowVectorXd P = bezier_eval(bezCurve, u(i-first));
  121. const double dist = (P-d.row(i)).squaredNorm();
  122. E(i-(first+1)) = dist;
  123. if (dist >= maxDist)
  124. {
  125. maxDist = dist;
  126. // Worst offender
  127. splitPoint = i;
  128. }
  129. }
  130. //const double half_total = E.array().sum()/2;
  131. //double run = 0;
  132. //for (int i = first + 1; i < last; i++)
  133. //{
  134. // run += E(i-(first+1));
  135. // if(run>half_total)
  136. // {
  137. // // When accumulated ½ the error --> more symmetric, but requires more
  138. // // curves
  139. // splitPoint = i;
  140. // break;
  141. // }
  142. //}
  143. return maxDist;
  144. };
  145. const static auto Straight = [](
  146. const Eigen::MatrixXd & d,
  147. const int first,
  148. const int last,
  149. const Eigen::RowVectorXd & tHat1,
  150. const Eigen::RowVectorXd & tHat2,
  151. Eigen::MatrixXd & bezCurve)
  152. {
  153. bezCurve.resize(4,d.cols());
  154. const double dist = (d.row(last)-d.row(first)).norm()/3.0;
  155. bezCurve.row(0) = d.row(first);
  156. bezCurve.row(1) = d.row(first) + tHat1*dist;
  157. bezCurve.row(2) = d.row(last) + tHat2*dist;
  158. bezCurve.row(3) = d.row(last);
  159. };
  160. const static auto GenerateBezier = [](
  161. const Eigen::MatrixXd & d,
  162. const int first,
  163. const int last,
  164. const Eigen::VectorXd & uPrime,
  165. const Eigen::RowVectorXd & tHat1,
  166. const Eigen::RowVectorXd & tHat2,
  167. Eigen::MatrixXd & bezCurve)
  168. {
  169. bezCurve.resize(4,d.cols());
  170. const int nPts = last - first + 1;
  171. const static auto B0 = [](const double u)->double
  172. { double tmp = 1.0 - u; return (tmp * tmp * tmp);};
  173. const static auto B1 = [](const double u)->double
  174. { double tmp = 1.0 - u; return (3 * u * (tmp * tmp));};
  175. const static auto B2 = [](const double u)->double
  176. { double tmp = 1.0 - u; return (3 * u * u * tmp); };
  177. const static auto B3 = [](const double u)->double
  178. { return (u * u * u); };
  179. /* Compute the A's */
  180. std::vector<std::vector<Eigen::RowVectorXd> > A(nPts);
  181. for (int i = 0; i < nPts; i++)
  182. {
  183. Eigen::RowVectorXd v1 = tHat1*B1(uPrime(i));
  184. Eigen::RowVectorXd v2 = tHat2*B2(uPrime(i));
  185. A[i] = {v1,v2};
  186. }
  187. /* Create the C and X matrices */
  188. Eigen::MatrixXd C(2,2);
  189. Eigen::VectorXd X(2);
  190. C(0,0) = 0.0;
  191. C(0,1) = 0.0;
  192. C(1,0) = 0.0;
  193. C(1,1) = 0.0;
  194. X(0) = 0.0;
  195. X(1) = 0.0;
  196. for( int i = 0; i < nPts; i++)
  197. {
  198. C(0,0) += A[i][0].dot(A[i][0]);
  199. C(0,1) += A[i][0].dot(A[i][1]);
  200. C(1,0) = C(0,1);
  201. C(1,1) += A[i][1].dot(A[i][1]);
  202. const Eigen::RowVectorXd tmp =
  203. d.row(first+i)-(
  204. d.row(first)*B0(uPrime(i))+
  205. d.row(first)*B1(uPrime(i))+
  206. d.row(last)*B2(uPrime(i))+
  207. d.row(last)*B3(uPrime(i)));
  208. X(0) += A[i][0].dot(tmp);
  209. X(1) += A[i][1].dot(tmp);
  210. }
  211. /* Compute the determinants of C and X */
  212. double det_C0_C1 = C(0,0) * C(1,1) - C(1,0) * C(0,1);
  213. const double det_C0_X = C(0,0) * X(1) - C(0,1) * X(0);
  214. const double det_X_C1 = X(0) * C(1,1) - X(1) * C(0,1);
  215. /* Finally, derive alpha values */
  216. if (det_C0_C1 == 0.0)
  217. {
  218. det_C0_C1 = (C(0,0) * C(1,1)) * 10e-12;
  219. }
  220. const double alpha_l = det_X_C1 / det_C0_C1;
  221. const double alpha_r = det_C0_X / det_C0_C1;
  222. /* If alpha negative, use the Wu/Barsky heuristic (see text) */
  223. /* (if alpha is 0, you get coincident control points that lead to
  224. * divide by zero in any subsequent NewtonRaphsonRootFind() call. */
  225. if (alpha_l < 1.0e-6 || alpha_r < 1.0e-6)
  226. {
  227. return Straight(d,first,last,tHat1,tHat2,bezCurve);
  228. }
  229. bezCurve.row(0) = d.row(first);
  230. bezCurve.row(1) = d.row(first) + tHat1*alpha_l;
  231. bezCurve.row(2) = d.row(last) + tHat2*alpha_r;
  232. bezCurve.row(3) = d.row(last);
  233. };
  234. const int maxIterations = 4;
  235. // This is a bad idea if error<1 ...
  236. //const double iterationError = error * error;
  237. const double iterationError = 100 * error;
  238. const int nPts = last - first + 1;
  239. /* Use heuristic if region only has two points in it */
  240. if(nPts == 2)
  241. {
  242. Eigen::MatrixXd bezCurve;
  243. Straight(d,first,last,tHat1,tHat2,bezCurve);
  244. cubics.push_back(bezCurve);
  245. return;
  246. }
  247. // ChordLengthParameterize
  248. Eigen::VectorXd u(last-first+1);
  249. u(0) = 0;
  250. for (int i = first+1; i <= last; i++)
  251. {
  252. u(i-first) = u(i-first-1) + (d.row(i)-d.row(i-1)).norm();
  253. }
  254. for (int i = first + 1; i <= last; i++)
  255. {
  256. u(i-first) = u(i-first) / u(last-first);
  257. }
  258. Eigen::MatrixXd bezCurve;
  259. GenerateBezier(d, first, last, u, tHat1, tHat2, bezCurve);
  260. int splitPoint;
  261. double maxError = ComputeMaxError(d, first, last, bezCurve, u, splitPoint);
  262. if (!force_split && maxError < error)
  263. {
  264. cubics.push_back(bezCurve);
  265. return;
  266. }
  267. /* If error not too large, try some reparameterization */
  268. /* and iteration */
  269. if (maxError < iterationError)
  270. {
  271. for (int i = 0; i < maxIterations; i++)
  272. {
  273. Eigen::VectorXd uPrime;
  274. // Reparameterize
  275. uPrime.resize(last-first+1);
  276. for (int i = first; i <= last; i++)
  277. {
  278. uPrime(i-first) = NewtonRaphsonRootFind(bezCurve, d.row(i), u(i- first));
  279. }
  280. GenerateBezier(d, first, last, uPrime, tHat1, tHat2, bezCurve);
  281. maxError = ComputeMaxError(d, first, last, bezCurve, uPrime, splitPoint);
  282. if (!force_split && maxError < error) {
  283. cubics.push_back(bezCurve);
  284. return;
  285. }
  286. u = uPrime;
  287. }
  288. }
  289. /* Fitting failed -- split at max error point and fit recursively */
  290. const Eigen::RowVectorXd tHatCenter =
  291. (d.row(splitPoint-1)-d.row(splitPoint+1)).normalized();
  292. //foobar
  293. fit_cubic_bezier_substring(
  294. d,first,splitPoint,tHat1,tHatCenter,error,false,cubics);
  295. fit_cubic_bezier_substring(
  296. d,splitPoint,last,(-tHatCenter).eval(),tHat2,error,false,cubics);
  297. }
  298. #ifdef IGL_STATIC_LIBRARY
  299. // Explicit template instantiation
  300. #endif