treeDataPrimitivePatch.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-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2020 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
30#include "indexedOctree.H"
31#include "triangleFuncs.H"
32#include "triSurfaceTools.H"
33#include "triFace.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class PatchType>
39{
40 if (cacheBb_)
41 {
42 bbs_.setSize(patch_.size());
43
44 forAll(patch_, i)
45 {
46 bbs_[i] = treeBoundBox(patch_.points(), patch_[i]);
47 }
48 }
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54template<class PatchType>
56(
57 const bool cacheBb,
58 const PatchType& patch,
59 const scalar planarTol
60)
61:
62 patch_(patch),
63 cacheBb_(cacheBb),
64 planarTol_(planarTol)
65{
66 update();
67}
68
69
70template<class PatchType>
72(
74)
75:
76 tree_(tree)
77{}
78
79
80template<class PatchType>
82(
84)
85:
86 tree_(tree)
87{}
88
89
90template<class PatchType>
92(
94 DynamicList<label>& shapeMask
95)
96:
97 tree_(tree),
98 shapeMask_(shapeMask)
99{}
100
101
102template<class PatchType>
105(
107 const label edgeID
108)
109:
110 tree_(tree),
111 edgeID_(edgeID)
112{}
113
114
115// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116
117template<class PatchType>
119{
120 pointField cc(patch_.size());
121
122 forAll(patch_, i)
123 {
124 cc[i] = patch_[i].centre(patch_.points());
125 }
126
127 return cc;
128}
129
130
131template<class PatchType>
133(
135 const point& sample
136) const
137{
138 // Need to determine whether sample is 'inside' or 'outside'
139 // Done by finding nearest face. This gives back a face which is
140 // guaranteed to contain nearest point. This point can be
141 // - in interior of face: compare to face normal
142 // - on edge of face: compare to edge normal
143 // - on point of face: compare to point normal
144 // Unfortunately the octree does not give us back the intersection point
145 // or where on the face it has hit so we have to recreate all that
146 // information.
147
148
149 // Find nearest face to sample
150 pointIndexHit info = oc.findNearest(sample, sqr(GREAT));
151
152 if (info.index() == -1)
153 {
155 << "Could not find " << sample << " in octree."
156 << abort(FatalError);
157 }
158
159 // Get actual intersection point on face
160 label facei = info.index();
161
162 if (debug & 2)
163 {
164 Pout<< "getSampleType : sample:" << sample
165 << " nearest face:" << facei;
166 }
167
168 const typename PatchType::face_type& localF = patch_.localFaces()[facei];
169 const typename PatchType::face_type& f = patch_[facei];
170 const pointField& points = patch_.points();
171 const labelList& mp = patch_.meshPoints();
172
173 // Retest to classify where on face info is. Note: could be improved. We
174 // already have point.
175
176 pointHit curHit = f.nearestPoint(sample, points);
177 const vector area = f.areaNormal(points);
178 const point& curPt = curHit.rawPoint();
179
180 //
181 // 1] Check whether sample is above face
182 //
183
184 if (curHit.hit())
185 {
186 // Nearest point inside face. Compare to face normal.
187
188 if (debug & 2)
189 {
190 Pout<< " -> face hit:" << curPt
191 << " comparing to face normal " << area << endl;
192 }
194 (
195 area,
196 sample - curPt
197 );
198 }
199
200 if (debug & 2)
201 {
202 Pout<< " -> face miss:" << curPt;
203 }
204
205 //
206 // 2] Check whether intersection is on one of the face vertices or
207 // face centre
208 //
209
210 const scalar typDimSqr = mag(area) + VSMALL;
211
212
213 forAll(f, fp)
214 {
215 if ((magSqr(points[f[fp]] - curPt)/typDimSqr) < planarTol_)
216 {
217 // Face intersection point equals face vertex fp
218
219 // Calculate point normal (wrong: uses face normals instead of
220 // triangle normals)
221
223 (
224 patch_.pointNormals()[localF[fp]],
225 sample - curPt
226 );
227 }
228 }
229
230 const point fc(f.centre(points));
231
232 if ((magSqr(fc - curPt)/typDimSqr) < planarTol_)
233 {
234 // Face intersection point equals face centre. Normal at face centre
235 // is already average of face normals
236
237 if (debug & 2)
238 {
239 Pout<< " -> centre hit:" << fc
240 << " distance:" << magSqr(fc - curPt)/typDimSqr << endl;
241 }
242
244 (
245 area,
246 sample - curPt
247 );
248 }
249
250
251
252 //
253 // 3] Get the 'real' edge the face intersection is on
254 //
255
256 const labelList& fEdges = patch_.faceEdges()[facei];
257
258 forAll(fEdges, fEdgeI)
259 {
260 label edgeI = fEdges[fEdgeI];
261 const edge& e = patch_.edges()[edgeI];
262 const linePointRef ln(points[mp[e.start()]], points[mp[e.end()]]);
263 pointHit edgeHit = ln.nearestDist(sample);
264
265 if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
266 {
267 // Face intersection point lies on edge e
268
269 // Calculate edge normal (wrong: uses face normals instead of
270 // triangle normals)
271 const labelList& eFaces = patch_.edgeFaces()[edgeI];
272
273 vector edgeNormal(Zero);
274
275 forAll(eFaces, i)
276 {
277 edgeNormal += patch_.faceNormals()[eFaces[i]];
278 }
279
280 if (debug & 2)
281 {
282 Pout<< " -> real edge hit point:" << edgeHit.rawPoint()
283 << " comparing to edge normal:" << edgeNormal
284 << endl;
285 }
286
287 // Found face intersection point on this edge. Compare to edge
288 // normal
290 (
291 edgeNormal,
292 sample - curPt
293 );
294 }
295 }
296
297
298 //
299 // 4] Get the internal edge the face intersection is on
300 //
301
302 forAll(f, fp)
303 {
304 pointHit edgeHit = linePointRef(points[f[fp]], fc).nearestDist(sample);
305
306 if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
307 {
308 // Face intersection point lies on edge between two face triangles
309
310 // Calculate edge normal as average of the two triangle normals
311 vector e = points[f[fp]] - fc;
312 vector ePrev = points[f[f.rcIndex(fp)]] - fc;
313 vector eNext = points[f[f.fcIndex(fp)]] - fc;
314
315 vector nLeft = normalised(ePrev ^ e);
316 vector nRight = normalised(e ^ eNext);
317
318 if (debug & 2)
319 {
320 Pout<< " -> internal edge hit point:" << edgeHit.rawPoint()
321 << " comparing to edge normal "
322 << 0.5*(nLeft + nRight)
323 << endl;
324 }
325
326 // Found face intersection point on this edge. Compare to edge
327 // normal
329 (
330 0.5*(nLeft + nRight),
331 sample - curPt
332 );
333 }
334 }
335
336 if (debug & 2)
337 {
338 Pout<< "Did not find sample " << sample
339 << " anywhere related to nearest face " << facei << endl
340 << "Face:";
341
342 forAll(f, fp)
343 {
344 Pout<< " vertex:" << f[fp]
345 << " coord:" << points[f[fp]]
346 << endl;
347 }
348 }
349
350 // Can't determine status of sample with respect to nearest face.
351 // Either
352 // - tolerances are wrong. (if e.g. face has zero area)
353 // - or (more likely) surface is not closed.
354
355 return volumeType::UNKNOWN;
356}
357
358
359// Check if any point on shape is inside cubeBb.
360template<class PatchType>
362(
363 const label index,
364 const treeBoundBox& cubeBb
365) const
366{
367 // 1. Quick rejection: bb does not intersect face bb at all
368 if
369 (
370 cacheBb_
371 ? !cubeBb.overlaps(bbs_[index])
372 : !cubeBb.overlaps(treeBoundBox(patch_.points(), patch_[index]))
373 )
374 {
375 return false;
376 }
377
378
379 // 2. Check if one or more face points inside
380
381 const pointField& points = patch_.points();
382 const typename PatchType::face_type& f = patch_[index];
383
384 if (cubeBb.containsAny(points, f))
385 {
386 return true;
387 }
388
389 // 3. Difficult case: all points are outside but connecting edges might
390 // go through cube. Use triangle-bounding box intersection.
391 const point fc = f.centre(points);
392
393 if (f.size() == 3)
394 {
396 (
397 points[f[0]],
398 points[f[1]],
399 points[f[2]],
400 cubeBb
401 );
402 }
403 else
404 {
405 forAll(f, fp)
406 {
407 bool triIntersects = triangleFuncs::intersectBb
408 (
409 points[f[fp]],
410 points[f[f.fcIndex(fp)]],
411 fc,
412 cubeBb
413 );
414
415 if (triIntersects)
416 {
417 return true;
418 }
419 }
420 }
421
422 return false;
423}
424
425
426// Check if any point on shape is inside sphere.
427template<class PatchType>
429(
430 const label index,
431 const point& centre,
432 const scalar radiusSqr
433) const
434{
435 // 1. Quick rejection: sphere does not intersect face bb at all
436 if
437 (
438 cacheBb_
439 ? !bbs_[index].overlaps(centre, radiusSqr)
440 : !treeBoundBox(patch_.points(),patch_[index]).overlaps(centre, radiusSqr)
441 )
442 {
443 return false;
444 }
445
446 const pointField& points = patch_.points();
447 const face& f = patch_[index];
448
449 pointHit nearHit = f.nearestPoint(centre, points);
450
451 // If the distance to the nearest point on the face from the sphere centres
452 // is within the radius, then the sphere touches the face.
453 if (sqr(nearHit.distance()) < radiusSqr)
454 {
455 return true;
456 }
457
458 return false;
459}
460
461
462template<class PatchType>
464(
465 const labelUList& indices,
466 const point& sample,
467
468 scalar& nearestDistSqr,
469 label& minIndex,
470 point& nearestPoint
471) const
472{
473 const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
474 const PatchType& patch = shape.patch();
475
476 const pointField& points = patch.points();
477
478 for (const label index : indices)
479 {
480 const typename PatchType::face_type& f = patch[index];
481
482 const pointHit nearHit = f.nearestPoint(sample, points);
483 const scalar distSqr = sqr(nearHit.distance());
484
485 if (distSqr < nearestDistSqr)
486 {
487 nearestDistSqr = distSqr;
488 minIndex = index;
489 nearestPoint = nearHit.rawPoint();
490 }
491 }
492}
493
494
495template<class PatchType>
497(
498 const labelUList& indices,
499 const linePointRef& ln,
500
501 treeBoundBox& tightest,
502 label& minIndex,
503 point& linePoint,
504 point& nearestPoint
505) const
506{
508}
509
510
511template<class PatchType>
513(
514 const label index,
515 const point& start,
516 const point& end,
517 point& intersectionPoint
518) const
519{
520 return findIntersection(tree_, index, start, end, intersectionPoint);
521}
522
523
524template<class PatchType>
526(
527 const label index,
528 const point& start,
529 const point& end,
530 point& intersectionPoint
531) const
532{
533 if (shapeMask_.found(index))
534 {
535 return false;
536 }
537
538 return findIntersection(tree_, index, start, end, intersectionPoint);
539}
540
541
542template<class PatchType>
544(
545 const label index,
546 const point& start,
547 const point& end,
548 point& intersectionPoint
549) const
550{
551 if (edgeID_ == -1)
552 {
554 << "EdgeID not set. Please set edgeID to the index of"
555 << " the edge you are testing"
556 << exit(FatalError);
557 }
558
559 const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
560 const PatchType& patch = shape.patch();
561
562 const typename PatchType::face_type& f = patch.localFaces()[index];
563 const edge& e = patch.edges()[edgeID_];
564
565 if (!f.found(e[0]) && !f.found(e[1]))
566 {
567 return findIntersection(tree_, index, start, end, intersectionPoint);
568 }
569
570 return false;
571}
572
573
574template<class PatchType>
576(
578 const label index,
579 const point& start,
580 const point& end,
581 point& intersectionPoint
582)
583{
584 const treeDataPrimitivePatch<PatchType>& shape = tree.shapes();
585 const PatchType& patch = shape.patch();
586
587 const pointField& points = patch.points();
588 const typename PatchType::face_type& f = patch[index];
589
590 // Do quick rejection test
591 if (shape.cacheBb_)
592 {
593 const treeBoundBox& faceBb = shape.bbs_[index];
594
595 if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
596 {
597 // start and end in same block outside of faceBb.
598 return false;
599 }
600 }
601
602 const vector dir(end - start);
603 pointHit inter;
604
605 if (f.size() == 3)
606 {
607 inter = triPointRef
608 (
609 points[f[0]],
610 points[f[1]],
611 points[f[2]]
612 ).intersection(start, dir, intersection::HALF_RAY, shape.planarTol_);
613 }
614 else
615 {
616 const pointField& faceCentres = patch.faceCentres();
617
618 inter = f.intersection
619 (
620 start,
621 dir,
622 faceCentres[index],
623 points,
625 shape.planarTol_
626 );
627 }
628
629 if (inter.hit() && inter.distance() <= 1)
630 {
631 // Note: no extra test on whether intersection is in front of us
632 // since using half_ray
633 intersectionPoint = inter.hitPoint();
634
635 return true;
636 }
637
638 return false;
639}
640
641
642// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
Minimal example by using system/controlDict.functions:
int overlaps
Flag to control which overlap calculations are performed.
Definition: PDRparams.H:97
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
const point_type & rawPoint() const noexcept
The point, no checks.
Definition: PointHit.H:172
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition: PointHit.H:145
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
label index() const noexcept
Return the hit index.
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
label rcIndex(const label i) const noexcept
Definition: UListI.H:67
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:249
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
A line primitive.
Definition: line.H:68
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
bool overlaps(const boundBox &bb) const
Overlaps other bounding box?
Definition: boundBoxI.H:221
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:393
Encapsulation of data needed to search on PrimitivePatches.
volumeType getVolumeType(const indexedOctree< treeDataPrimitivePatch< PatchType > > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does shape at index overlap bb.
static bool findIntersection(const indexedOctree< treeDataPrimitivePatch< PatchType > > &tree, const label index, const point &start, const point &end, point &intersectionPoint)
Helper: find intersection of line with shapes.
const PatchType & patch() const
Return access to the underlying patch.
static bool intersectBb(const point &p0, const point &p1, const point &p2, const treeBoundBox &cubeBb)
Intersect triangle with bounding box.
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:61
@ UNKNOWN
Unknown state.
Definition: volumeType.H:67
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:933
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333