treeDataPoint.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) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "treeDataPoint.H"
30 #include "treeBoundBox.H"
31 #include "indexedOctree.H"
32 #include "triangleFuncs.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(treeDataPoint, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 :
46  points_(points),
47  useSubset_(false)
48 {}
49 
50 
52 (
53  const pointField& points,
54  const labelUList& pointLabels,
55  const bool useSubsetPoints
56 )
57 :
58  points_(points),
59  pointLabels_(pointLabels),
60  useSubset_(useSubsetPoints)
61 {}
62 
63 
65 (
66  const pointField& points,
68  const bool useSubsetPoints
69 )
70 :
71  points_(points),
72  pointLabels_(std::move(pointLabels)),
73  useSubset_(useSubsetPoints)
74 {}
75 
76 
78 (
80 )
81 :
82  tree_(tree)
83 {}
84 
85 
87 (
89 )
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (useSubset_)
98  {
99  return pointField(points_, pointLabels_);
100  }
101 
102  return points_;
103 }
104 
105 
107 (
109  const point& sample
110 ) const
111 {
112  return volumeType::UNKNOWN;
113 }
114 
115 
117 (
118  const label index,
119  const treeBoundBox& cubeBb
120 ) const
121 {
122  return cubeBb.contains(shapePoint(index));
123 }
124 
125 
127 (
128  const label index,
129  const point& centre,
130  const scalar radiusSqr
131 ) const
132 {
133  return (magSqr(shapePoint(index) - centre) <= radiusSqr);
134 }
135 
136 
137 void Foam::treeDataPoint::findNearestOp::operator()
138 (
139  const labelUList& indices,
140  const point& sample,
141 
142  scalar& nearestDistSqr,
143  label& minIndex,
144  point& nearestPoint
145 ) const
146 {
147  const treeDataPoint& shape = tree_.shapes();
148 
149  for (const label index : indices)
150  {
151  const point& pt = shape.shapePoint(index);
152 
153  const scalar distSqr = magSqr(pt - sample);
154 
155  if (distSqr < nearestDistSqr)
156  {
157  nearestDistSqr = distSqr;
158  minIndex = index;
159  nearestPoint = pt;
160  }
161  }
162 }
163 
164 
165 void Foam::treeDataPoint::findNearestOp::operator()
166 (
167  const labelUList& indices,
168  const linePointRef& ln,
169 
170  treeBoundBox& tightest,
171  label& minIndex,
172  point& linePoint,
173  point& nearestPoint
174 ) const
175 {
176  const treeDataPoint& shape = tree_.shapes();
177 
178  // Best so far
179  scalar nearestDistSqr = GREAT;
180  if (minIndex >= 0)
181  {
182  nearestDistSqr = magSqr(linePoint - nearestPoint);
183  }
184 
185  for (const label index : indices)
186  {
187  const point& shapePt = shape.shapePoint(index);
188 
189  if (tightest.contains(shapePt))
190  {
191  // Nearest point on line
192  pointHit pHit = ln.nearestDist(shapePt);
193  const scalar distSqr = sqr(pHit.distance());
194 
195  if (distSqr < nearestDistSqr)
196  {
197  nearestDistSqr = distSqr;
198  minIndex = index;
199  linePoint = pHit.rawPoint();
200  nearestPoint = shapePt;
201 
202  {
203  point& minPt = tightest.min();
204  minPt = min(ln.start(), ln.end());
205  minPt.x() -= pHit.distance();
206  minPt.y() -= pHit.distance();
207  minPt.z() -= pHit.distance();
208  }
209  {
210  point& maxPt = tightest.max();
211  maxPt = max(ln.start(), ln.end());
212  maxPt.x() += pHit.distance();
213  maxPt.y() += pHit.distance();
214  maxPt.z() += pHit.distance();
215  }
216  }
217  }
218  }
219 }
220 
221 
222 bool Foam::treeDataPoint::findIntersectOp::operator()
223 (
224  const label index,
225  const point& start,
226  const point& end,
227  point& result
228 ) const
229 {
231  return false;
232 }
233 
234 
235 // ************************************************************************* //
Foam::treeDataPoint::shapePoint
const point & shapePoint(const label index) const
Point at specified index.
Definition: treeDataPoint.H:210
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::treeDataPoint::findIntersectOp::findIntersectOp
findIntersectOp(const indexedOctree< treeDataPoint > &tree)
Definition: treeDataPoint.C:87
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
indexedOctree.H
Foam::PointHit::distance
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
Foam::treeDataPoint::treeDataPoint
treeDataPoint(const pointField &points)
Construct from pointField.
Definition: treeDataPoint.C:44
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::treeDataPoint
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
Definition: treeDataPoint.H:62
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Foam::VectorSpace::min
static const Form min
Definition: VectorSpace.H:118
Foam::volumeType
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:60
Foam::Field< vector >
treeDataPoint.H
treeBoundBox.H
Foam::volumeType::UNKNOWN
Unknown state.
Definition: volumeType.H:67
Foam::treeDataPoint::shapePoints
pointField shapePoints() const
Definition: treeDataPoint.C:95
Foam::PointHit::rawPoint
const point_type & rawPoint() const noexcept
The point, no checks.
Definition: PointHit.H:172
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
triangleFuncs.H
Foam::treeDataPoint::findNearestOp::findNearestOp
findNearestOp(const indexedOctree< treeDataPoint > &tree)
Definition: treeDataPoint.C:78
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::treeBoundBox::contains
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:320
Foam::VectorSpace::max
static const Form max
Definition: VectorSpace.H:117
Foam::Vector< scalar >
Foam::List< label >
Foam::UList< label >
Foam::treeDataPoint::overlaps
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does (bb of) shape at index overlap bb.
Definition: treeDataPoint.C:117
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::line
A line primitive.
Definition: line.H:53
Foam::treeDataPoint::getVolumeType
volumeType getVolumeType(const indexedOctree< treeDataPoint > &os, const point &sample) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataPoint.C:107
Foam::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:925
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointLabels
labelList pointLabels(nPoints, -1)
sample
Minimal example by using system/controlDict.functions: