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-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
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 }
129
130 return true;
131}
132
133
134// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135
137:
138 surface_(surface),
139 tolerance_(indexedOctree<treeDataTriSurface>::perturbTol()),
140 maxTreeDepth_(10),
141 treePtr_(nullptr)
142{}
143
144
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
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 * * * * * * * * * * * * * * * //
190
192{
193 clearOut();
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 = bb.extend(rndGen, 1e-4);
234 bb.min() -= point::uniform(ROOTVSMALL);
235 bb.max() += point::uniform(ROOTVSMALL);
236 }
237
238 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
240
241 treePtr_.reset
242 (
244 (
245 treeDataTriSurface(false, surface_, tolerance_),
246 bb,
247 maxTreeDepth_, // maxLevel
248 10, // leafsize
249 3.0 // duplicity
250 )
251 );
252
254 }
255
256 return *treePtr_;
257}
258
259
260// Determine inside/outside for samples
262(
263 const pointField& samples
264) const
265{
266 boolList inside(samples.size());
267
268 forAll(samples, sampleI)
269 {
270 const point& sample = samples[sampleI];
271
272 if (!tree().bb().contains(sample))
273 {
274 inside[sampleI] = false;
275 }
276 else if (tree().getVolumeType(sample) == volumeType::INSIDE)
277 {
278 inside[sampleI] = true;
279 }
280 else
281 {
282 inside[sampleI] = false;
283 }
284 }
285 return inside;
286}
287
288
290(
291 const pointField& samples,
292 const scalarField& nearestDistSqr,
294) const
295{
296 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
298
299 const indexedOctree<treeDataTriSurface>& octree = tree();
300
301 const treeDataTriSurface::findNearestOp fOp(octree);
302
303 info.setSize(samples.size());
304
305 forAll(samples, i)
306 {
307 info[i] = octree.findNearest
308 (
309 samples[i],
310 nearestDistSqr[i],
311 fOp
312 );
313 }
314
316}
317
318
320(
321 const point& pt,
322 const vector& span
323)
324const
325{
326 const scalar nearestDistSqr = 0.25*magSqr(span);
327
328 return tree().findNearest(pt, nearestDistSqr);
329}
330
331
333(
334 const pointField& start,
335 const pointField& end,
337) const
338{
339 const indexedOctree<treeDataTriSurface>& octree = tree();
340
341 info.setSize(start.size());
342
343 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
345
346 forAll(start, i)
347 {
348 info[i] = octree.findLine(start[i], end[i]);
349 }
350
352}
353
354
356(
357 const pointField& start,
358 const pointField& end,
360) const
361{
362 const indexedOctree<treeDataTriSurface>& octree = tree();
363
364 info.setSize(start.size());
365
366 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
368
369 forAll(start, i)
370 {
371 info[i] = octree.findLineAny(start[i], end[i]);
372 }
373
375}
376
377
379(
380 const pointField& start,
381 const pointField& end,
383) const
384{
385 const indexedOctree<treeDataTriSurface>& octree = tree();
386
387 info.setSize(start.size());
388
389 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
391
392 // Work array
394
395 DynamicList<label> shapeMask;
396
397 treeDataTriSurface::findAllIntersectOp allIntersectOp(octree, shapeMask);
398
399 forAll(start, pointi)
400 {
401 hits.clear();
402 shapeMask.clear();
403
404 while (true)
405 {
406 // See if any intersection between pt and end
407 pointIndexHit inter = octree.findLine
408 (
409 start[pointi],
410 end[pointi],
411 allIntersectOp
412 );
413
414 if (inter.hit())
415 {
416 const vector lineVec = normalised(end[pointi] - start[pointi]);
417
418 if
419 (
420 checkUniqueHit
421 (
422 inter,
423 hits,
424 lineVec
425 )
426 )
427 {
428 hits.append(inter);
429 }
430
431 shapeMask.append(inter.index());
432 }
433 else
434 {
435 break;
436 }
437 }
438
439 info[pointi].transfer(hits);
440 }
441
443}
444
445
446// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Minimal example by using system/controlDict.functions:
void setSize(const label n)
Alias for resize()
Definition: List.H:218
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 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:60
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
static scalar & perturbTol()
Get the perturbation tolerance.
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
label nearest() const
The nearest control point, or -1 if invalid.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Encapsulation of data needed to search on PrimitivePatches.
Helper class to search on triSurface.
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.
const triSurface & surface() const
Return reference to the surface.
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
~triSurfaceSearch()
Destructor.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void clearOut()
Clear storage.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
Definition: triSurface.H:79
@ POINT
Close to point.
Definition: triangle.H:96
@ EDGE
Close to edge.
Definition: triangle.H:97
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
@ INSIDE
A location inside the volume.
Definition: volumeType.H:68
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:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
scalarField samples(nIntervals, Zero)
Random rndGen
Definition: createFields.H:23