All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
treeDataEdge.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "treeDataEdge.H"
27 #include "indexedOctree.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(treeDataEdge, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 Foam::treeBoundBox Foam::treeDataEdge::calcBb(const label edgeI) const
40 {
41  const edge& e = edges_[edgeI];
42  const point& p0 = points_[e[0]];
43  const point& p1 = points_[e[1]];
44 
45  return treeBoundBox(min(p0, p1), max(p0, p1));
46 }
47 
48 
49 void Foam::treeDataEdge::update()
50 {
51  if (cacheBb_)
52  {
53  bbs_.setSize(edgeLabels_.size());
54 
55  forAll(edgeLabels_, i)
56  {
57  bbs_[i] = calcBb(edgeLabels_[i]);
58  }
59  }
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const bool cacheBb,
68  const edgeList& edges,
69  const pointField& points,
70  const labelUList& edgeLabels
71 )
72 :
73  edges_(edges),
74  points_(points),
75  edgeLabels_(edgeLabels),
76  cacheBb_(cacheBb)
77 {
78  update();
79 }
80 
81 
83 (
84  const bool cacheBb,
85  const edgeList& edges,
86  const pointField& points,
87  labelList&& edgeLabels
88 )
89 :
90  edges_(edges),
91  points_(points),
92  edgeLabels_(move(edgeLabels)),
93  cacheBb_(cacheBb)
94 {
95  update();
96 }
97 
98 
100 (
101  const indexedOctree<treeDataEdge>& tree
102 )
103 :
104  tree_(tree)
105 {}
106 
107 
109 (
110  const indexedOctree<treeDataEdge>& tree
111 )
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  pointField eMids(edgeLabels_.size());
120 
121  forAll(edgeLabels_, i)
122  {
123  const edge& e = edges_[edgeLabels_[i]];
124 
125  eMids[i] = e.centre(points_);
126  }
127  return eMids;
128 }
129 
130 
132 (
133  const indexedOctree<treeDataEdge>& oc,
134  const point& sample
135 ) const
136 {
137  return volumeType::unknown;
138 }
139 
140 
142 (
143  const label index,
144  const treeBoundBox& cubeBb
145 ) const
146 {
147  const edge& e = edges_[edgeLabels_[index]];
148 
149  const point& start = points_[e.start()];
150  const point& end = points_[e.end()];
151 
152  point intersect;
153 
154  return cubeBb.intersects(start, end, intersect);
155 }
156 
157 
159 (
160  const label index,
161  const point& centre,
162  const scalar radiusSqr
163 ) const
164 {
165  const edge& e = edges_[edgeLabels_[index]];
166 
167  const pointHit nearHit = e.line(points_).nearestDist(centre);
168 
169  const scalar distSqr = sqr(nearHit.distance());
170 
171  if (distSqr <= radiusSqr)
172  {
173  return true;
174  }
175 
176  return false;
177 }
178 
179 
180 void Foam::treeDataEdge::findNearestOp::operator()
181 (
182  const labelUList& indices,
183  const point& sample,
184 
185  scalar& nearestDistSqr,
186  label& minIndex,
187  point& nearestPoint
188 ) const
189 {
190  const treeDataEdge& shape = tree_.shapes();
191 
192  forAll(indices, i)
193  {
194  const label index = indices[i];
195 
196  const edge& e = shape.edges()[shape.edgeLabels()[index]];
197 
198  pointHit nearHit = e.line(shape.points()).nearestDist(sample);
199 
200  scalar distSqr = sqr(nearHit.distance());
201 
202  if (distSqr < nearestDistSqr)
203  {
204  nearestDistSqr = distSqr;
205  minIndex = index;
206  nearestPoint = nearHit.rawPoint();
207  }
208  }
209 }
210 
211 
212 void Foam::treeDataEdge::findNearestOp::operator()
213 (
214  const labelUList& indices,
215  const linePointRef& ln,
216 
217  treeBoundBox& tightest,
218  label& minIndex,
219  point& linePoint,
220  point& nearestPoint
221 ) const
222 {
223  const treeDataEdge& shape = tree_.shapes();
224 
225  // Best so far
226  scalar nearestDistSqr = magSqr(linePoint - nearestPoint);
227 
228  forAll(indices, i)
229  {
230  const label index = indices[i];
231 
232  const edge& e = shape.edges()[shape.edgeLabels()[index]];
233 
234  // Note: could do bb test ? Worthwhile?
235 
236  // Nearest point on line
237  point ePoint, lnPt;
238  scalar dist = e.line(shape.points()).nearestDist(ln, ePoint, lnPt);
239  scalar distSqr = sqr(dist);
240 
241  if (distSqr < nearestDistSqr)
242  {
243  nearestDistSqr = distSqr;
244  minIndex = index;
245  linePoint = lnPt;
246  nearestPoint = ePoint;
247 
248  {
249  point& minPt = tightest.min();
250  minPt = min(ln.start(), ln.end());
251  minPt.x() -= dist;
252  minPt.y() -= dist;
253  minPt.z() -= dist;
254  }
255  {
256  point& maxPt = tightest.max();
257  maxPt = max(ln.start(), ln.end());
258  maxPt.x() += dist;
259  maxPt.y() += dist;
260  maxPt.z() += dist;
261  }
262  }
263  }
264 }
265 
266 
267 bool Foam::treeDataEdge::findIntersectOp::operator()
268 (
269  const label index,
270  const point& start,
271  const point& end,
272  point& result
273 ) const
274 {
276  return false;
277 }
278 
279 
280 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A line primitive.
Definition: line.H:56
const pointField & points() const
Definition: treeDataEdge.H:174
static const Form max
Definition: VectorSpace.H:115
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
point centre(const pointField &) const
Return centre (centroid)
Definition: edgeI.H:169
static const Form min
Definition: VectorSpace.H:116
findNearestOp(const indexedOctree< treeDataEdge > &tree)
Definition: treeDataEdge.C:100
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:187
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:236
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
findIntersectOp(const indexedOctree< treeDataEdge > &tree)
Definition: treeDataEdge.C:109
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does (bb of) shape at index overlap bb.
Definition: treeDataEdge.C:142
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const Cmpt & x() const
Definition: VectorI.H:75
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:912
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
treeDataEdge(const bool cacheBb, const edgeList &edges, const pointField &points, const labelUList &edgeLabels)
Construct from selected edges. !Holds references to edges and points.
Definition: treeDataEdge.C:66
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const edgeList & edges() const
Definition: treeDataEdge.H:169
vector point
Point is a vector.
Definition: point.H:41
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
pointField shapePoints() const
Get representative point cloud for all shapes inside.
Definition: treeDataEdge.C:117
label end() const
Return end vertex label.
Definition: edgeI.H:92
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
const labelList & edgeLabels() const
Definition: treeDataEdge.H:179
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:95
label start() const
Return start vertex label.
Definition: edgeI.H:81
Namespace for OpenFOAM.
volumeType getVolumeType(const indexedOctree< treeDataEdge > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataEdge.C:132