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