Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Friends | List of all members
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< label >
typedef SubList< labelsubList
 Declare type of subList. More...
 
- Public Types inherited from UList< T >
typedef T value_type
 Type of values the UList contains. More...
 
typedef Treference
 Type that can be used for storing into. More...
 
typedef const Tconst_reference
 Type that can be used for storing into. More...
 
typedef label difference_type
 The type that can represent the difference between any two. More...
 
typedef label size_type
 The type that can represent the size of a UList. More...
 
typedef Titerator
 Random access iterator for traversing UList. More...
 
typedef const Tconst_iterator
 Random access iterator for traversing UList. More...
 
typedef Treverse_iterator
 Reverse iterator for reverse traversal of UList. More...
 
typedef const Tconst_reverse_iterator
 Reverse iterator for reverse traversal of constant UList. 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< label >
label size () const
 Return the number of elements in the UList. More...
 
 List ()
 Null constructor. More...
 
 List (const label)
 Construct with given size. More...
 
 List (const label, const label &)
 Construct with given size and value for all elements. More...
 
 List (const label, const zero)
 Construct with given size initializing all elements to zero. More...
 
 List (const List< label > &)
 Copy constructor. More...
 
 List (const Xfer< List< label >> &)
 Construct by transferring the parameter contents. More...
 
 List (List< label > &, bool reuse)
 Construct as copy or re-use as specified. More...
 
 List (const UList< label > &, const labelUList &mapAddressing)
 Construct as subset. More...
 
 List (const FixedList< label, Size > &)
 Construct as copy of FixedList<T, Size> More...
 
 List (const PtrList< label > &)
 Construct as copy of PtrList<T> More...
 
 List (const SLList< label > &)
 Construct as copy of SLList<T> More...
 
 List (const UIndirectList< label > &)
 Construct as copy of UIndirectList<T> More...
 
 List (const BiIndirectList< label > &)
 Construct as copy of BiIndirectList<T> More...
 
 List (Istream &)
 Construct from Istream. More...
 
autoPtr< List< label > > clone () const
 Clone. More...
 
 ~List ()
 Destructor. More...
 
void resize (const label)
 Alias for setSize(const label) More...
 
void resize (const label, const label &)
 Alias for setSize(const label, const T&) More...
 
void setSize (const label)
 Reset size of List. More...
 
void setSize (const label, const label &)
 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 label &)
 Append an element at the end of the list. More...
 
void append (const UList< label > &)
 Append a List at the end of this list. More...
 
void append (const UIndirectList< label > &)
 Append a UIndirectList at the end of this list. More...
 
void transfer (List< label > &)
 Transfer the contents of the argument List into this list. More...
 
void transfer (DynamicList< label, SizeInc, SizeMult, SizeDiv > &)
 Transfer the contents of the argument List into this list. More...
 
void transfer (SortableList< label > &)
 Transfer the contents of the argument List into this list. More...
 
Xfer< List< label > > xfer ()
 Transfer contents to the Xfer container. More...
 
labelnewElmt (const label)
 Return subscript-checked element of UList. More...
 
void shallowCopy (const UList< label > &)=delete
 Disallow implicit shallowCopy. More...
 
void operator= (const UList< label > &)
 Assignment from UList operator. Takes linear time. More...
 
void operator= (const List< label > &)
 Assignment operator. Takes linear time. More...
 
void operator= (const SLList< label > &)
 Assignment from SLList operator. Takes linear time. More...
 
void operator= (const UIndirectList< label > &)
 Assignment from UIndirectList operator. Takes linear time. More...
 
void operator= (const BiIndirectList< label > &)
 Assignment from BiIndirectList operator. Takes linear time. More...
 
void operator= (const label &)
 Assignment of all entries to the given value. More...
 
void operator= (const zero)
 Assignment of all entries to zero. More...
 
- Public Member Functions inherited from UList< T >
 UList ()
 Null constructor. More...
 
 UList (T *__restrict__ v, label size)
 Construct from components. More...
 
label fcIndex (const label i) const
 Return the forward circular index, i.e. the next index. More...
 
label rcIndex (const label i) const
 Return the reverse circular index, i.e. the previous index. More...
 
std::streamsize byteSize () const
 Return the binary size in number of characters of the UList. More...
 
const Tcdata () const
 Return a const pointer to the first data element,. More...
 
Tdata ()
 Return a pointer to the first data element,. More...
 
Tfirst ()
 Return the first element of the list. More...
 
const Tfirst () const
 Return first element of the list. More...
 
Tlast ()
 Return the last element of the list. More...
 
const Tlast () const
 Return the last element of the list. More...
 
void checkStart (const label start) const
 Check start is within valid range (0 ... size-1) More...
 
void checkSize (const label size) const
 Check size is within valid range (0 ... size) More...
 
void checkIndex (const label i) const
 Check index i is within valid range (0 ... size-1) More...
 
void shallowCopy (const UList< T > &)
 Copy the pointer held by the given UList. More...
 
void deepCopy (const UList< T > &)
 Copy elements of the given UList. More...
 
void writeEntry (Ostream &) const
 Write the UList as a dictionary entry. More...
 
void writeEntry (const word &keyword, Ostream &) const
 Write the UList as a dictionary entry with keyword. More...
 
Toperator[] (const label)
 Return element of UList. More...
 
const Toperator[] (const label) const
 Return element of constant UList. More...
 
 operator const Foam::List< T > & () const
 Allow cast to a const List<T>&. More...
 
void operator= (const T &)
 Assignment of all entries to the given value. More...
 
void operator= (const zero)
 Assignment of all entries to zero. More...
 
iterator begin ()
 Return an iterator to begin traversing the UList. More...
 
iterator end ()
 Return an iterator to end traversing the UList. More...
 
const_iterator cbegin () const
 Return const_iterator to begin traversing the constant UList. More...
 
const_iterator cend () const
 Return const_iterator to end traversing the constant UList. More...
 
const_iterator begin () const
 Return const_iterator to begin traversing the constant UList. More...
 
const_iterator end () const
 Return const_iterator to end traversing the constant UList. More...
 
reverse_iterator rbegin ()
 Return reverse_iterator to begin reverse traversing the UList. More...
 
reverse_iterator rend ()
 Return reverse_iterator to end reverse traversing the UList. More...
 
const_reverse_iterator crbegin () const
 Return const_reverse_iterator to begin reverse traversing the UList. More...
 
const_reverse_iterator crend () const
 Return const_reverse_iterator to end reverse traversing the UList. More...
 
const_reverse_iterator rbegin () const
 Return const_reverse_iterator to begin reverse traversing the UList. More...
 
const_reverse_iterator rend () const
 Return const_reverse_iterator to end reverse traversing the UList. More...
 
label size () const
 Return the number of elements in the UList. More...
 
label max_size () const
 Return size of the largest possible UList. More...
 
bool empty () const
 Return true if the UList is empty (ie, size() is zero) More...
 
void swap (UList< T > &)
 Swap two ULists of the same type in constant time. More...
 
bool operator== (const UList< T > &) const
 Equality operation on ULists of the same type. More...
 
bool operator!= (const UList< T > &) const
 The opposite of the equality operation. Takes linear time. More...
 
bool operator< (const UList< T > &) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator> (const UList< T > &) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator<= (const UList< T > &) const
 Return true if !(a > b). Takes linear time. More...
 
bool operator>= (const UList< T > &) const
 Return true if !(a < b). Takes linear time. More...
 
template<>
const bool & operator[] (const label i) const
 

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< label >
static const List< label > & null ()
 Return a null List. More...
 
- Static Public Member Functions inherited from UList< T >
static const UList< T > & null ()
 Return a null UList. More...
 

Static Public Attributes

static const char *const typeName = "face"
 

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< label >
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

enum proxType

Return types for classify.

Enumerator
NONE 
POINT 
EDGE 

Definition at line 134 of file face.H.

Constructor & Destructor Documentation

face ( )
inline

Construct null.

Definition at line 44 of file faceI.H.

Referenced by face::reverseFace().

Here is the caller graph for this function:

face ( label  s)
inlineexplicit

Construct given size.

Definition at line 48 of file faceI.H.

face ( const labelUList lst)
inlineexplicit

Construct from list of labels.

Definition at line 54 of file faceI.H.

face ( const labelList lst)
inlineexplicit

Construct from list of labels.

Definition at line 60 of file faceI.H.

face ( const Xfer< labelList > &  lst)
inlineexplicit

Construct by transferring the parameter contents.

Definition at line 66 of file faceI.H.

face ( const triFace f)

Copy construct from triFace.

Definition at line 290 of file face.C.

face ( Istream is)
inline

Construct from Istream.

Definition at line 72 of file faceI.H.

Member Function Documentation

Foam::label collapse ( )

Collapse face by removing duplicate point labels.

return the collapsed size

Definition at line 449 of file face.C.

References UList< T >::operator[](), List< label >::setSize(), and List< label >::size().

Referenced by STARCD::readCells().

Here is the call graph for this function:

Here is the caller graph for this function:

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, List< label >::size(), and Foam::Swap().

Referenced by hexRef8::faceLevel().

Here is the call graph for this function:

Here is the caller graph for this function:

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, UList< T >::operator[](), p, and List< label >::size().

Referenced by AMIMethod< SourcePatch, TargetPatch >::writeIntersectionOBJ().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::point centre ( const pointField points) const

Centre point of face.

Definition at line 488 of file face.C.

References Foam::mag(), nPoints, UList< T >::operator[](), List< label >::size(), and Foam::Zero.

Referenced by searchableBox::boundingSpheres(), edgeCollapser::checkMeshQuality(), hexRef8::faceLevel(), AMIMethod< SourcePatch, TargetPatch >::findTargetFace(), face::inertia(), cellCuts::nonAnchorPoints(), and face::sweptVol().

Here is the call graph for this function:

Here is the caller graph for this function:

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 Foam::mag(), nPoints, and Foam::Zero.

Referenced by face::triangles().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Magnitude of face area.

Definition at line 97 of file faceI.H.

References Foam::mag(), and face::normal().

Referenced by face::areaInContact(), edgeCollapser::checkMeshQuality(), faceAreaWeightAMI< SourcePatch, TargetPatch >::interArea(), and listPlusEqOp< T >::operator()().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::vector normal ( const pointField p) const

Vector normal; magnitude is equal to area of face.

Definition at line 552 of file face.C.

References n, nPoints, UList< T >::operator[](), List< label >::size(), and Foam::Zero.

Referenced by hexRef8::faceLevel(), faceTriangulation::faceTriangulation(), face::mag(), cellCuts::nonAnchorPoints(), and snappySnapDriver::repatchToSurface().

Here is the call graph for this function:

Here is the caller graph for this function:

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(), face::face(), List< T >::size(), List< label >::size(), and Foam::xferMove().

Referenced by NamedEnum< Enum, nEnum >::names(), ifEqEqOp< value >::operator()(), hexRef8::setInstance(), removeCells::setRefinement(), and createShellMesh::setRefinement().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::label which ( const label  globalIndex) const

Navigation through face vertices.

Which vertex on face (face index given a global index)

returns -1 if not found

Definition at line 630 of file face.C.

References f(), forAll, and face::sweptVol().

Here is the call graph for this function:

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

Previous vertex on face.

Definition at line 124 of file faceI.H.

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

Referenced by edgeCollapser::checkMeshQuality(), and tetDecomposer::setRefinement().

Here is the call graph for this function:

Here is the caller graph for this function:

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 face::centre(), face::inertia(), nPoints, Foam::constant::mathematical::pi(), and List< label >::size().

Referenced by face::which().

Here is the call graph for this function:

Here is the caller graph for this function:

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 face::centre(), UList< T >::fcIndex(), forAll, List< label >::size(), and Foam::Zero.

Referenced by edgeCollapser::checkMeshQuality(), and face::sweptVol().

Here is the call graph for this function:

Here is the caller graph for this function:

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 34 of file faceIntersection.C.

References Foam::average(), PointHit< Point >::distance(), PointHit< Point >::eligibleMiss(), f(), PointHit< Point >::hit(), PointHit< Point >::hitPoint(), face::intersection(), Foam::mag(), PointHit< Point >::missPoint(), nPoints, points, PointHit< Point >::setDistance(), PointHit< Point >::setHit(), PointHit< Point >::setMiss(), and PointHit< Point >::setPoint().

Referenced by treeDataCell::findIntersectOp::operator()(), AMIInterpolation< SourcePatch, TargetPatch >::srcPointFace(), and AMIInterpolation< SourcePatch, TargetPatch >::tgtPointFace().

Here is the call graph for this function:

Here is the caller graph for this function:

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< Point >::distance(), f(), forAll, PointHit< Point >::hit(), PointHit< Point >::hitPoint(), Foam::mag(), face::nearestPoint(), PointHit< Point >::setDistance(), PointHit< Point >::setHit(), and PointHit< Point >::setPoint().

Referenced by face::ray().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return nearest point to face.

Definition at line 199 of file faceIntersection.C.

References face::nearestPointClassify().

Referenced by meshSearch::findNearestBoundaryFace(), sampledSet::findNearFace(), treeDataFace::getVolumeType(), face::intersection(), nearestEqOp::operator()(), treeDataFace::findNearestOp::operator()(), and treeDataPrimitivePatch< PatchType >::overlaps().

Here is the call graph for this function:

Here is the caller graph for this function:

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< Point >::distance(), f(), PointHit< Point >::hit(), PointHit< Point >::hitPoint(), Foam::mag(), PointHit< Point >::missPoint(), triangle< Point, PointRef >::nearestPointClassify(), nPoints, PointHit< Point >::setDistance(), PointHit< Point >::setHit(), PointHit< Point >::setMiss(), PointHit< Point >::setPoint(), and List< T >::size().

Referenced by face::nearestPoint().

Here is the call graph for this function:

Here is the caller graph for this function:

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(), and p.

Here is the call graph for this function:

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

Return area in contact, given the displacement in vertices.

Definition at line 33 of file faceAreaInContact.C.

References forAll, face::mag(), Foam::mag(), List< T >::setSize(), and List< T >::size().

Here is the call graph for this function:

Foam::label nEdges ( ) const
inline

Return number of edges.

Definition at line 103 of file faceI.H.

References List< label >::size().

Here is the call graph for this function:

Foam::edgeList edges ( ) const

Return edges in face point ordering,.

i.e. edges()[0] is edge between [0] and [1]

Definition at line 761 of file face.C.

References Foam::e, UList< T >::last(), and List< T >::size().

Referenced by edgeCollapser::checkMeshQuality(), Foam::longestEdge(), and Foam::polyMeshZipUpCells().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return n-th face edge.

Definition at line 110 of file faceI.H.

References UList< T >::fcIndex().

Here is the call graph for this function:

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 edge::end(), UList< T >::fcIndex(), forAll, UList< T >::rcIndex(), and edge::start().

Referenced by particle< Type >::crossEdgeConnectedFace(), patchFaceOrientation::updateEdge(), and patchFaceOrientation::updateFace().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::label nTriangles ( ) const
inline

Number of triangles after splitting.

Definition at line 130 of file faceI.H.

References List< label >::size().

Referenced by boundaryMesh::getNTris(), face::nTriangles(), STARCDsurfaceFormat< Face >::read(), boundaryMesh::triangulate(), triSurfaceTools::triangulate(), and boundaryMesh::triangulateLocal().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::label nTriangles ( const pointField points) const

Number of triangles after splitting.

Definition at line 822 of file face.C.

References face::nTriangles(), and face::triangles().

Here is the call graph for this function:

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 829 of file face.C.

References face::nTrianglesQuads().

Referenced by faceAreaIntersect::calc(), NamedEnum< Enum, nEnum >::names(), face::nTriangles(), STARCDsurfaceFormat< Face >::read(), boundaryMesh::triangulate(), MeshedSurface< Face >::triangulate(), triSurfaceTools::triangulate(), boundaryMesh::triangulateLocal(), and patchInjectionBase::updateMesh().

Here is the call graph for this function:

Here is the caller graph for this function:

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.

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 843 of file face.C.

References face::trianglesQuads().

Referenced by face::triangles().

Here is the call graph for this function:

Here is the caller graph for this function:

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 857 of file face.C.

Referenced by face::nTrianglesQuads().

Here is the caller graph for this function:

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, ConstCirculator< ContainerType >::circulate(), CirculatorBase::CLOCKWISE, ConstCirculator< ContainerType >::setFulcrumToIterator(), ConstCirculator< ContainerType >::setIteratorToFulcrum(), and List< T >::size().

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

Here is the call graph for this function:

Here is the caller graph for this function:

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 forAll, and List< T >::size().

Here is the call graph for this function:

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

Definition at line 33 of file faceTemplates.C.

References face::average(), DynamicList< T, SizeInc, SizeMult, SizeDiv >::setSize(), and List< T >::size().

Here is the call graph for this function:

Friends And Related Function Documentation

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

Member Data Documentation

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: