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 PointField< Type > &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 64 of file pointMVCWeight.H.
pointMVCWeight | ( | const polyMesh & | mesh, |
const vector & | position, | ||
const label | celli, | ||
const label | facei = -1 |
||
) |
Construct from components.
Definition at line 245 of file pointMVCWeight.C.
References pointMVCWeight::calcWeights(), pointMVCWeight::cellIndex_, primitiveMesh::cellPoints(), f(), polyMesh::faces(), forAll, HashTable< T, Key, Hash >::insert(), Foam::mag(), Foam::pid(), polyMesh::points(), List< T >::setSize(), List< T >::size(), pointMVCWeight::tol, and pointMVCWeight::weights_.
|
protected |
Calculate weights from single face's vertices only.
Definition at line 41 of file pointMVCWeight.C.
References Foam::asin(), f(), forAll, Foam::mag(), Foam::pid(), List< T >::setSize(), HashTable< T, Key, Hash >::size(), Foam::tan(), pointMVCWeight::tol, and pointMVCWeight::weights().
Referenced by pointMVCWeight::pointMVCWeight().
|
protected |
Calculate weights from all cell's vertices.
Definition at line 82 of file pointMVCWeight.C.
References alpha(), Foam::asin(), primitiveMesh::cells(), f(), polyMesh::faces(), UList< T >::fcIndex(), forAll, Foam::mag(), Foam::min(), Foam::pid(), UList< T >::rcIndex(), Foam::sin(), Foam::sum(), Foam::tan(), and VectorSpace< Form, Cmpt, Ncmpts >::zero.
ClassName | ( | "pointMVCWeight" | ) |
Type information.
|
inline |
Interpolation weights (in order of cellPoints)
Definition at line 132 of file pointMVCWeight.H.
References pointMVCWeight::weights_.
Referenced by pointMVCWeight::calcWeights().
|
inline |
Interpolate field.
Definition at line 31 of file pointMVCWeightI.H.
References pointMVCWeight::cellIndex_, forAll, DimensionedField< Type, GeoMesh >::mesh(), Foam::vertices(), pointMVCWeight::weights_, and Foam::Zero.
Referenced by interpolationPointMVC< Type >::interpolate().
|
protected |
Cell index.
Definition at line 71 of file pointMVCWeight.H.
Referenced by pointMVCWeight::cell(), pointMVCWeight::interpolate(), and pointMVCWeight::pointMVCWeight().
|
protected |
Weights applied to cell vertices.
Definition at line 74 of file pointMVCWeight.H.
Referenced by pointMVCWeight::interpolate(), pointMVCWeight::pointMVCWeight(), and pointMVCWeight::weights().
|
static |
Tolerance used in calculating barycentric co-ordinates.
(applied to normalised values)
Definition at line 108 of file pointMVCWeight.H.
Referenced by pointMVCWeight::calcWeights(), and pointMVCWeight::pointMVCWeight().