Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coordinates. More...
Public Member Functions | |
ClassName ("pointMVCWeight") | |
Type information. More... | |
pointMVCWeight (const polyMesh &mesh, const vector &position, const label celli, const label facei=-1) | |
Construct from components. More... | |
label | cell () const |
Cell index. More... | |
const scalarField & | weights () const |
Interpolation weights (in order of cellPoints) More... | |
template<class Type > | |
Type | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &psip) const |
Interpolate field. More... | |
Static Public Attributes | |
static scalar | tol |
Tolerance used in calculating barycentric co-ordinates. More... | |
Protected Member Functions | |
void | calcWeights (const Map< label > &toLocal, const face &f, const DynamicList< point > &u, const scalarField &dist, scalarField &weights) const |
Calculate weights from single face's vertices only. More... | |
void | calcWeights (const polyMesh &mesh, const labelList &toGlobal, const Map< label > &toLocal, const vector &position, const vectorField &uVec, const scalarField &dist, scalarField &weights) const |
Calculate weights from all cell's vertices. More... | |
Protected Attributes | |
const label | cellIndex_ |
Cell index. More... | |
scalarField | weights_ |
Weights applied to cell vertices. More... | |
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coordinates.
Based on (VTK's vtkMeanValueCoordinatesInterpolator's) implementation of "Spherical Barycentric Coordinates" 2006 paper Eurographics Symposium on Geometry Processing by Torsten Langer, Alexander Belyaev and Hans-Peter Seide
Definition at line 67 of file pointMVCWeight.H.
pointMVCWeight | ( | const polyMesh & | mesh, |
const vector & | position, | ||
const label | celli, | ||
const label | facei = -1 |
||
) |
Construct from components.
Definition at line 246 of file pointMVCWeight.C.
References primitiveMesh::cellPoints(), polyMesh::faces(), forAll, HashTable< T, Key, Hash >::insert(), Foam::mag(), Foam::pid(), polyMesh::points(), and List< T >::size().
Referenced by pointMVCWeight::calcWeights().
|
protected |
Calculate weights from single face's vertices only.
Definition at line 42 of file pointMVCWeight.C.
References Foam::asin(), UList< T >::fcIndex(), forAll, Foam::mag(), Foam::pid(), List< T >::setSize(), List< T >::size(), HashTable< T, Key, Hash >::size(), and Foam::tan().
|
protected |
Calculate weights from all cell's vertices.
Definition at line 83 of file pointMVCWeight.C.
References alpha(), Foam::asin(), primitiveMesh::cells(), polyMesh::faces(), UList< T >::fcIndex(), forAll, Foam::mag(), Foam::min(), Foam::pid(), pointMVCWeight::pointMVCWeight(), UList< T >::rcIndex(), Foam::sin(), Foam::sum(), Foam::tan(), and VectorSpace< Vector< scalar >, scalar, 3 >::zero.
ClassName | ( | "pointMVCWeight" | ) |
Type information.
|
inline |
Interpolation weights (in order of cellPoints)
Definition at line 135 of file pointMVCWeight.H.
References pointMVCWeight::interpolate(), and pointMVCWeight::weights_.
|
inline |
Interpolate field.
Definition at line 32 of file pointMVCWeightI.H.
References forAll, DimensionedField< Type, GeoMesh >::mesh(), Foam::vertices(), and Foam::Zero.
Referenced by interpolationPointMVC< Type >::interpolate(), and pointMVCWeight::weights().
|
protected |
|
protected |
Weights applied to cell vertices.
Definition at line 77 of file pointMVCWeight.H.
Referenced by pointMVCWeight::weights().
|
static |
Tolerance used in calculating barycentric co-ordinates.
(applied to normalised values)
Definition at line 111 of file pointMVCWeight.H.