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