plane.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-2016 OpenFOAM Foundation
9  Copyright (C) 2017-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 Class
28  Foam::plane
29 
30 Description
31  Geometric class that creates a 3D plane and can return the intersection
32  point between a line and the plane.
33 
34  Construction from a dictionary is driven by the \c planeType
35 
36  For \c planeType as \c pointAndNormal :
37  \verbatim
38  pointAndNormalDict
39  {
40  point <point>; // or basePoint
41  normal <vector>; // or normalVector
42  }
43  \endverbatim
44 
45  For \c planeType as \c embeddedPoints :
46  \verbatim
47  embeddedPointsDict
48  {
49  point1 <point>;
50  point2 <point>;
51  point3 <point>;
52  }
53  \endverbatim
54 
55  For \c planeType with \c planeEquation coefficients
56  \f$ ax + by + cz + d = 0 \f$ :
57 
58  \verbatim
59  planeEquationDict
60  {
61  a <scalar>;
62  b <scalar>;
63  c <scalar>;
64  d <scalar>;
65  }
66  \endverbatim
67 
68 SourceFiles
69  plane.C
70 
71 \*---------------------------------------------------------------------------*/
72 
73 #ifndef plane_H
74 #define plane_H
75 
76 #include "point.H"
77 #include "scalarList.H"
78 #include "dictionary.H"
79 #include "line.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 namespace Foam
84 {
85 
86 /*---------------------------------------------------------------------------*\
87  Class plane Declaration
88 \*---------------------------------------------------------------------------*/
89 
90 class plane
91 {
92 public:
93 
94  //- Side of the plane
95  enum side
96  {
97  FRONT = 1,
98  BACK = -1,
99  NORMAL = 1,
100  FLIP = -1
101  };
102 
103 
104  //- A reference point and direction
105  class ray
106  {
107  point pt_;
108  vector dir_;
109 
110  public:
111 
112  ray(const point& pt, const vector& dir)
113  :
114  pt_(pt),
115  dir_(dir)
116  {}
117 
118  const point& refPoint() const
119  {
120  return pt_;
121  }
122 
123  const vector& dir() const
124  {
125  return dir_;
126  }
127  };
128 
129 
130 private:
131 
132  // Private Data
133 
134  //- The unit normal of the plane
135  vector normal_;
136 
137  //- The origin or reference point for the plane
138  point origin_;
139 
140 
141  // Private Member Functions
142 
143  //- Normalise normal_ and emit error if its mag is less than VSMALL
144  // Optionally pass as test only, without normalisation.
145  void makeUnitNormal
146  (
147  const char * const caller,
148  const bool notTest = true
149  );
150 
151  //- Calculates point and normal given plane coefficients.
152  void calcFromCoeffs
153  (
154  const scalar a,
155  const scalar b,
156  const scalar c,
157  const scalar d,
158  const char * const caller
159  );
160 
161  //- Calculates point and normal vector given three points
162  //- Normal vector determined using right hand rule
163  void calcFromEmbeddedPoints
164  (
165  const point& point1,
166  const point& point2,
167  const point& point3,
168  const char * const caller
169  );
170 
171 
172 public:
173 
174  // Constructors
175 
176  //- Construct zero-initialised.
177  inline plane();
178 
179  //- Construct from normal vector through the origin.
180  // The vector is normalised to a unit vector on input.
181  explicit plane(const vector& normalVector);
182 
183  //- Construct from normal vector and point in plane.
184  // By default, the vector is normalised to a unit vector on input.
185  plane
186  (
187  const point& originPoint,
188  const vector& normalVector,
189  const bool doNormalise = true
190  );
191 
192  //- Construct from three points in plane
193  plane
194  (
195  const point& point1,
196  const point& point2,
197  const point& point3
198  );
199 
200  //- Construct from coefficients for the plane equation:
201  //- ax + by + cz + d = 0
202  explicit plane(const scalarList& coeffs);
203 
204  //- Construct from coefficients for the plane equation:
205  //- ax + by + cz + d = 0
206  explicit plane(const FixedList<scalar,4>& coeffs);
207 
208  //- Construct from dictionary
209  explicit plane(const dictionary& dict);
210 
211  //- Construct from Istream. Assumes (normal) (origin) input.
212  explicit plane(Istream& is);
213 
214 
215  // Member Functions
216 
217  //- The plane unit normal
218  inline const vector& normal() const;
219 
220  //- The plane base point
221  inline const point& origin() const;
222 
223  //- The plane base point, for modification
224  inline point& origin();
225 
226  //- The plane base point (same as origin)
227  inline const point& refPoint() const;
228 
229  //- Flip the plane by reversing the normal
230  inline void flip();
231 
232  //- Return coefficients for the plane equation:
233  //- ax + by + cz + d = 0
235 
236  //- Return nearest point in the plane for the given point
237  inline point nearestPoint(const point& p) const;
238 
239  //- Return distance (magnitude) from the given point to the plane
240  inline scalar distance(const point& p) const;
241 
242  //- Return distance from the given point to the plane
243  inline scalar signedDistance(const point& p) const;
244 
245  //- Return cut coefficient for plane and line defined by
246  //- origin and direction
247  scalar normalIntersect(const point& pnt0, const vector& dir) const;
248 
249  //- Return cut coefficient for plane and ray
250  scalar normalIntersect(const ray& r) const
251  {
252  return normalIntersect(r.refPoint(), r.dir());
253  }
254 
255  //- Return the cutting point between the plane and
256  // a line passing through the supplied points
257  template<class Point, class PointRef>
258  scalar lineIntersect(const line<Point, PointRef>& l) const
259  {
260  return normalIntersect(l.start(), l.vec());
261  }
262 
263  //- Return the cutting line between this plane and another.
264  // Returned as direction vector and point line goes through.
265  ray planeIntersect(const plane& plane2) const;
266 
267  //- Return the cutting point between this plane and two other planes
269  (
270  const plane& plane2,
271  const plane& plane3
272  )
273  const;
274 
275  //- Return a point somewhere on the plane, a distance from the base
276  point somePointInPlane(const scalar dist = 1e-3) const;
277 
278  //- Return the side of the plane that the point is on.
279  // If the point is on the plane, then returns NORMAL.
280  inline side sideOfPlane(const point& p) const;
281 
282  //- The sign for the side of the plane that the point is on.
283  // Uses the supplied tolerance for rounding around zero.
284  // \return
285  // - 0: on plane
286  // - +1: front-side
287  // - -1: back-side
288  inline int sign(const point& p, const scalar tol = SMALL) const;
289 
290  //- Mirror the supplied point in the plane. Return the mirrored point.
291  point mirror(const point& p) const;
292 
293  //- Write to dictionary
294  void writeDict(Ostream& os) const;
295 };
296 
297 
298 // IOstream Operators
299 
300 //- Write plane normal, origin
301 Ostream& operator<<(Ostream& os, const plane& pln);
302 
303 
304 // Global Operators
305 
306 //- Test for equality of origin and normal
307 inline bool operator==(const plane& a, const plane& b);
308 
309 //- Test for inequality of origin or normal
310 inline bool operator!=(const plane& a, const plane& b);
311 
312 //- Compare origin
313 inline bool operator<(const plane& a, const plane& b);
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 } // End namespace Foam
319 
320 #include "planeI.H"
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 #endif
325 
326 // ************************************************************************* //
Foam::plane::distance
scalar distance(const point &p) const
Return distance (magnitude) from the given point to the plane.
Definition: planeI.H:75
Foam::plane::planePlaneIntersect
point planePlaneIntersect(const plane &plane2, const plane &plane3) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:370
Foam::plane::normal
const vector & normal() const
The plane unit normal.
Definition: planeI.H:39
Foam::plane::writeDict
void writeDict(Ostream &os) const
Write to dictionary.
Definition: plane.C:435
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::plane::nearestPoint
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: planeI.H:69
Foam::plane::FLIP
Same as BACK.
Definition: plane.H:99
Foam::plane::ray::dir
const vector & dir() const
Definition: plane.H:122
point.H
Foam::line::start
PointRef start() const noexcept
Return first point.
Definition: lineI.H:85
Foam::plane::ray::ray
ray(const point &pt, const vector &dir)
Definition: plane.H:111
line.H
planeI.H
scalarList.H
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Foam::plane::NORMAL
Same as FRONT.
Definition: plane.H:98
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::plane::side
side
Side of the plane.
Definition: plane.H:94
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::plane::normalIntersect
scalar normalIntersect(const ray &r) const
Return cut coefficient for plane and ray.
Definition: plane.H:249
Foam::plane::sign
int sign(const point &p, const scalar tol=SMALL) const
The sign for the side of the plane that the point is on.
Definition: planeI.H:95
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::plane::signedDistance
scalar signedDistance(const point &p) const
Return distance from the given point to the plane.
Definition: planeI.H:81
Foam::plane::flip
void flip()
Flip the plane by reversing the normal.
Definition: planeI.H:63
Foam::plane::planeCoeffs
FixedList< scalar, 4 > planeCoeffs() const
Definition: plane.C:239
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:123
Foam::plane::lineIntersect
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:257
os
OBJstream os(runTime.globalPath()/outputName)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::plane::ray
A reference point and direction.
Definition: plane.H:104
Foam::plane::origin
const point & origin() const
The plane base point.
Definition: planeI.H:45
Foam::plane::planeIntersect
ray planeIntersect(const plane &plane2) const
Return the cutting line between this plane and another.
Definition: plane.C:301
Foam::plane::BACK
The back (negative normal) side of the plane.
Definition: plane.H:97
Foam::Vector< scalar >
Foam::List< scalar >
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::plane::refPoint
const point & refPoint() const
The plane base point (same as origin)
Definition: planeI.H:57
dictionary.H
Foam::plane::ray::refPoint
const point & refPoint() const
Definition: plane.H:117
Foam::line
A line primitive.
Definition: line.H:53
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::plane::FRONT
The front (positive normal) side of the plane.
Definition: plane.H:96
Foam::operator<
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
Definition: IOstreamOption.H:398
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::plane::sideOfPlane
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: planeI.H:87
Foam::plane::normalIntersect
scalar normalIntersect(const point &pnt0, const vector &dir) const
Definition: plane.C:290
Foam::plane::plane
plane()
Construct zero-initialised.
Definition: planeI.H:30
Foam::line::vec
Point vec() const
Return start-to-end vector.
Definition: lineI.H:112
Foam::plane::somePointInPlane
point somePointInPlane(const scalar dist=1e-3) const
Return a point somewhere on the plane, a distance from the base.
Definition: plane.C:392
Foam::plane::mirror
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:420