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-2020 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  {
69  FLIP
70  };
71 
72  //- A direction and a reference point
73  class ray
74  {
75  point refPoint_;
76  vector dir_;
77 
78  public:
79 
80  ray(const point& refPoint, const vector& dir)
81  :
82  refPoint_(refPoint),
83  dir_(dir)
84  {}
85 
86  const point& refPoint() const
87  {
88  return refPoint_;
89  }
90 
91  const vector& dir() const
92  {
93  return dir_;
94  }
95  };
96 
97 
98 private:
99 
100  // Private Data
101 
102  //- Normal
103  vector normal_;
104 
105  //- Reference point
106  point point_;
107 
108 
109  // Private Member Functions
110 
111  //- Calculates basePoint and normal vector given plane coefficients
112  void calcPntAndVec(const scalarList& C);
113 
114  //- Calculates basePoint and normal vector given three points
115  //- Normal vector determined using right hand rule
116  void calcPntAndVec
117  (
118  const point& point1,
119  const point& point2,
120  const point& point3
121  );
122 
123 
124 public:
125 
126  // Constructors
127 
128  //- Construct from normal vector through the origin
129  explicit plane(const vector& normalVector);
130 
131  //- Construct from normal vector and point in plane
132  plane(const point& basePoint, const vector& normalVector);
133 
134  //- Construct from three points in plane
135  plane
136  (
137  const point& point1,
138  const point& point2,
139  const point& point3
140  );
141 
142  //- Construct from coefficients for the
143  // plane equation: ax + by + cz + d = 0
144  explicit plane(const scalarList& C);
145 
146  //- Construct from dictionary
147  explicit plane(const dictionary& planeDict);
148 
149  //- Construct from Istream. Assumes the base + normal notation.
150  explicit plane(Istream& is);
151 
152 
153  // Member Functions
154 
155  //- Return plane normal
156  const vector& normal() const;
157 
158  //- Return or return plane base point
159  const point& refPoint() const;
160 
161  //- Return coefficients for the
162  // plane equation: ax + by + cz + d = 0
164 
165  //- Return a point on the plane
166  point aPoint() const;
167 
168  //- Return nearest point in the plane for the given point
169  point nearestPoint(const point& p) const;
170 
171  //- Return distance from the given point to the plane
172  scalar distance(const point& p) const;
173 
174  //- Return cut coefficient for plane and line defined by
175  // origin and direction
176  scalar normalIntersect(const point& pnt0, const vector& dir) const;
177 
178  //- Return cut coefficient for plane and ray
179  scalar normalIntersect(const ray& r) const
180  {
181  return normalIntersect(r.refPoint(), r.dir());
182  }
183 
184  //- Return the cutting point between the plane and
185  // a line passing through the supplied points
186  template<class Point, class PointRef>
187  scalar lineIntersect(const line<Point, PointRef>& l) const
188  {
189  return normalIntersect(l.start(), l.vec());
190  }
191 
192  //- Return the cutting line between this plane and another.
193  // Returned as direction vector and point line goes through.
194  ray planeIntersect(const plane&) const;
195 
196  //- Return the cutting point between this plane and two other planes
197  point planePlaneIntersect(const plane&, const plane&) const;
198 
199  //- Return the side of the plane that the point is on.
200  // If the point is on the plane, then returns NORMAL.
201  side sideOfPlane(const point& p) const;
202 
203  //- Mirror the supplied point in the plane. Return the mirrored point.
204  point mirror(const point& p) const;
205 
206  //- Write to dictionary
207  void writeDict(Ostream&) const;
208 
209 
210  // friend Operators
211 
212  friend bool operator==(const plane&, const plane&);
213  friend bool operator!=(const plane&, const plane&);
214 
215 
216  // IOstream Operators
217 
218  //- Write plane properties
219  friend Ostream& operator<<(Ostream&, const plane&);
220 };
221 
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 } // End namespace Foam
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:48
A line primitive.
Definition: line.H:56
A direction and a reference point.
Definition: plane.H:72
friend Ostream & operator<<(Ostream &, const plane &)
Write plane properties.
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const vector & dir() const
Definition: plane.H:90
friend bool operator!=(const plane &, const plane &)
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: plane.C:361
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
ray planeIntersect(const plane &) const
Return the cutting line between this plane and another.
Definition: plane.C:386
point planePlaneIntersect(const plane &, const plane &) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:454
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
plane(const vector &normalVector)
Construct from normal vector through the origin.
Definition: plane.C:97
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: plane.C:476
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by.
Definition: plane.C:374
point aPoint() const
Return a point on the plane.
Definition: plane.C:294
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
ray(const point &refPoint, const vector &dir)
Definition: plane.H:79
scalar distance(const point &p) const
Return distance from the given point to the plane.
Definition: plane.C:367
PointRef start() const
Return first vertex.
Definition: lineI.H:60
const point & refPoint() const
Definition: plane.H:85
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:484
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:186
Point vec() const
Return start-end vector.
Definition: lineI.H:87
friend bool operator==(const plane &, const plane &)
side
Side of the plane.
Definition: plane.H:65
const vector & normal() const
Return plane normal.
Definition: plane.C:235
Ostream & operator<<(Ostream &, const ensightPart &)
void writeDict(Ostream &) const
Write to dictionary.
Definition: plane.C:499
FixedList< scalar, 4 > planeCoeffs() const
Return coefficients for the.
Definition: plane.C:247
volScalarField & p
bool operator!=(const particle &, const particle &)
Definition: particle.C:1205
Namespace for OpenFOAM.