You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

648 lines
19 KiB

  1. #include "delaunator.hpp"
  2. #include <iostream>
  3. #include <algorithm>
  4. #include <assert.h>
  5. #include <cmath>
  6. #include <numeric>
  7. #include <limits>
  8. #include <stdexcept>
  9. #include <tuple>
  10. #include <vector>
  11. namespace delaunator {
  12. //@see https://stackoverflow.com/questions/33333363/built-in-mod-vs-custom-mod-function-improve-the-performance-of-modulus-op/33333636#33333636
  13. inline size_t fast_mod(const size_t i, const size_t c) {
  14. return i >= c ? i % c : i;
  15. }
  16. // Kahan and Babuska summation, Neumaier variant; accumulates less FP error
  17. inline double sum(const std::vector<double>& x) {
  18. double sum = x[0];
  19. double err = 0.0;
  20. for (size_t i = 1; i < x.size(); i++) {
  21. const double k = x[i];
  22. const double m = sum + k;
  23. err += std::fabs(sum) >= std::fabs(k) ? sum - m + k : k - m + sum;
  24. sum = m;
  25. }
  26. return sum + err;
  27. }
  28. inline double dist(
  29. const double ax,
  30. const double ay,
  31. const double bx,
  32. const double by) {
  33. const double dx = ax - bx;
  34. const double dy = ay - by;
  35. return dx * dx + dy * dy;
  36. }
  37. inline double circumradius(const Point& p1, const Point& p2, const Point& p3)
  38. {
  39. Point d = Point::vector(p1, p2);
  40. Point e = Point::vector(p1, p3);
  41. const double bl = d.magnitude2();
  42. const double cl = e.magnitude2();
  43. const double det = Point::determinant(d, e);
  44. Point radius((e.y() * bl - d.y() * cl) * 0.5 / det,
  45. (d.x() * cl - e.x() * bl) * 0.5 / det);
  46. if ((bl > 0.0 || bl < 0.0) &&
  47. (cl > 0.0 || cl < 0.0) &&
  48. (det > 0.0 || det < 0.0))
  49. return radius.magnitude2();
  50. return (std::numeric_limits<double>::max)();
  51. }
  52. inline double circumradius(
  53. const double ax,
  54. const double ay,
  55. const double bx,
  56. const double by,
  57. const double cx,
  58. const double cy) {
  59. const double dx = bx - ax;
  60. const double dy = by - ay;
  61. const double ex = cx - ax;
  62. const double ey = cy - ay;
  63. const double bl = dx * dx + dy * dy;
  64. const double cl = ex * ex + ey * ey;
  65. const double d = dx * ey - dy * ex;
  66. const double x = (ey * bl - dy * cl) * 0.5 / d;
  67. const double y = (dx * cl - ex * bl) * 0.5 / d;
  68. if ((bl > 0.0 || bl < 0.0) && (cl > 0.0 || cl < 0.0) && (d > 0.0 || d < 0.0)) {
  69. return x * x + y * y;
  70. } else {
  71. return (std::numeric_limits<double>::max)();
  72. }
  73. }
  74. inline bool clockwise(const Point& p0, const Point& p1, const Point& p2)
  75. {
  76. Point v0 = Point::vector(p0, p1);
  77. Point v1 = Point::vector(p0, p2);
  78. double det = Point::determinant(v0, v1);
  79. double dist = v0.magnitude2() + v1.magnitude2();
  80. double dist2 = Point::dist2(v0, v1);
  81. if (det == 0)
  82. {
  83. return false;
  84. }
  85. double reldet = std::abs(dist / det);
  86. if (reldet > 1e14)
  87. return false;
  88. return det < 0;
  89. }
  90. inline bool clockwise(double px, double py, double qx, double qy,
  91. double rx, double ry)
  92. {
  93. Point p0(px, py);
  94. Point p1(qx, qy);
  95. Point p2(rx, ry);
  96. return clockwise(p0, p1, p2);
  97. }
  98. inline bool counterclockwise(const Point& p0, const Point& p1, const Point& p2)
  99. {
  100. Point v0 = Point::vector(p0, p1);
  101. Point v1 = Point::vector(p0, p2);
  102. double det = Point::determinant(v0, v1);
  103. double dist = v0.magnitude2() + v1.magnitude2();
  104. double dist2 = Point::dist2(v0, v1);
  105. if (det == 0)
  106. return false;
  107. double reldet = std::abs(dist / det);
  108. if (reldet > 1e14)
  109. return false;
  110. return det > 0;
  111. }
  112. inline bool counterclockwise(double px, double py, double qx, double qy,
  113. double rx, double ry)
  114. {
  115. Point p0(px, py);
  116. Point p1(qx, qy);
  117. Point p2(rx, ry);
  118. return counterclockwise(p0, p1, p2);
  119. }
  120. inline Point circumcenter(
  121. const double ax,
  122. const double ay,
  123. const double bx,
  124. const double by,
  125. const double cx,
  126. const double cy) {
  127. const double dx = bx - ax;
  128. const double dy = by - ay;
  129. const double ex = cx - ax;
  130. const double ey = cy - ay;
  131. const double bl = dx * dx + dy * dy;
  132. const double cl = ex * ex + ey * ey;
  133. //ABELL - This is suspect for div-by-0.
  134. const double d = dx * ey - dy * ex;
  135. const double x = ax + (ey * bl - dy * cl) * 0.5 / d;
  136. const double y = ay + (dx * cl - ex * bl) * 0.5 / d;
  137. return Point(x, y);
  138. }
  139. inline bool in_circle(
  140. const double ax,
  141. const double ay,
  142. const double bx,
  143. const double by,
  144. const double cx,
  145. const double cy,
  146. const double px,
  147. const double py) {
  148. const double dx = ax - px;
  149. const double dy = ay - py;
  150. const double ex = bx - px;
  151. const double ey = by - py;
  152. const double fx = cx - px;
  153. const double fy = cy - py;
  154. const double ap = dx * dx + dy * dy;
  155. const double bp = ex * ex + ey * ey;
  156. const double cp = fx * fx + fy * fy;
  157. return (dx * (ey * cp - bp * fy) -
  158. dy * (ex * cp - bp * fx) +
  159. ap * (ex * fy - ey * fx)) < 0.0;
  160. }
  161. constexpr double EPSILON = std::numeric_limits<double>::epsilon();
  162. inline bool check_pts_equal(double x1, double y1, double x2, double y2) {
  163. return std::fabs(x1 - x2) <= EPSILON &&
  164. std::fabs(y1 - y2) <= EPSILON;
  165. }
  166. // monotonically increases with real angle, but doesn't need expensive trigonometry
  167. inline double pseudo_angle(const double dx, const double dy) {
  168. const double p = dx / (std::abs(dx) + std::abs(dy));
  169. return (dy > 0.0 ? 3.0 - p : 1.0 + p) / 4.0; // [0..1)
  170. }
  171. Delaunator::Delaunator(std::vector<double> const& in_coords)
  172. : coords(in_coords), m_points(in_coords)
  173. {
  174. std::size_t n = coords.size() >> 1;
  175. std::vector<std::size_t> ids(n);
  176. std::iota(ids.begin(), ids.end(), 0);
  177. double max_x = std::numeric_limits<double>::lowest();
  178. double max_y = std::numeric_limits<double>::lowest();
  179. double min_x = (std::numeric_limits<double>::max)();
  180. double min_y = (std::numeric_limits<double>::max)();
  181. for (const Point& p : m_points)
  182. {
  183. min_x = std::min(p.x(), min_x);
  184. min_y = std::min(p.y(), min_y);
  185. max_x = std::max(p.x(), max_x);
  186. max_y = std::max(p.y(), max_y);
  187. }
  188. double width = max_x - min_x;
  189. double height = max_y - min_y;
  190. double span = width * width + height * height; // Everything is square dist.
  191. Point center((min_x + max_x) / 2, (min_y + max_y) / 2);
  192. std::size_t i0 = INVALID_INDEX;
  193. std::size_t i1 = INVALID_INDEX;
  194. std::size_t i2 = INVALID_INDEX;
  195. // pick a seed point close to the centroid
  196. double min_dist = (std::numeric_limits<double>::max)();
  197. for (size_t i = 0; i < m_points.size(); ++i)
  198. {
  199. const Point& p = m_points[i];
  200. const double d = Point::dist2(center, p);
  201. if (d < min_dist) {
  202. i0 = i;
  203. min_dist = d;
  204. }
  205. }
  206. const Point& p0 = m_points[i0];
  207. min_dist = (std::numeric_limits<double>::max)();
  208. // find the point closest to the seed
  209. for (std::size_t i = 0; i < n; i++) {
  210. if (i == i0) continue;
  211. const double d = Point::dist2(p0, m_points[i]);
  212. if (d < min_dist && d > 0.0) {
  213. i1 = i;
  214. min_dist = d;
  215. }
  216. }
  217. const Point& p1 = m_points[i1];
  218. double min_radius = (std::numeric_limits<double>::max)();
  219. // find the third point which forms the smallest circumcircle
  220. // with the first two
  221. for (std::size_t i = 0; i < n; i++) {
  222. if (i == i0 || i == i1) continue;
  223. const double r = circumradius(p0, p1, m_points[i]);
  224. if (r < min_radius) {
  225. i2 = i;
  226. min_radius = r;
  227. }
  228. }
  229. if (!(min_radius < (std::numeric_limits<double>::max)())) {
  230. throw std::runtime_error("not triangulation");
  231. }
  232. const Point& p2 = m_points[i2];
  233. if (counterclockwise(p0, p1, p2))
  234. std::swap(i1, i2);
  235. double i0x = p0.x();
  236. double i0y = p0.y();
  237. double i1x = m_points[i1].x();
  238. double i1y = m_points[i1].y();
  239. double i2x = m_points[i2].x();
  240. double i2y = m_points[i2].y();
  241. m_center = circumcenter(i0x, i0y, i1x, i1y, i2x, i2y);
  242. // Calculate the distances from the center once to avoid having to
  243. // calculate for each compare. This used to be done in the comparator,
  244. // but GCC 7.5+ would copy the comparator to iterators used in the
  245. // sort, and this was excruciatingly slow when there were many points
  246. // because you had to copy the vector of distances.
  247. std::vector<double> dists;
  248. dists.reserve(m_points.size());
  249. for (const Point& p : m_points)
  250. dists.push_back(dist(p.x(), p.y(), m_center.x(), m_center.y()));
  251. // sort the points by distance from the seed triangle circumcenter
  252. std::sort(ids.begin(), ids.end(),
  253. [&dists](std::size_t i, std::size_t j)
  254. { return dists[i] < dists[j]; });
  255. // initialize a hash table for storing edges of the advancing convex hull
  256. m_hash_size = static_cast<std::size_t>(std::ceil(std::sqrt(n)));
  257. m_hash.resize(m_hash_size);
  258. std::fill(m_hash.begin(), m_hash.end(), INVALID_INDEX);
  259. // initialize arrays for tracking the edges of the advancing convex hull
  260. hull_prev.resize(n);
  261. hull_next.resize(n);
  262. hull_tri.resize(n);
  263. hull_start = i0;
  264. size_t hull_size = 3;
  265. hull_next[i0] = hull_prev[i2] = i1;
  266. hull_next[i1] = hull_prev[i0] = i2;
  267. hull_next[i2] = hull_prev[i1] = i0;
  268. hull_tri[i0] = 0;
  269. hull_tri[i1] = 1;
  270. hull_tri[i2] = 2;
  271. m_hash[hash_key(i0x, i0y)] = i0;
  272. m_hash[hash_key(i1x, i1y)] = i1;
  273. m_hash[hash_key(i2x, i2y)] = i2;
  274. // ABELL - Why are we doing this is n < 3? There is no triangulation if
  275. // there is no triangle.
  276. std::size_t max_triangles = n < 3 ? 1 : 2 * n - 5;
  277. triangles.reserve(max_triangles * 3);
  278. halfedges.reserve(max_triangles * 3);
  279. add_triangle(i0, i1, i2, INVALID_INDEX, INVALID_INDEX, INVALID_INDEX);
  280. double xp = std::numeric_limits<double>::quiet_NaN();
  281. double yp = std::numeric_limits<double>::quiet_NaN();
  282. // Go through points based on distance from the center.
  283. for (std::size_t k = 0; k < n; k++) {
  284. const std::size_t i = ids[k];
  285. const double x = coords[2 * i];
  286. const double y = coords[2 * i + 1];
  287. // skip near-duplicate points
  288. if (k > 0 && check_pts_equal(x, y, xp, yp))
  289. continue;
  290. xp = x;
  291. yp = y;
  292. //ABELL - This is dumb. We have the indices. Use them.
  293. // skip seed triangle points
  294. if (check_pts_equal(x, y, i0x, i0y) ||
  295. check_pts_equal(x, y, i1x, i1y) ||
  296. check_pts_equal(x, y, i2x, i2y)) continue;
  297. // find a visible edge on the convex hull using edge hash
  298. std::size_t start = 0;
  299. size_t key = hash_key(x, y);
  300. for (size_t j = 0; j < m_hash_size; j++) {
  301. start = m_hash[fast_mod(key + j, m_hash_size)];
  302. // ABELL - Not sure how hull_next[start] could ever equal start
  303. // I *think* hull_next is just a representation of the hull in one
  304. // direction.
  305. if (start != INVALID_INDEX && start != hull_next[start])
  306. break;
  307. }
  308. //ABELL
  309. // Make sure what we found is on the hull.
  310. assert(hull_prev[start] != start);
  311. assert(hull_prev[start] != INVALID_INDEX);
  312. start = hull_prev[start];
  313. size_t e = start;
  314. size_t q;
  315. // Advance until we find a place in the hull where our current point
  316. // can be added.
  317. while (true)
  318. {
  319. q = hull_next[e];
  320. if (Point::equal(m_points[i], m_points[e], span) ||
  321. Point::equal(m_points[i], m_points[q], span))
  322. {
  323. e = INVALID_INDEX;
  324. break;
  325. }
  326. if (counterclockwise(x, y, coords[2 * e], coords[2 * e + 1],
  327. coords[2 * q], coords[2 * q + 1]))
  328. break;
  329. e = q;
  330. if (e == start) {
  331. e = INVALID_INDEX;
  332. break;
  333. }
  334. }
  335. // ABELL
  336. // This seems wrong. Perhaps we should check what's going on?
  337. if (e == INVALID_INDEX) // likely a near-duplicate point; skip it
  338. continue;
  339. // add the first triangle from the point
  340. std::size_t t = add_triangle(
  341. e,
  342. i,
  343. hull_next[e],
  344. INVALID_INDEX,
  345. INVALID_INDEX,
  346. hull_tri[e]);
  347. hull_tri[i] = legalize(t + 2); // Legalize the triangle we just added.
  348. hull_tri[e] = t;
  349. hull_size++;
  350. // walk forward through the hull, adding more triangles and
  351. // flipping recursively
  352. std::size_t next = hull_next[e];
  353. while (true)
  354. {
  355. q = hull_next[next];
  356. if (!counterclockwise(x, y, coords[2 * next], coords[2 * next + 1],
  357. coords[2 * q], coords[2 * q + 1]))
  358. break;
  359. t = add_triangle(next, i, q,
  360. hull_tri[i], INVALID_INDEX, hull_tri[next]);
  361. hull_tri[i] = legalize(t + 2);
  362. hull_next[next] = next; // mark as removed
  363. hull_size--;
  364. next = q;
  365. }
  366. // walk backward from the other side, adding more triangles and flipping
  367. if (e == start) {
  368. while (true)
  369. {
  370. q = hull_prev[e];
  371. if (!counterclockwise(x, y, coords[2 * q], coords[2 * q + 1],
  372. coords[2 * e], coords[2 * e + 1]))
  373. break;
  374. t = add_triangle(q, i, e,
  375. INVALID_INDEX, hull_tri[e], hull_tri[q]);
  376. legalize(t + 2);
  377. hull_tri[q] = t;
  378. hull_next[e] = e; // mark as removed
  379. hull_size--;
  380. e = q;
  381. }
  382. }
  383. // update the hull indices
  384. hull_prev[i] = e;
  385. hull_start = e;
  386. hull_prev[next] = i;
  387. hull_next[e] = i;
  388. hull_next[i] = next;
  389. m_hash[hash_key(x, y)] = i;
  390. m_hash[hash_key(coords[2 * e], coords[2 * e + 1])] = e;
  391. }
  392. }
  393. double Delaunator::get_hull_area()
  394. {
  395. std::vector<double> hull_area;
  396. size_t e = hull_start;
  397. size_t cnt = 1;
  398. do {
  399. hull_area.push_back((coords[2 * e] - coords[2 * hull_prev[e]]) *
  400. (coords[2 * e + 1] + coords[2 * hull_prev[e] + 1]));
  401. cnt++;
  402. e = hull_next[e];
  403. } while (e != hull_start);
  404. return sum(hull_area);
  405. }
  406. double Delaunator::get_triangle_area()
  407. {
  408. std::vector<double> vals;
  409. for (size_t i = 0; i < triangles.size(); i += 3)
  410. {
  411. const double ax = coords[2 * triangles[i]];
  412. const double ay = coords[2 * triangles[i] + 1];
  413. const double bx = coords[2 * triangles[i + 1]];
  414. const double by = coords[2 * triangles[i + 1] + 1];
  415. const double cx = coords[2 * triangles[i + 2]];
  416. const double cy = coords[2 * triangles[i + 2] + 1];
  417. double val = std::fabs((by - ay) * (cx - bx) - (bx - ax) * (cy - by));
  418. vals.push_back(val);
  419. }
  420. return sum(vals);
  421. }
  422. std::size_t Delaunator::legalize(std::size_t a) {
  423. std::size_t i = 0;
  424. std::size_t ar = 0;
  425. m_edge_stack.clear();
  426. // recursion eliminated with a fixed-size stack
  427. while (true) {
  428. const size_t b = halfedges[a];
  429. /* if the pair of triangles doesn't satisfy the Delaunay condition
  430. * (p1 is inside the circumcircle of [p0, pl, pr]), flip them,
  431. * then do the same check/flip recursively for the new pair of triangles
  432. *
  433. * pl pl
  434. * /||\ / \
  435. * al/ || \bl al/ \a
  436. * / || \ / \
  437. * / a||b \ flip /___ar___\
  438. * p0\ || /p1 => p0\---bl---/p1
  439. * \ || / \ /
  440. * ar\ || /br b\ /br
  441. * \||/ \ /
  442. * pr pr
  443. */
  444. const size_t a0 = 3 * (a / 3);
  445. ar = a0 + (a + 2) % 3;
  446. if (b == INVALID_INDEX) {
  447. if (i > 0) {
  448. i--;
  449. a = m_edge_stack[i];
  450. continue;
  451. } else {
  452. //i = INVALID_INDEX;
  453. break;
  454. }
  455. }
  456. const size_t b0 = 3 * (b / 3);
  457. const size_t al = a0 + (a + 1) % 3;
  458. const size_t bl = b0 + (b + 2) % 3;
  459. const std::size_t p0 = triangles[ar];
  460. const std::size_t pr = triangles[a];
  461. const std::size_t pl = triangles[al];
  462. const std::size_t p1 = triangles[bl];
  463. const bool illegal = in_circle(
  464. coords[2 * p0],
  465. coords[2 * p0 + 1],
  466. coords[2 * pr],
  467. coords[2 * pr + 1],
  468. coords[2 * pl],
  469. coords[2 * pl + 1],
  470. coords[2 * p1],
  471. coords[2 * p1 + 1]);
  472. if (illegal) {
  473. triangles[a] = p1;
  474. triangles[b] = p0;
  475. auto hbl = halfedges[bl];
  476. // Edge swapped on the other side of the hull (rare).
  477. // Fix the halfedge reference
  478. if (hbl == INVALID_INDEX) {
  479. std::size_t e = hull_start;
  480. do {
  481. if (hull_tri[e] == bl) {
  482. hull_tri[e] = a;
  483. break;
  484. }
  485. e = hull_prev[e];
  486. } while (e != hull_start);
  487. }
  488. link(a, hbl);
  489. link(b, halfedges[ar]);
  490. link(ar, bl);
  491. std::size_t br = b0 + (b + 1) % 3;
  492. if (i < m_edge_stack.size()) {
  493. m_edge_stack[i] = br;
  494. } else {
  495. m_edge_stack.push_back(br);
  496. }
  497. i++;
  498. } else {
  499. if (i > 0) {
  500. i--;
  501. a = m_edge_stack[i];
  502. continue;
  503. } else {
  504. break;
  505. }
  506. }
  507. }
  508. return ar;
  509. }
  510. std::size_t Delaunator::hash_key(const double x, const double y) const {
  511. const double dx = x - m_center.x();
  512. const double dy = y - m_center.y();
  513. return fast_mod(
  514. static_cast<std::size_t>(std::llround(std::floor(pseudo_angle(dx, dy) * static_cast<double>(m_hash_size)))),
  515. m_hash_size);
  516. }
  517. std::size_t Delaunator::add_triangle(
  518. std::size_t i0,
  519. std::size_t i1,
  520. std::size_t i2,
  521. std::size_t a,
  522. std::size_t b,
  523. std::size_t c) {
  524. std::size_t t = triangles.size();
  525. triangles.push_back(i0);
  526. triangles.push_back(i1);
  527. triangles.push_back(i2);
  528. link(t, a);
  529. link(t + 1, b);
  530. link(t + 2, c);
  531. return t;
  532. }
  533. void Delaunator::link(const std::size_t a, const std::size_t b) {
  534. std::size_t s = halfedges.size();
  535. if (a == s) {
  536. halfedges.push_back(b);
  537. } else if (a < s) {
  538. halfedges[a] = b;
  539. } else {
  540. throw std::runtime_error("Cannot link edge");
  541. }
  542. if (b != INVALID_INDEX) {
  543. std::size_t s2 = halfedges.size();
  544. if (b == s2) {
  545. halfedges.push_back(a);
  546. } else if (b < s2) {
  547. halfedges[b] = a;
  548. } else {
  549. throw std::runtime_error("Cannot link edge");
  550. }
  551. }
  552. }
  553. } //namespace delaunator