plane Class Reference

Geometric class that creates a 3D plane and can return the intersection point between a line and the plane. More...

Inheritance diagram for plane:

Classes

class  ray
 A reference point and direction. More...
 

Public Types

enum  side { FRONT = 1, BACK = -1, NORMAL = FRONT, FLIP = BACK }
 Side of the plane. More...
 

Public Member Functions

 plane ()
 Construct zero-initialised. More...
 
 plane (const vector &normalVector)
 Construct from normal vector through the origin. More...
 
 plane (const point &originPoint, const vector &normalVector, const bool doNormalise=true)
 Construct from normal vector and point in plane. More...
 
 plane (const point &point1, const point &point2, const point &point3)
 Construct from three points in plane. More...
 
 plane (const UList< scalar > &coeffs)
 Construct from coefficients for the plane equation: ax + by + cz + d = 0. More...
 
 plane (const FixedList< scalar, 4 > &coeffs)
 Construct from coefficients for the plane equation: ax + by + cz + d = 0. More...
 
 plane (const dictionary &dict)
 Construct from dictionary. More...
 
 plane (Istream &is)
 Construct from Istream. Assumes (normal) (origin) input. More...
 
const vectornormal () const noexcept
 The plane unit normal. More...
 
const pointorigin () const noexcept
 The plane base point. More...
 
pointorigin () noexcept
 The plane base point, for modification. More...
 
void flip ()
 Flip the plane by reversing the normal. More...
 
FixedList< scalar, 4 > planeCoeffs () const
 Return coefficients for the plane equation: ax + by + cz + d = 0. More...
 
point nearestPoint (const point &p) const
 Return nearest point in the plane for the given point. More...
 
scalar distance (const point &p) const
 Return distance (magnitude) from the given point to the plane. More...
 
scalar signedDistance (const point &p) const
 Return distance from the given point to the plane. More...
 
scalar normalIntersect (const point &pnt0, const vector &dir) const
 Return cut coefficient for plane and line defined by origin and direction. More...
 
scalar normalIntersect (const ray &r) const
 Return cut coefficient for plane and ray. More...
 
template<class PointType , class PointRef >
scalar lineIntersect (const line< PointType, PointRef > &l) const
 Return the cutting point between the plane and a line passing through the supplied points. More...
 
ray planeIntersect (const plane &plane2) const
 Return the cutting line between this plane and another. More...
 
point planePlaneIntersect (const plane &plane2, const plane &plane3) const
 Return the cutting point between this plane and two other planes. More...
 
point somePointInPlane (const scalar dist=1e-3) const
 Return a point somewhere on the plane, a distance from the base. More...
 
side whichSide (const point &p) const
 Return the side of the plane that the point is on. More...
 
int sign (const point &p, const scalar tol=SMALL) const
 The sign for the side of the plane that the point is on. More...
 
point mirror (const point &p) const
 Mirror the supplied point in the plane. Return the mirrored point. More...
 
void writeDict (Ostream &os) const
 Write to dictionary. More...
 
const pointrefPoint () const noexcept
 The plane base point (same as origin) More...
 
side sideOfPlane (const point &p) const
 Same as whichSide() More...
 

Detailed Description

Geometric class that creates a 3D plane and can return the intersection point between a line and the plane.

Construction from a dictionary is driven by the planeType. If planeType is missing, pointAndNormal is used and the pointAndNormalDict becomes optional.

For planeType as pointAndNormal :

pointAndNormalDict
{
    point   <point>;   // basePoint (1612 and earlier)
    normal  <vector>;  // normalVector (1612 and earlier)
}

For planeType as embeddedPoints :

embeddedPointsDict
{
    point1  <point>;
    point2  <point>;
    point3  <point>;
}

For planeType with planeEquation coefficients $ ax + by + cz + d = 0 $ :

planeEquationDict
{
   a   <scalar>;
   b   <scalar>;
   c   <scalar>;
   d   <scalar>;
}
Source files

Definition at line 90 of file plane.H.

Member Enumeration Documentation

◆ side

enum side

Side of the plane.

Enumerator
FRONT 

The front (positive normal) side of the plane.

BACK 

The back (negative normal) side of the plane.

NORMAL 

Alias for FRONT.

FLIP 

Alias for BACK.

Definition at line 99 of file plane.H.

Constructor & Destructor Documentation

◆ plane() [1/8]

plane ( )
inline

Construct zero-initialised.

Definition at line 23 of file planeI.H.

◆ plane() [2/8]

plane ( const vector normalVector)
explicit

Construct from normal vector through the origin.

The vector is normalised to a unit vector on input.

Definition at line 115 of file plane.C.

References FUNCTION_NAME.

◆ plane() [3/8]

plane ( const point originPoint,
const vector normalVector,
const bool  doNormalise = true 
)

Construct from normal vector and point in plane.

By default, the vector is normalised to a unit vector on input.

Definition at line 125 of file plane.C.

References FUNCTION_NAME.

◆ plane() [4/8]

plane ( const point point1,
const point point2,
const point point3 
)

Construct from three points in plane.

Definition at line 164 of file plane.C.

References Foam::constant::physicoChemical::b, Foam::constant::universal::c, and FUNCTION_NAME.

◆ plane() [5/8]

plane ( const UList< scalar > &  coeffs)
explicit

Construct from coefficients for the plane equation: ax + by + cz + d = 0.

Definition at line 138 of file plane.C.

References FUNCTION_NAME.

◆ plane() [6/8]

plane ( const FixedList< scalar, 4 > &  coeffs)
explicit

Construct from coefficients for the plane equation: ax + by + cz + d = 0.

Definition at line 151 of file plane.C.

References FUNCTION_NAME.

◆ plane() [7/8]

plane ( const dictionary dict)
explicit

Construct from dictionary.

Definition at line 170 of file plane.C.

References Foam::abort(), dict, Foam::FatalIOError, FatalIOErrorInFunction, dictionary::get(), dictionary::getCompat(), and Foam::nl.

Here is the call graph for this function:

◆ plane() [8/8]

plane ( Istream is)
explicit

Construct from Istream. Assumes (normal) (origin) input.

Definition at line 231 of file plane.C.

References FUNCTION_NAME.

Member Function Documentation

◆ normal()

const Foam::vector & normal ( ) const
inlinenoexcept

The plane unit normal.

Definition at line 32 of file planeI.H.

Referenced by geomCellLooper::cut(), Foam::operator<<(), Foam::operator==(), plane::planeIntersect(), and sampledPlane::print().

Here is the caller graph for this function:

◆ origin() [1/2]

const Foam::point & origin ( ) const
inlinenoexcept

The plane base point.

Definition at line 38 of file planeI.H.

Referenced by searchablePlane::coordinates(), searchableDisk::coordinates(), Foam::operator<(), Foam::operator<<(), Foam::operator==(), plane::planeIntersect(), and sampledPlane::print().

Here is the caller graph for this function:

◆ origin() [2/2]

Foam::point & origin ( )
inlinenoexcept

The plane base point, for modification.

Definition at line 44 of file planeI.H.

◆ flip()

void flip ( )
inline

Flip the plane by reversing the normal.

Definition at line 50 of file planeI.H.

◆ planeCoeffs()

Foam::FixedList< Foam::scalar, 4 > planeCoeffs ( ) const

Return coefficients for the plane equation: ax + by + cz + d = 0.

Definition at line 242 of file plane.C.

References Foam::mag().

Referenced by plane::planePlaneIntersect().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nearestPoint()

Foam::point nearestPoint ( const point p) const
inline

Return nearest point in the plane for the given point.

Definition at line 56 of file planeI.H.

References p.

◆ distance()

Foam::scalar distance ( const point p) const
inline

Return distance (magnitude) from the given point to the plane.

Definition at line 62 of file planeI.H.

References Foam::mag(), and p.

Referenced by geomCellLooper::cut().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ signedDistance()

Foam::scalar signedDistance ( const point p) const
inline

Return distance from the given point to the plane.

Definition at line 68 of file planeI.H.

References p.

◆ normalIntersect() [1/2]

Foam::scalar normalIntersect ( const point pnt0,
const vector dir 
) const

Return cut coefficient for plane and line defined by origin and direction.

Definition at line 293 of file plane.C.

References Foam::stabilise().

Referenced by plane::lineIntersect(), and plane::normalIntersect().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ normalIntersect() [2/2]

scalar normalIntersect ( const ray r) const
inline

Return cut coefficient for plane and ray.

Definition at line 307 of file plane.H.

References plane::normalIntersect().

Here is the call graph for this function:

◆ lineIntersect()

scalar lineIntersect ( const line< PointType, PointRef > &  l) const
inline

Return the cutting point between the plane and a line passing through the supplied points.

Definition at line 317 of file plane.H.

References plane::normalIntersect(), line< Point, PointRef >::start(), and line< Point, PointRef >::vec().

Referenced by slidingInterface::modifyMotionPoints(), cuttingPlane::performCut(), and surfaceFeatures::subsetPlane().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ planeIntersect()

Foam::plane::ray planeIntersect ( const plane plane2) const

Return the cutting line between this plane and another.

Returned as direction vector and point line goes through.

Definition at line 304 of file plane.C.

References Foam::mag(), plane::normal(), and plane::origin().

Referenced by searchableSurfacesQueries::findNearest().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ planePlaneIntersect()

Foam::point planePlaneIntersect ( const plane plane2,
const plane plane3 
) const

Return the cutting point between this plane and two other planes.

Definition at line 373 of file plane.C.

References Foam::constant::physicoChemical::b, Foam::inv(), and plane::planeCoeffs().

Here is the call graph for this function:

◆ somePointInPlane()

Foam::point somePointInPlane ( const scalar  dist = 1e-3) const

Return a point somewhere on the plane, a distance from the base.

Definition at line 395 of file plane.C.

References Foam::mag(), and p.

Here is the call graph for this function:

◆ whichSide()

Foam::plane::side whichSide ( const point p) const
inline

Return the side of the plane that the point is on.

If the point is on the plane, then returns FRONT.

Returns
  • +1 (FRONT): above or on plane (front-side)
  • -1 (BACK): below plane (back-side)

Definition at line 74 of file planeI.H.

References p.

Referenced by plane::sideOfPlane().

Here is the caller graph for this function:

◆ sign()

int sign ( const point p,
const scalar  tol = SMALL 
) const
inline

The sign for the side of the plane that the point is on.

Uses the supplied tolerance for rounding around zero.

Returns
  • 0: on plane
  • +1 (FRONT): above plane (front-side)
  • -1 (BACK): below plane (back-side)

Definition at line 82 of file planeI.H.

References p.

◆ mirror()

Foam::point mirror ( const point p) const

Mirror the supplied point in the plane. Return the mirrored point.

Definition at line 423 of file plane.C.

References Foam::distance(), and p.

Here is the call graph for this function:

◆ writeDict()

void writeDict ( Ostream os) const

Write to dictionary.

Definition at line 438 of file plane.C.

References Ostream::beginBlock(), Ostream::endBlock(), os(), and Ostream::writeEntry().

Here is the call graph for this function:

◆ refPoint()

const point& refPoint ( ) const
inlinenoexcept

The plane base point (same as origin)

Definition at line 381 of file plane.H.

◆ sideOfPlane()

side sideOfPlane ( const point p) const
inline

Same as whichSide()

Definition at line 386 of file plane.H.

References p, and plane::whichSide().

Here is the call graph for this function:

The documentation for this class was generated from the following files: