FreeFOAM The Cross-Platform CFD Toolkit
plane Class Reference

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

#include <OpenFOAM/plane.H>


Detailed Description

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

Source files

Definition at line 61 of file plane.H.

+ Inheritance diagram for plane:

List of all members.

Classes

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

Public Member Functions

 plane (const vector &normalVector)
 Construct from normal vector through the origin.
 plane (const point &basePoint, const vector &normalVector)
 Construct from normal vector and point in plane.
 plane (const point &point1, const point &point2, const point &point3)
 Construct from three points in plane.
 plane (const scalarList &C)
 Construct from coefficients for the.
 plane (const dictionary &planeDict)
 Construct from dictionary.
 plane (Istream &is)
 Construct from Istream. Assumes the base + normal notation.
const vectornormal () const
 Return plane normal.
const pointrefPoint () const
 Return or return plane base point.
FixedList< scalar, 4 > planeCoeffs () const
 Return coefficients for the.
point nearestPoint (const point &p) const
 Return nearest point in the plane for the given point.
scalar distance (const point &p) const
 Return distance from the given point to the plane.
scalar normalIntersect (const point &pnt0, const vector &dir) const
 Return cut coefficient for plane and line defined by.
scalar normalIntersect (const ray &r) const
 Return cut coefficient for plane and ray.
template<class Point , class PointRef >
scalar lineIntersect (const line< Point, PointRef > &l) const
 Return the cutting point between the plane and.
ray planeIntersect (const plane &) const
 Return the cutting line between this plane and another.
point planePlaneIntersect (const plane &, const plane &) const
 Return the cutting point between this plane and two other planes.
void writeDict (Ostream &) const
 Write to dictionary.

Friends

bool operator== (const plane &, const plane &)
bool operator!= (const plane &, const plane &)
Ostreamoperator<< (Ostream &, const plane &)
 Write plane properties.

Constructor & Destructor Documentation

plane ( const vector normalVector)

Construct from normal vector through the origin.

Definition at line 126 of file plane.C.

plane ( const point basePoint,
const vector normalVector 
)

Construct from normal vector and point in plane.

Definition at line 134 of file plane.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and Foam::mag().

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

Construct from three points in plane.

Definition at line 163 of file plane.C.

plane ( const scalarList C)

Construct from coefficients for the.

plane equation: ax + by + cz + d = 0

Definition at line 155 of file plane.C.

plane ( const dictionary planeDict)

Construct from dictionary.

Definition at line 174 of file plane.C.

References Foam::abort(), Foam::FatalIOError, FatalIOErrorIn, dictionary::lookup(), Foam::mag(), Foam::readScalar(), and dictionary::subDict().

plane ( Istream is)

Construct from Istream. Assumes the base + normal notation.

Definition at line 226 of file plane.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and Foam::mag().


Member Function Documentation

const Foam::vector & normal ( ) const

Return plane normal.

Definition at line 249 of file plane.C.

Referenced by geomCellLooper::cut(), and plane::planeIntersect().

const Foam::point & refPoint ( ) const

Return or return plane base point.

Definition at line 256 of file plane.C.

Referenced by searchablePlane::coordinates(), and plane::planeIntersect().

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

Return coefficients for the.

plane equation: ax + by + cz + d = 0

Definition at line 263 of file plane.C.

References Foam::mag().

Referenced by plane::planePlaneIntersect().

Foam::point nearestPoint ( const point p) const

Return nearest point in the plane for the given point.

Definition at line 311 of file plane.C.

Foam::scalar distance ( const point p) const

Return distance from the given point to the plane.

Definition at line 318 of file plane.C.

References Foam::mag().

Referenced by geomCellLooper::cut().

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 326 of file plane.C.

References Foam::stabilise().

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

scalar normalIntersect ( const ray r) const
inline

Return cut coefficient for plane and ray.

Definition at line 165 of file plane.H.

References plane::ray::dir(), plane::normalIntersect(), and plane::ray::refPoint().

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

Return the cutting point between the plane and.

a line passing through the supplied points

Definition at line 173 of file plane.H.

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

Referenced by slidingInterface::modifyMotionPoints().

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 338 of file plane.C.

References Foam::mag(), plane::normal(), plane::refPoint(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

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

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

Definition at line 407 of file plane.C.

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

void writeDict ( Ostream os) const

Friends And Related Function Documentation

bool operator== ( const plane ,
const plane  
)
friend
bool operator!= ( const plane ,
const plane  
)
friend
Ostream& operator<< ( Ostream ,
const plane  
)
friend

Write plane properties.


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