53 #ifndef searchableSphere_H 54 #define searchableSphere_H 68 class searchableSphere
70 public searchableSurface
109 const struct componentOrder order_;
118 inline static componentOrder getOrdering(
const vector&
radii);
126 (
p.x() - origin_.x()) / radii_.x(),
127 (
p.y() - origin_.y()) / radii_.y(),
128 (
p.z() - origin_.z()) / radii_.z()
137 p.x() * radii_.x() + origin_.x(),
138 p.y() * radii_.y() + origin_.y(),
139 p.z() * radii_.z() + origin_.z()
152 const scalar nearestDistSqr
166 searchableSphere(
const searchableSphere&) =
delete;
169 void operator=(
const searchableSphere&) =
delete;
221 return radii_[order_.major];
264 virtual label
size()
const 296 virtual void findNearest
318 virtual void findLineAll
virtual volumeType outsideVolumeType() const
What is type of points outside bounds.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual tmp< pointField > points() const
Get the points that define the surface.
vector surfaceNormal(const scalar theta, const scalar phi) const
Surface normal on the sphere at given location.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
virtual ~searchableSphere()=default
Destructor.
A bounding box defined in terms of min/max extrema points.
const vector & radii() const noexcept
The radii of the spheroid.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
Oblate (major = mezzo > minor)
const point & centre() const noexcept
The centre (origin) of the sphere.
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
virtual void getRegion(const List< pointIndexHit > &, labelList ®ion) const
From a set of points and indices get the region.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
point surfacePoint(const scalar theta, const scalar phi) const
A point on the sphere at given location.
Sphere (all components equal)
enum shapeType shape() const noexcept
The type of shape.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
A location outside the volume.
shapeType
The type of shape.
bool writeData(Ostream &) const
Pure virtual writeData function.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
scalar radius() const noexcept
The radius of the sphere, or major radius of the spheroid.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
TypeName("searchableSphere")
Runtime type information.
List< word > wordList
List of word.
virtual bool hasVolumeType() const
Whether supports volume type (below)
vector point
Point is a vector.
virtual const wordList & regions() const
Names of regions.
virtual void boundingSpheres(pointField ¢res, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual label size() const
Range of local indices that can be returned.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Prolate (major > mezzo = minor)