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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "dictionary.H"
30 #include "plane.H"
31 #include "tensor.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::plane::makeUnitNormal
36 (
37  const char * const caller,
38  const bool notTest
39 )
40 {
41  const scalar magNormal(Foam::mag(normal_));
42 
43  if (magNormal < VSMALL)
44  {
46  << "Plane normal has zero length.\nCalled from " << caller
47  << abort(FatalError);
48  }
49 
50  if (notTest)
51  {
52  normal_ /= magNormal;
53  }
54 }
55 
56 
57 void Foam::plane::calcFromCoeffs
58 (
59  const scalar a,
60  const scalar b,
61  const scalar c,
62  const scalar d,
63  const char * const caller
64 )
65 {
66  if (mag(a) > VSMALL)
67  {
68  origin_ = vector((-d/a), 0, 0);
69  }
70  else if (mag(b) > VSMALL)
71  {
72  origin_ = vector(0, (-d/b), 0);
73  }
74  else if (mag(c) > VSMALL)
75  {
76  origin_ = vector(0, 0, (-d/c));
77  }
78  else
79  {
81  << "At least one plane coefficient must have a value"
82  << abort(FatalError);
83  }
84 
85  normal_ = vector(a, b, c);
86  makeUnitNormal(caller);
87 }
88 
89 
90 void Foam::plane::calcFromEmbeddedPoints
91 (
92  const point& point1,
93  const point& point2,
94  const point& point3,
95  const char * const caller
96 )
97 {
98  origin_ = (point1 + point2 + point3)/3;
99  const vector line12 = point1 - point2;
100  const vector line23 = point2 - point3;
101 
102  if
103  (
104  mag(line12) < VSMALL
105  || mag(line23) < VSMALL
106  || mag(point3-point1) < VSMALL
107  )
108  {
110  << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
111  << abort(FatalError);
112  }
113 
114  normal_ = line12 ^ line23;
116  makeUnitNormal(caller);
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 Foam::plane::plane(const vector& normalVector)
123 :
124  normal_(normalVector),
125  origin_(Zero)
126 {
127  makeUnitNormal(FUNCTION_NAME);
128 }
129 
130 
132 (
133  const point& originPoint,
134  const vector& normalVector,
135  const bool doNormalise
136 )
137 :
138  normal_(normalVector),
139  origin_(originPoint)
140 {
141  makeUnitNormal(FUNCTION_NAME, doNormalise);
142 }
143 
144 
145 Foam::plane::plane(const UList<scalar>& coeffs)
146 {
147  calcFromCoeffs
148  (
149  coeffs[0],
150  coeffs[1],
151  coeffs[2],
152  coeffs[3],
154  );
155 }
156 
157 
159 {
160  calcFromCoeffs
161  (
162  coeffs[0],
163  coeffs[1],
164  coeffs[2],
165  coeffs[3],
167  );
168 }
169 
171 Foam::plane::plane(const point& a, const point& b, const point& c)
172 {
173  calcFromEmbeddedPoints(a, b, c, FUNCTION_NAME);
174 }
175 
176 
178 :
179  normal_(Zero),
180  origin_(Zero)
181 {
182  word planeType;
183  dict.readIfPresent("planeType", planeType);
184 
185  if (planeType.empty())
186  {
187  const dictionary& coeffs = dict.optionalSubDict("pointAndNormalDict");
188 
189  origin_ = coeffs.get<point>("point");
190  normal_ = coeffs.get<point>("normal");
191 
192  makeUnitNormal("point/normal");
193  }
194  else if (planeType == "pointAndNormal")
195  {
196  const dictionary& coeffs = dict.subDict("pointAndNormalDict");
197 
198  origin_ = coeffs.getCompat<point>("point", {{"basePoint", 1612}});
199  normal_ = coeffs.getCompat<point>("normal", {{"normalVector", 1612}});
200 
201  makeUnitNormal("point/normal");
202  }
203  else if (planeType == "planeEquation")
204  {
205  const dictionary& subDict = dict.subDict("planeEquationDict");
206 
207  calcFromCoeffs
208  (
209  subDict.get<scalar>("a"),
210  subDict.get<scalar>("b"),
211  subDict.get<scalar>("c"),
212  subDict.get<scalar>("d"),
213  "planeEquation" // caller name for makeUnitNormal
214  );
215  }
216  else if (planeType == "embeddedPoints")
217  {
218  const dictionary& subDict = dict.subDict("embeddedPointsDict");
219 
220  calcFromEmbeddedPoints
221  (
222  subDict.get<point>("point1"),
223  subDict.get<point>("point2"),
224  subDict.get<point>("point3"),
225  "embeddedPoints" // caller name for makeUnitNormal
226  );
227  }
228  else
229  {
231  << "Invalid plane type: " << planeType << nl
232  << "Valid options: (planeEquation embeddedPoints pointAndNormal)"
233  << abort(FatalIOError);
234  }
235 }
236 
237 
238 Foam::plane::plane(Istream& is)
239 :
240  normal_(is),
241  origin_(is)
242 {
243  makeUnitNormal(FUNCTION_NAME);
244 }
245 
246 
247 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248 
250 {
251  FixedList<scalar, 4> coeffs;
252 
253  const scalar magX = mag(normal_.x());
254  const scalar magY = mag(normal_.y());
255  const scalar magZ = mag(normal_.z());
256 
257  if (magX > magY)
258  {
259  if (magX > magZ)
260  {
261  coeffs[0] = 1;
262  coeffs[1] = normal_.y()/normal_.x();
263  coeffs[2] = normal_.z()/normal_.x();
264  }
265  else
266  {
267  coeffs[0] = normal_.x()/normal_.z();
268  coeffs[1] = normal_.y()/normal_.z();
269  coeffs[2] = 1;
270  }
271  }
272  else
273  {
274  if (magY > magZ)
275  {
276  coeffs[0] = normal_.x()/normal_.y();
277  coeffs[1] = 1;
278  coeffs[2] = normal_.z()/normal_.y();
279  }
280  else
281  {
282  coeffs[0] = normal_.x()/normal_.z();
283  coeffs[1] = normal_.y()/normal_.z();
284  coeffs[2] = 1;
285  }
286  }
287 
288  coeffs[3] =
289  (
290  - coeffs[0] * origin_.x()
291  - coeffs[1] * origin_.y()
292  - coeffs[2] * origin_.z()
293  );
294 
295  return coeffs;
296 }
297 
298 
299 Foam::scalar Foam::plane::normalIntersect
300 (
301  const point& pnt0,
302  const vector& dir
303 ) const
304 {
305  const scalar denom = stabilise((dir & normal_), VSMALL);
306 
307  return ((origin_ - pnt0) & normal_)/denom;
308 }
309 
310 
312 {
313  // Mathworld plane-plane intersection. Assume there is a point on the
314  // intersection line with z=0 and solve the two plane equations
315  // for that (now 2x2 equation in x and y)
316  // Better: use either z=0 or x=0 or y=0.
317 
318  const vector& n1 = this->normal();
319  const vector& n2 = plane2.normal();
320 
321  const point& p1 = this->origin();
322  const point& p2 = plane2.origin();
323 
324  const scalar n1p1 = n1 & p1;
325  const scalar n2p2 = n2 & p2;
326 
327  const vector dir = n1 ^ n2;
328 
329  // Determine zeroed out direction (can be x,y or z) by looking at which
330  // has the largest component in dir.
331  const scalar magX = mag(dir.x());
332  const scalar magY = mag(dir.y());
333  const scalar magZ = mag(dir.z());
334 
335  direction iZero, i1, i2;
336 
337  if (magX > magY)
338  {
339  if (magX > magZ)
340  {
341  iZero = 0;
342  i1 = 1;
343  i2 = 2;
344  }
345  else
346  {
347  iZero = 2;
348  i1 = 0;
349  i2 = 1;
350  }
351  }
352  else
353  {
354  if (magY > magZ)
355  {
356  iZero = 1;
357  i1 = 2;
358  i2 = 0;
359  }
360  else
361  {
362  iZero = 2;
363  i1 = 0;
364  i2 = 1;
365  }
366  }
367 
368 
369  vector pt;
370 
371  pt[iZero] = 0;
372  pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
373  pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);
374 
375  return ray(pt, dir);
376 }
377 
378 
380 (
381  const plane& plane2,
382  const plane& plane3
383 ) const
384 {
385  FixedList<scalar, 4> coeffs1(planeCoeffs());
386  FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
387  FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
388 
389  tensor a
390  (
391  coeffs1[0],coeffs1[1],coeffs1[2],
392  coeffs2[0],coeffs2[1],coeffs2[2],
393  coeffs3[0],coeffs3[1],coeffs3[2]
394  );
396  vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
397 
398  return (inv(a) & (-b));
399 }
400 
401 
402 Foam::point Foam::plane::somePointInPlane(const scalar dist) const
403 {
404  // ax + by + cz + d = 0
405  const FixedList<scalar, 4> coeff(planeCoeffs());
406 
407  // Perturb the base-point
408  point p = origin_ + point::uniform(dist);
409 
410  if (Foam::mag(coeff[2]) < SMALL)
411  {
412  if (Foam::mag(coeff[1]) < SMALL)
413  {
414  p[0] = -1.0*(coeff[1]*p[1] + coeff[2]*p[2] + coeff[3])/coeff[0];
415  }
416  else
417  {
418  p[1] = -1.0*(coeff[0]*p[0] + coeff[2]*p[2] + coeff[3])/coeff[1];
419  }
420  }
421  else
422  {
423  p[2] = -1.0*(coeff[0]*p[0] + coeff[1]*p[1] + coeff[3])/coeff[2];
424  }
425 
426  return p;
427 }
428 
429 
431 {
432  const vector mirroredPtDir = p - nearestPoint(p);
433 
434  if ((normal() & mirroredPtDir) > 0)
435  {
436  return p - 2.0*distance(p)*normal();
437  }
438  else
439  {
440  return p + 2.0*distance(p)*normal();
441  }
442 }
443 
444 
445 void Foam::plane::writeDict(Ostream& os) const
446 {
447  os.writeEntry("planeType", "pointAndNormal");
448 
449  os.beginBlock("pointAndNormalDict");
450 
451  os.writeEntry("point", origin_);
452  os.writeEntry("normal", normal_);
454  os.endBlock();
455 }
456 
457 
458 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
459 
460 Foam::Ostream& Foam::operator<<(Ostream& os, const plane& pln)
461 {
462  os << pln.normal() << token::SPACE << pln.origin();
463  return os;
464 }
465 
466 
467 // ************************************************************************* //
dictionary dict
uint8_t direction
Definition: direction.H:46
A reference point and direction.
Definition: plane.H:113
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const point & origin() const noexcept
The plane base point.
Definition: planeI.H:38
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by origin and direction.
Definition: plane.C:293
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:108
Space [isspace].
Definition: token.H:131
void writeDict(Ostream &os) const
Write to dictionary.
Definition: plane.C:438
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
plane()
Construct zero-initialised.
Definition: planeI.H:23
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:423
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
vector point
Point is a vector.
Definition: point.H:37
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
const dimensionedScalar c
Speed of light in a vacuum.
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:90
ray planeIntersect(const plane &plane2) const
Return the cutting line between this plane and another.
Definition: plane.C:304
point somePointInPlane(const scalar dist=1e-3) const
Return a point somewhere on the plane, a distance from the base.
Definition: plane.C:395
FixedList< scalar, 4 > planeCoeffs() const
Return coefficients for the plane equation: ax + by + cz + d = 0.
Definition: plane.C:242
volScalarField & p
point planePlaneIntersect(const plane &plane2, const plane &plane3) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:373
Tensor of scalars, i.e. Tensor<scalar>.
const vector & normal() const noexcept
The plane unit normal.
Definition: planeI.H:32
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127