Loading...
Searching...
No Matches
treeDataFace.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) 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 "treeDataFace.H"
30#include "polyMesh.H"
31#include "triangle.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
39 scalar treeDataFace::tolSqr = sqr(1e-6);
40}
41
42
43// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48// Bound boxes corresponding to specified faces
49template<class ElementIds>
51(
52 const primitiveMesh& mesh,
53 const ElementIds& elemIds
54)
55{
56 treeBoundBoxList bbs(elemIds.size());
57
58 std::transform
59 (
60 elemIds.cbegin(),
61 elemIds.cend(),
62 bbs.begin(),
63 [&](label facei)
64 {
65 return treeBoundBox(mesh.points(), mesh.faces()[facei]);
66 }
67 );
69 return bbs;
70}
71
72
73// Overall bound box for specified faces
74template<class ElementIds>
75static treeBoundBox boundsImpl
76(
77 const primitiveMesh& mesh,
78 const ElementIds& elemIds
79)
80{
81 treeBoundBox bb;
82
83 for (const label facei : elemIds)
84 {
85 bb.add(mesh.points(), mesh.faces()[facei]);
86 }
87
88 return bb;
89}
91} // End namespace Foam
92
93
94// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
95
99 // All faces
100 return boxesImpl(mesh, labelRange(mesh.nFaces()));
101}
102
103
106(
107 const primitiveMesh& mesh,
108 const labelRange& range
117(
118 const primitiveMesh& mesh,
119 const labelUList& faceIds
121{
122 return boxesImpl(mesh, faceIds);
123}
124
125
128(
129 const primitiveMesh& mesh,
130 const labelRange& range
132{
133 return boundsImpl(mesh, range);
134}
135
136
139(
140 const primitiveMesh& mesh,
141 const labelUList& faceIds
142)
143{
144 return boundsImpl(mesh, faceIds);
145}
146
147
148// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
149
150inline bool
151Foam::treeDataFace::usesFace(const label facei) const
152{
153 return (!useSubset_ || isTreeFace_.test(facei));
154}
155
156
157inline Foam::treeBoundBox
158Foam::treeDataFace::getBounds(const label index) const
159{
160 const label facei = objectIndex(index);
161 return treeBoundBox(mesh_.points(), mesh_.faces()[facei]);
162}
163
164
165void Foam::treeDataFace::update()
166{
167 if (cacheBb_)
168 {
169 if (useSubset_)
170 {
171 bbs_ = treeDataFace::boxes(mesh_, faceLabels_);
172 }
173 else
174 {
175 bbs_ = treeDataFace::boxes(mesh_);
177 }
178}
179
180
181// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
182
184(
185 const bool cacheBb,
186 const primitiveMesh& mesh
187)
188:
189 mesh_(mesh),
190 faceLabels_(),
191 isTreeFace_(),
192 useSubset_(false),
193 cacheBb_(cacheBb)
194{
195 update();
196}
197
198
200(
201 const bool cacheBb,
202 const primitiveMesh& mesh,
203 const labelRange& range
204)
205:
206 mesh_(mesh),
207 faceLabels_(identity(range)),
208 isTreeFace_(range),
209 useSubset_(true),
210 cacheBb_(cacheBb)
211{
212 update();
213}
214
215
217(
218 const bool cacheBb,
219 const polyPatch& patch
220)
221:
222 treeDataFace(cacheBb, patch.boundaryMesh().mesh(), patch.range())
223{}
224
225
227(
228 const bool cacheBb,
229 const primitiveMesh& mesh,
231)
232:
233 mesh_(mesh),
234 faceLabels_(faceLabels),
235 isTreeFace_(mesh_.nFaces(), faceLabels_),
236 useSubset_(true),
237 cacheBb_(cacheBb)
238{
239 update();
240}
241
242
244(
245 const bool cacheBb,
246 const primitiveMesh& mesh,
248)
249:
250 mesh_(mesh),
251 faceLabels_(std::move(faceLabels)),
252 isTreeFace_(mesh_.nFaces(), faceLabels_),
253 useSubset_(true),
254 cacheBb_(cacheBb)
255{
257}
258
259
260// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
261
263Foam::treeDataFace::bounds(const labelUList& indices) const
264{
265 if (useSubset_)
266 {
267 treeBoundBox bb;
268
269 for (const label index : indices)
270 {
271 const label facei = faceLabels_[index];
272
273 bb.add(mesh_.points(), mesh_.faces()[facei]);
274 }
275
276 return bb;
277 }
278
279 return treeDataFace::bounds(mesh_, indices);
280}
281
282
284{
285 if (useSubset_)
286 {
287 return tmp<pointField>::New(mesh_.faceCentres(), faceLabels_);
288 }
289
290 return mesh_.faceCentres();
291}
292
293
295(
296 const indexedOctree<treeDataFace>& oc,
297 const point& sample
298) const
299{
300 // Need to determine whether sample is 'inside' or 'outside'
301 // Done by finding nearest face. This gives back a face which is
302 // guaranteed to contain nearest point. This point can be
303 // - in interior of face: compare to face normal
304 // - on edge of face: compare to edge normal
305 // - on point of face: compare to point normal
306 // Unfortunately the octree does not give us back the intersection point
307 // or where on the face it has hit so we have to recreate all that
308 // information.
309
310
311 // Find nearest face to sample
312 pointIndexHit info = oc.findNearest(sample, sqr(GREAT));
313
314 const label index = info.index();
315
316 if (index == -1)
317 {
319 << "Could not find " << sample << " in octree."
320 << abort(FatalError);
321 }
322
323
324 // Get actual intersection point on face
325 const label facei = objectIndex(index);
326
327 if (debug & 2)
328 {
329 Pout<< "getSampleType : sample:" << sample
330 << " nearest face:" << facei;
331 }
332
333 const pointField& points = mesh_.points();
334
335 // Retest to classify where on face info is. Note: could be improved. We
336 // already have point.
337
338 const face& f = mesh_.faces()[facei];
339 const vector& area = mesh_.faceAreas()[facei];
340 const point& fc = mesh_.faceCentres()[facei];
341
342 pointHit curHit = f.nearestPoint(sample, points);
343 const point& curPt = curHit.point();
344
345 //
346 // 1] Check whether sample is above face
347 //
348
349 if (curHit.hit())
350 {
351 // Nearest point inside face. Compare to face normal.
352
353 if (debug & 2)
354 {
355 Pout<< " -> face hit:" << curPt
356 << " comparing to face normal " << area << endl;
357 }
358 return indexedOctree<treeDataFace>::getSide(area, sample - curPt);
359 }
360
361 if (debug & 2)
362 {
363 Pout<< " -> face miss:" << curPt;
364 }
365
366 //
367 // 2] Check whether intersection is on one of the face vertices or
368 // face centre
369 //
370
371 const scalar typDimSqr = mag(area) + VSMALL;
372
373 for (const label fp : f)
374 {
375 const scalar relDistSqr = (magSqr(points[fp] - curPt)/typDimSqr);
376
377 if (relDistSqr < tolSqr)
378 {
379 // Face intersection point equals face vertex
380
381 // Calculate point normal (wrong: uses face normals instead of
382 // triangle normals)
383
384 vector pointNormal(Zero);
385
386 for (const label ptFacei : mesh_.pointFaces()[fp])
387 {
388 if (usesFace(ptFacei))
389 {
390 pointNormal += normalised(mesh_.faceAreas()[ptFacei]);
391 }
392 }
393
394 if (debug & 2)
395 {
396 Pout<< " -> face point hit :" << points[fp]
397 << " point normal:" << pointNormal
398 << " distance:" << relDistSqr << endl;
399 }
401 (
402 pointNormal,
403 sample - curPt
404 );
405 }
406 }
407
408 const scalar relDistSqr = (magSqr(fc - curPt)/typDimSqr);
409
410 if (relDistSqr < tolSqr)
411 {
412 // Face intersection point equals face centre. Normal at face centre
413 // is already average of face normals
414
415 if (debug & 2)
416 {
417 Pout<< " -> centre hit:" << fc
418 << " distance:" << relDistSqr << endl;
419 }
420
421 return indexedOctree<treeDataFace>::getSide(area, sample - curPt);
422 }
423
424
425 //
426 // 3] Get the 'real' edge the face intersection is on
427 //
428
429 for (const label edgei : mesh_.faceEdges()[facei])
430 {
431 pointHit edgeHit =
432 mesh_.edges()[edgei].line(points).nearestDist(sample);
433
434 const scalar relDistSqr = edgeHit.point().distSqr(curPt)/typDimSqr;
435
436 if (relDistSqr < tolSqr)
437 {
438 // Face intersection point lies on edge e
439
440 // Calculate edge normal (wrong: uses face normals instead of
441 // triangle normals)
442
443 vector edgeNormal(Zero);
444
445 for (const label eFacei : mesh_.edgeFaces()[edgei])
446 {
447 if (usesFace(eFacei))
448 {
449 edgeNormal += normalised(mesh_.faceAreas()[eFacei]);
450 }
451 }
452
453 if (debug & 2)
454 {
455 Pout<< " -> real edge hit point:" << edgeHit.point()
456 << " comparing to edge normal:" << edgeNormal
457 << endl;
458 }
459
460 // Found face intersection point on this edge. Compare to edge
461 // normal
463 (
464 edgeNormal,
465 sample - curPt
466 );
467 }
468 }
469
470
471 //
472 // 4] Get the internal edge the face intersection is on
473 //
474
475 forAll(f, fp)
476 {
477 pointHit edgeHit = linePointRef(points[f[fp]], fc).nearestDist(sample);
478
479 const scalar relDistSqr = edgeHit.point().distSqr(curPt)/typDimSqr;
480
481 if (relDistSqr < tolSqr)
482 {
483 // Face intersection point lies on edge between two face triangles
484
485 // Calculate edge normal as average of the two triangle normals
486 vector e = points[f[fp]] - fc;
487 vector ePrev = points[f[f.rcIndex(fp)]] - fc;
488 vector eNext = points[f[f.fcIndex(fp)]] - fc;
489
490 vector nLeft = normalised(ePrev ^ e);
491 vector nRight = normalised(e ^ eNext);
492
493 if (debug & 2)
494 {
495 Pout<< " -> internal edge hit point:" << edgeHit.point()
496 << " comparing to edge normal "
497 << 0.5*(nLeft + nRight)
498 << endl;
499 }
500
501 // Found face intersection point on this edge. Compare to edge
502 // normal
504 (
505 0.5*(nLeft + nRight),
506 sample - curPt
507 );
508 }
509 }
510
511 if (debug & 2)
512 {
513 Pout<< "Did not find sample " << sample
514 << " anywhere related to nearest face " << facei << endl
515 << "Face:";
516
517 forAll(f, fp)
518 {
519 Pout<< " vertex:" << f[fp] << " coord:" << points[f[fp]]
520 << endl;
521 }
522 }
523
524 // Can't determine status of sample with respect to nearest face.
525 // Either
526 // - tolerances are wrong. (if e.g. face has zero area)
527 // - or (more likely) surface is not closed.
528
529 return volumeType::UNKNOWN;
530}
531
532
533// Check if any point on shape is inside searchBox.
535(
536 const label index,
537 const treeBoundBox& searchBox
538) const
539{
540 // 1. Quick rejection: bb does not intersect face bb at all
541 if
542 (
543 cacheBb_
544 ? !searchBox.overlaps(bbs_[index])
545 : !searchBox.overlaps(getBounds(index))
546 )
547 {
548 return false;
549 }
550
551 const pointField& points = mesh_.points();
552
553 // 2. Check if one or more face points inside
554 const label facei = objectIndex(index);
555
556 const face& f = mesh_.faces()[facei];
557
558 if (f.size() == 3)
559 {
560 const triPointRef tri(points, f[0], f[1], f[2]);
561
562 return searchBox.intersects(tri);
563 }
564
565 if (searchBox.containsAny(points, f))
566 {
567 return true;
568 }
569
570 // 3. Difficult case: all points are outside but connecting edges might
571 // go through cube. Use triangle-bounding box intersection.
572
573 const point& fc = mesh_.faceCentres()[facei];
574
575 forAll(f, fp)
576 {
577 const triPointRef tri
578 (
579 points[f.thisLabel(fp)], points[f.nextLabel(fp)], fc
580 );
581
582 if (searchBox.intersects(tri))
583 {
584 return true;
585 }
587
588 return false;
589}
590
591
592// Check if any point on shape is inside sphere.
594(
595 const label index,
596 const point& centre,
597 const scalar radiusSqr
598) const
599{
600 // 1. Quick rejection: sphere does not intersect face bb at all
601 if
602 (
603 cacheBb_
604 ? !bbs_[index].overlaps(centre, radiusSqr)
605 : !getBounds(index).overlaps(centre, radiusSqr)
606 )
607 {
608 return false;
609 }
610
611 const label facei = objectIndex(index);
612
613 const face& f = mesh().faces()[facei];
614
615 pointHit nearHit = f.nearestPoint(centre, mesh().points());
616
617 if (sqr(nearHit.distance()) < radiusSqr)
618 {
619 return true;
620 }
622 return false;
623}
624
625
626// * * * * * * * * * * * * * * * * Searching * * * * * * * * * * * * * * * * //
627
629(
631)
632:
633 tree_(tree)
634{}
635
636
638(
640)
641:
642 tree_(tree)
643{}
644
645
647(
648 const labelUList& indices,
649 const point& sample,
650
651 scalar& nearestDistSqr,
652 label& minIndex,
653 point& nearestPoint
654) const
655{
656 for (const label index : indices)
657 {
658 const label facei = objectIndex(index);
659
660 const face& f = mesh().faces()[facei];
661
662 pointHit nearHit = f.nearestPoint(sample, mesh().points());
663 scalar distSqr = sqr(nearHit.distance());
664
665 if (distSqr < nearestDistSqr)
666 {
667 nearestDistSqr = distSqr;
668 minIndex = index;
669 nearestPoint = nearHit.point();
670 }
671 }
672}
673
674
675void Foam::treeDataFace::findNearestOp::operator()
676(
677 const labelUList& indices,
678 const point& sample,
679
680 scalar& nearestDistSqr,
681 label& minIndex,
682 point& nearestPoint
683) const
684{
685 tree_.shapes().findNearest
686 (
687 indices,
688 sample,
689 nearestDistSqr,
690 minIndex,
691 nearestPoint
692 );
693}
694
695
696void Foam::treeDataFace::findNearestOp::operator()
697(
698 const labelUList& indices,
699 const linePointRef& ln,
700
701 treeBoundBox& tightest,
702 label& minIndex,
703 point& linePoint,
704 point& nearestPoint
705) const
706{
708}
709
710
711bool Foam::treeDataFace::findIntersectOp::operator()
712(
713 const label index,
714 const point& start,
715 const point& end,
716 point& intersectionPoint
717) const
718{
719 const treeDataFace& shape = tree_.shapes();
720
721 // Do quick rejection test
722 if (shape.cacheBb_)
723 {
724 const treeBoundBox& faceBb = shape.bbs_[index];
725
726 if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
727 {
728 // start and end in same block outside of faceBb.
729 return false;
730 }
731 }
732
733 const label facei = shape.objectIndex(index);
734
735 const vector dir(end - start);
736
737 pointHit inter = shape.mesh_.faces()[facei].intersection
738 (
739 start,
740 dir,
741 shape.mesh_.faceCentres()[facei],
742 shape.mesh_.points(),
744 );
745
746 if (inter.hit() && inter.distance() <= 1)
747 {
748 // Note: no extra test on whether intersection is in front of us
749 // since using half_ray
750 intersectionPoint = inter.point();
751 return true;
752 }
753
754 return false;
755}
756
757
758// ************************************************************************* //
scalar range
labelList faceLabels(nFaceLabels)
Minimal example by using system/controlDict.functions:
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
bool hit() const noexcept
Is there a hit.
Definition pointHit.H:145
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
label index() const noexcept
Return the hit index.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition UListI.H:106
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge).
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Non-pointer based hierarchical recursive searching.
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition lineI.H:180
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Cell-face mesh analysis engine.
const vectorField & faceCentres() const
virtual const faceList & faces() const =0
Return faces.
virtual const pointField & points() const =0
Return mesh points.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Standard boundBox with extra functionality for use in octree.
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
Definition boundBoxI.H:439
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,.
direction posBits(const point &pt) const
Position of point relative to bounding box.
findIntersectOp(const indexedOctree< treeDataFace > &tree)
findNearestOp(const indexedOctree< treeDataFace > &tree)
Encapsulation of data for searching on faces.
tmp< pointField > centres() const
Representative point cloud for contained shapes. One point per shape, corresponding to the face centr...
const labelList & faceLabels() const noexcept
The subset of face ids to use.
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
static treeBoundBox bounds(const primitiveMesh &mesh, const labelRange &range)
Return bounding box of specified range of faces.
const point & centre(const label index) const
Representative point (face centre) at shape index.
label objectIndex(const label index) const
Map from shape index to original (non-subset) face label.
treeDataFace(const bool cacheBb, const primitiveMesh &mesh)
Construct from mesh, using all faces.
volumeType getVolumeType(const indexedOctree< treeDataFace > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
static treeBoundBoxList boxes(const primitiveMesh &mesh)
Calculate and return bounding boxes for all mesh faces.
const primitiveMesh & mesh() const noexcept
Reference to the supporting mesh.
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does (bb of) shape at index overlap searchBox.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition volumeType.H:56
@ UNKNOWN
Unknown state.
Definition volumeType.H:64
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
mesh update()
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
Namespace for OpenFOAM.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
static treeBoundBox boundsImpl(const primitiveMesh &mesh, const ElementIds &elemIds)
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
static treeBoundBoxList boxesImpl(const primitiveMesh &mesh, const ElementIds &elemIds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition POSIX.C:1239
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299