Loading...
Searching...
No Matches
triangleFuncs.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "triangleFuncs.H"
30#include "triangle.H"
31#include "treeBoundBox.H"
32#include "SortableList.H"
33#include "boolList.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37void Foam::triangleFuncs::setIntersection
38(
39 const point& oppositeSidePt,
40 const scalar oppositeSign,
41
42 const point& thisSidePt,
43 const scalar thisSign,
44
45 const scalar tol,
46
47 point& pt
48)
49{
50 const scalar denom = oppositeSign - thisSign;
51
52 if (mag(denom) < tol)
53 {
54 // If almost does not cut choose one which certainly cuts.
55 pt = oppositeSidePt;
56 }
57 else
58 {
59 pt = oppositeSidePt + oppositeSign/denom*(thisSidePt - oppositeSidePt);
60 }
61}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
67(
68 const point& V0,
69 const point& V10,
70 const point& V20,
71 const label i0,
72 const pointField& origin,
73 const scalar maxLength,
74 point& pInter
75)
76{
77 // Based on Graphics Gems - Fast Ray Triangle intersection.
78 // Since direction is coordinate axis there is no need to do projection,
79 // we can directly check u,v components for inclusion in triangle.
80
81 // Get other components
82 label i1 = (i0 + 1) % 3;
83 label i2 = (i1 + 1) % 3;
84
85 scalar u1 = V10[i1];
86 scalar v1 = V10[i2];
87
88 scalar u2 = V20[i1];
89 scalar v2 = V20[i2];
90
91 scalar localScale = mag(u1)+mag(v1)+mag(u2)+mag(v2);
92
93 scalar det = v2*u1 - u2*v1;
94
95 // Fix for V0:(-31.71428 0 -15.10714)
96 // V10:(-1.285715 8.99165e-16 -1.142858)
97 // V20:(0 0 -1.678573)
98 // i0:0
99 if (localScale < VSMALL || Foam::mag(det)/localScale < SMALL)
100 {
101 // Triangle parallel to dir
102 return false;
103 }
104
105 forAll(origin, originI)
106 {
107 const point& P = origin[originI];
108
109 scalar u0 = P[i1] - V0[i1];
110 scalar v0 = P[i2] - V0[i2];
111
112 scalar alpha = 0;
113 scalar beta = 0;
114 bool inter = false;
115
116 if (Foam::mag(u1) < ROOTVSMALL)
117 {
118 beta = u0/u2;
119 if ((beta >= 0) && (beta <= 1))
120 {
121 alpha = (v0 - beta*v2)/v1;
122 inter = ((alpha >= 0) && ((alpha + beta) <= 1));
123 }
124 }
125 else
126 {
127 beta = (v0*u1 - u0*v1)/det;
128 if ((beta >= 0) && (beta <= 1))
129 {
130 alpha = (u0 - beta*u2)/u1;
131 inter = ((alpha >= 0) && ((alpha + beta) <= 1));
132 }
133 }
134
135 if (inter)
136 {
137 pInter = V0 + alpha*V10 + beta*V20;
138 scalar s = (pInter - origin[originI])[i0];
139
140 if ((s >= 0) && (s <= maxLength))
141 {
142 return true;
144 }
145 }
146 return false;
147}
148
149
151(
152 const triPointRef& tri,
153 const treeBoundBox& cubeBb
154)
155{
156 // Slow (edge by edge) bounding box intersection. TBD: replace with call
157 // to above intersectAxesBundle. However this function is not fully
158 // correct and misses intersection between some triangles.
159 {
160 const pointField points(cubeBb.points());
161
162 for (const edge& e : treeBoundBox::edges)
163 {
164 const point& start = points[e[0]];
165 const point& end = points[e[1]];
166
167 pointHit inter = tri.intersection
168 (
169 start,
170 end-start,
172 );
173
174 if (inter.hit() && inter.distance() <= 1)
175 {
176 return true;
177 }
178 }
179 }
180
181
182 // Intersect triangle edges with bounding box
183 point pInter;
184
185 return
186 (
187 cubeBb.intersects(tri.a(), tri.b(), pInter)
188 || cubeBb.intersects(tri.b(), tri.c(), pInter)
189 || cubeBb.intersects(tri.c(), tri.a(), pInter)
190 );
191}
192
193
195(
196 const point& p0,
197 const point& p1,
198 const point& p2,
199 const treeBoundBox& cubeBb
200)
202 const triPointRef tri(p0, p1, p2);
203
204 return intersectBb(tri, cubeBb);
205}
206
207
209(
210 const point& va0,
211 const point& va10,
212 const point& va20,
213
214 const point& base,
215 const point& normal,
216
217 point& pInter0,
218 point& pInter1
219)
220{
221 // Get triangle normal
222 vector na = va10 ^ va20;
223 scalar magArea = mag(na);
224 na/magArea;
225
226 if (mag(na & normal) > (1 - SMALL))
227 {
228 // Parallel
229 return false;
230 }
231
232 const point va1 = va0 + va10;
233 const point va2 = va0 + va20;
234
235 // Find the triangle point on the other side.
236 scalar sign0 = (va0 - base) & normal;
237 scalar sign1 = (va1 - base) & normal;
238 scalar sign2 = (va2 - base) & normal;
239
240 label oppositeVertex = -1;
241
242 if (sign0 < 0)
243 {
244 if (sign1 < 0)
245 {
246 if (sign2 < 0)
247 {
248 // All on same side of plane
249 return false;
250 }
251 else // sign2 >= 0
252 {
253 // 2 on opposite side.
254 oppositeVertex = 2;
255 }
256 }
257 else // sign1 >= 0
258 {
259 if (sign2 < 0)
260 {
261 // 1 on opposite side.
262 oppositeVertex = 1;
263 }
264 else
265 {
266 // 0 on opposite side.
267 oppositeVertex = 0;
268 }
269 }
270 }
271 else // sign0 >= 0
272 {
273 if (sign1 < 0)
274 {
275 if (sign2 < 0)
276 {
277 // 0 on opposite side.
278 oppositeVertex = 0;
279 }
280 else // sign2 >= 0
281 {
282 // 1 on opposite side.
283 oppositeVertex = 1;
284 }
285 }
286 else // sign1 >= 0
287 {
288 if (sign2 < 0)
289 {
290 // 2 on opposite side.
291 oppositeVertex = 2;
292 }
293 else // sign2 >= 0
294 {
295 // All on same side of plane
296 return false;
297 }
298 }
299 }
300
301 scalar tol = SMALL*Foam::sqrt(magArea);
302
303 if (oppositeVertex == 0)
304 {
305 // 0 on opposite side. Cut edges 01 and 02
306 setIntersection(va0, sign0, va1, sign1, tol, pInter0);
307 setIntersection(va0, sign0, va2, sign2, tol, pInter1);
308 }
309 else if (oppositeVertex == 1)
310 {
311 // 1 on opposite side. Cut edges 10 and 12
312 setIntersection(va1, sign1, va0, sign0, tol, pInter0);
313 setIntersection(va1, sign1, va2, sign2, tol, pInter1);
314 }
315 else // oppositeVertex == 2
316 {
317 // 2 on opposite side. Cut edges 20 and 21
318 setIntersection(va2, sign2, va0, sign0, tol, pInter0);
319 setIntersection(va2, sign2, va1, sign1, tol, pInter1);
320 }
321
322 return true;
323}
324
325
327(
328 const point& va0,
329 const point& va10,
330 const point& va20,
331
332 const point& vb0,
333 const point& vb10,
334 const point& vb20,
335
336 point& pInter0,
337 point& pInter1
338)
339{
340 // Get triangle normals
341 vector na = va10 ^ va20;
342 na/mag(na);
343
344 vector nb = vb10 ^ vb20;
345 nb/mag(nb);
346
347 // Calculate intersection of triangle a with plane of b
348 point planeB0;
349 point planeB1;
350 if (!intersect(va0, va10, va20, vb0, nb, planeB0, planeB1))
351 {
352 return false;
353 }
354
355 // ,, triangle b with plane of a
356 point planeA0;
357 point planeA1;
358 if (!intersect(vb0, vb10, vb20, va0, na, planeA0, planeA1))
359 {
360 return false;
361 }
362
363 // Now check if intersections overlap (w.r.t. intersection of the two
364 // planes)
365
366 vector intersection(na ^ nb);
367
368 scalar coordB0 = planeB0 & intersection;
369 scalar coordB1 = planeB1 & intersection;
370
371 scalar coordA0 = planeA0 & intersection;
372 scalar coordA1 = planeA1 & intersection;
373
374 // Put information in indexable form for sorting.
375 List<point*> pts(4);
376 boolList isFromB(4);
377 SortableList<scalar> sortCoords(4);
378
379 pts[0] = &planeB0;
380 isFromB[0] = true;
381 sortCoords[0] = coordB0;
382
383 pts[1] = &planeB1;
384 isFromB[1] = true;
385 sortCoords[1] = coordB1;
386
387 pts[2] = &planeA0;
388 isFromB[2] = false;
389 sortCoords[2] = coordA0;
390
391 pts[3] = &planeA1;
392 isFromB[3] = false;
393 sortCoords[3] = coordA1;
394
395 sortCoords.sort();
396
397 const labelList& indices = sortCoords.indices();
398
399 if (isFromB[indices[0]] == isFromB[indices[1]])
400 {
401 // Entry 0 and 1 are of same region (both a or both b). Hence that
402 // region does not overlap.
403 return false;
404 }
405 else
406 {
407 // Different type. Start of overlap at indices[1], end at indices[2]
408 pInter0 = *pts[indices[1]];
409 pInter1 = *pts[indices[2]];
410
411 return true;
412 }
413}
414
415
416// ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
bool hit() const noexcept
Is there a hit.
Definition pointHit.H:145
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const noexcept
Return the list of sorted indices. Updated every sort.
void sort()
Forward (stable) sort the list (if changed after construction).
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Foam::intersection.
Standard boundBox with extra functionality for use in octree.
static const edgeList edges
Edge to point addressing, using octant corner points.
tmp< pointField > points() const
Vertex coordinates. In octant coding.
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
static bool intersectBb(const triPointRef &tri, const treeBoundBox &cubeBb)
Intersect triangle with bounding box.
static bool intersect(const point &va0, const point &va10, const point &va20, const point &basePoint, const vector &normal, point &pInter0, point &pInter1)
Intersect triangle with plane.
static bool intersectAxesBundle(const point &V0, const point &V10, const point &V20, const label i0, const pointField &origin, const scalar maxLength, point &pInter)
Intersect triangle with parallel edges aligned with axis i0.
const Point & b() const noexcept
The second vertex.
Definition triangle.H:403
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition triangleI.H:672
const Point & c() const noexcept
The third vertex.
Definition triangle.H:408
const Point & a() const noexcept
The first vertex.
Definition triangle.H:398
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & p0
Definition EEqn.H:36
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const char * end
Definition SVGTools.H:223
dimensionedScalar det(const dimensionedSphericalTensor &dt)
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
scalar u0
volScalarField & alpha
const pointField & pts
volScalarField & e
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299