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