Loading...
Searching...
No Matches
treeDataEdge.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) 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
30#include "indexedOctree.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37}
38
39
40// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
41
42namespace Foam
43{
44
45// Bound boxes corresponding to specified edges
46template<class ElementIds>
48(
49 const edgeList& edges,
50 const pointField& points,
51 const ElementIds& elemIds
52)
53{
54 treeBoundBoxList bbs(elemIds.size());
55
56 std::transform
57 (
58 elemIds.cbegin(),
59 elemIds.cend(),
60 bbs.begin(),
61 [&](label edgei)
62 {
63 return treeBoundBox(edges[edgei].box(points));
64 }
65 );
67 return bbs;
68}
69
70
71// Overall bound box for specified edges
72template<class ElementIds>
73static treeBoundBox boundsImpl
74(
75 const edgeList& edges,
76 const pointField& points,
77 const ElementIds& elemIds
78)
79{
80 treeBoundBox bb;
81
82 for (const label edgei : elemIds)
83 {
84 const edge& e = edges[edgei];
85
86 bb.add(points[e.first()], points[e.second()]);
87 }
88
89 return bb;
90}
92} // End namespace Foam
93
94
95// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
96
99(
100 const edgeList& edges,
101 const pointField& points
102)
103{
105
106 std::transform
107 (
108 edges.cbegin(),
109 edges.cend(),
110 bbs.begin(),
111 [&](const edge& e) { return treeBoundBox(e.box(points)); }
112 );
113
114 return bbs;
115}
116
117
120(
121 const edgeList& edges,
122 const pointField& points,
123 const labelRange& range
125{
126 return boxesImpl(edges, points, range);
127}
128
129
132(
133 const edgeList& edges,
134 const pointField& points,
135 const labelUList& edgeIds
137{
138 return boxesImpl(edges, points, edgeIds);
139}
140
141
144(
145 const edgeList& edges,
146 const pointField& points
147)
148{
149 treeBoundBox bb;
150
151 for (const edge& e : edges)
152 {
153 bb.add(points[e.first()], points[e.second()]);
155
156 return bb;
157}
158
159
162(
163 const edgeList& edges,
164 const pointField& points,
165 const labelRange& range
167{
168 return boundsImpl(edges, points, range);
169}
170
171
174(
175 const edgeList& edges,
176 const pointField& points,
177 const labelUList& edgeIds
178)
179{
180 return boundsImpl(edges, points, edgeIds);
181}
182
183
184// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
185
186void Foam::treeDataEdge::update()
187{
188 if (cacheBb_)
189 {
190 if (useSubset_)
191 {
192 bbs_ = treeDataEdge::boxes(edges_, points_, edgeLabels_);
193 }
194 else
195 {
196 bbs_ = treeDataEdge::boxes(edges_, points_);
198 }
199}
200
201
202// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
203
205(
206 const bool cacheBb,
207 const edgeList& edges,
208 const pointField& points
209)
210:
211 points_(points),
212 edges_(edges),
213 edgeLabels_(),
214 useSubset_(false),
215 cacheBb_(cacheBb)
216{
217 update();
218}
219
220
222(
223 const bool cacheBb,
224 const edgeList& edges,
225 const pointField& points,
226 const labelRange& range
227)
228:
229 points_(points),
230 edges_(edges),
231 edgeLabels_(identity(range)),
232 useSubset_(true),
233 cacheBb_(cacheBb)
234{
235 update();
236}
237
238
240(
241 const bool cacheBb,
242 const edgeList& edges,
243 const pointField& points,
244 const labelUList& edgeLabels
245)
246:
247 points_(points),
248 edges_(edges),
249 edgeLabels_(edgeLabels),
250 useSubset_(true),
251 cacheBb_(cacheBb)
252{
253 update();
254}
255
256
258(
259 const bool cacheBb,
260 const edgeList& edges,
261 const pointField& points,
262 labelList&& edgeLabels
263)
264:
265 points_(points),
266 edges_(edges),
267 edgeLabels_(std::move(edgeLabels)),
268 useSubset_(true),
269 cacheBb_(cacheBb)
271 update();
272}
273
274
275// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
276
278{
279 if (useSubset_)
280 {
281 treeBoundBox bb;
282
283 for (const label index : indices)
284 {
285 const edge& e = edges_[edgeLabels_[index]];
286
287 bb.add(points_[e.first()], points_[e.second()]);
288 }
289
290 return bb;
291 }
292
293 return treeDataEdge::bounds(edges_, points_, indices);
294}
295
296
298{
299 tmp<pointField> tpts;
300
301 if (useSubset_)
302 {
303 tpts = tmp<pointField>::New(edgeLabels_.size());
304
305 std::transform
306 (
307 edgeLabels_.cbegin(),
308 edgeLabels_.cend(),
309 tpts.ref().begin(),
310 [&](label edgei) { return edges_[edgei].centre(points_); }
311 );
312 }
313 else
314 {
315 tpts = tmp<pointField>::New(edges_.size());
316
317 std::transform
318 (
319 edges_.cbegin(),
320 edges_.cend(),
321 tpts.ref().begin(),
322 [&](const edge& e) { return e.centre(points_); }
323 );
324 }
325
326 return tpts;
327}
328
329
331(
332 const indexedOctree<treeDataEdge>& oc,
334) const
335{
336 return volumeType::UNKNOWN;
337}
338
339
341(
342 const label index,
343 const treeBoundBox& searchBox
344) const
345{
346 point intersect;
347 return searchBox.intersects(this->line(index), intersect);
348}
349
350
352(
353 const label index,
354 const point& centre,
355 const scalar radiusSqr
356) const
357{
358 const pointHit nearHit = this->line(index).nearestDist(centre);
360 return (sqr(nearHit.distance()) <= radiusSqr);
361}
362
363
364// * * * * * * * * * * * * * * * * Searching * * * * * * * * * * * * * * * * //
365
367(
369)
370:
371 tree_(tree)
372{}
373
374
376(
378)
379{}
380
381
383(
384 const labelUList& indices,
385 const point& sample,
386
387 scalar& nearestDistSqr,
388 label& minIndex,
389 point& nearestPoint
390) const
391{
392 for (const label index : indices)
393 {
394 pointHit nearHit = this->line(index).nearestDist(sample);
395
396 const scalar distSqr = sqr(nearHit.distance());
397
398 if (distSqr < nearestDistSqr)
399 {
400 nearestDistSqr = distSqr;
401 minIndex = index;
402 nearestPoint = nearHit.point();
403 }
404 }
405}
406
407
408void Foam::treeDataEdge::findNearestOp::operator()
409(
410 const labelUList& indices,
411 const point& sample,
412
413 scalar& nearestDistSqr,
414 label& minIndex,
415 point& nearestPoint
416) const
417{
418 tree_.shapes().findNearest
419 (
420 indices,
421 sample,
422 nearestDistSqr,
423 minIndex,
424 nearestPoint
425 );
426}
427
428
429void Foam::treeDataEdge::findNearestOp::operator()
430(
431 const labelUList& indices,
432 const linePointRef& ln,
433
434 treeBoundBox& tightest,
435 label& minIndex,
436 point& linePoint,
437 point& nearestPoint
438) const
439{
440 const treeDataEdge& shape = tree_.shapes();
441
442 const treeBoundBox lnBb(ln.box());
443
444 // Best so far
445 scalar nearestDistSqr = linePoint.distSqr(nearestPoint);
446
447 for (const label index : indices)
448 {
449 // Note: could do bb test ? Worthwhile?
450
451 // Nearest point on line
452 point ePoint, lnPt;
453 const scalar dist = shape.line(index).nearestDist(ln, ePoint, lnPt);
454
455 const scalar distSqr = sqr(dist);
456
457 if (distSqr < nearestDistSqr)
458 {
459 nearestDistSqr = distSqr;
460 minIndex = index;
461 linePoint = lnPt;
462 nearestPoint = ePoint;
463
464 tightest = lnBb;
465 tightest.grow(dist);
466 }
467 }
468}
469
470
471bool Foam::treeDataEdge::findIntersectOp::operator()
472(
473 const label index,
474 const point& start,
475 const point& end,
476 point& result
477) const
478{
480 return false;
481}
482
483
484// ************************************************************************* //
scalar range
Minimal example by using system/controlDict.functions:
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
Definition UListI.H:468
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition UListI.H:424
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition boundBoxI.H:367
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
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
Non-pointer based hierarchical recursive searching.
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
A line primitive.
Definition line.H:180
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition lineI.H:180
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Standard boundBox with extra functionality for use in octree.
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,.
findIntersectOp(const indexedOctree< treeDataEdge > &tree)
findNearestOp(const indexedOctree< treeDataEdge > &tree)
Holds data for octree to work on an edges subset.
point centre(const label index) const
Representative point (edge centre) at shape index.
static treeBoundBoxList boxes(const edgeList &edges, const pointField &points)
Calculate and return bounding boxes for all edges.
tmp< pointField > centres() const
Representative point cloud for contained shapes. One point per shape, corresponding to the edge centr...
treeDataEdge(const bool cacheBb, const edgeList &edges, const pointField &points)
Construct from all edges.
static treeBoundBox bounds(const edgeList &edges, const pointField &points)
Return bounding box of all edges.
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
const labelList & edgeLabels() const noexcept
The subset of edge ids to use.
const edgeList & edges() const noexcept
The original list of edges.
const linePointRef line(const label index) const
Geometric line for edge at specified shape index. Frequently used.
const pointField & points() const noexcept
The reference point field.
volumeType getVolumeType(const indexedOctree< treeDataEdge > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does (bb of) shape at index searchBox.
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 defineTypeName(Type)
Define the typeName.
Definition className.H:113
mesh update()
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
const pointField & points
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
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)
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)
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...
vector point
Point is a vector.
Definition point.H:37
vectorField pointField
pointField is a vectorField.
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
Tree tree(triangles.begin(), triangles.end())
volScalarField & e