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-2018 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 "plane.H"
30 #include "tensor.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 void Foam::plane::makeUnitNormal
35 (
36  const char * const caller,
37  const bool notTest
38 )
39 {
40  const scalar magNormal(Foam::mag(normal_));
41 
42  if (magNormal < VSMALL)
43  {
45  << "Plane normal has zero length.\nCalled from " << caller
46  << abort(FatalError);
47  }
48 
49  if (notTest)
50  {
51  normal_ /= magNormal;
52  }
53 }
54 
55 
56 void Foam::plane::calcFromCoeffs
57 (
58  const scalar a,
59  const scalar b,
60  const scalar c,
61  const scalar d,
62  const char * const caller
63 )
64 {
65  if (mag(a) > VSMALL)
66  {
67  origin_ = vector((-d/a), 0, 0);
68  }
69  else if (mag(b) > VSMALL)
70  {
71  origin_ = vector(0, (-d/b), 0);
72  }
73  else if (mag(c) > VSMALL)
74  {
75  origin_ = vector(0, 0, (-d/c));
76  }
77  else
78  {
80  << "At least one plane coefficient must have a value"
81  << abort(FatalError);
82  }
83 
84  normal_ = vector(a, b, c);
85  makeUnitNormal(caller);
86 }
87 
88 
89 void Foam::plane::calcFromEmbeddedPoints
90 (
91  const point& point1,
92  const point& point2,
93  const point& point3,
94  const char * const caller
95 )
96 {
97  origin_ = (point1 + point2 + point3)/3;
98  const vector line12 = point1 - point2;
99  const vector line23 = point2 - point3;
100 
101  if
102  (
103  mag(line12) < VSMALL
104  || mag(line23) < VSMALL
105  || mag(point3-point1) < VSMALL
106  )
107  {
109  << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
110  << abort(FatalError);
111  }
112 
113  normal_ = line12 ^ line23;
114 
115  makeUnitNormal(caller);
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
121 Foam::plane::plane(const vector& normalVector)
122 :
123  normal_(normalVector),
124  origin_(Zero)
125 {
126  makeUnitNormal(FUNCTION_NAME);
127 }
128 
129 
131 (
132  const point& originPoint,
133  const vector& normalVector,
134  const bool doNormalise
135 )
136 :
137  normal_(normalVector),
138  origin_(originPoint)
139 {
140  makeUnitNormal(FUNCTION_NAME, doNormalise);
141 }
142 
143 
145 {
146  calcFromCoeffs
147  (
148  coeffs[0],
149  coeffs[1],
150  coeffs[2],
151  coeffs[3],
153  );
154 }
155 
156 
158 {
159  calcFromCoeffs
160  (
161  coeffs[0],
162  coeffs[1],
163  coeffs[2],
164  coeffs[3],
166  );
167 }
168 
169 
170 Foam::plane::plane(const point& a, const point& b, const point& c)
171 {
172  calcFromEmbeddedPoints(a, b, c, FUNCTION_NAME);
173 }
174 
175 
177 :
178  normal_(Zero),
179  origin_(Zero)
180 {
181  const word planeType(dict.get<word>("planeType"));
182 
183  if (planeType == "planeEquation")
184  {
185  const dictionary& subDict = dict.subDict("planeEquationDict");
186 
187  calcFromCoeffs
188  (
189  subDict.get<scalar>("a"),
190  subDict.get<scalar>("b"),
191  subDict.get<scalar>("c"),
192  subDict.get<scalar>("d"),
193  "planeEquationDict" // caller name for makeUnitNormal
194  );
195  }
196  else if (planeType == "embeddedPoints")
197  {
198  const dictionary& subDict = dict.subDict("embeddedPointsDict");
199 
200  calcFromEmbeddedPoints
201  (
202  subDict.get<point>("point1"),
203  subDict.get<point>("point2"),
204  subDict.get<point>("point3"),
205  "embeddedPointsDict" // caller name for makeUnitNormal
206  );
207 
208  }
209  else if (planeType == "pointAndNormal")
210  {
211  const dictionary& subDict = dict.subDict("pointAndNormalDict");
212 
213  origin_ = subDict.getCompat<point>("point", {{"basePoint", 1612}});
214  normal_ = subDict.getCompat<point>("normal", {{"normalVector", 1612}});
215 
216  makeUnitNormal("pointAndNormalDict");
217  }
218  else
219  {
221  << "Invalid plane type: " << planeType << nl
222  << "Valid options: (planeEquation embeddedPoints pointAndNormal)"
223  << abort(FatalIOError);
224  }
225 }
226 
227 
229 :
230  normal_(is),
231  origin_(is)
232 {
233  makeUnitNormal(FUNCTION_NAME);
234 }
235 
236 
237 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
238 
240 {
241  FixedList<scalar, 4> coeffs(4);
242 
243  const scalar magX = mag(normal_.x());
244  const scalar magY = mag(normal_.y());
245  const scalar magZ = mag(normal_.z());
246 
247  if (magX > magY)
248  {
249  if (magX > magZ)
250  {
251  coeffs[0] = 1;
252  coeffs[1] = normal_.y()/normal_.x();
253  coeffs[2] = normal_.z()/normal_.x();
254  }
255  else
256  {
257  coeffs[0] = normal_.x()/normal_.z();
258  coeffs[1] = normal_.y()/normal_.z();
259  coeffs[2] = 1;
260  }
261  }
262  else
263  {
264  if (magY > magZ)
265  {
266  coeffs[0] = normal_.x()/normal_.y();
267  coeffs[1] = 1;
268  coeffs[2] = normal_.z()/normal_.y();
269  }
270  else
271  {
272  coeffs[0] = normal_.x()/normal_.z();
273  coeffs[1] = normal_.y()/normal_.z();
274  coeffs[2] = 1;
275  }
276  }
277 
278  coeffs[3] =
279  (
280  - coeffs[0] * origin_.x()
281  - coeffs[1] * origin_.y()
282  - coeffs[2] * origin_.z()
283  );
284 
285  return coeffs;
286 }
287 
288 
289 Foam::scalar Foam::plane::normalIntersect
290 (
291  const point& pnt0,
292  const vector& dir
293 ) const
294 {
295  const scalar denom = stabilise((dir & normal_), VSMALL);
296 
297  return ((origin_ - pnt0) & normal_)/denom;
298 }
299 
300 
302 {
303  // Mathworld plane-plane intersection. Assume there is a point on the
304  // intersection line with z=0 and solve the two plane equations
305  // for that (now 2x2 equation in x and y)
306  // Better: use either z=0 or x=0 or y=0.
307 
308  const vector& n1 = this->normal();
309  const vector& n2 = plane2.normal();
310 
311  const point& p1 = this->origin();
312  const point& p2 = plane2.origin();
313 
314  const scalar n1p1 = n1 & p1;
315  const scalar n2p2 = n2 & p2;
316 
317  const vector dir = n1 ^ n2;
318 
319  // Determine zeroed out direction (can be x,y or z) by looking at which
320  // has the largest component in dir.
321  const scalar magX = mag(dir.x());
322  const scalar magY = mag(dir.y());
323  const scalar magZ = mag(dir.z());
324 
325  direction iZero, i1, i2;
326 
327  if (magX > magY)
328  {
329  if (magX > magZ)
330  {
331  iZero = 0;
332  i1 = 1;
333  i2 = 2;
334  }
335  else
336  {
337  iZero = 2;
338  i1 = 0;
339  i2 = 1;
340  }
341  }
342  else
343  {
344  if (magY > magZ)
345  {
346  iZero = 1;
347  i1 = 2;
348  i2 = 0;
349  }
350  else
351  {
352  iZero = 2;
353  i1 = 0;
354  i2 = 1;
355  }
356  }
357 
358 
359  vector pt;
360 
361  pt[iZero] = 0;
362  pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
363  pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);
364 
365  return ray(pt, dir);
366 }
367 
368 
370 (
371  const plane& plane2,
372  const plane& plane3
373 ) const
374 {
375  FixedList<scalar, 4> coeffs1(planeCoeffs());
376  FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
377  FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
378 
379  tensor a
380  (
381  coeffs1[0],coeffs1[1],coeffs1[2],
382  coeffs2[0],coeffs2[1],coeffs2[2],
383  coeffs3[0],coeffs3[1],coeffs3[2]
384  );
385 
386  vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
387 
388  return (inv(a) & (-b));
389 }
390 
391 
393 {
394  // ax + by + cz + d = 0
395  const FixedList<scalar, 4> coeff(planeCoeffs());
396 
397  // Perturb the base-point
398  point p = origin_ + point::uniform(dist);
399 
400  if (Foam::mag(coeff[2]) < SMALL)
401  {
402  if (Foam::mag(coeff[1]) < SMALL)
403  {
404  p[0] = -1.0*(coeff[1]*p[1] + coeff[2]*p[2] + coeff[3])/coeff[0];
405  }
406  else
407  {
408  p[1] = -1.0*(coeff[0]*p[0] + coeff[2]*p[2] + coeff[3])/coeff[1];
409  }
410  }
411  else
412  {
413  p[2] = -1.0*(coeff[0]*p[0] + coeff[1]*p[1] + coeff[3])/coeff[2];
414  }
415 
416  return p;
417 }
418 
419 
421 {
422  const vector mirroredPtDir = p - nearestPoint(p);
423 
424  if ((normal() & mirroredPtDir) > 0)
425  {
426  return p - 2.0*distance(p)*normal();
427  }
428  else
429  {
430  return p + 2.0*distance(p)*normal();
431  }
432 }
433 
434 
436 {
437  os.writeEntry("planeType", "pointAndNormal");
438 
439  os.beginBlock("pointAndNormalDict");
440 
441  os.writeEntry("point", origin_);
442  os.writeEntry("normal", normal_);
443 
444  os.endBlock();
445 }
446 
447 
448 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
449 
451 {
452  os << pln.normal() << token::SPACE << pln.origin();
453  return os;
454 }
455 
456 
457 // ************************************************************************* //
Foam::plane::planePlaneIntersect
point planePlaneIntersect(const plane &plane2, const plane &plane3) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:370
Foam::Tensor< scalar >
Foam::plane::normal
const vector & normal() const
The plane unit normal.
Definition: planeI.H:39
Foam::plane::writeDict
void writeDict(Ostream &os) const
Write to dictionary.
Definition: plane.C:435
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::FatalIOError
IOerror FatalIOError
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::dictionary::getCompat
T getCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, enum keyType::option=keyType::REGEX) const
Definition: dictionaryTemplates.C:134
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
tensor.H
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
plane.H
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:43
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::plane::planeCoeffs
FixedList< scalar, 4 > planeCoeffs() const
Definition: plane.C:239
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::plane::ray
A reference point and direction.
Definition: plane.H:104
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::plane::origin
const point & origin() const
The plane base point.
Definition: planeI.H:45
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::plane::planeIntersect
ray planeIntersect(const plane &plane2) const
Return the cutting line between this plane and another.
Definition: plane.C:301
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::List< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::token::SPACE
Space [isspace].
Definition: token.H:125
Foam::direction
uint8_t direction
Definition: direction.H:52
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::plane::normalIntersect
scalar normalIntersect(const point &pnt0, const vector &dir) const
Definition: plane.C:290
Foam::plane::plane
plane()
Construct zero-initialised.
Definition: planeI.H:30
Foam::plane::somePointInPlane
point somePointInPlane(const scalar dist=1e-3) const
Return a point somewhere on the plane, a distance from the base.
Definition: plane.C:392
Foam::plane::mirror
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:420