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