Loading...
Searching...
No Matches
triSurfaceSearch.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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 "triSurfaceSearch.H"
30#include "triSurface.H"
31#include "PatchTools.H"
32#include "volumeType.H"
33
34// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35
36bool Foam::triSurfaceSearch::checkUniqueHit
37(
38 const pointIndexHit& currHit,
39 const UList<pointIndexHit>& hits,
40 const vector& lineVec
41) const
42{
43 // Classify the hit
44 label nearType = -1;
45 label nearLabel = -1;
46
47 const labelledTri& f = surface()[currHit.index()];
48
49 f.nearestPointClassify
50 (
51 currHit.hitPoint(),
52 surface().points(),
53 nearType,
54 nearLabel
55 );
56
57 if (nearType == triPointRef::POINT)
58 {
59 // near point
60
61 const label nearPointi = f[nearLabel];
62
63 const labelList& pointFaces =
64 surface().pointFaces()[surface().meshPointMap()[nearPointi]];
65
66 forAll(pointFaces, pI)
67 {
68 const label pointFacei = pointFaces[pI];
69
70 if (pointFacei != currHit.index())
71 {
72 forAll(hits, hI)
73 {
74 const pointIndexHit& hit = hits[hI];
75
76 if (hit.index() == pointFacei)
77 {
78 return false;
79 }
80 }
81 }
82 }
83 }
84 else if (nearType == triPointRef::EDGE)
85 {
86 // near edge
87 // check if the other face of the edge is already hit
88
89 const labelList& fEdges = surface().faceEdges()[currHit.index()];
90
91 const label edgeI = fEdges[nearLabel];
92
93 const labelList& edgeFaces = surface().edgeFaces()[edgeI];
94
95 forAll(edgeFaces, fI)
96 {
97 const label edgeFacei = edgeFaces[fI];
98
99 if (edgeFacei != currHit.index())
100 {
101 forAll(hits, hI)
102 {
103 const pointIndexHit& hit = hits[hI];
104
105 if (hit.index() == edgeFacei)
106 {
107 // Check normals
108 const vector currHitNormal =
109 surface().faceNormals()[currHit.index()];
110
111 const vector existingHitNormal =
112 surface().faceNormals()[edgeFacei];
113
114 const label signCurrHit =
115 pos0(currHitNormal & lineVec);
116
117 const label signExistingHit =
118 pos0(existingHitNormal & lineVec);
119
120 if (signCurrHit == signExistingHit)
121 {
122 return false;
123 }
124 }
125 }
126 }
127 }
128 }
130 return true;
131}
132
133
134// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135
136Foam::triSurfaceSearch::triSurfaceSearch(const triSurface& surface)
137:
138 surface_(surface),
139 tolerance_(indexedOctree<treeDataTriSurface>::perturbTol()),
140 maxTreeDepth_(10),
141 treePtr_(nullptr)
142{}
143
144
145Foam::triSurfaceSearch::triSurfaceSearch
146(
147 const triSurface& surface,
148 const dictionary& dict
149)
150:
151 surface_(surface),
152 tolerance_(indexedOctree<treeDataTriSurface>::perturbTol()),
153 maxTreeDepth_(10),
154 treePtr_(nullptr)
155{
156 // Have optional non-standard search tolerance for gappy surfaces.
157 if (dict.readIfPresent("tolerance", tolerance_) && tolerance_ > 0)
158 {
159 Info<< " using intersection tolerance " << tolerance_ << endl;
160 }
161
162 // Have optional non-standard tree-depth to limit storage.
163 if (dict.readIfPresent("maxTreeDepth", maxTreeDepth_) && maxTreeDepth_ > 0)
164 {
165 Info<< " using maximum tree depth " << maxTreeDepth_ << endl;
166 }
167}
168
169
170Foam::triSurfaceSearch::triSurfaceSearch
171(
172 const triSurface& surface,
173 const scalar tolerance,
174 const label maxTreeDepth
175)
176:
177 surface_(surface),
178 tolerance_(tolerance),
179 maxTreeDepth_(maxTreeDepth),
180 treePtr_(nullptr)
181{
182 if (tolerance_ < 0)
183 {
185 }
186}
187
188
189// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
194}
195
196
198{
199 treePtr_.clear();
200}
201
202
203// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204
207{
208 if (!treePtr_)
209 {
210 // Calculate bb without constructing local point numbering.
212
213 if (surface().size())
214 {
215 label nPoints;
216 PatchTools::calcBounds(surface(), bb, nPoints);
217
218 if (nPoints != surface().points().size())
219 {
221 << "Surface does not have compact point numbering."
222 << " Of " << surface().points().size()
223 << " only " << nPoints
224 << " are used. This might give problems in some routines."
225 << endl;
226 }
227
228 // Random number generator. Bit dodgy since not exactly random ;-)
229 Random rndGen(65431);
230
231 // Slightly extended bb. Slightly off-centred just so on symmetric
232 // geometry there are less face/edge aligned items.
233 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
234 }
235
236 const scalar oldTol =
238
239 treePtr_.reset
240 (
242 (
243 treeDataTriSurface(surface_, tolerance_),
244 bb,
245 maxTreeDepth_, // maxLevel
246 10, // leafsize
247 3.0 // duplicity
248 )
249 );
250
252 }
253
254 return *treePtr_;
255}
256
257
259{
260 if (treePtr_)
261 {
262 PackedList<2>& nodeTypes = treePtr_->nodeTypes();
263 forAll(nodeTypes, i)
264 {
265 if (nodeTypes[i] == volumeType::INSIDE)
266 {
267 nodeTypes[i] = volumeType::OUTSIDE;
268 }
269 else if (nodeTypes[i] == volumeType::OUTSIDE)
270 {
271 nodeTypes[i] = volumeType::INSIDE;
273 }
274 }
275}
276
277
278// Determine inside/outside for samples
280(
281 const pointField& samples
282) const
283{
284 boolList inside(samples.size());
285
286 forAll(samples, sampleI)
287 {
288 const point& sample = samples[sampleI];
289
290 if (!tree().bb().contains(sample))
291 {
292 inside[sampleI] = false;
293 }
294 else if (tree().getVolumeType(sample) == volumeType::INSIDE)
295 {
296 inside[sampleI] = true;
297 }
298 else
299 {
300 inside[sampleI] = false;
301 }
302 }
303 return inside;
304}
305
306
308(
309 const pointField& samples,
310 const scalarField& nearestDistSqr,
312) const
313{
314 const scalar oldTol =
316
317 const indexedOctree<treeDataTriSurface>& octree = tree();
318
319 const treeDataTriSurface::findNearestOp fOp(octree);
320
321 info.setSize(samples.size());
322
323 forAll(samples, i)
324 {
325 info[i] = octree.findNearest
326 (
327 samples[i],
328 nearestDistSqr[i],
329 fOp
330 );
331 }
332
334}
335
336
338(
339 const point& pt,
340 const vector& span
341)
342const
344 const scalar nearestDistSqr = 0.25*magSqr(span);
345
346 return tree().findNearest(pt, nearestDistSqr);
347}
348
349
351(
352 const pointField& start,
353 const pointField& end,
355) const
356{
357 const indexedOctree<treeDataTriSurface>& octree = tree();
358
359 info.setSize(start.size());
360
361 const scalar oldTol =
363
364 forAll(start, i)
365 {
366 info[i] = octree.findLine(start[i], end[i]);
367 }
368
370}
371
372
374(
375 const pointField& start,
376 const pointField& end,
377 List<pointIndexHit>& info
378) const
379{
380 const indexedOctree<treeDataTriSurface>& octree = tree();
381
382 info.setSize(start.size());
383
384 const scalar oldTol =
386
387 forAll(start, i)
388 {
389 info[i] = octree.findLineAny(start[i], end[i]);
390 }
391
393}
394
395
397(
398 const pointField& start,
399 const pointField& end,
400 List<List<pointIndexHit>>& info
401) const
402{
403 const indexedOctree<treeDataTriSurface>& octree = tree();
404
405 info.setSize(start.size());
406
407 const scalar oldTol =
409
410 // Work array
411 DynamicList<pointIndexHit> hits;
412
413 DynamicList<label> shapeMask;
414
415 treeDataTriSurface::findAllIntersectOp allIntersectOp(octree, shapeMask);
416
417 forAll(start, pointi)
418 {
419 hits.clear();
420 shapeMask.clear();
421
422 while (true)
423 {
424 // See if any intersection between pt and end
425 pointIndexHit inter = octree.findLine
426 (
427 start[pointi],
428 end[pointi],
429 allIntersectOp
430 );
431
432 if (inter.hit())
433 {
434 const vector lineVec = normalised(end[pointi] - start[pointi]);
435
436 if
437 (
438 checkUniqueHit
439 (
440 inter,
441 hits,
442 lineVec
443 )
444 )
445 {
446 hits.append(inter);
447 }
448
449 shapeMask.append(inter.index());
450 }
451 else
452 {
453 break;
454 }
455 }
456
457 info[pointi].transfer(hits);
458 }
459
461}
462
463
464// ************************************************************************* //
if(patchID !=-1)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
Minimal example by using system/controlDict.functions:
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
void setSize(label n)
Alias for resize().
Definition List.H:536
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition PackedList.H:146
static void calcBounds(const PrimitivePatch< FaceList, PointField > &p, boundBox &bb, label &nPoints)
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
const Map< label > & meshPointMap() const
Mesh point map.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const labelListList & pointFaces() const
Return point-face addressing.
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
Random number generator.
Definition Random.H:56
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static scalar & perturbTol() noexcept
Get the perturbation tolerance.
Non-pointer based hierarchical recursive searching.
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Standard boundBox with extra functionality for use in octree.
void flip()
Flip orientation (if cached on octree).
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &info) const
Calculate all intersections from start to end.
boolList calcInside(const pointField &searchPoints) const
Calculate for each searchPoint inside/outside status.
scalar tolerance() const
Return tolerance to use in searches.
const triSurface & surface() const
Return reference to the surface.
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
label maxTreeDepth() const
Return max tree depth of octree.
void clearOut()
Clear storage.
pointIndexHit nearest(const point &pt, const vector &span) const
Calculate nearest point on surface for single searchPoint. Returns.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
Definition triSurface.H:74
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
@ OUTSIDE
A location outside the volume.
Definition volumeType.H:66
@ INSIDE
A location inside the volume.
Definition volumeType.H:65
const pointField & points
label nPoints
#define WarningInFunction
Report a warning using Foam::Warning.
treeDataPrimitivePatch< triSurface > treeDataTriSurface
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
dictionary dict
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
scalarField samples(nIntervals, Zero)
Random rndGen