Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Friends
face Class Reference

A face is a list of labels corresponding to mesh vertices. More...

Inheritance diagram for face:
Inheritance graph
[legend]
Collaboration diagram for face:
Collaboration graph
[legend]

Public Types

enum  proxType { NONE, POINT, EDGE }
 Return types for classify. More...
 
- Public Types inherited from List
typedef SubList< TsubList
 Declare type of subList. More...
 

Public Member Functions

 face ()
 Construct null. More...
 
 face (label)
 Construct given size. More...
 
 face (const labelUList &)
 Construct from list of labels. More...
 
 face (const labelList &)
 Construct from list of labels. More...
 
 face (const Xfer< labelList > &)
 Construct by transferring the parameter contents. More...
 
 face (const triFace &)
 Copy construct from triFace. More...
 
 face (Istream &)
 Construct from Istream. More...
 
label collapse ()
 Collapse face by removing duplicate point labels. More...
 
void flip ()
 Flip the face in-place. More...
 
pointField points (const pointField &) const
 Return the points corresponding to this face. More...
 
point centre (const pointField &) const
 Centre point of face. More...
 
template<class Type >
Type average (const pointField &, const Field< Type > &) const
 Calculate average value at centroid of face. More...
 
scalar mag (const pointField &) const
 Magnitude of face area. More...
 
vector normal (const pointField &) const
 Vector normal; magnitude is equal to area of face. More...
 
face reverseFace () const
 Return face with reverse direction. More...
 
label which (const label globalIndex) const
 Navigation through face vertices. More...
 
label nextLabel (const label i) const
 Next vertex on face. More...
 
label prevLabel (const label i) const
 Previous vertex on face. More...
 
scalar sweptVol (const pointField &oldPoints, const pointField &newPoints) const
 Return the volume swept out by the face when its points move. More...
 
tensor inertia (const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
 Return the inertia tensor, with optional reference. More...
 
pointHit ray (const point &p, const vector &n, const pointField &, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
 Return potential intersection with face with a ray starting. More...
 
pointHit intersection (const point &p, const vector &q, const point &ctr, const pointField &, const intersection::algorithm alg, const scalar tol=0.0) const
 Fast intersection with a ray. More...
 
pointHit nearestPoint (const point &p, const pointField &) const
 Return nearest point to face. More...
 
pointHit nearestPointClassify (const point &p, const pointField &, label &nearType, label &nearLabel) const
 Return nearest point to face and classify it: More...
 
scalar contactSphereDiameter (const point &p, const vector &n, const pointField &) const
 Return contact sphere diameter. More...
 
scalar areaInContact (const pointField &, const scalarField &v) const
 Return area in contact, given the displacement in vertices. More...
 
label nEdges () const
 Return number of edges. More...
 
edgeList edges () const
 Return edges in face point ordering,. More...
 
edge faceEdge (const label n) const
 Return n-th face edge. More...
 
int edgeDirection (const edge &) const
 Return the edge direction on the face. More...
 
label nTriangles () const
 Number of triangles after splitting. More...
 
label nTriangles (const pointField &points) const
 Number of triangles after splitting. More...
 
label triangles (const pointField &points, label &triI, faceList &triFaces) const
 Split into triangles using existing points. More...
 
template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
label triangles (const pointField &points, DynamicList< face, SizeInc, SizeMult, SizeDiv > &triFaces) const
 Split into triangles using existing points. More...
 
label nTrianglesQuads (const pointField &points, label &nTris, label &nQuads) const
 Number of triangles and quads after splitting. More...
 
label trianglesQuads (const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
 Split into triangles and quads. More...
 
template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
Foam::label triangles (const pointField &points, DynamicList< face, SizeInc, SizeMult, SizeDiv > &triFaces) const
 
- Public Member Functions inherited from List
 List ()
 Null constructor. More...
 
 List (const label)
 Construct with given size. More...
 
 List (const label, const T &)
 Construct with given size and value for all elements. More...
 
 List (const List< T > &)
 Copy constructor. More...
 
 List (const Xfer< List< T > > &)
 Construct by transferring the parameter contents. More...
 
 List (List< T > &, bool reUse)
 Construct as copy or re-use as specified. More...
 
 List (const UList< T > &, const labelUList &mapAddressing)
 Construct as subset. More...
 
template<class InputIterator >
 List (InputIterator first, InputIterator last)
 Construct given start and end iterators. More...
 
template<unsigned Size>
 List (const FixedList< T, Size > &)
 Construct as copy of FixedList<T, Size> More...
 
 List (const PtrList< T > &)
 Construct as copy of PtrList<T> More...
 
 List (const SLList< T > &)
 Construct as copy of SLList<T> More...
 
 List (const UIndirectList< T > &)
 Construct as copy of UIndirectList<T> More...
 
 List (const BiIndirectList< T > &)
 Construct as copy of BiIndirectList<T> More...
 
 List (Istream &)
 Construct from Istream. More...
 
autoPtr< List< T > > clone () const
 Clone. More...
 
 ~List ()
 Destructor. More...
 
label size () const
 Return the number of elements in the UList. More...
 
void resize (const label)
 Alias for setSize(const label) More...
 
void resize (const label, const T &)
 Alias for setSize(const label, const T&) More...
 
void setSize (const label)
 Reset size of List. More...
 
void setSize (const label, const T &)
 Reset size of List and value for new elements. More...
 
void clear ()
 Clear the list, i.e. set size to zero. More...
 
void append (const T &)
 Append an element at the end of the list. More...
 
void append (const UList< T > &)
 Append a List at the end of this list. More...
 
void append (const UIndirectList< T > &)
 Append a UIndirectList at the end of this list. More...
 
void transfer (List< T > &)
 Transfer the contents of the argument List into this list. More...
 
template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
void transfer (DynamicList< T, SizeInc, SizeMult, SizeDiv > &)
 Transfer the contents of the argument List into this list. More...
 
void transfer (SortableList< T > &)
 Transfer the contents of the argument List into this list. More...
 
Xfer< List< T > > xfer ()
 Transfer contents to the Xfer container. More...
 
TnewElmt (const label)
 Return subscript-checked element of UList. More...
 
void operator= (const UList< T > &)
 Assignment from UList operator. Takes linear time. More...
 
void operator= (const List< T > &)
 Assignment operator. Takes linear time. More...
 
void operator= (const SLList< T > &)
 Assignment from SLList operator. Takes linear time. More...
 
void operator= (const UIndirectList< T > &)
 Assignment from UIndirectList operator. Takes linear time. More...
 
void operator= (const BiIndirectList< T > &)
 Assignment from BiIndirectList operator. Takes linear time. More...
 
void operator= (const T &)
 Assignment of all entries to the given value. More...
 
template<class T >
 List (const label s)
 
template<class T >
 List (const label s, const T &a)
 
template<class T >
 List (const List< T > &a)
 
template<class T >
 List (const Xfer< List< T > > &lst)
 
template<class T >
 List (List< T > &a, bool reUse)
 
template<class T >
 List (const UList< T > &a, const labelUList &map)
 
template<class T >
 List (const PtrList< T > &lst)
 
template<class T >
 List (const SLList< T > &lst)
 
template<class T >
 List (const UIndirectList< T > &lst)
 
template<class T >
 List (const BiIndirectList< T > &lst)
 
template<class T >
void transfer (List< T > &a)
 
template<class T >
void transfer (SortableList< T > &a)
 
template<class T >
 List ()
 
template<class T >
 List (Istream &is)
 

Static Public Member Functions

static int compare (const face &, const face &)
 Compare faces. More...
 
static bool sameVertices (const face &, const face &)
 Return true if the faces have the same vertices. More...
 
- Static Public Member Functions inherited from List
static const List< T > & null ()
 Return a null List. More...
 

Static Public Attributes

static const char *const typeName = "face"
 

Private Types

enum  splitMode { COUNTTRIANGLE, COUNTQUAD, SPLITTRIANGLE, SPLITQUAD }
 Enumeration listing the modes for split() More...
 

Private Member Functions

label right (const label i) const
 Edge to the right of face vertex i. More...
 
label left (const label i) const
 Edge to the left of face vertex i. More...
 
tmp< vectorFieldcalcEdges (const pointField &points) const
 Construct list of edge vectors for face. More...
 
scalar edgeCos (const vectorField &edges, const label index) const
 Cos between neighbouring edges. More...
 
label mostConcaveAngle (const pointField &points, const vectorField &edges, scalar &edgeCos) const
 Find index of largest internal angle on face. More...
 
label split (const splitMode mode, const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
 Split face into triangles or triangles&quads. More...
 

Friends

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

Additional Inherited Members

- Protected Member Functions inherited from List
void size (const label)
 Override size to be inconsistent with allocated storage. More...
 

Detailed Description

A face is a list of labels corresponding to mesh vertices.

See also
Foam::triFace
Source files

Definition at line 75 of file face.H.

Member Enumeration Documentation

◆ splitMode

enum splitMode
private

Enumeration listing the modes for split()

Enumerator
COUNTTRIANGLE 
COUNTQUAD 
SPLITTRIANGLE 
SPLITQUAD 

Definition at line 109 of file face.H.

◆ proxType

enum proxType

Return types for classify.

Enumerator
NONE 
POINT 
EDGE 

Definition at line 134 of file face.H.

Constructor & Destructor Documentation

◆ face() [1/7]

face ( )
inline

Construct null.

Definition at line 44 of file faceI.H.

◆ face() [2/7]

face ( label  s)
inlineexplicit

Construct given size.

Definition at line 48 of file faceI.H.

◆ face() [3/7]

face ( const labelUList lst)
inlineexplicit

Construct from list of labels.

Definition at line 54 of file faceI.H.

◆ face() [4/7]

face ( const labelList lst)
inlineexplicit

Construct from list of labels.

Definition at line 60 of file faceI.H.

◆ face() [5/7]

face ( const Xfer< labelList > &  lst)
inlineexplicit

Construct by transferring the parameter contents.

Definition at line 66 of file faceI.H.

◆ face() [6/7]

face ( const triFace f)

Copy construct from triFace.

Definition at line 290 of file face.C.

◆ face() [7/7]

face ( Istream is)
inline

Construct from Istream.

Definition at line 72 of file faceI.H.

Member Function Documentation

◆ right()

Foam::label right ( const label  i) const
inlineprivate

Edge to the right of face vertex i.

Definition at line 29 of file faceI.H.

◆ left()

Foam::label left ( const label  i) const
inlineprivate

Edge to the left of face vertex i.

Definition at line 36 of file faceI.H.

◆ calcEdges()

Foam::tmp< Foam::vectorField > calcEdges ( const pointField points) const
private

Construct list of edge vectors for face.

Definition at line 41 of file face.C.

References face::edges(), forAll, Foam::mag(), face::points(), and List::size().

Here is the call graph for this function:

◆ edgeCos()

Foam::scalar edgeCos ( const vectorField edges,
const label  index 
) const
private

Cos between neighbouring edges.

Definition at line 64 of file face.C.

◆ mostConcaveAngle()

Foam::label mostConcaveAngle ( const pointField points,
const vectorField edges,
scalar &  edgeCos 
) const
private

Find index of largest internal angle on face.

Definition at line 78 of file face.C.

References Foam::acos(), forAll, Foam::max(), Foam::min(), n, Foam::constant::mathematical::pi(), and points.

Here is the call graph for this function:

◆ split()

Foam::label split ( const splitMode  mode,
const pointField points,
label triI,
label quadI,
faceList triFaces,
faceList quadFaces 
) const
private

Split face into triangles or triangles&quads.

Stores results quadFaces[quadI], triFaces[triI] Returns number of new faces created

Definition at line 125 of file face.C.

References Foam::abort(), Foam::acos(), Foam::diff(), Foam::FatalError, FatalErrorInFunction, Foam::mag(), Foam::max(), Foam::min(), Foam::mode(), Foam::constant::mathematical::pi(), points, face::split(), and triFace().

Referenced by face::split().

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

◆ collapse()

Foam::label collapse ( )

Collapse face by removing duplicate point labels.

return the collapsed size

Definition at line 449 of file face.C.

References setSize().

Here is the call graph for this function:

◆ flip()

void flip ( )

Flip the face in-place.

The starting points of the original and flipped face are identical.

Definition at line 474 of file face.C.

References n, and Foam::Swap().

Referenced by slidingInterface::decoupleInterface(), main(), layerAdditionRemoval::removeCellLayer(), layerAdditionRemoval::setLayerPairing(), and hexRef8::storeMidPointInfo().

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

◆ points()

Foam::pointField points ( const pointField meshPoints) const
inline

Return the points corresponding to this face.

Definition at line 80 of file faceI.H.

References forAll, and p.

Referenced by directAMI< SourcePatch, TargetPatch >::appendToDirectSeeds(), enrichedPatch::calcCutFaces(), face::calcEdges(), faceCoupleInfo::checkMatch(), faceCoupleInfo::growCutFaces(), and AMIMethod< SourcePatch, TargetPatch >::writeIntersectionOBJ().

Here is the caller graph for this function:

◆ centre()

Foam::point centre ( const pointField points) const

◆ average()

Type average ( const pointField meshPoints,
const Field< Type > &  fld 
) const

Calculate average value at centroid of face.

Definition at line 51 of file faceTemplates.C.

References fld(), Foam::mag(), and nPoints.

Here is the call graph for this function:

◆ mag()

Foam::scalar mag ( const pointField p) const
inline

Magnitude of face area.

Definition at line 97 of file faceI.H.

References Foam::mag(), and p.

Referenced by face::areaInContact(), hexCellLooper::cut(), faceAreaWeightAMI< SourcePatch, TargetPatch >::interArea(), and autoSnapDriver::preventFaceSqueeze().

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

◆ normal()

Foam::vector normal ( const pointField p) const

◆ reverseFace()

Foam::face reverseFace ( ) const

Return face with reverse direction.

The starting points of the original and reverse face are identical.

Definition at line 611 of file face.C.

References f(), List::size(), and Foam::xferMove().

Referenced by hexRef8::addFace(), meshCutAndRemove::addFace(), meshCutter::addFace(), slidingInterface::coupleInterface(), meshOctreeAddressing::createOctreeFaces(), decomposeFaces::decomposeMeshFaces(), extrudedMesh::extrudedFaces(), meshCutAndRemove::modFace(), hexRef8::modFace(), meshCutter::modFace(), surfaceMorpherCells::morphInternalFaces(), perfectInterface::setRefinement(), and createShellMesh::setRefinement().

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

◆ which()

Foam::label which ( const label  globalIndex) const

◆ nextLabel()

Foam::label nextLabel ( const label  i) const
inline

◆ prevLabel()

Foam::label prevLabel ( const label  i) const
inline

◆ sweptVol()

Foam::scalar sweptVol ( const pointField oldPoints,
const pointField newPoints 
) const

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

Definition at line 647 of file face.C.

References nPoints, and Foam::constant::mathematical::pi().

Here is the call graph for this function:

◆ inertia()

Foam::tensor inertia ( const pointField p,
const point refPt = vector::zero,
scalar  density = 1.0 
) const

Return the inertia tensor, with optional reference.

point and density specification

Definition at line 726 of file face.C.

References forAll, p, and Tensor::zero.

◆ ray()

Foam::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-centre decomposition and returns triangle intersection point closest to p. Face-centre is calculated from point average. 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

Definition at line 43 of file faceIntersection.C.

References Foam::average(), PointHit::distance(), PointHit::eligibleMiss(), f(), PointHit::hit(), PointHit::hitPoint(), Foam::mag(), PointHit::missPoint(), n, nPoints, p, points, PointHit::setDistance(), PointHit::setHit(), PointHit::setMiss(), and PointHit::setPoint().

Referenced by directAMI< SourcePatch, TargetPatch >::appendToDirectSeeds().

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

◆ intersection()

Foam::pointHit intersection ( const point p,
const vector q,
const point ctr,
const pointField meshPoints,
const intersection::algorithm  alg,
const scalar  tol = 0.0 
) const

Fast intersection with a ray.

Does face-centre decomposition and returns triangle intersection point closest to p. See triangle::intersection for details.

Definition at line 141 of file faceIntersection.C.

References PointHit::distance(), f(), forAll, PointHit::hit(), PointHit::hitPoint(), Foam::mag(), p, PointHit::setDistance(), PointHit::setHit(), and PointHit::setPoint().

Here is the call graph for this function:

◆ nearestPoint()

Foam::pointHit nearestPoint ( const point p,
const pointField meshPoints 
) const

Return nearest point to face.

Definition at line 199 of file faceIntersection.C.

References p.

Referenced by faceCoupleInfo::maxDistance(), and PairCollision< CloudType >::wallInteraction().

Here is the caller graph for this function:

◆ nearestPointClassify()

Foam::pointHit nearestPointClassify ( const point p,
const pointField meshPoints,
label nearType,
label nearLabel 
) const

Return nearest point to face and classify it:

+ near point (nearType=POINT, nearLabel=0, 1, 2) + near edge (nearType=EDGE, nearLabel=0, 1, 2) Note: edges are counted from starting vertex so e.g. edge n is from f[n] to f[0], where the face has n + 1 points

Definition at line 213 of file faceIntersection.C.

References PointHit::distance(), f(), PointHit::hit(), PointHit::hitPoint(), Foam::mag(), PointHit::missPoint(), triangle::nearestPointClassify(), nPoints, p, PointHit::setDistance(), PointHit::setHit(), PointHit::setMiss(), and PointHit::setPoint().

Here is the call graph for this function:

◆ contactSphereDiameter()

Foam::scalar contactSphereDiameter ( const point p,
const vector n,
const pointField meshPoints 
) const

Return contact sphere diameter.

Definition at line 34 of file faceContactSphere.C.

References Foam::mag(), n, and p.

Here is the call graph for this function:

◆ areaInContact()

Foam::scalar areaInContact ( const pointField meshPoints,
const scalarField v 
) const

Return area in contact, given the displacement in vertices.

Definition at line 36 of file faceAreaInContact.C.

References forAll, face::mag(), Foam::mag(), and List::size().

Here is the call graph for this function:

◆ nEdges()

Foam::label nEdges ( ) const
inline

Return number of edges.

Definition at line 103 of file faceI.H.

◆ edges()

Foam::edgeList edges ( ) const

◆ faceEdge()

Foam::edge faceEdge ( const label  n) const
inline

◆ edgeDirection()

int edgeDirection ( const edge e) const

Return the edge direction on the face.

Returns:

  • 0: edge not found on the face
  • +1: forward (counter-clockwise) on the face
  • -1: reverse (clockwise) on the face

Definition at line 779 of file face.C.

References Foam::e, and forAll.

◆ nTriangles() [1/2]

Foam::label nTriangles ( ) const
inline

Number of triangles after splitting.

Definition at line 130 of file faceI.H.

◆ nTriangles() [2/2]

Foam::label nTriangles ( const pointField points) const

Number of triangles after splitting.

Definition at line 823 of file face.C.

◆ triangles() [1/3]

Foam::label 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. Returns number of faces created.

Definition at line 830 of file face.C.

References points.

Referenced by faceAreaIntersect::calc(), and thresholdCellFaces::calculate().

Here is the caller graph for this function:

◆ triangles() [2/3]

label triangles ( const pointField points,
DynamicList< face, SizeInc, SizeMult, SizeDiv > &  triFaces 
) const

Split into triangles using existing points.

Append to DynamicList. Returns number of faces created.

◆ nTrianglesQuads()

Foam::label nTrianglesQuads ( const pointField points,
label nTris,
label nQuads 
) const

Number of triangles and quads after splitting.

Returns the sum of both

Definition at line 844 of file face.C.

References points.

◆ trianglesQuads()

Foam::label trianglesQuads ( const pointField points,
label triI,
label quadI,
faceList triFaces,
faceList quadFaces 
) const

Split into triangles and quads.

Results in triFaces (starting at triI) and quadFaces (starting at quadI). Returns number of new faces created.

Definition at line 858 of file face.C.

References points.

◆ compare()

int compare ( const face a,
const face b 
)
static

Compare faces.

0: different +1: identical -1: same face, but different orientation

Definition at line 298 of file face.C.

References CirculatorBase::ANTICLOCKWISE, Foam::constant::physicoChemical::b, ConstCirculator< ContainerType >::circulate(), CirculatorBase::CLOCKWISE, ConstCirculator< ContainerType >::setFulcrumToIterator(), ConstCirculator< ContainerType >::setIteratorToFulcrum(), and List::size().

Referenced by meshReader::createPolyCells(), main(), Foam::operator!=(), and Foam::operator==().

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

◆ sameVertices()

bool sameVertices ( const face a,
const face b 
)
static

Return true if the faces have the same vertices.

Definition at line 400 of file face.C.

References Foam::constant::physicoChemical::b, forAll, and List::size().

Here is the call graph for this function:

◆ triangles() [3/3]

Foam::label triangles ( const pointField points,
DynamicList< face, SizeInc, SizeMult, SizeDiv > &  triFaces 
) const

Definition at line 33 of file faceTemplates.C.

References points, and DynamicList::setSize().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator==

bool operator== ( const face a,
const face b 
)
friend

◆ operator!=

bool operator!= ( const face a,
const face b 
)
friend

◆ operator>>

Istream& operator>> ( Istream ,
face  
)
friend

Field Documentation

◆ typeName

const char *const typeName = "face"
static

Definition at line 143 of file face.H.


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