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  basePoint_ = vector((-C[3]/C[0]), 0, 0);
36  }
37  else
38  {
39  if (mag(C[1]) > VSMALL)
40  {
41  basePoint_ = vector(0, (-C[3]/C[1]), 0);
42  }
43  else
44  {
45  if (mag(C[2]) > VSMALL)
46  {
47  basePoint_ = 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  unitVector_ = vector(C[0], C[1], C[2]);
59  scalar magUnitVector(mag(unitVector_));
60 
61  if (magUnitVector < VSMALL)
62  {
64  << "Plane normal defined with zero length"
65  << abort(FatalError);
66  }
67 
68  unitVector_ /= 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  basePoint_ = (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  unitVector_ = line12 ^ line23;
96  scalar magUnitVector(mag(unitVector_));
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  unitVector_ /= magUnitVector;
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
112 Foam::plane::plane(const vector& normalVector)
113 :
114  unitVector_(normalVector),
115  basePoint_(Zero)
116 {
117  scalar magUnitVector(mag(unitVector_));
118 
119  if (magUnitVector > VSMALL)
120  {
121  unitVector_ /= magUnitVector;
122  }
123  else
124  {
126  << "plane normal has zero length. basePoint:" << basePoint_
127  << abort(FatalError);
128  }
129 }
130 
131 
132 Foam::plane::plane(const point& basePoint, const vector& normalVector)
133 :
134  unitVector_(normalVector),
135  basePoint_(basePoint)
136 {
137  scalar magUnitVector(mag(unitVector_));
138 
139  if (magUnitVector > VSMALL)
140  {
141  unitVector_ /= magUnitVector;
142  }
143  else
144  {
146  << "plane normal has zero length. basePoint:" << basePoint_
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  unitVector_(Zero),
172  basePoint_(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  basePoint_ = subDict.lookup("basePoint");
204  unitVector_ = subDict.lookup("normalVector");
205  unitVector_ /= mag(unitVector_);
206  }
207  else
208  {
210  << "Invalid plane type: " << planeType << nl
211  << "Valid options include: planeEquation, embeddedPoints and "
212  << "pointAndNormal"
213  << abort(FatalIOError);
214  }
215 }
216 
217 
219 :
220  unitVector_(is),
221  basePoint_(is)
222 {
223  scalar magUnitVector(mag(unitVector_));
224 
225  if (magUnitVector > VSMALL)
226  {
227  unitVector_ /= magUnitVector;
228  }
229  else
230  {
232  << "plane normal has zero length. basePoint:" << basePoint_
233  << abort(FatalError);
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239 
241 {
242  return unitVector_;
243 }
244 
245 
247 {
248  return basePoint_;
249 }
250 
251 
253 {
255 
256  scalar magX = mag(unitVector_.x());
257  scalar magY = mag(unitVector_.y());
258  scalar magZ = mag(unitVector_.z());
259 
260  if (magX > magY)
261  {
262  if (magX > magZ)
263  {
264  C[0] = 1;
265  C[1] = unitVector_.y()/unitVector_.x();
266  C[2] = unitVector_.z()/unitVector_.x();
267  }
268  else
269  {
270  C[0] = unitVector_.x()/unitVector_.z();
271  C[1] = unitVector_.y()/unitVector_.z();
272  C[2] = 1;
273  }
274  }
275  else
276  {
277  if (magY > magZ)
278  {
279  C[0] = unitVector_.x()/unitVector_.y();
280  C[1] = 1;
281  C[2] = unitVector_.z()/unitVector_.y();
282  }
283  else
284  {
285  C[0] = unitVector_.x()/unitVector_.z();
286  C[1] = unitVector_.y()/unitVector_.z();
287  C[2] = 1;
288  }
289  }
290 
291  C[3] = - C[0] * basePoint_.x()
292  - C[1] * basePoint_.y()
293  - C[2] * basePoint_.z();
294 
295  return C;
296 }
297 
298 
300 {
301  return p - unitVector_*((p - basePoint_) & unitVector_);
302 }
303 
304 
305 Foam::scalar Foam::plane::distance(const point& p) const
306 {
307  return mag((p - basePoint_) & unitVector_);
308 }
309 
310 
311 Foam::scalar Foam::plane::normalIntersect
312 (
313  const point& pnt0,
314  const vector& dir
315 ) const
316 {
317  scalar denom = stabilise((dir & unitVector_), VSMALL);
318 
319  return ((basePoint_ - pnt0) & unitVector_)/denom;
320 }
321 
322 
324 {
325  // Mathworld plane-plane intersection. Assume there is a point on the
326  // intersection line with z=0 and solve the two plane equations
327  // for that (now 2x2 equation in x and y)
328  // Better: use either z=0 or x=0 or y=0.
329 
330  const vector& n1 = normal();
331  const vector& n2 = plane2.normal();
332 
333  const point& p1 = refPoint();
334  const point& p2 = plane2.refPoint();
335 
336  scalar n1p1 = n1&p1;
337  scalar n2p2 = n2&p2;
338 
339  vector dir = n1 ^ n2;
340 
341  // Determine zeroed out direction (can be x,y or z) by looking at which
342  // has the largest component in dir.
343  scalar magX = mag(dir.x());
344  scalar magY = mag(dir.y());
345  scalar magZ = mag(dir.z());
346 
347  direction iZero, i1, i2;
348 
349  if (magX > magY)
350  {
351  if (magX > magZ)
352  {
353  iZero = 0;
354  i1 = 1;
355  i2 = 2;
356  }
357  else
358  {
359  iZero = 2;
360  i1 = 0;
361  i2 = 1;
362  }
363  }
364  else
365  {
366  if (magY > magZ)
367  {
368  iZero = 1;
369  i1 = 2;
370  i2 = 0;
371  }
372  else
373  {
374  iZero = 2;
375  i1 = 0;
376  i2 = 1;
377  }
378  }
379 
380  vector pt;
381 
382  pt[iZero] = 0;
383  pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
384  pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);
385 
386  return ray(pt, dir);
387 }
388 
389 
391 (
392  const plane& plane2,
393  const plane& plane3
394 ) const
395 {
397  FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
398  FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
399 
400  tensor a
401  (
402  coeffs1[0],coeffs1[1],coeffs1[2],
403  coeffs2[0],coeffs2[1],coeffs2[2],
404  coeffs3[0],coeffs3[1],coeffs3[2]
405  );
406 
407  vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
408 
409  return (inv(a) & (-b));
410 }
411 
412 
414 {
415  const scalar angle((p - basePoint_) & unitVector_);
416 
417  return (angle < 0 ? FLIP : NORMAL);
418 }
419 
420 
422 {
423  const vector mirroredPtDir = p - nearestPoint(p);
424 
425  if ((normal() & mirroredPtDir) > 0)
426  {
427  return p - 2.0*distance(p)*normal();
428  }
429  else
430  {
431  return p + 2.0*distance(p)*normal();
432  }
433 }
434 
435 
437 {
438  os.writeKeyword("planeType") << "pointAndNormal"
439  << token::END_STATEMENT << nl;
440  os << indent << "pointAndNormalDict" << nl
442  os.writeKeyword("basePoint") << basePoint_ << token::END_STATEMENT << nl;
443  os.writeKeyword("normalVector") << unitVector_ << token::END_STATEMENT
444  << nl;
445  os << decrIndent << indent << token::END_BLOCK << endl;
446 }
447 
448 
449 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
450 
451 bool Foam::operator==(const plane& a, const plane& b)
452 {
453  if (a.basePoint_ == b.basePoint_ && a.unitVector_ == b.unitVector_)
454  {
455  return true;
456  }
457  else
458  {
459  return false;
460  }
461 }
462 
463 bool Foam::operator!=(const plane& a, const plane& b)
464 {
465  return !(a == b);
466 }
467 
468 
469 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
470 
472 {
473  os << a.unitVector_ << token::SPACE << a.basePoint_;
474 
475  return os;
476 }
477 
478 
479 // ************************************************************************* //
const Cmpt & z() const
Definition: VectorI.H:87
dictionary dict
uint8_t direction
Definition: direction.H:46
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
A direction and a reference point.
Definition: plane.H:73
const Cmpt & x() const
Definition: VectorI.H:75
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
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
const point & refPoint() const
Return or return plane base point.
Definition: plane.C:246
side sideOfPlane(const point &p) const
Return the side of the plane that the point is on.
Definition: plane.C:413
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: plane.C:299
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
const vector & normal() const
Return plane normal.
Definition: plane.C:240
plane(const vector &normalVector)
Construct from normal vector through the origin.
Definition: plane.C:112
point mirror(const point &p) const
Mirror the supplied point in the plane. Return the mirrored point.
Definition: plane.C:421
ray planeIntersect(const plane &) const
Return the cutting line between this plane and another.
Definition: plane.C:323
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
scalar distance(const point &p) const
Return distance from the given point to the plane.
Definition: plane.C:305
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
const Cmpt & y() const
Definition: VectorI.H:81
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
point planePlaneIntersect(const plane &, const plane &) const
Return the cutting point between this plane and two other planes.
Definition: plane.C:391
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalar normalIntersect(const point &pnt0, const vector &dir) const
Return cut coefficient for plane and line defined by.
Definition: plane.C:312
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
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
side
Side of the plane.
Definition: plane.H:65
vector point
Point is a vector.
Definition: point.H:41
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
FixedList< scalar, 4 > planeCoeffs() const
Return coefficients for the.
Definition: plane.C:252
Ostream & operator<<(Ostream &, const ensightPart &)
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
bool operator!=(const particle &, const particle &)
Definition: particle.C:145
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230
void writeDict(Ostream &) const
Write to dictionary.
Definition: plane.C:436
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
IOerror FatalIOError