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