Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coordinates. More...
Public Member Functions | |
pointMVCWeight (const polyMesh &mesh, const vector &position, const label celli, const label facei=-1) | |
Construct from components. More... | |
label | cell () const noexcept |
Cell index. More... | |
const scalarField & | weights () const noexcept |
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 int | debug |
Debug switch. More... | |
static scalar | tol |
Tolerance used in calculating barycentric coordinates. 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 239 of file pointMVCWeight.C.
References f(), forAll, Foam::invertToMap(), Foam::mag(), mesh, Foam::pid(), and UList< T >::size().
|
protected |
Calculate weights from single face's vertices only.
Definition at line 37 of file pointMVCWeight.C.
References Foam::asin(), f(), forAll, Foam::mag(), Foam::pid(), List< T >::setSize(), HashTable< T, Key, Hash >::size(), and Foam::tan().
|
protected |
Calculate weights from all cell's vertices.
Definition at line 78 of file pointMVCWeight.C.
References Foam::constant::atomic::alpha, Foam::asin(), f(), forAll, Foam::mag(), mesh, Foam::min(), Foam::normalised(), Foam::pid(), Foam::sin(), Foam::sum(), Foam::tan(), and Foam::Zero.
|
inlinenoexcept |
|
inlinenoexcept |
Interpolation weights (in order of cellPoints)
Definition at line 151 of file pointMVCWeight.H.
References pointMVCWeight::weights_.
|
inline |
Interpolate field.
Definition at line 27 of file pointMVCWeightI.H.
References forAll, DimensionedField< Type, GeoMesh >::mesh(), Foam::vertices(), and Foam::Zero.
Referenced by interpolationPointMVC< Type >::interpolate(), and surfaceAlignedSBRStressFvMotionSolver::solve().
|
protected |
|
protected |
Weights applied to cell vertices.
Definition at line 78 of file pointMVCWeight.H.
Referenced by pointMVCWeight::weights().
|
static |
Debug switch.
Definition at line 114 of file pointMVCWeight.H.
|
static |
Tolerance used in calculating barycentric coordinates.
(applied to normalised values)
Definition at line 121 of file pointMVCWeight.H.