searchableCone.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) 2015-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::searchableCone
28 
29 Description
30  Searching on (optionally hollow) cone.
31 
32  \heading Dictionary parameters
33  \table
34  Property | Description | Required | Default
35  type | cone | selector |
36  point1 | coordinate of endpoint | yes |
37  radius1 | radius at point1 | yes |
38  innerRadius1| inner radius at point1 | no | 0
39  point2 | coordinate of endpoint | yes |
40  radius2 | radius at point2 | yes |
41  innerRadius2| inner radius at point2 | no | 0
42  \endtable
43 
44 Note
45  Initial implementation, might suffer from robustness (e.g. radius1==radius2)
46 
47 Note
48  Longer type name : \c searchableCone
49 
50 SourceFiles
51  searchableCone.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef searchableCone_H
56 #define searchableCone_H
57 
58 #include "searchableSurface.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 
65 /*---------------------------------------------------------------------------*\
66  Class searchableCone Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 class searchableCone
70 :
71  public searchableSurface
72 {
73  // Private Member Data
74 
75  //- 'Left' point
76  const point point1_;
77 
78  //- Outer radius at point1
79  const scalar radius1_;
80 
81  //- Inner radius at point1
82  const scalar innerRadius1_;
83 
84 
85  //- 'Right' point
86  const point point2_;
87 
88  //- Outer radius at point2
89  const scalar radius2_;
90 
91  //- Inner radius at point2
92  const scalar innerRadius2_;
93 
94 
95  //- Length of vector point2-point1
96  const scalar magDir_;
97 
98  //- Normalised vector point2-point1
99  const vector unitDir_;
100 
101 
102  //- Names of regions
103  mutable wordList regions_;
105 
106  // Private Member Functions
107 
108  //- Find nearest point on cylinder.
109  void findNearestAndNormal
110  (
111  const point& sample,
112  const scalar nearestDistSqr,
113  pointIndexHit & nearInfo,
114  vector& normal
115  ) const;
116 
117  //- Determine radial coordinate (squared)
118  static scalar radius2(const searchableCone& cone, const point& pt);
119 
120  //- Find both intersections with cone. innerRadii supplied externally.
121  void findLineAll
122  (
123  const searchableCone& cone,
124  const scalar innerRadius1,
125  const scalar innerRadius2,
126  const point& start,
127  const point& end,
128  pointIndexHit& near,
129  pointIndexHit& far
130  ) const;
131 
132  //- Insert a hit if it differs (by a tolerance) from the existing ones
133  void insertHit
134  (
135  const point& start,
136  const point& end,
138  const pointIndexHit& hit
139  ) const;
140 
141  //- Return the boundBox of the cylinder
142  boundBox calcBounds() const;
143 
144  //- No copy construct
145  searchableCone(const searchableCone&) = delete;
146 
147  //- No copy assignment
148  void operator=(const searchableCone&) = delete;
149 
150 
151 public:
152 
153  //- Runtime type information
154  TypeName("searchableCone");
155 
156 
157  // Constructors
158 
159  //- Construct from components
161  (
162  const IOobject& io,
163  const point& point1,
164  const scalar radius1,
165  const scalar innerRadius1,
166  const point& point2,
167  const scalar radius2,
168  const scalar innerRadius2
169  );
170 
171  //- Construct from dictionary (used by searchableSurface)
173  (
174  const IOobject& io,
175  const dictionary& dict
176  );
177 
178 
179  //- Destructor
180  virtual ~searchableCone() = default;
181 
182 
183  // Member Functions
184 
185  //- Names of regions
186  virtual const wordList& regions() const;
187 
188  //- Whether supports volume type (below)
189  virtual bool hasVolumeType() const
190  {
191  return true;
192  }
193 
194  //- What is type of points outside bounds
195  virtual volumeType outsideVolumeType() const
196  {
197  return volumeType::OUTSIDE;
198  }
199 
200  //- Range of local indices that can be returned.
201  virtual label size() const
202  {
203  return 1;
204  }
205 
206  //- Get representative set of element coordinates
207  // Usually the element centres (should be of length size()).
208  virtual tmp<pointField> coordinates() const;
209 
210  //- Get bounding spheres (centre and radius squared), one per element.
211  // Any point on element is guaranteed to be inside.
212  virtual void boundingSpheres
213  (
214  pointField& centres,
215  scalarField& radiusSqr
216  ) const;
217 
218  //- Get the points that define the surface.
219  virtual tmp<pointField> points() const;
220 
221  //- Does any part of the surface overlap the supplied bound box?
222  virtual bool overlaps(const boundBox& bb) const
223  {
225  return false;
226  }
227 
228 
229  // Multiple point queries.
230 
231  //- Find nearest point on cylinder
232  virtual void findNearest
233  (
234  const pointField& sample,
235  const scalarField& nearestDistSqr,
236  List<pointIndexHit>&
237  ) const;
238 
239  //- Find nearest intersection on line from start to end
240  virtual void findLine
241  (
242  const pointField& start,
243  const pointField& end,
244  List<pointIndexHit>&
245  ) const;
246 
247  //- Find all intersections in order from start to end
248  virtual void findLineAll
249  (
250  const pointField& start,
251  const pointField& end,
252  List<List<pointIndexHit>>&
253  ) const;
254 
255  //- Find any intersection on line from start to end
256  virtual void findLineAny
257  (
258  const pointField& start,
259  const pointField& end,
260  List<pointIndexHit>&
261  ) const;
262 
263  //- From a set of points and indices get the region
264  virtual void getRegion
265  (
266  const List<pointIndexHit>&,
267  labelList& region
268  ) const;
269 
270  //- From a set of points and indices get the normal
271  virtual void getNormal
272  (
273  const List<pointIndexHit>&,
274  vectorField& normal
275  ) const;
277  //- Determine type (inside/outside/mixed) for point.
278  // Unknown if cannot be determined (e.g. non-manifold surface)
279  virtual void getVolumeType
280  (
281  const pointField& points,
282  List<volumeType>& volType
283  ) const;
285 
286  // regIOobject implementation
287 
288  virtual bool writeData(Ostream&) const
289  {
291  return false;
292  }
293 };
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace Foam
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 #endif
303 
304 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< pointField > points() const
Get the points that define the surface.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find nearest intersection on line from start to end.
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
virtual bool hasVolumeType() const
Whether supports volume type (below)
virtual const wordList & regions() const
Names of regions.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) 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.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
virtual ~searchableCone()=default
Destructor.
virtual bool writeData(Ostream &) const
Pure virtual writeData function.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
Find nearest point on cylinder.
Searching on (optionally hollow) cone.
Vector< scalar > vector
Definition: vector.H:57
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
A location outside the volume.
Definition: volumeType.H:66
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition: IOobject.H:955
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
TypeName("searchableCone")
Runtime type information.
List< word > wordList
List of word.
Definition: fileName.H:59
vector point
Point is a vector.
Definition: point.H:37
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find any intersection on line from start to end.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual volumeType outsideVolumeType() const
What is type of points outside bounds.
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
virtual label size() const
Range of local indices that can be returned.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.