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