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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "treeDataEdge.H"
29#include "indexedOctree.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
36}
37
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41Foam::treeBoundBox Foam::treeDataEdge::calcBb(const label edgeI) const
42{
43 const edge& e = edges_[edgeI];
44 const point& p0 = points_[e[0]];
45 const point& p1 = points_[e[1]];
46
47 return treeBoundBox(min(p0, p1), max(p0, p1));
48}
49
50
51void Foam::treeDataEdge::update()
52{
53 if (cacheBb_)
54 {
55 bbs_.setSize(edgeLabels_.size());
56
57 forAll(edgeLabels_, i)
58 {
59 bbs_[i] = calcBb(edgeLabels_[i]);
60 }
61 }
62}
63
64
65// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66
68(
69 const bool cacheBb,
70 const edgeList& edges,
71 const pointField& points,
72 const labelUList& edgeLabels
73)
74:
75 edges_(edges),
76 points_(points),
77 edgeLabels_(edgeLabels),
78 cacheBb_(cacheBb)
79{
80 update();
81}
82
83
85(
86 const bool cacheBb,
87 const edgeList& edges,
88 const pointField& points,
89 labelList&& edgeLabels
90)
91:
92 edges_(edges),
93 points_(points),
94 edgeLabels_(std::move(edgeLabels)),
95 cacheBb_(cacheBb)
96{
97 update();
98}
99
100
102(
104)
105:
106 tree_(tree)
107{}
108
109
111(
113)
114{}
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
120{
121 pointField eMids(edgeLabels_.size());
122
123 forAll(edgeLabels_, i)
124 {
125 const edge& e = edges_[edgeLabels_[i]];
126
127 eMids[i] = e.centre(points_);
128 }
129 return eMids;
130}
131
132
134(
136 const point& sample
137) const
138{
139 return volumeType::UNKNOWN;
140}
141
142
144(
145 const label index,
146 const treeBoundBox& cubeBb
147) const
148{
149 const edge& e = edges_[edgeLabels_[index]];
150
151 const point& start = points_[e.start()];
152 const point& end = points_[e.end()];
153
154 point intersect;
155
156 return cubeBb.intersects(start, end, intersect);
157}
158
159
161(
162 const label index,
163 const point& centre,
164 const scalar radiusSqr
165) const
166{
167 const edge& e = edges_[edgeLabels_[index]];
168
169 const pointHit nearHit = e.line(points_).nearestDist(centre);
170
171 const scalar distSqr = sqr(nearHit.distance());
172
173 return (distSqr <= radiusSqr);
174}
175
176
178(
179 const labelUList& indices,
180 const point& sample,
181
182 scalar& nearestDistSqr,
183 label& minIndex,
184 point& nearestPoint
185) const
186{
187 const treeDataEdge& shape = tree_.shapes();
188
189 for (const label index : indices)
190 {
191 const edge& e = shape.edges()[shape.edgeLabels()[index]];
192
193 pointHit nearHit = e.line(shape.points()).nearestDist(sample);
194
195 const scalar distSqr = sqr(nearHit.distance());
196
197 if (distSqr < nearestDistSqr)
198 {
199 nearestDistSqr = distSqr;
200 minIndex = index;
201 nearestPoint = nearHit.rawPoint();
202 }
203 }
204}
205
206
208(
209 const labelUList& indices,
210 const linePointRef& ln,
211
212 treeBoundBox& tightest,
213 label& minIndex,
214 point& linePoint,
215 point& nearestPoint
216) const
217{
218 const treeDataEdge& shape = tree_.shapes();
219
220 // Best so far
221 scalar nearestDistSqr = magSqr(linePoint - nearestPoint);
222
223 for (const label index : indices)
224 {
225 const edge& e = shape.edges()[shape.edgeLabels()[index]];
226
227 // Note: could do bb test ? Worthwhile?
228
229 // Nearest point on line
230 point ePoint, lnPt;
231 scalar dist = e.line(shape.points()).nearestDist(ln, ePoint, lnPt);
232 scalar distSqr = sqr(dist);
233
234 if (distSqr < nearestDistSqr)
235 {
236 nearestDistSqr = distSqr;
237 minIndex = index;
238 linePoint = lnPt;
239 nearestPoint = ePoint;
240
241 {
242 point& minPt = tightest.min();
243 minPt = min(ln.start(), ln.end());
244 minPt.x() -= dist;
245 minPt.y() -= dist;
246 minPt.z() -= dist;
247 }
248 {
249 point& maxPt = tightest.max();
250 maxPt = max(ln.start(), ln.end());
251 maxPt.x() += dist;
252 maxPt.y() += dist;
253 maxPt.z() += dist;
254 }
255 }
256 }
257}
258
259
261(
262 const label index,
263 const point& start,
264 const point& end,
265 point& result
266) const
267{
269 return false;
270}
271
272
273// ************************************************************************* //
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
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static const Form max
Definition: VectorSpace.H:117
static const Form min
Definition: VectorSpace.H:118
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
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 > &)
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
A line primitive.
Definition: line.H:68
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
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,.
Definition: treeBoundBox.C:162
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:57
const edgeList & edges() const
Definition: treeDataEdge.H:173
const labelList & edgeLabels() const
Definition: treeDataEdge.H:183
volumeType getVolumeType(const indexedOctree< treeDataEdge > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataEdge.C:134
pointField shapePoints() const
Definition: treeDataEdge.C:119
const pointField & points() const
Definition: treeDataEdge.H:178
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 defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & p0
Definition: EEqn.H:36
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
const pointField & points
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
vector point
Point is a vector.
Definition: point.H:43
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
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)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333