Foam::cellPointWeight. More...


Public Member Functions | |
| cellPointWeight (const polyMesh &mesh, const vector &position, const label celli, const label facei=-1) | |
| Construct from components. More... | |
| label | cell () const |
| Cell index. More... | |
| const List< scalar > & | weights () const |
| Interpolation weights. More... | |
| const List< label > & | faceVertices () const |
| Interpolation addressing for points on face. More... | |
Static Public Attributes | |
| static int | debug |
| Debug switch. More... | |
| static scalar | tol |
| Tolerance used in calculating barycentric co-ordinates. More... | |
Protected Member Functions | |
| void | findTetrahedron (const polyMesh &mesh, const vector &position, const label celli) |
| void | findTriangle (const polyMesh &mesh, const vector &position, const label facei) |
Protected Attributes | |
| const label | celli_ |
| Cell index. More... | |
| List< scalar > | weights_ |
| Weights applied to tet vertices. More... | |
| List< label > | faceVertices_ |
| Face vertex indices. More... | |
Definition at line 50 of file cellPointWeight.H.
| cellPointWeight | ( | const polyMesh & | mesh, |
| const vector & | position, | ||
| const label | celli, | ||
| const label | facei = -1 |
||
| ) |
Construct from components.
Definition at line 256 of file cellPointWeight.C.
Referenced by cellPointWeight::findTriangle().

|
protected |
Definition at line 39 of file cellPointWeight.C.
References tetrahedron< Point, PointRef >::barycentric(), primitiveMesh::cellVolumes(), Foam::det(), PointHit< Point >::distance(), Foam::endl(), f(), tetIndices::face(), tetIndices::faceBasePt(), tetIndices::facePtA(), tetIndices::facePtB(), polyMesh::faces(), cellPointWeight::findTriangle(), forAll, Foam::mag(), tetrahedron< Point, PointRef >::nearestPoint(), Foam::nl, pFaces, Foam::Pout, and tetIndices::tet().

Definition at line 141 of file cellPointWeight.C.
References triangle< Point, PointRef >::barycentric(), cellPointWeight::cellPointWeight(), Foam::det(), PointHit< Point >::distance(), Foam::endl(), f(), primitiveMesh::faceAreas(), tetIndices::faceBasePt(), polyMesh::faceOwner(), tetIndices::facePtA(), tetIndices::facePtB(), polyMesh::faces(), tetIndices::faceTri(), forAll, Foam::mag(), Foam::magSqr(), triangle< Point, PointRef >::nearestPoint(), Foam::nl, and Foam::Pout.
Referenced by cellPointWeight::findTetrahedron().


Cell index.
Definition at line 108 of file cellPointWeight.H.
References cellPointWeight::celli_.
Referenced by interpolationCellPointWallModified< Type >::interpolate(), and interpolationCellPoint< Type >::interpolate().

|
inline |
Interpolation weights.
Definition at line 114 of file cellPointWeight.H.
References cellPointWeight::weights_.
Referenced by interpolationCellPointWallModified< Type >::interpolate(), and interpolationCellPoint< Type >::interpolate().

Interpolation addressing for points on face.
Definition at line 120 of file cellPointWeight.H.
References cellPointWeight::faceVertices_.
Referenced by interpolationCellPointWallModified< Type >::interpolate(), and interpolationCellPoint< Type >::interpolate().

|
protected |
|
protected |
Weights applied to tet vertices.
Definition at line 60 of file cellPointWeight.H.
Referenced by cellPointWeight::weights().
Face vertex indices.
Definition at line 63 of file cellPointWeight.H.
Referenced by cellPointWeight::faceVertices().
|
static |
Debug switch.
Definition at line 86 of file cellPointWeight.H.
|
static |
Tolerance used in calculating barycentric co-ordinates.
(applied to normalised values)
Definition at line 90 of file cellPointWeight.H.
1.8.11