searchableSphere.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) 2018 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 "searchableSphere.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(searchableSphere, 0);
38  (
39  searchableSurface,
40  searchableSphere,
41  dict
42  );
44  (
45  searchableSurface,
46  searchableSphere,
47  dict,
48  sphere
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 Foam::pointIndexHit Foam::searchableSphere::findNearest
56 (
57  const point& sample,
58  const scalar nearestDistSqr
59 ) const
60 {
61  pointIndexHit info(false, sample, -1);
62 
63  const vector n(sample - origin_);
64  scalar magN = mag(n);
65 
66  if (nearestDistSqr >= sqr(magN - radius_))
67  {
68  if (magN < ROOTVSMALL)
69  {
70  info.rawPoint() = origin_ + vector(1,0,0)*radius_;
71  }
72  else
73  {
74  info.rawPoint() = origin_ + n/magN*radius_;
75  }
76  info.setHit();
77  info.setIndex(0);
78  }
79 
80  return info;
81 }
82 
83 
84 // From Graphics Gems - intersection of sphere with ray
85 void Foam::searchableSphere::findLineAll
86 (
87  const point& start,
88  const point& end,
89  pointIndexHit& near,
90  pointIndexHit& far
91 ) const
92 {
93  near.setMiss();
94  far.setMiss();
95 
96  vector dir(end-start);
97  scalar magSqrDir = magSqr(dir);
98 
99  if (magSqrDir > ROOTVSMALL)
100  {
101  const vector toCentre(origin_ - start);
102  scalar magSqrToCentre = magSqr(toCentre);
103 
104  dir /= Foam::sqrt(magSqrDir);
105 
106  scalar v = (toCentre & dir);
107 
108  scalar disc = sqr(radius_) - (magSqrToCentre - sqr(v));
109 
110  if (disc >= 0)
111  {
112  scalar d = Foam::sqrt(disc);
113 
114  scalar nearParam = v-d;
115 
116  if (nearParam >= 0 && sqr(nearParam) <= magSqrDir)
117  {
118  near.setHit();
119  near.setPoint(start + nearParam*dir);
120  near.setIndex(0);
121  }
122 
123  scalar farParam = v+d;
124 
125  if (farParam >= 0 && sqr(farParam) <= magSqrDir)
126  {
127  far.setHit();
128  far.setPoint(start + farParam*dir);
129  far.setIndex(0);
130  }
131  }
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137 
138 Foam::searchableSphere::searchableSphere
139 (
140  const IOobject& io,
141  const point& origin,
142  const scalar radius
143 )
144 :
145  searchableSurface(io),
146  origin_(origin),
147  radius_(radius)
148 {
149  bounds() = boundBox
150  (
151  origin_ - radius_*vector::one,
152  origin_ + radius_*vector::one
153  );
154 }
155 
156 
157 Foam::searchableSphere::searchableSphere
158 (
159  const IOobject& io,
160  const dictionary& dict
161 )
162 :
164  (
165  io,
166  dict.getCompat<vector>("origin", {{"centre", -1806}}),
167  dict.get<scalar>("radius")
168  )
169 {}
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
175 {
176  return bb.overlaps(origin_, sqr(radius_));
177 }
178 
179 
181 {
182  if (regions_.empty())
183  {
184  regions_.resize(1);
185  regions_.first() = "region0";
186  }
187  return regions_;
188 }
189 
190 
191 
193 (
194  pointField& centres,
195  scalarField& radiusSqr
196 ) const
197 {
198  centres.resize(1);
199  centres[0] = origin_;
200 
201  radiusSqr.resize(1);
202  radiusSqr[0] = Foam::sqr(radius_);
203 
204  // Add a bit to make sure all points are tested inside
205  radiusSqr += Foam::sqr(SMALL);
206 }
207 
208 
209 void Foam::searchableSphere::findNearest
210 (
211  const pointField& samples,
212  const scalarField& nearestDistSqr,
213  List<pointIndexHit>& info
214 ) const
215 {
216  info.setSize(samples.size());
217 
218  forAll(samples, i)
219  {
220  info[i] = findNearest(samples[i], nearestDistSqr[i]);
221  }
222 }
223 
224 
226 (
227  const pointField& start,
228  const pointField& end,
229  List<pointIndexHit>& info
230 ) const
231 {
232  info.setSize(start.size());
233 
235 
236  forAll(start, i)
237  {
238  // Pick nearest intersection. If none intersected take second one.
239  findLineAll(start[i], end[i], info[i], b);
240  if (!info[i].hit() && b.hit())
241  {
242  info[i] = b;
243  }
244  }
245 }
246 
247 
249 (
250  const pointField& start,
251  const pointField& end,
252  List<pointIndexHit>& info
253 ) const
254 {
255  info.setSize(start.size());
256 
258 
259  forAll(start, i)
260  {
261  // Discard far intersection
262  findLineAll(start[i], end[i], info[i], b);
263  if (!info[i].hit() && b.hit())
264  {
265  info[i] = b;
266  }
267  }
268 }
269 
270 
271 void Foam::searchableSphere::findLineAll
272 (
273  const pointField& start,
274  const pointField& end,
276 ) const
277 {
278  info.setSize(start.size());
279 
280  forAll(start, i)
281  {
282  pointIndexHit near, far;
283  findLineAll(start[i], end[i], near, far);
284 
285  if (near.hit())
286  {
287  if (far.hit())
288  {
289  info[i].setSize(2);
290  info[i][0] = near;
291  info[i][1] = far;
292  }
293  else
294  {
295  info[i].setSize(1);
296  info[i][0] = near;
297  }
298  }
299  else
300  {
301  if (far.hit())
302  {
303  info[i].setSize(1);
304  info[i][0] = far;
305  }
306  else
307  {
308  info[i].clear();
309  }
310  }
311  }
312 }
313 
314 
316 (
317  const List<pointIndexHit>& info,
318  labelList& region
319 ) const
320 {
321  region.setSize(info.size());
322  region = 0;
323 }
324 
325 
327 (
328  const List<pointIndexHit>& info,
329  vectorField& normal
330 ) const
331 {
332  normal.setSize(info.size());
333  normal = Zero;
334 
335  forAll(info, i)
336  {
337  if (info[i].hit())
338  {
339  normal[i] = normalised(info[i].hitPoint() - origin_);
340  }
341  else
342  {
343  // Set to what?
344  }
345  }
346 }
347 
348 
350 (
351  const pointField& points,
352  List<volumeType>& volType
353 ) const
354 {
355  volType.setSize(points.size());
356 
357  const scalar rad2 = sqr(radius_);
358 
359  forAll(points, pointi)
360  {
361  const point& pt = points[pointi];
362 
363  volType[pointi] =
364  (
365  (magSqr(pt - origin_) <= rad2)
367  );
368  }
369 }
370 
371 
372 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::one
static const Vector< scalar > one
Definition: VectorSpace.H:116
Foam::searchableSphere::findLine
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
Definition: searchableSphere.C:226
searchableSphere.H
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::searchableSphere::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: searchableSphere.C:316
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:107
Foam::searchableSphere::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchableSphere.C:180
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::getCompat
T getCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, enum keyType::option=keyType::REGEX) const
Definition: dictionaryTemplates.C:108
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:519
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
sphere
Specialization of rigidBody to construct a sphere given the mass and radius.
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:55
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::searchableSphere::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: searchableSphere.C:327
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
samples
scalarField samples(nIntervals, Zero)
Foam::searchableSphere::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchableSphere.C:193
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::searchableSphere::findLineAny
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
Definition: searchableSphere.C:249
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::searchableSphere::getVolumeType
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
Definition: searchableSphere.C:350
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:115
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:496
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::Vector< scalar >
Foam::List< word >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::volumeType::INSIDE
A location inside the volume.
Definition: volumeType.H:68
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::searchableSphere::overlaps
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Definition: searchableSphere.C:174
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::volumeType::OUTSIDE
A location outside the volume.
Definition: volumeType.H:69
Foam::searchableSphere
Searching on sphere.
Definition: searchableSphere.H:86