Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coordinates. More...
#include <pointMVCWeight.H>

Public Member Functions | |
| pointMVCWeight (const polyMesh &mesh, const vector &position, const label celli, const label facei=-1) | |
| Construct from components. | |
| label | cell () const noexcept |
| Cell index. | |
| const scalarField & | weights () const noexcept |
| Interpolation weights (in order of cellPoints). | |
| template<class Type> | |
| Type | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &psip) const |
| Interpolate field. | |
Static Public Attributes | |
| static int | debug |
| Debug switch. | |
| static scalar | tol |
| Tolerance used in calculating barycentric coordinates. | |
Protected Member Functions | |
| void | calcWeights (const Map< label > &toLocal, const face &f, const UList< point > &u, const scalarField &dist, scalarField &weights) const |
| Calculate weights from single face's vertices only. | |
| 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. | |
Protected Attributes | |
| const label | cellIndex_ |
| Cell index. | |
| scalarField | weights_ |
| Weights applied to cell vertices. | |
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 238 of file pointMVCWeight.C.
References calcWeights(), cellIndex_, f(), forAll, Foam::invertToMap(), Foam::mag(), mesh, Foam::pid(), UList< T >::size(), tol, and weights_.

|
protected |
Calculate weights from single face's vertices only.
Definition at line 36 of file pointMVCWeight.C.
References Foam::asin(), f(), forAll, Foam::mag(), Foam::pid(), HashTable< T, Key, Hash >::size(), Foam::tan(), tol, and weights().
Referenced by calcWeights(), and pointMVCWeight().


|
protected |
Calculate weights from all cell's vertices.
Definition at line 77 of file pointMVCWeight.C.
References alpha, Foam::asin(), calcWeights(), cellIndex_, f(), forAll, Foam::mag(), mesh, Foam::min(), Foam::normalised(), Foam::pid(), Foam::sin(), Foam::sum(), Foam::tan(), tol, weights(), and Foam::Zero.

|
inlinenoexcept |
Cell index.
Definition at line 143 of file pointMVCWeight.H.
References cellIndex_, and Foam::noexcept.
|
inlinenoexcept |
Interpolation weights (in order of cellPoints).
Definition at line 151 of file pointMVCWeight.H.
References Foam::noexcept, and weights_.
Referenced by calcWeights(), and calcWeights().

|
inline |
Interpolate field.
Definition at line 26 of file pointMVCWeightI.H.
References cellIndex_, forAll, DimensionedField< Type, GeoMesh >::mesh(), Foam::vertices(), weights_, and Foam::Zero.
Referenced by interpolationPointMVC< Type >::interpolate(), and surfaceAlignedSBRStressFvMotionSolver::solve().


|
protected |
Cell index.
Definition at line 73 of file pointMVCWeight.H.
Referenced by calcWeights(), cell(), interpolate(), and pointMVCWeight().
|
protected |
Weights applied to cell vertices.
Definition at line 78 of file pointMVCWeight.H.
Referenced by interpolate(), pointMVCWeight(), and 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.
Referenced by calcWeights(), calcWeights(), and pointMVCWeight().