plane.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::plane
26 
27 Description
28  Geometric class that creates a 2D plane and can return the intersection
29  point between a line and the plane.
30 
31 SourceFiles
32  plane.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef plane_H
37 #define plane_H
38 
39 #include "point.H"
40 #include "scalarList.H"
41 #include "dictionary.H"
42 #include "line.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // Forward declaration of friend functions and operators
50 
51 class plane;
52 bool operator==(const plane&, const plane&);
53 bool operator!=(const plane&, const plane&);
54 Ostream& operator<<(Ostream&, const plane&);
55 
56 
57 /*---------------------------------------------------------------------------*\
58  Class plane Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class plane
62 {
63 public:
64 
65  //- Side of the plane
66  enum side
67  {
68  NORMAL,
70  };
71 
72  //- A direction and a reference point
73  class ray
74  {
75  point refPoint_;
76 
77  vector dir_;
78 
79  public:
80 
81  ray(const point& refPoint, const vector& dir)
82  :
83  refPoint_(refPoint),
84  dir_(dir)
85  {}
86 
87  const point& refPoint() const
88  {
89  return refPoint_;
90  }
91 
92  const vector& dir() const
93  {
94  return dir_;
95  }
96  };
97 
98 
99 private:
100 
101  // Private Data
102 
103  //- Normal
104  vector normal_;
105 
106  //- Reference point
107  point point_;
108 
109 
110  // Private Member Functions
111 
112  //- Calculates basePoint and normal vector given plane coefficients
113  void calcPntAndVec
114  (
115  const scalar C1,
116  const scalar C2,
117  const scalar C3,
118  const scalar C4
119  );
120 
121  //- Calculates basePoint and normal vector given three points
122  void calcPntAndVec
123  (
124  const point& point1,
125  const point& point2,
126  const point& point3
127  );
128 
129 
130 public:
131 
132  // Constructors
133 
134  //- Construct from normal vector through the origin
135  explicit plane(const vector& normalVector);
136 
137  //- Construct from normal vector and point in plane
138  plane(const point& basePoint, const vector& normalVector);
139 
140  //- Construct from three points in plane
141  plane
142  (
143  const point& point1,
144  const point& point2,
145  const point& point3
146  );
147 
148  //- Construct from coefficients for the plane equation,
149  // ax + by + cz + d = 0
150  explicit plane
151  (
152  const scalar a,
153  const scalar b,
154  const scalar c,
155  const scalar d
156  );
157 
158  //- Construct from dictionary
159  explicit plane(const dictionary& planeDict);
160 
161  //- Construct from Istream
162  explicit plane(Istream& is);
163 
164 
165  // Member Functions
166 
167  //- Return whether or not this plane is valid
168  inline bool valid() const
169  {
170  return normal_ != vector::zero && point_ != point::max;
171  }
172 
173  //- Return the plane normal
174  inline const vector& normal() const
175  {
176  return normal_;
177  }
178 
179  //- Return the plane base point
180  inline const point& refPoint() const
181  {
182  return point_;
183  }
184 
185  //- Return coefficients for the plane equation, ax + by + cz + d = 0
187 
188  //- Return a point on the plane
189  point aPoint() const;
190 
191  //- Return nearest point in the plane for the given point
192  point nearestPoint(const point& p) const;
193 
194  //- Return distance from the given point to the plane
195  scalar distance(const point& p) const;
196 
197  //- Return signed distance from the given point to the plane
198  scalar signedDistance(const point& p) const;
199 
200  //- Return cut coefficient for plane and line defined by
201  // origin and direction
202  scalar normalIntersect(const point& pnt0, const vector& dir) const;
203 
204  //- Return cut coefficient for plane and ray
205  inline scalar normalIntersect(const ray& r) const
206  {
207  return normalIntersect(r.refPoint(), r.dir());
208  }
209 
210  //- Return the cutting point between the plane and
211  // a line passing through the supplied points
212  template<class Point, class PointRef>
213  inline scalar lineIntersect(const line<Point, PointRef>& l) const
214  {
215  return normalIntersect(l.start(), l.vec());
216  }
217 
218  //- Return the cutting line between this plane and another.
219  // Returned as direction vector and point line goes through.
220  ray planeIntersect(const plane&) const;
221 
222  //- Return the cutting point between this plane and two other planes
223  point planePlaneIntersect(const plane&, const plane&) const;
224 
225  //- Return the side of the plane that the point is on.
226  // If the point is on the plane, then returns NORMAL.
227  side sideOfPlane(const point& p) const;
228 
229  //- Mirror the supplied point in the plane. Return the mirrored point.
230  point mirror(const point& p) const;
231 
232  //- Write to dictionary
233  void writeDict(Ostream&) const;
234 
235 
236  // friend Operators
237 
238  friend bool operator==(const plane&, const plane&);
239  friend bool operator!=(const plane&, const plane&);
240 
241 
242  // IOstream Operators
243 
244  //- Write plane properties
245  friend Ostream& operator<<(Ostream&, const plane&);
246 };
247 
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 } // End namespace Foam
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 #endif
256 
257 // ************************************************************************* //
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
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:57
static const Form zero
Definition: VectorSpace.H:113
static const Form max
Definition: VectorSpace.H:115
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A line primitive.
Definition: line.H:71
PointRef start() const
Return first vertex.
Definition: lineI.H:60
Point vec() const
Return start-end vector.
Definition: lineI.H:87
A direction and a reference point.
Definition: plane.H:73
const point & refPoint() const
Definition: plane.H:86
ray(const point &refPoint, const vector &dir)
Definition: plane.H:80
const vector & dir() const
Definition: plane.H:91
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:61
const point & refPoint() const
Return the plane base point.
Definition: plane.H:179
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:439
FixedList< scalar, 4 > planeCoeffs() const
Return coefficients for the plane equation, ax + by + cz + d = 0.
Definition: plane.C:196
friend bool operator!=(const plane &, const plane &)
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by.
Definition: plane.C:329
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: plane.C:431
bool valid() const
Return whether or not this plane is valid.
Definition: plane.H:167
plane(const vector &normalVector)
Construct from normal vector through the origin.
Definition: plane.C:80
friend Ostream & operator<<(Ostream &, const plane &)
Write plane properties.
scalar distance(const point &p) const
Return distance from the given point to the plane.
Definition: plane.C:316
const vector & normal() const
Return the plane normal.
Definition: plane.H:173
point aPoint() const
Return a point on the plane.
Definition: plane.C:243
void writeDict(Ostream &) const
Write to dictionary.
Definition: plane.C:454
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:212
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: plane.C:310
ray planeIntersect(const plane &) const
Return the cutting line between this plane and another.
Definition: plane.C:341
friend bool operator==(const plane &, const plane &)
point planePlaneIntersect(const plane &, const plane &) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:409
scalar signedDistance(const point &p) const
Return signed distance from the given point to the plane.
Definition: plane.C:322
side
Side of the plane.
Definition: plane.H:66
@ NORMAL
Definition: plane.H:67
@ FLIP
Definition: plane.H:68
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1257
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & operator<<(Ostream &, const ensightPart &)
volScalarField & p