PolyFEM
Loading...
Searching...
No Matches
GeometryUtils.cpp
Go to the documentation of this file.
1#include "GeometryUtils.hpp"
2
4
5#include <Eigen/Geometry>
6
7#include <igl/barycentric_coordinates.h>
8
9#include <ipc/distance/point_edge.hpp>
10#include <ipc/distance/point_triangle.hpp>
11
12namespace polyfem::utils
13{
15 const Eigen::Vector2d &a,
16 const Eigen::Vector2d &b,
17 const Eigen::Vector2d &c)
18 {
19 return ((b.x() - a.x()) * (c.y() - a.y()) - (c.x() - a.x()) * (b.y() - a.y())) / 2.0;
20 }
21
23 const Eigen::Vector3d &a,
24 const Eigen::Vector3d &b,
25 const Eigen::Vector3d &c)
26 {
27 return (b - a).cross(c - a).norm() / 2.0;
28 }
29
30 double triangle_area(const Eigen::MatrixXd V)
31 {
32 assert(V.rows() == 3);
33 if (V.cols() == 2)
34 {
35 return triangle_area_2D(V.row(0), V.row(1), V.row(2));
36 }
37 else
38 {
39 return triangle_area_3D(V.row(0), V.row(1), V.row(2));
40 }
41 }
42
44 const Eigen::Vector3d &a,
45 const Eigen::Vector3d &b,
46 const Eigen::Vector3d &c,
47 const Eigen::Vector3d &d)
48 {
49 return (b - a).cross(c - a).dot(d - a) / 6.0;
50 }
51
52 double tetrahedron_volume(const Eigen::MatrixXd V)
53 {
54 assert(V.rows() == 4 && V.cols() == 3);
55 return tetrahedron_volume(V.row(0), V.row(1), V.row(2), V.row(3));
56 }
57
58 Eigen::MatrixXd triangle_to_clockwise_order(const Eigen::MatrixXd &triangle)
59 {
60 // Only works for 2D triangles.
61 assert(triangle.rows() == 3 && triangle.cols() == 2);
62
63 if (triangle_area(triangle) <= 0)
64 return triangle; // triangle aleady in clockwise order.
65
66 Eigen::Matrix<double, 3, 2> triangle_clockwise;
67 triangle_clockwise.row(0) = triangle.row(2);
68 triangle_clockwise.row(1) = triangle.row(1);
69 triangle_clockwise.row(2) = triangle.row(0);
70
71 return triangle_clockwise;
72 }
73
74 std::vector<Eigen::MatrixXd> triangle_fan(const Eigen::MatrixXd &convex_polygon)
75 {
76 assert(convex_polygon.rows() >= 3);
77 std::vector<Eigen::MatrixXd> triangles;
78 for (int i = 1; i < convex_polygon.rows() - 1; ++i)
79 {
80 triangles.emplace_back(3, convex_polygon.cols());
81 triangles.back().row(0) = convex_polygon.row(0);
82 triangles.back().row(1) = convex_polygon.row(i);
83 triangles.back().row(2) = convex_polygon.row(i + 1);
84 }
85 return triangles;
86 }
87
88 Eigen::VectorXd barycentric_coordinates(
89 const Eigen::VectorXd &p,
90 const Eigen::MatrixXd &V)
91 {
92 Eigen::MatrixXd A(p.size() + 1, p.size() + 1);
93 A.topRows(p.size()) = V.transpose();
94 A.bottomRows(1).setOnes();
95
96 Eigen::VectorXd rhs(p.size() + 1);
97 rhs.topRows(p.size()) = p;
98 rhs.bottomRows(1).setOnes();
99
100 // TODO: Can we use better than LU?
101 const Eigen::VectorXd bc = A.partialPivLu().solve(rhs);
102 assert((A * bc - rhs).norm() / rhs.norm() < 1e-12);
103
104 return bc;
105 }
106
108 const Eigen::Vector2d &t0,
109 const Eigen::Vector2d &t1,
110 const Eigen::Vector2d &t2,
111 const Eigen::Vector2d &center,
112 const double radius)
113 {
114 Eigen::RowVector3d bc;
115 igl::barycentric_coordinates(
116 center.transpose(), t0.transpose(), t1.transpose(), t2.transpose(),
117 bc);
118
119 if (bc.minCoeff() >= 0) // point is inside triangle
120 {
121 assert(bc.maxCoeff() <= 1 + 1e-12); // bc.sum() == 1
122 return true;
123 }
124
125 const double radius_sqr = radius * radius;
126
127 if (ipc::point_edge_distance(center, t0, t1) <= radius_sqr
128 || ipc::point_edge_distance(center, t1, t2) <= radius_sqr
129 || ipc::point_edge_distance(center, t2, t0) <= radius_sqr)
130 {
131 return true;
132 }
133
134 return false;
135 }
136
138 const Eigen::Vector3d &t0,
139 const Eigen::Vector3d &t1,
140 const Eigen::Vector3d &t2,
141 const Eigen::Vector3d &t3,
142 const Eigen::Vector3d &center,
143 const double radius)
144 {
145 Eigen::RowVector4d bc;
146 igl::barycentric_coordinates(
147 center.transpose(), t0.transpose(), t1.transpose(),
148 t2.transpose(), t3.transpose(), bc);
149
150 if (bc.minCoeff() >= 0)
151 {
152 assert(bc.maxCoeff() <= 1 + 1e-12); // bc.sum() == 1
153 return true;
154 }
155
156 const double radius_sqr = radius * radius;
157
158 if (ipc::point_triangle_distance(center, t0, t1, t2) <= radius_sqr
159 || ipc::point_triangle_distance(center, t0, t1, t3) <= radius_sqr
160 || ipc::point_triangle_distance(center, t0, t2, t3) <= radius_sqr
161 || ipc::point_triangle_distance(center, t1, t2, t3) <= radius_sqr)
162 {
163 return true;
164 }
165
166 return false;
167 }
168
170 const VectorNd &ea0,
171 const VectorNd &ea1,
172 const VectorNd &eb0,
173 const VectorNd &eb1,
174 const double tol)
175 {
176 assert((ea0 - ea1).norm() != 0 && (eb1 - eb0).norm() != 0);
177 const VectorNd ea = (ea1 - ea0).normalized();
178 const VectorNd eb = (eb1 - eb0).normalized();
179 return abs(ea.dot(eb)) > 1 - tol;
180 }
181
183 const Eigen::Vector3d &t00,
184 const Eigen::Vector3d &t01,
185 const Eigen::Vector3d &t02,
186 const Eigen::Vector3d &t10,
187 const Eigen::Vector3d &t11,
188 const Eigen::Vector3d &t12,
189 const double tol)
190 {
191 const Eigen::Vector3d n0 = (t01 - t00).cross(t02 - t00).normalized();
192 const Eigen::Vector3d n1 = (t11 - t10).cross(t12 - t10).normalized();
193 return abs(n0.dot(n1)) > 1 - tol;
194 }
195
197 const VectorNd &aabb0_min,
198 const VectorNd &aabb0_max,
199 const VectorNd &aabb1_min,
200 const VectorNd &aabb1_max)
201 {
202 return (aabb0_min.array() <= aabb1_max.array()).all() && (aabb1_min.array() <= aabb0_max.array()).all();
203 }
204
205 // =========================================================================
206 // The following automatically generated using SymPy
207 // =========================================================================
208
209 void triangle_area_2D_gradient(double ax, double ay, double bx, double by, double cx, double cy, double g[6])
210 {
211 const auto t0 = -cy;
212 const auto t1 = -cx;
213 g[0] = (1.0 / 2.0) * by + (1.0 / 2.0) * t0;
214 g[1] = -1.0 / 2.0 * bx - 1.0 / 2.0 * t1;
215 g[2] = -1.0 / 2.0 * ay - 1.0 / 2.0 * t0;
216 g[3] = (1.0 / 2.0) * ax + (1.0 / 2.0) * t1;
217 g[4] = (1.0 / 2.0) * (ay - by);
218 g[5] = (1.0 / 2.0) * (-ax + bx);
219 }
220
221 void triangle_area_2D_hessian(double ax, double ay, double bx, double by, double cx, double cy, double H[36])
222 {
223 H[0] = 0;
224 H[1] = 0;
225 H[2] = 0;
226 H[3] = 1.0 / 2.0;
227 H[4] = 0;
228 H[5] = -1.0 / 2.0;
229 H[6] = 0;
230 H[7] = 0;
231 H[8] = -1.0 / 2.0;
232 H[9] = 0;
233 H[10] = 1.0 / 2.0;
234 H[11] = 0;
235 H[12] = 0;
236 H[13] = -1.0 / 2.0;
237 H[14] = 0;
238 H[15] = 0;
239 H[16] = 0;
240 H[17] = 1.0 / 2.0;
241 H[18] = 1.0 / 2.0;
242 H[19] = 0;
243 H[20] = 0;
244 H[21] = 0;
245 H[22] = -1.0 / 2.0;
246 H[23] = 0;
247 H[24] = 0;
248 H[25] = 1.0 / 2.0;
249 H[26] = 0;
250 H[27] = -1.0 / 2.0;
251 H[28] = 0;
252 H[29] = 0;
253 H[30] = -1.0 / 2.0;
254 H[31] = 0;
255 H[32] = 1.0 / 2.0;
256 H[33] = 0;
257 H[34] = 0;
258 H[35] = 0;
259 }
260
261 void tetrahedron_volume_gradient(double ax, double ay, double az, double bx, double by, double bz, double cx, double cy, double cz, double dx, double dy, double dz, double g[12])
262 {
263 const auto t0 = az - dz;
264 const auto t1 = -cy;
265 const auto t2 = by + t1;
266 const auto t3 = ay - dy;
267 const auto t4 = -cz;
268 const auto t5 = bz + t4;
269 const auto t6 = ay - by;
270 const auto t7 = az + t4;
271 const auto t8 = ay + t1;
272 const auto t9 = az - bz;
273 const auto t10 = t6 * t7 - t8 * t9;
274 const auto t11 = -cx;
275 const auto t12 = bx + t11;
276 const auto t13 = ax - dx;
277 const auto t14 = ax - bx;
278 const auto t15 = ax + t11;
279 const auto t16 = t14 * t7 - t15 * t9;
280 const auto t17 = t14 * t8 - t15 * t6;
281 g[0] = -1.0 / 6.0 * t0 * t2 - 1.0 / 6.0 * t10 + (1.0 / 6.0) * t3 * t5;
282 g[1] = (1.0 / 6.0) * t0 * t12 - 1.0 / 6.0 * t13 * t5 + (1.0 / 6.0) * t16;
283 g[2] = -1.0 / 6.0 * t12 * t3 + (1.0 / 6.0) * t13 * t2 - 1.0 / 6.0 * t17;
284 g[3] = (1.0 / 6.0) * t0 * t8 - 1.0 / 6.0 * t3 * t7;
285 g[4] = -1.0 / 6.0 * t0 * t15 + (1.0 / 6.0) * t13 * t7;
286 g[5] = -1.0 / 6.0 * t13 * t8 + (1.0 / 6.0) * t15 * t3;
287 g[6] = -1.0 / 6.0 * t0 * t6 + (1.0 / 6.0) * t3 * t9;
288 g[7] = (1.0 / 6.0) * t0 * t14 - 1.0 / 6.0 * t13 * t9;
289 g[8] = (1.0 / 6.0) * t13 * t6 - 1.0 / 6.0 * t14 * t3;
290 g[9] = (1.0 / 6.0) * t10;
291 g[10] = -1.0 / 6.0 * t16;
292 g[11] = (1.0 / 6.0) * t17;
293 }
294
295 void tetrahedron_volume_hessian(double ax, double ay, double az, double bx, double by, double bz, double cx, double cy, double cz, double dx, double dy, double dz, double H[144])
296 {
297 const auto t0 = -dz;
298 const auto t1 = cz + t0;
299 const auto t2 = -1.0 / 6.0 * t1;
300 const auto t3 = -dy;
301 const auto t4 = cy + t3;
302 const auto t5 = (1.0 / 6.0) * t4;
303 const auto t6 = bz + t0;
304 const auto t7 = (1.0 / 6.0) * t6;
305 const auto t8 = by + t3;
306 const auto t9 = -1.0 / 6.0 * t8;
307 const auto t10 = -cz;
308 const auto t11 = bz + t10;
309 const auto t12 = -1.0 / 6.0 * t11;
310 const auto t13 = -cy;
311 const auto t14 = by + t13;
312 const auto t15 = (1.0 / 6.0) * t14;
313 const auto t16 = (1.0 / 6.0) * t1;
314 const auto t17 = -dx;
315 const auto t18 = cx + t17;
316 const auto t19 = -1.0 / 6.0 * t18;
317 const auto t20 = -1.0 / 6.0 * t6;
318 const auto t21 = bx + t17;
319 const auto t22 = (1.0 / 6.0) * t21;
320 const auto t23 = (1.0 / 6.0) * t11;
321 const auto t24 = -cx;
322 const auto t25 = bx + t24;
323 const auto t26 = -1.0 / 6.0 * t25;
324 const auto t27 = -1.0 / 6.0 * t4;
325 const auto t28 = (1.0 / 6.0) * t18;
326 const auto t29 = (1.0 / 6.0) * t8;
327 const auto t30 = -1.0 / 6.0 * t21;
328 const auto t31 = -1.0 / 6.0 * t14;
329 const auto t32 = (1.0 / 6.0) * t25;
330 const auto t33 = az + t0;
331 const auto t34 = -1.0 / 6.0 * t33;
332 const auto t35 = ay + t3;
333 const auto t36 = (1.0 / 6.0) * t35;
334 const auto t37 = az + t10;
335 const auto t38 = (1.0 / 6.0) * t37;
336 const auto t39 = ay + t13;
337 const auto t40 = -1.0 / 6.0 * t39;
338 const auto t41 = (1.0 / 6.0) * t33;
339 const auto t42 = ax + t17;
340 const auto t43 = -1.0 / 6.0 * t42;
341 const auto t44 = -1.0 / 6.0 * t37;
342 const auto t45 = ax + t24;
343 const auto t46 = (1.0 / 6.0) * t45;
344 const auto t47 = -1.0 / 6.0 * t35;
345 const auto t48 = (1.0 / 6.0) * t42;
346 const auto t49 = (1.0 / 6.0) * t39;
347 const auto t50 = -1.0 / 6.0 * t45;
348 const auto t51 = az - bz;
349 const auto t52 = -1.0 / 6.0 * t51;
350 const auto t53 = ay - by;
351 const auto t54 = (1.0 / 6.0) * t53;
352 const auto t55 = (1.0 / 6.0) * t51;
353 const auto t56 = ax - bx;
354 const auto t57 = -1.0 / 6.0 * t56;
355 const auto t58 = -1.0 / 6.0 * t53;
356 const auto t59 = (1.0 / 6.0) * t56;
357 H[0] = 0;
358 H[1] = 0;
359 H[2] = 0;
360 H[3] = 0;
361 H[4] = t2;
362 H[5] = t5;
363 H[6] = 0;
364 H[7] = t7;
365 H[8] = t9;
366 H[9] = 0;
367 H[10] = t12;
368 H[11] = t15;
369 H[12] = 0;
370 H[13] = 0;
371 H[14] = 0;
372 H[15] = t16;
373 H[16] = 0;
374 H[17] = t19;
375 H[18] = t20;
376 H[19] = 0;
377 H[20] = t22;
378 H[21] = t23;
379 H[22] = 0;
380 H[23] = t26;
381 H[24] = 0;
382 H[25] = 0;
383 H[26] = 0;
384 H[27] = t27;
385 H[28] = t28;
386 H[29] = 0;
387 H[30] = t29;
388 H[31] = t30;
389 H[32] = 0;
390 H[33] = t31;
391 H[34] = t32;
392 H[35] = 0;
393 H[36] = 0;
394 H[37] = t16;
395 H[38] = t27;
396 H[39] = 0;
397 H[40] = 0;
398 H[41] = 0;
399 H[42] = 0;
400 H[43] = t34;
401 H[44] = t36;
402 H[45] = 0;
403 H[46] = t38;
404 H[47] = t40;
405 H[48] = t2;
406 H[49] = 0;
407 H[50] = t28;
408 H[51] = 0;
409 H[52] = 0;
410 H[53] = 0;
411 H[54] = t41;
412 H[55] = 0;
413 H[56] = t43;
414 H[57] = t44;
415 H[58] = 0;
416 H[59] = t46;
417 H[60] = t5;
418 H[61] = t19;
419 H[62] = 0;
420 H[63] = 0;
421 H[64] = 0;
422 H[65] = 0;
423 H[66] = t47;
424 H[67] = t48;
425 H[68] = 0;
426 H[69] = t49;
427 H[70] = t50;
428 H[71] = 0;
429 H[72] = 0;
430 H[73] = t20;
431 H[74] = t29;
432 H[75] = 0;
433 H[76] = t41;
434 H[77] = t47;
435 H[78] = 0;
436 H[79] = 0;
437 H[80] = 0;
438 H[81] = 0;
439 H[82] = t52;
440 H[83] = t54;
441 H[84] = t7;
442 H[85] = 0;
443 H[86] = t30;
444 H[87] = t34;
445 H[88] = 0;
446 H[89] = t48;
447 H[90] = 0;
448 H[91] = 0;
449 H[92] = 0;
450 H[93] = t55;
451 H[94] = 0;
452 H[95] = t57;
453 H[96] = t9;
454 H[97] = t22;
455 H[98] = 0;
456 H[99] = t36;
457 H[100] = t43;
458 H[101] = 0;
459 H[102] = 0;
460 H[103] = 0;
461 H[104] = 0;
462 H[105] = t58;
463 H[106] = t59;
464 H[107] = 0;
465 H[108] = 0;
466 H[109] = t23;
467 H[110] = t31;
468 H[111] = 0;
469 H[112] = t44;
470 H[113] = t49;
471 H[114] = 0;
472 H[115] = t55;
473 H[116] = t58;
474 H[117] = 0;
475 H[118] = 0;
476 H[119] = 0;
477 H[120] = t12;
478 H[121] = 0;
479 H[122] = t32;
480 H[123] = t38;
481 H[124] = 0;
482 H[125] = t50;
483 H[126] = t52;
484 H[127] = 0;
485 H[128] = t59;
486 H[129] = 0;
487 H[130] = 0;
488 H[131] = 0;
489 H[132] = t15;
490 H[133] = t26;
491 H[134] = 0;
492 H[135] = t40;
493 H[136] = t46;
494 H[137] = 0;
495 H[138] = t54;
496 H[139] = t57;
497 H[140] = 0;
498 H[141] = 0;
499 H[142] = 0;
500 H[143] = 0;
501 }
502} // namespace polyfem::utils
int V
Eigen::Matrix< double, dim, 1 > cross(const Eigen::Matrix< double, dim, 1 > &x, const Eigen::Matrix< double, dim, 1 > &y)
bool tetrahedron_intersects_ball(const Eigen::Vector3d &t0, const Eigen::Vector3d &t1, const Eigen::Vector3d &t2, const Eigen::Vector3d &t3, const Eigen::Vector3d &center, const double radius)
Determine if a 3D tetrahedron intersects a 3D ball.
bool are_aabbs_intersecting(const VectorNd &aabb0_min, const VectorNd &aabb0_max, const VectorNd &aabb1_min, const VectorNd &aabb1_max)
Determine if two axis-aligned bounding boxes intersect.
void triangle_area_2D_hessian(double ax, double ay, double bx, double by, double cx, double cy, double H[36])
Compute the Hessian of the signed area of a 2D triangle defined by three points.
void tetrahedron_volume_gradient(double ax, double ay, double az, double bx, double by, double bz, double cx, double cy, double cz, double dx, double dy, double dz, double g[12])
Compute the gradient of the signed volume of a tetrahedron defined by four points.
void tetrahedron_volume_hessian(double ax, double ay, double az, double bx, double by, double bz, double cx, double cy, double cz, double dx, double dy, double dz, double H[144])
Compute the gradient of the signed area of a 2D triangle defined by three points.
double tetrahedron_volume(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const Eigen::Vector3d &d)
Compute the signed volume of a tetrahedron defined by four points.
void triangle_area_2D_gradient(double ax, double ay, double bx, double by, double cx, double cy, double g[6])
Compute the gradient of the signed area of a 2D triangle defined by three points.
bool are_edges_collinear(const VectorNd &ea0, const VectorNd &ea1, const VectorNd &eb0, const VectorNd &eb1, const double tol)
Determine if two edges are collinear.
double triangle_area_3D(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c)
Compute the area of a 3D triangle defined by three points.
Eigen::MatrixXd triangle_to_clockwise_order(const Eigen::MatrixXd &triangle)
Reorder the vertices of a triangle so they are in clockwise order.
std::vector< Eigen::MatrixXd > triangle_fan(const Eigen::MatrixXd &convex_polygon)
Convert a convex polygon to a list of triangles fanned around the first vertex.
bool are_triangles_coplanar(const Eigen::Vector3d &t00, const Eigen::Vector3d &t01, const Eigen::Vector3d &t02, const Eigen::Vector3d &t10, const Eigen::Vector3d &t11, const Eigen::Vector3d &t12, const double tol)
Determine if two triangles are coplanar.
Eigen::VectorXd barycentric_coordinates(const Eigen::VectorXd &p, const Eigen::MatrixXd &V)
Compute barycentric coordinates for point p with respect to a simplex.
double triangle_area(const Eigen::MatrixXd V)
Compute the signed area of a triangle defined by three points.
double triangle_area_2D(const Eigen::Vector2d &a, const Eigen::Vector2d &b, const Eigen::Vector2d &c)
Compute the signed area of a 2D triangle defined by three points.
bool triangle_intersects_disk(const Eigen::Vector2d &t0, const Eigen::Vector2d &t1, const Eigen::Vector2d &t2, const Eigen::Vector2d &center, const double radius)
Determine if a 2D triangle intersects a 2D disk.
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11