plane.C
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-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "plane.H"
27 #include "tensor.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Calculate base point and unit normal vector from plane equation
33 {
34  if (mag(C[0]) > VSMALL)
35  {
36  basePoint_ = vector((-C[3]/C[0]), 0, 0);
37  }
38  else
39  {
40  if (mag(C[1]) > VSMALL)
41  {
42  basePoint_ = vector(0, (-C[3]/C[1]), 0);
43  }
44  else
45  {
46  if (mag(C[2]) > VSMALL)
47  {
48  basePoint_ = vector(0, 0, (-C[3]/C[2]));
49  }
50  else
51  {
53  << "At least one plane coefficient must have a value"
54  << abort(FatalError);
55  }
56  }
57  }
58 
59  unitVector_ = vector(C[0], C[1], C[2]);
60  scalar magUnitVector(mag(unitVector_));
61 
62  if (magUnitVector < VSMALL)
63  {
65  << "Plane normal defined with zero length"
66  << abort(FatalError);
67  }
68 
69  unitVector_ /= magUnitVector;
70 }
71 
72 
74 (
75  const point& point1,
76  const point& point2,
77  const point& point3
78 )
79 {
80  basePoint_ = (point1 + point2 + point3)/3;
81  vector line12 = point1 - point2;
82  vector line23 = point2 - point3;
83 
84  if
85  (
86  mag(line12) < VSMALL
87  || mag(line23) < VSMALL
88  || mag(point3-point1) < VSMALL
89  )
90  {
92  << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
93  << abort(FatalError);
94  }
95 
96  unitVector_ = line12 ^ line23;
97  scalar magUnitVector(mag(unitVector_));
98 
99  if (magUnitVector < VSMALL)
100  {
102  << "Plane normal defined with zero length" << nl
103  << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
104  << abort(FatalError);
105  }
106 
107  unitVector_ /= magUnitVector;
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
113 // Construct from normal vector through the origin
114 Foam::plane::plane(const vector& normalVector)
115 :
116  unitVector_(normalVector),
117  basePoint_(vector::zero)
118 {
119  scalar magUnitVector(mag(unitVector_));
120 
121  if (magUnitVector > VSMALL)
122  {
123  unitVector_ /= magUnitVector;
124  }
125  else
126  {
128  << "plane normal has zero length. basePoint:" << basePoint_
129  << abort(FatalError);
130  }
131 }
132 
133 
134 // Construct from point and normal vector
135 Foam::plane::plane(const point& basePoint, const vector& normalVector)
136 :
137  unitVector_(normalVector),
138  basePoint_(basePoint)
139 {
140  scalar magUnitVector(mag(unitVector_));
141 
142  if (magUnitVector > VSMALL)
143  {
144  unitVector_ /= magUnitVector;
145  }
146  else
147  {
149  << "plane normal has zero length. basePoint:" << basePoint_
150  << abort(FatalError);
151  }
152 }
153 
154 
155 // Construct from plane equation
157 {
158  calcPntAndVec(C);
159 }
160 
161 
162 // Construct from three points
164 (
165  const point& a,
166  const point& b,
167  const point& c
168 )
169 {
170  calcPntAndVec(a, b, c);
171 }
172 
173 
174 // Construct from dictionary
176 :
177  unitVector_(vector::zero),
178  basePoint_(point::zero)
179 {
180  const word planeType(dict.lookup("planeType"));
181 
182  if (planeType == "planeEquation")
183  {
184  const dictionary& subDict = dict.subDict("planeEquationDict");
185  scalarList C(4);
186 
187  C[0] = readScalar(subDict.lookup("a"));
188  C[1] = readScalar(subDict.lookup("b"));
189  C[2] = readScalar(subDict.lookup("c"));
190  C[3] = readScalar(subDict.lookup("d"));
191 
192  calcPntAndVec(C);
193 
194  }
195  else if (planeType == "embeddedPoints")
196  {
197  const dictionary& subDict = dict.subDict("embeddedPointsDict");
198 
199  point point1(subDict.lookup("point1"));
200  point point2(subDict.lookup("point2"));
201  point point3(subDict.lookup("point3"));
202 
203  calcPntAndVec(point1, point2, point3);
204  }
205  else if (planeType == "pointAndNormal")
206  {
207  const dictionary& subDict = dict.subDict("pointAndNormalDict");
208 
209  basePoint_ = subDict.lookup("basePoint");
210  unitVector_ = subDict.lookup("normalVector");
212  }
213  else
214  {
216  << "Invalid plane type: " << planeType << nl
217  << "Valid options include: planeEquation, embeddedPoints and "
218  << "pointAndNormal"
219  << abort(FatalIOError);
220  }
221 }
222 
223 
224 // Construct from Istream. Assumes point and normal vector.
226 :
227  unitVector_(is),
228  basePoint_(is)
229 {
230  scalar magUnitVector(mag(unitVector_));
231 
232  if (magUnitVector > VSMALL)
233  {
234  unitVector_ /= magUnitVector;
235  }
236  else
237  {
239  << "plane normal has zero length. basePoint:" << basePoint_
240  << abort(FatalError);
241  }
242 }
243 
244 
245 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
246 
247 // Return plane normal vector
249 {
250  return unitVector_;
251 }
252 
253 
254 // Return plane base point
256 {
257  return basePoint_;
258 }
259 
260 
261 // Return coefficients for plane equation: ax + by + cz + d = 0
263 {
265 
266  scalar magX = mag(unitVector_.x());
267  scalar magY = mag(unitVector_.y());
268  scalar magZ = mag(unitVector_.z());
269 
270  if (magX > magY)
271  {
272  if (magX > magZ)
273  {
274  C[0] = 1;
275  C[1] = unitVector_.y()/unitVector_.x();
276  C[2] = unitVector_.z()/unitVector_.x();
277  }
278  else
279  {
280  C[0] = unitVector_.x()/unitVector_.z();
281  C[1] = unitVector_.y()/unitVector_.z();
282  C[2] = 1;
283  }
284  }
285  else
286  {
287  if (magY > magZ)
288  {
289  C[0] = unitVector_.x()/unitVector_.y();
290  C[1] = 1;
291  C[2] = unitVector_.z()/unitVector_.y();
292  }
293  else
294  {
295  C[0] = unitVector_.x()/unitVector_.z();
296  C[1] = unitVector_.y()/unitVector_.z();
297  C[2] = 1;
298  }
299  }
300 
301  C[3] = - C[0] * basePoint_.x()
302  - C[1] * basePoint_.y()
303  - C[2] * basePoint_.z();
304 
305  return C;
306 }
307 
308 
309 // Return nearest point in the plane for the given point
311 {
312  return p - unitVector_*((p - basePoint_) & unitVector_);
313 }
314 
315 
316 // Return distance from the given point to the plane
317 Foam::scalar Foam::plane::distance(const point& p) const
318 {
319  return mag((p - basePoint_) & unitVector_);
320 }
321 
322 
323 // Cutting point for plane and line defined by origin and direction
324 Foam::scalar Foam::plane::normalIntersect
325 (
326  const point& pnt0,
327  const vector& dir
328 ) const
329 {
330  scalar denom = stabilise((dir & unitVector_), VSMALL);
331 
332  return ((basePoint_ - pnt0) & unitVector_)/denom;
333 }
334 
335 
336 // Cutting line of two planes
338 {
339  // Mathworld plane-plane intersection. Assume there is a point on the
340  // intersection line with z=0 and solve the two plane equations
341  // for that (now 2x2 equation in x and y)
342  // Better: use either z=0 or x=0 or y=0.
343 
344  const vector& n1 = normal();
345  const vector& n2 = plane2.normal();
346 
347  const point& p1 = refPoint();
348  const point& p2 = plane2.refPoint();
349 
350  scalar n1p1 = n1&p1;
351  scalar n2p2 = n2&p2;
352 
353  vector dir = n1 ^ n2;
354 
355  // Determine zeroed out direction (can be x,y or z) by looking at which
356  // has the largest component in dir.
357  scalar magX = mag(dir.x());
358  scalar magY = mag(dir.y());
359  scalar magZ = mag(dir.z());
360 
361  direction iZero, i1, i2;
362 
363  if (magX > magY)
364  {
365  if (magX > magZ)
366  {
367  iZero = 0;
368  i1 = 1;
369  i2 = 2;
370  }
371  else
372  {
373  iZero = 2;
374  i1 = 0;
375  i2 = 1;
376  }
377  }
378  else
379  {
380  if (magY > magZ)
381  {
382  iZero = 1;
383  i1 = 2;
384  i2 = 0;
385  }
386  else
387  {
388  iZero = 2;
389  i1 = 0;
390  i2 = 1;
391  }
392  }
393 
394  vector pt;
395 
396  pt[iZero] = 0;
397  pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
398  pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);
399 
400  return ray(pt, dir);
401 }
402 
403 
404 // Cutting point of three planes
406 (
407  const plane& plane2,
408  const plane& plane3
409 ) const
410 {
411  FixedList<scalar, 4> coeffs1(planeCoeffs());
412  FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
413  FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
414 
415  tensor a
416  (
417  coeffs1[0],coeffs1[1],coeffs1[2],
418  coeffs2[0],coeffs2[1],coeffs2[2],
419  coeffs3[0],coeffs3[1],coeffs3[2]
420  );
421 
422  vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
423 
424  return (inv(a) & (-b));
425 }
426 
427 
429 {
430  const scalar angle((p - basePoint_) & unitVector_);
431 
432  return (angle < 0 ? FLIP : NORMAL);
433 }
434 
435 
437 {
438  const vector mirroredPtDir = p - nearestPoint(p);
439 
440  if ((normal() & mirroredPtDir) > 0)
441  {
442  return p - 2.0*distance(p)*normal();
443  }
444  else
445  {
446  return p + 2.0*distance(p)*normal();
447  }
448 }
449 
450 
452 {
453  os.writeKeyword("planeType") << "pointAndNormal"
454  << token::END_STATEMENT << nl;
455  os << indent << "pointAndNormalDict" << nl
457  os.writeKeyword("basePoint") << basePoint_ << token::END_STATEMENT << nl;
458  os.writeKeyword("normalVector") << unitVector_ << token::END_STATEMENT
459  << nl;
460  os << decrIndent << indent << token::END_BLOCK << endl;
461 }
462 
463 
464 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
465 
466 bool Foam::operator==(const plane& a, const plane& b)
467 {
468  if (a.basePoint_ == b.basePoint_ && a.unitVector_ == b.unitVector_)
469  {
470  return true;
471  }
472  else
473  {
474  return false;
475  }
476 }
477 
478 bool Foam::operator!=(const plane& a, const plane& b)
479 {
480  return !(a == b);
481 }
482 
483 
484 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
485 
487 {
488  os << a.unitVector_ << token::SPACE << a.basePoint_;
489 
490  return os;
491 }
492 
493 
494 // ************************************************************************* //
Foam::plane::distance
scalar distance(const point &p) const
Return distance from the given point to the plane.
Definition: plane.C:317
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::plane::normal
const vector & normal() const
Return plane normal.
Definition: plane.C:248
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::plane::calcPntAndVec
void calcPntAndVec(const scalarList &C)
Calculates basePoint and normal vector given plane coefficients.
Definition: plane.C:32
Foam::plane::nearestPoint
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: plane.C:310
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::plane::unitVector_
vector unitVector_
Plane normal.
Definition: plane.H:104
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::plane::plane
plane(const vector &normalVector)
Construct from normal vector through the origin.
Definition: plane.C:114
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
tensor.H
C
volScalarField & C
Definition: readThermalProperties.H:104
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Foam::plane::side
side
Side of the plane.
Definition: plane.H:65
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
plane.H
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:40
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:71
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::plane::writeDict
void writeDict(Ostream &) const
Write to dictionary.
Definition: plane.C:451
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
Foam::plane::planeCoeffs
FixedList< scalar, 4 > planeCoeffs() const
Return coefficients for the.
Definition: plane.C:262
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::operator!=
bool operator!=(const particle &, const particle &)
Definition: particle.C:147
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::plane::basePoint_
point basePoint_
Base point.
Definition: plane.H:107
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::plane::ray
A direction and a reference point.
Definition: plane.H:73
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::plane::planeIntersect
ray planeIntersect(const plane &) const
Return the cutting line between this plane and another.
Definition: plane.C:337
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::token::END_BLOCK
@ END_BLOCK
Definition: token.H:105
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::plane::refPoint
const point & refPoint() const
Return or return plane base point.
Definition: plane.C:255
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::token::BEGIN_BLOCK
@ BEGIN_BLOCK
Definition: token.H:104
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::plane::sideOfPlane
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: plane.C:428
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::plane::normalIntersect
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by.
Definition: plane.C:325
Foam::plane::planePlaneIntersect
point planePlaneIntersect(const plane &, const plane &) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:406
Foam::plane::mirror
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:436
normal
A normal distribution model.
Foam::token::SPACE
@ SPACE
Definition: token.H:95
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47