pointEdgePoint Class Reference

Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave) To be used in wall distance calculation. More...

Inheritance diagram for pointEdgePoint:

Public Member Functions

 pointEdgePoint ()
 Default construct. Max point. More...
 
 pointEdgePoint (const point &origin, const scalar distSqr)
 Construct from origin, distance. More...
 
const pointorigin () const noexcept
 
pointorigin () noexcept
 
scalar distSqr () const noexcept
 
scalar & distSqr () noexcept
 
bool operator== (const pointEdgePoint &) const
 Test for equality. More...
 
bool operator!= (const pointEdgePoint &) const
 Test for inequality. More...
 
template<class TrackingData >
bool valid (TrackingData &td) const
 Changed or contains original (invalid) value. More...
 
template<class TrackingData >
bool sameGeometry (const pointEdgePoint &, const scalar tol, TrackingData &td) const
 Check for identical geometrical data (eg, cyclics checking) More...
 
template<class TrackingData >
void leaveDomain (const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
 Convert origin to relative vector to leaving point. More...
 
template<class TrackingData >
void enterDomain (const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
 Convert relative origin to absolute by adding entering point. More...
 
template<class TrackingData >
void transform (const tensor &rotTensor, TrackingData &td)
 Apply rotation matrix to origin. More...
 
template<class TrackingData >
bool updatePoint (const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
 Influence of edge on point. More...
 
template<class TrackingData >
bool updatePoint (const polyMesh &mesh, const label pointi, const pointEdgePoint &newPointInfo, const scalar tol, TrackingData &td)
 Influence of different value on same point. More...
 
template<class TrackingData >
bool updatePoint (const pointEdgePoint &newPointInfo, const scalar tol, TrackingData &td)
 Influence of different value on same point. More...
 
template<class TrackingData >
bool updateEdge (const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
 Influence of point on edge. More...
 
template<class TrackingData >
bool equal (const pointEdgePoint &, TrackingData &td) const
 Test for equality, with TrackingData. More...
 

Friends

Ostreamoperator<< (Ostream &, const pointEdgePoint &)
 
Istreamoperator>> (Istream &, pointEdgePoint &)
 

Detailed Description

Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave) To be used in wall distance calculation.

Source files

Definition at line 61 of file pointEdgePoint.H.

Constructor & Destructor Documentation

◆ pointEdgePoint() [1/2]

pointEdgePoint ( )
inline

Default construct. Max point.

Definition at line 116 of file pointEdgePointI.H.

◆ pointEdgePoint() [2/2]

pointEdgePoint ( const point origin,
const scalar  distSqr 
)
inline

Construct from origin, distance.

Definition at line 124 of file pointEdgePointI.H.

Member Function Documentation

◆ origin() [1/2]

const point& origin ( ) const
inlinenoexcept

Definition at line 126 of file pointEdgePoint.H.

◆ origin() [2/2]

point& origin ( )
inlinenoexcept

Definition at line 128 of file pointEdgePoint.H.

◆ distSqr() [1/2]

scalar distSqr ( ) const
inlinenoexcept

Definition at line 130 of file pointEdgePoint.H.

◆ distSqr() [2/2]

scalar& distSqr ( )
inlinenoexcept

Definition at line 132 of file pointEdgePoint.H.

◆ operator==()

bool operator== ( const pointEdgePoint rhs) const
inline

Test for equality.

Definition at line 287 of file pointEdgePointI.H.

◆ operator!=()

bool operator!= ( const pointEdgePoint rhs) const
inline

Test for inequality.

Definition at line 296 of file pointEdgePointI.H.

◆ valid()

bool valid ( TrackingData &  td) const
inline

Changed or contains original (invalid) value.

Definition at line 137 of file pointEdgePointI.H.

References Foam::max().

Here is the call graph for this function:

◆ sameGeometry()

bool sameGeometry ( const pointEdgePoint w2,
const scalar  tol,
TrackingData &  td 
) const
inline

Check for identical geometrical data (eg, cyclics checking)

Definition at line 146 of file pointEdgePointI.H.

References Foam::diff(), Foam::mag(), and w2.

Here is the call graph for this function:

◆ leaveDomain()

void leaveDomain ( const polyPatch patch,
const label  patchPointi,
const point pos,
TrackingData &  td 
)
inline

Convert origin to relative vector to leaving point.

(= point coordinate)

Definition at line 174 of file pointEdgePointI.H.

◆ enterDomain()

void enterDomain ( const polyPatch patch,
const label  patchPointi,
const point pos,
TrackingData &  td 
)
inline

Convert relative origin to absolute by adding entering point.

Definition at line 200 of file pointEdgePointI.H.

◆ transform()

void transform ( const tensor rotTensor,
TrackingData &  td 
)
inline

Apply rotation matrix to origin.

Definition at line 187 of file pointEdgePointI.H.

References Foam::transform().

Here is the call graph for this function:

◆ updatePoint() [1/3]

bool updatePoint ( const polyMesh mesh,
const label  pointi,
const label  edgeI,
const pointEdgePoint edgeInfo,
const scalar  tol,
TrackingData &  td 
)
inline

Influence of edge on point.

Definition at line 215 of file pointEdgePointI.H.

References mesh, polyMesh::points(), and update().

Here is the call graph for this function:

◆ updatePoint() [2/3]

bool updatePoint ( const polyMesh mesh,
const label  pointi,
const pointEdgePoint newPointInfo,
const scalar  tol,
TrackingData &  td 
)
inline

Influence of different value on same point.

Merge new and old info.

Definition at line 231 of file pointEdgePointI.H.

References mesh, polyMesh::points(), and update().

Here is the call graph for this function:

◆ updatePoint() [3/3]

bool updatePoint ( const pointEdgePoint newPointInfo,
const scalar  tol,
TrackingData &  td 
)
inline

Influence of different value on same point.

No information about current position whatsoever.

Definition at line 246 of file pointEdgePointI.H.

References update().

Here is the call graph for this function:

◆ updateEdge()

bool updateEdge ( const polyMesh mesh,
const label  edgeI,
const label  pointi,
const pointEdgePoint pointInfo,
const scalar  tol,
TrackingData &  td 
)
inline

Influence of point on edge.

Definition at line 259 of file pointEdgePointI.H.

References Foam::constant::electromagnetic::e, primitiveMesh::edges(), mesh, polyMesh::points(), and update().

Here is the call graph for this function:

◆ equal()

bool equal ( const pointEdgePoint rhs,
TrackingData &  td 
) const
inline

Test for equality, with TrackingData.

Definition at line 275 of file pointEdgePointI.H.

References Foam::operator==().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator<<

Ostream& operator<< ( Ostream ,
const pointEdgePoint  
)
friend

◆ operator>>

Istream& operator>> ( Istream ,
pointEdgePoint  
)
friend

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