OpenFOAM logo
Open Source CFD Toolkit

face Class Reference

Inheritance diagram for face:

Inheritance graph
[legend]
Collaboration diagram for face:

Collaboration graph
[legend]
List of all members.

Public Member Functions

 face ()
 Construct null.
 face (label)
 Construct given size.
 face (const labelList &)
 Construct from components.
 face (Istream &)
 Construct from Istream.
void collapse ()
 Collapse face by removing duplicate point labels.
pointField points (const pointField &meshPoints) const
 Return the points corresponding to this face.
point centre (const pointField &) const
 Centre point of face.
scalar mag (const pointField &) const
 Scalar magnitude.
vector normal (const pointField &) const
 Vector normal; magnitude is equal to area of face.
face reverseFace () const
 Return face with reverse direction.
label which (const label globalIndex) const
 Navigation through face vertices Which vertex on face (face index given a global index).
label nextLabel (const label i) const
 Next vertex on face.
label prevLabel (const label i) const
 Previous vertex on face.
scalar sweptVol (const pointField &oldPoints, const pointField &newPoints) const
 Return the volume swept out by the face when its points move.
pointHit ray (const point &p, const vector &n, const pointField &meshPoints, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
 Return potential intersection with face with a ray starting.
pointHit nearestPoint (const point &p, const pointField &meshPoints) const
 Return nearest point to face.
scalar contactSphereDiameter (const point &p, const vector &n, const pointField &meshPoints) const
 Return contact sphere diameter.
scalar areaInContact (const pointField &points, const scalarField &v) const
 Return area in contact, given the displacement in vertices.
label nEdges () const
 Return number of edges.
edgeList edges () const
 Return edges in face point ordering, i.e. edges()[0] is edge.
edge faceEdge (const label n) const
 Return n-th face edge.
label nTriangles (const pointField &points) const
 Number of triangles after splitting.
void triangles (const pointField &points, label &triI, faceList &triFaces) const
 Split into triangles using existing points. Result in.
void nTrianglesQuads (const pointField &points, label &nTris, label &nQuads) const
 Number of triangles and quads after splitting.
void trianglesQuads (const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
 Split into triangles and quads. Result in triFaces (starting at.

Friends

bool operator== (const face &, const face &)
bool operator!= (const face &, const face &)
Istreamoperator>> (Istream &, face &)

Constructor & Destructor Documentation

face  )  [inline]
 

Construct null.

Definition at line 61 of file faceI.H.

face label   )  [inline, explicit]
 

Construct given size.

Definition at line 66 of file faceI.H.

References Foam::labelList.

face const labelList  )  [inline, explicit]
 

Construct from components.

Definition at line 73 of file faceI.H.

face Istream  )  [inline]
 

Construct from Istream.

Definition at line 80 of file faceI.H.

References forAll, p, Foam::pointField, and UList::size().

Here is the call graph for this function:


Member Function Documentation

scalar areaInContact const pointField points,
const scalarField v
const
 

Return area in contact, given the displacement in vertices.

point centre const pointField  )  const
 

Centre point of face.

void collapse  ) 
 

Collapse face by removing duplicate point labels.

scalar contactSphereDiameter const point p,
const vector n,
const pointField meshPoints
const
 

Return contact sphere diameter.

edgeList edges  )  const
 

Return edges in face point ordering, i.e. edges()[0] is edge.

between [0] and [1]

edge faceEdge const label  n  )  const [inline]
 

Return n-th face edge.

Definition at line 118 of file faceI.H.

References UList::fcIndex(), and UList::operator[]().

Here is the call graph for this function:

scalar mag const pointField  )  const [inline]
 

Scalar magnitude.

Definition at line 105 of file faceI.H.

pointHit nearestPoint const point p,
const pointField meshPoints
const
 

Return nearest point to face.

label nEdges  )  const [inline]
 

Return number of edges.

Definition at line 111 of file faceI.H.

References UList::fcIndex().

Here is the call graph for this function:

label nextLabel const label  i  )  const [inline]
 

Next vertex on face.

Definition at line 125 of file faceI.H.

References UList::operator[](), and UList::rcIndex().

Here is the call graph for this function:

vector normal const pointField  )  const
 

Vector normal; magnitude is equal to area of face.

label nTriangles const pointField points  )  const
 

Number of triangles after splitting.

void nTrianglesQuads const pointField points,
label nTris,
label nQuads
const
 

Number of triangles and quads after splitting.

pointField points const pointField meshPoints  )  const [inline]
 

Return the points corresponding to this face.

Definition at line 88 of file faceI.H.

References UList::operator[](), and p.

Here is the call graph for this function:

label prevLabel const label  i  )  const [inline]
 

Previous vertex on face.

Definition at line 132 of file faceI.H.

References Istream::readBegin(), and IOstream::version().

Here is the call graph for this function:

pointHit ray const point p,
const vector n,
const pointField meshPoints,
const intersection::algorithm  alg = intersection::FULL_RAY,
const intersection::direction  dir = intersection::VECTOR
const
 

Return potential intersection with face with a ray starting.

at p, direction n (does not need to be normalized) Does face-center decomposition and returns triangle intersection point closest to p. For a hit, the distance is signed. Positive number represents the point in front of triangle In case of miss the point is the nearest point on the face and the distance is the distance between the intersection point and the original point. The half-ray or full-ray intersection and the contact sphere adjustment of the projection vector is set by the intersection parameters

face reverseFace  )  const
 

Return face with reverse direction.

scalar sweptVol const pointField oldPoints,
const pointField newPoints
const
 

Return the volume swept out by the face when its points move.

void triangles const pointField points,
label triI,
faceList triFaces
const
 

Split into triangles using existing points. Result in.

triFaces[triI..triI+nTri]. Splits intelligently to maximize triangle quality.

void trianglesQuads const pointField points,
label triI,
label quadI,
faceList triFaces,
faceList quadFaces
const
 

Split into triangles and quads. Result in triFaces (starting at.

triI and quadFaces (starting at quadI)

label which const label  globalIndex  )  const
 

Navigation through face vertices Which vertex on face (face index given a global index).


Friends And Related Function Documentation

bool operator!= const face ,
const face
[friend]
 

bool operator== const face ,
const face
[friend]
 

Istream& operator>> Istream is,
face f
[friend]
 

Definition at line 140 of file faceI.H.


The documentation for this class was generated from the following files:
For further information go to www.openfoam.org