searchableSphere.H
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-2020 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
27Class
28 Foam::searchableSphere
29
30Description
31 Searching on general spheroid.
32
33 \heading Dictionary parameters
34 \table
35 Property | Description | Required | Default
36 type | sphere | selector |
37 origin | The origin (centre) of the sphere | yes |
38 radius | The (outside) radius/radiii of sphere | yes |
39 centre | Alternative name for 'origin' | no |
40 \endtable
41
42Note
43 The \c radius can be specified as a single \em scalar (for a sphere)
44 or a \em vector of three values (for a general spheroid).
45
46 Longer type name : \c searchableSphere
47
48SourceFiles
49 searchableSphere.C
50
51\*---------------------------------------------------------------------------*/
52
53#ifndef searchableSphere_H
54#define searchableSphere_H
55
56#include "treeBoundBox.H"
57#include "searchableSurface.H"
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61namespace Foam
62{
63
64/*---------------------------------------------------------------------------*\
65 Class searchableSphere Declaration
66\*---------------------------------------------------------------------------*/
67
68class searchableSphere
69:
70 public searchableSurface
71{
72public:
73
74 // Public Types
75
76 //- The type of shape
77 enum shapeType : uint8_t
78 {
79 SPHERE = 0,
80 OBLATE = 1,
81 PROLATE = 2,
82 GENERAL = 3
83 };
84
85
86private:
87
88 // Private Types
89
90 //- Component order (largest to smallest)
91 struct componentOrder
92 {
93 uint8_t major;
94 uint8_t mezzo;
95 uint8_t minor;
96 shapeType shape;
97 };
98
99
100 // Private Data
102 //- Centre point of the sphere
103 const point origin_;
105 //- The outer radii of the spheroid
106 const vector radii_;
108 //- The canonical (sorted) order and shape
109 const struct componentOrder order_;
110
111 //- Names of regions
112 mutable wordList regions_;
113
114
115 // Private Member Functions
116
117 //- Determine sorted order and classify the shape
118 inline static componentOrder getOrdering(const vector& radii);
119
120 //- Shift point relative to origin
121 //- and scale relative to spheroid dimensions
122 inline point scalePoint(const point& p) const
123 {
124 return point
125 (
126 (p.x() - origin_.x()) / radii_.x(),
127 (p.y() - origin_.y()) / radii_.y(),
128 (p.z() - origin_.z()) / radii_.z()
129 );
130 }
131
132 //- Undo scalePoint(): unscale point and unshift relative to origin
133 inline point unscalePoint(const point& p) const
134 {
135 return point
136 (
137 p.x() * radii_.x() + origin_.x(),
138 p.y() * radii_.y() + origin_.y(),
139 p.z() * radii_.z() + origin_.z()
140 );
141 }
142
143
144 //- Inherit findNearest from searchableSurface
146
147 //- Find nearest point on general spheroid.
148 // With some optimization for special shapes
149 pointIndexHit findNearest
150 (
151 const point& sample,
152 const scalar nearestDistSqr
153 ) const;
154
155 //- Find intersection with general spheroid
156 void findLineAll
157 (
158 const point& start,
159 const point& end,
160 pointIndexHit& near,
161 pointIndexHit& far
162 ) const;
163
164
165 //- No copy construct
166 searchableSphere(const searchableSphere&) = delete;
167
168 //- No copy assignment
169 void operator=(const searchableSphere&) = delete;
170
171
172public:
173
174 //- Runtime type information
175 TypeName("searchableSphere");
176
177
178 // Constructors
179
180 //- Construct a sphere from components
181 searchableSphere
182 (
183 const IOobject& io,
184 const point& origin,
185 const scalar radius
186 );
187
188 //- Construct a spheroid from components
189 searchableSphere
190 (
191 const IOobject& io,
192 const point& centre,
193 const vector& radii
194 );
195
196 //- Construct from dictionary (used by searchableSurface)
197 searchableSphere
198 (
199 const IOobject& io,
200 const dictionary& dict
201 );
202
203
204 //- Destructor
205 virtual ~searchableSphere() = default;
206
207
208 // Member Functions
209
210 // Geometric
211
212 //- The centre (origin) of the sphere
213 const point& centre() const noexcept
214 {
215 return origin_;
216 }
217
218 //- The radius of the sphere, or major radius of the spheroid
219 scalar radius() const noexcept
220 {
221 return radii_[order_.major];
222 }
223
224 //- The radii of the spheroid
225 const vector& radii() const noexcept
226 {
227 return radii_;
228 }
230 //- The type of shape
231 enum shapeType shape() const noexcept
232 {
233 return order_.shape;
234 }
235
236 //- A point on the sphere at given location
237 // theta [-pi,pi], phi [0,pi]
238 point surfacePoint(const scalar theta, const scalar phi) const;
239
240 //- Surface normal on the sphere at given location
241 // theta [-pi,pi], phi [0,pi]
242 vector surfaceNormal(const scalar theta, const scalar phi) const;
244
245
246 // Searching
247
248 //- Names of regions
249 virtual const wordList& regions() const;
250
251 //- Whether supports volume type (below)
252 virtual bool hasVolumeType() const
253 {
254 return true;
256
257 //- What is type of points outside bounds
258 virtual volumeType outsideVolumeType() const
259 {
260 return volumeType::OUTSIDE;
261 }
262
263 //- Range of local indices that can be returned.
264 virtual label size() const
265 {
266 return 1;
267 }
268
269 //- Get representative set of element coordinates
270 // Usually the element centres (should be of length size()).
271 virtual tmp<pointField> coordinates() const
272 {
273 return tmp<pointField>::New(1, origin_);
274 }
275
276 //- Get bounding spheres (centre and radius squared), one per element.
277 // Any point on element is guaranteed to be inside.
278 virtual void boundingSpheres
279 (
280 pointField& centres,
281 scalarField& radiusSqr
282 ) const;
283
284 //- Get the points that define the surface.
285 virtual tmp<pointField> points() const
286 {
287 return coordinates();
289
290 //- Does any part of the surface overlap the supplied bound box?
291 virtual bool overlaps(const boundBox& bb) const;
292
293
294 // Multiple point queries.
296 virtual void findNearest
297 (
298 const pointField& sample,
299 const scalarField& nearestDistSqr,
301 ) const;
302
303 virtual void findLine
304 (
305 const pointField& start,
306 const pointField& end,
308 ) const;
310 virtual void findLineAny
311 (
312 const pointField& start,
313 const pointField& end,
315 ) const;
316
317 //- Get all intersections in order from start to end.
318 virtual void findLineAll
319 (
320 const pointField& start,
321 const pointField& end,
323 ) const;
324
325 //- From a set of points and indices get the region
326 virtual void getRegion
327 (
328 const List<pointIndexHit>&,
329 labelList& region
330 ) const;
331
332 //- From a set of points and indices get the normal
333 virtual void getNormal
334 (
335 const List<pointIndexHit>&,
336 vectorField& normal
337 ) const;
338
339 //- Determine type (inside/outside/mixed) for point.
340 virtual void getVolumeType
341 (
342 const pointField& points,
343 List<volumeType>& volType
344 ) const;
345
346
347 // Output
348
349 // Implementation for regIOobject. NotImplemented
350 bool writeData(Ostream&) const
351 {
353 return false;
354 }
355};
356
357
358// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
359
360} // End namespace Foam
361
362// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363
364#endif
365
366// ************************************************************************* //
#define minor(dev)
Definition: fileStat.C:40
#define major(dev)
Definition: fileStat.C:39
surfaceScalarField & phi
Minimal example by using system/controlDict.functions:
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
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
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Searching on general spheroid.
virtual label size() const
Range of local indices that can be returned.
TypeName("searchableSphere")
Runtime type information.
const vector & radii() const noexcept
The radii of the spheroid.
bool writeData(Ostream &) const
Pure virtual writeData function.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
vector surfaceNormal(const scalar theta, const scalar phi) const
Surface normal on the sphere at given location.
enum shapeType shape() const noexcept
The type of shape.
virtual ~searchableSphere()=default
Destructor.
const point & centre() const noexcept
The centre (origin) of the sphere.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
scalar radius() const noexcept
The radius of the sphere, or major radius of the spheroid.
virtual tmp< pointField > points() const
Get the points that define the surface.
virtual volumeType outsideVolumeType() const
What is type of points outside bounds.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
virtual const wordList & regions() const
Names of regions.
point surfacePoint(const scalar theta, const scalar phi) const
A point on the sphere at given location.
shapeType
The type of shape.
@ PROLATE
Prolate (major > mezzo = minor)
@ GENERAL
General spheroid.
@ SPHERE
Sphere (all components equal)
@ OBLATE
Oblate (major = mezzo > minor)
virtual bool hasVolumeType() const
Whether supports volume type (below)
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
A class for managing temporary objects.
Definition: tmp.H:65
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:61
@ OUTSIDE
A location outside the volume.
Definition: volumeType.H:69
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const direction noexcept
Definition: Scalar.H:223
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73