DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
gjk.h
1// Copyright (C) 2020 Chris Richardson
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <array>
10#include <concepts>
11#include <dolfinx/common/math.h>
12#include <limits>
13#include <numeric>
14#include <span>
15#include <utility>
16#include <vector>
17
18namespace dolfinx::geometry
19{
20
21namespace impl_gjk
22{
23
27template <std::floating_point T>
28std::pair<std::vector<T>, std::array<T, 3>>
29nearest_simplex(std::span<const T> s)
30{
31 assert(s.size() % 3 == 0);
32 const std::size_t s_rows = s.size() / 3;
33 switch (s_rows)
34 {
35 case 2:
36 {
37 // Compute lm = dot(s0, ds / |ds|)
38 auto s0 = s.template subspan<0, 3>();
39 auto s1 = s.template subspan<3, 3>();
40 std::array ds = {s1[0] - s0[0], s1[1] - s0[1], s1[2] - s0[2]};
41 const T lm = -(s0[0] * ds[0] + s0[1] * ds[1] + s0[2] * ds[2])
42 / (ds[0] * ds[0] + ds[1] * ds[1] + ds[2] * ds[2]);
43 if (lm >= 0.0 and lm <= 1.0)
44 {
45 // The origin is between A and B
46 // v = s0 + lm * (s1 - s0);
47 std::array v
48 = {s0[0] + lm * ds[0], s0[1] + lm * ds[1], s0[2] + lm * ds[2]};
49 return {std::vector<T>(s.begin(), s.end()), v};
50 }
51
52 if (lm < 0.0)
53 return {std::vector<T>(s0.begin(), s0.end()), {s0[0], s0[1], s0[2]}};
54 else
55 return {std::vector<T>(s1.begin(), s1.end()), {s1[0], s1[1], s1[2]}};
56 }
57 case 3:
58 {
59 auto a = s.template subspan<0, 3>();
60 auto b = s.template subspan<3, 3>();
61 auto c = s.template subspan<6, 3>();
62 auto length = [](auto& x, auto& y)
63 {
64 return std::transform_reduce(x.begin(), x.end(), y.begin(), 0.0,
65 std::plus{}, [](auto x, auto y)
66 { return (x - y) * (x - y); });
67 };
68 const T ab2 = length(a, b);
69 const T ac2 = length(a, c);
70 const T bc2 = length(b, c);
71
72 // Helper to compute dot(x, x - y)
73 auto helper = [](auto& x, auto& y)
74 {
75 return std::transform_reduce(x.begin(), x.end(), y.begin(), 0.0,
76 std::plus{},
77 [](auto x, auto y) { return x * (x - y); });
78 };
79 const std::array lm
80 = {helper(a, b) / ab2, helper(a, c) / ac2, helper(b, c) / bc2};
81
82 T caba = 0;
83 for (std::size_t i = 0; i < 3; ++i)
84 caba += (c[i] - a[i]) * (b[i] - a[i]);
85
86 // Calculate triangle ABC
87 const T c2 = 1 - caba * caba / (ab2 * ac2);
88 const T lbb = (lm[0] - lm[1] * caba / ab2) / c2;
89 const T lcc = (lm[1] - lm[0] * caba / ac2) / c2;
90
91 // Intersects triangle
92 if (lbb >= 0.0 and lcc >= 0.0 and (lbb + lcc) <= 1.0)
93 {
94 // Calculate intersection more accurately
95 // v = (c - a) x (b - a)
96 std::array dx0 = {c[0] - a[0], c[1] - a[1], c[2] - a[2]};
97 std::array dx1 = {b[0] - a[0], b[1] - a[1], b[2] - a[2]};
98 std::array v = math::cross(dx0, dx1);
99
100 // Barycentre of triangle
101 std::array p = {(a[0] + b[0] + c[0]) / 3, (a[1] + b[1] + c[1]) / 3,
102 (a[2] + b[2] + c[2]) / 3};
103
104 T sum = v[0] * p[0] + v[1] * p[1] + v[2] * p[2];
105 T vnorm2 = v[0] * v[0] + v[1] * v[1] + v[2] * v[2];
106 for (std::size_t i = 0; i < 3; ++i)
107 v[i] *= sum / vnorm2;
108
109 return {std::vector(s.begin(), s.end()), v};
110 }
111
112 // Get closest point
113 std::size_t pos = 0;
114 {
115 T norm0 = std::numeric_limits<T>::max();
116 for (std::size_t i = 0; i < s.size(); i += 3)
117 {
118 std::span<const T, 3> p(s.data() + i, 3);
119 T norm = p[0] * p[0] + p[1] * p[1] + p[2] * p[2];
120 if (norm < norm0)
121 {
122 pos = i / 3;
123 norm0 = norm;
124 }
125 }
126 }
127
128 std::array vmin = {s[3 * pos], s[3 * pos + 1], s[3 * pos + 2]};
129 T qmin = 0;
130 for (std::size_t k = 0; k < 3; ++k)
131 qmin += vmin[k] * vmin[k];
132
133 std::vector smin = {vmin[0], vmin[1], vmin[2]};
134
135 // Check if edges are closer
136 constexpr int f[3][2] = {{0, 1}, {0, 2}, {1, 2}};
137 for (std::size_t i = 0; i < s_rows; ++i)
138 {
139 auto s0 = s.subspan(3 * f[i][0], 3);
140 auto s1 = s.subspan(3 * f[i][1], 3);
141 if (lm[i] > 0 and lm[i] < 1)
142 {
143 std::array<T, 3> v;
144 for (std::size_t k = 0; k < 3; ++k)
145 v[k] = s0[k] + lm[i] * (s1[k] - s0[k]);
146 T qnorm = 0;
147 for (std::size_t k = 0; k < 3; ++k)
148 qnorm += v[k] * v[k];
149 if (qnorm < qmin)
150 {
151 std::copy(v.begin(), v.end(), vmin.begin());
152 qmin = qnorm;
153 smin.resize(2 * 3);
154 std::span<T, 3> smin0(smin.data(), 3);
155 std::copy(s0.begin(), s0.end(), smin0.begin());
156 std::span<T, 3> smin1(smin.data() + 3, 3);
157 std::copy(s1.begin(), s1.end(), smin1.begin());
158 }
159 }
160 }
161 return {std::move(smin), vmin};
162 }
163 case 4:
164 {
165 auto s0 = s.template subspan<0, 3>();
166 auto s1 = s.template subspan<3, 3>();
167 auto s2 = s.template subspan<6, 3>();
168 auto s3 = s.template subspan<9, 3>();
169 auto W1 = math::cross(s0, s1);
170 auto W2 = math::cross(s2, s3);
171
172 std::array<T, 4> B;
173 B[0] = std::transform_reduce(s2.begin(), s2.end(), W1.begin(), 0.0);
174 B[1] = -std::transform_reduce(s3.begin(), s3.end(), W1.begin(), 0.0);
175 B[2] = std::transform_reduce(s0.begin(), s0.end(), W2.begin(), 0.0);
176 B[3] = -std::transform_reduce(s1.begin(), s1.end(), W2.begin(), 0.0);
177
178 const bool signDetM = std::signbit(std::reduce(B.begin(), B.end(), 0.0));
179 std::array<bool, 4> f_inside;
180 for (int i = 0; i < 4; ++i)
181 f_inside[i] = (std::signbit(B[i]) == signDetM);
182
183 if (f_inside[1] and f_inside[2] and f_inside[3])
184 {
185 if (f_inside[0]) // The origin is inside the tetrahedron
186 return {std::vector<T>(s.begin(), s.end()), {0, 0, 0}};
187 else // The origin projection P faces BCD
188 return nearest_simplex<T>(s.template subspan<0, 3 * 3>());
189 }
190
191 // Test ACD, ABD and/or ABC
192 std::vector<T> smin;
193 std::array<T, 3> vmin = {0, 0, 0};
194 constexpr int facets[3][3] = {{0, 1, 3}, {0, 2, 3}, {1, 2, 3}};
195 T qmin = std::numeric_limits<T>::max();
196 std::vector<T> M(9);
197 for (int i = 0; i < 3; ++i)
198 {
199 if (f_inside[i + 1] == false)
200 {
201 std::copy_n(std::next(s.begin(), 3 * facets[i][0]), 3, M.begin());
202 std::copy_n(std::next(s.begin(), 3 * facets[i][1]), 3,
203 std::next(M.begin(), 3));
204 std::copy_n(std::next(s.begin(), 3 * facets[i][2]), 3,
205 std::next(M.begin(), 6));
206
207 const auto [snew, v] = nearest_simplex<T>(M);
208 T q = std::transform_reduce(v.begin(), v.end(), v.begin(), 0);
209 if (q < qmin)
210 {
211 qmin = q;
212 vmin = v;
213 smin = snew;
214 }
215 }
216 }
217 return {smin, vmin};
218 }
219 default:
220 throw std::runtime_error("Number of rows defining simplex not supported.");
221 }
222}
223
225template <std::floating_point T>
226std::array<T, 3> support(std::span<const T> bd, std::array<T, 3> v)
227{
228 int i = 0;
229 T qmax = bd[0] * v[0] + bd[1] * v[1] + bd[2] * v[2];
230 for (std::size_t m = 1; m < bd.size() / 3; ++m)
231 {
232 T q = bd[3 * m] * v[0] + bd[3 * m + 1] * v[1] + bd[3 * m + 2] * v[2];
233 if (q > qmax)
234 {
235 qmax = q;
236 i = m;
237 }
238 }
239
240 return {bd[3 * i], bd[3 * i + 1], bd[3 * i + 2]};
241}
242} // namespace impl_gjk
243
254template <std::floating_point T>
255std::array<T, 3> compute_distance_gjk(std::span<const T> p,
256 std::span<const T> q)
257{
258 assert(p.size() % 3 == 0);
259 assert(q.size() % 3 == 0);
260
261 constexpr int maxk = 15; // Maximum number of iterations of the GJK algorithm
262
263 // Tolerance
264 constexpr T eps = 1.0e4 * std::numeric_limits<T>::epsilon();
265
266 // Initialise vector and simplex
267 std::array<T, 3> v = {p[0] - q[0], p[1] - q[1], p[2] - q[2]};
268 std::vector<T> s = {v[0], v[1], v[2]};
269
270 // Begin GJK iteration
271 int k;
272 for (k = 0; k < maxk; ++k)
273 {
274 // Support function
275 std::array w1 = impl_gjk::support(p, {-v[0], -v[1], -v[2]});
276 std::array w0 = impl_gjk::support(q, {v[0], v[1], v[2]});
277 const std::array w = {w1[0] - w0[0], w1[1] - w0[1], w1[2] - w0[2]};
278
279 // Break if any existing points are the same as w
280 assert(s.size() % 3 == 0);
281 std::size_t m;
282 for (m = 0; m < s.size() / 3; ++m)
283 {
284 auto it = std::next(s.begin(), 3 * m);
285 if (std::equal(it, std::next(it, 3), w.begin(), w.end()))
286 break;
287 }
288
289 if (m != s.size() / 3)
290 break;
291
292 // 1st exit condition (v - w).v = 0
293 const T vnorm2 = v[0] * v[0] + v[1] * v[1] + v[2] * v[2];
294 const T vw = vnorm2 - (v[0] * w[0] + v[1] * w[1] + v[2] * w[2]);
295 if (vw < (eps * vnorm2) or vw < eps)
296 break;
297
298 // Add new vertex to simplex
299 s.insert(s.end(), w.begin(), w.end());
300
301 // Find nearest subset of simplex
302 auto [snew, vnew] = impl_gjk::nearest_simplex<T>(s);
303 s.assign(snew.data(), snew.data() + snew.size());
304 v = {vnew[0], vnew[1], vnew[2]};
305
306 // 2nd exit condition - intersecting or touching
307 if ((v[0] * v[0] + v[1] * v[1] + v[2] * v[2]) < eps * eps)
308 break;
309 }
310
311 if (k == maxk)
312 throw std::runtime_error("GJK error - max iteration limit reached");
313
314 return v;
315}
316
317} // namespace dolfinx::geometry
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:21
std::array< T, 3 > compute_distance_gjk(std::span< const T > p, std::span< const T > q)
Compute the distance between two convex bodies p and q, each defined by a set of points.
Definition gjk.h:255
auto norm(const V &x, Norm type=Norm::l2)
Definition Vector.h:267