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