planeI.H
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) 2018-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 inline Foam::plane::plane()
31 :
32  normal_(Zero),
33  origin_(Zero)
34 {}
35 
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 inline const Foam::vector& Foam::plane::normal() const noexcept
40 {
41  return normal_;
42 }
43 
44 
45 inline const Foam::point& Foam::plane::origin() const noexcept
46 {
47  return origin_;
48 }
49 
50 
52 {
53  return origin_;
54 }
55 
56 
57 inline void Foam::plane::flip()
58 {
59  normal_ = -normal_;
60 }
61 
62 
63 inline Foam::point Foam::plane::nearestPoint(const point& p) const
64 {
65  return p - normal_*((p - origin_) & normal_);
66 }
67 
68 
69 inline Foam::scalar Foam::plane::distance(const point& p) const
70 {
71  return mag(signedDistance(p));
72 }
73 
74 
75 inline Foam::scalar Foam::plane::signedDistance(const point& p) const
76 {
77  return ((p - origin_) & normal_);
78 }
79 
80 
82 {
83  const scalar dist = signedDistance(p);
84 
85  return (dist < 0 ? side::BACK : side::FRONT);
86 }
87 
88 
89 inline int Foam::plane::sign(const point& p, const scalar tol) const
90 {
91  const scalar dist = signedDistance(p);
92 
93  return ((dist < -tol) ? -1 : (dist > tol) ? +1 : 0);
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
98 
99 inline bool Foam::operator==(const plane& a, const plane& b)
100 {
101  return (a.origin() == b.origin() && a.normal() == b.normal());
102 }
103 
105 inline bool Foam::operator!=(const plane& a, const plane& b)
106 {
107  return !(a == b);
108 }
109 
110 
111 inline bool Foam::operator<(const plane& a, const plane& b)
112 {
113  return (a.origin() < b.origin());
114 }
115 
116 
117 // ************************************************************************* //
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar signedDistance(const point &p) const
Return distance from the given point to the plane.
Definition: planeI.H:68
void flip()
Flip the plane by reversing the normal.
Definition: planeI.H:50
point nearestPoint(const point &p) const
Return nearest point in the plane for the given point.
Definition: planeI.H:56
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
const point & origin() const noexcept
The plane base point.
Definition: planeI.H:38
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
scalar distance(const point &p) const
Return distance (magnitude) from the given point to the plane.
Definition: planeI.H:62
side whichSide(const point &p) const
Return the side of the plane that the point is on.
Definition: planeI.H:74
plane()
Construct zero-initialised.
Definition: planeI.H:23
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
side
Side of the plane.
Definition: plane.H:99
int sign(const point &p, const scalar tol=SMALL) const
The sign for the side of the plane that the point is on.
Definition: planeI.H:82
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const vector & normal() const noexcept
The plane unit normal.
Definition: planeI.H:32
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127