Go to the documentation of this file.
43 #ifndef pointMVCWeight_H
44 #define pointMVCWeight_H
60 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
123 const label faceI = -1
pointMVCWeight(const polyMesh &mesh, const vector &position, const label cellI, const label faceI=-1)
Construct from components.
static int debug
Debug switch.
scalarField weights_
Weights applied to cell vertices.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Type interpolate(const GeometricField< Type, pointPatchField, pointMesh > &psip) const
Interpolate field.
Abstract base class for point-mesh patch fields.
Mesh consisting of general polyhedral cells.
static scalar tol
Tolerance used in calculating barycentric co-ordinates.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coor...
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Mesh representing a set of points created from polyMesh.
const label cellIndex_
Cell index.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
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.
label cell() const
Cell index.
const scalarField & weights() const
Interpolation weights (in order of cellPoints)
A face is a list of labels corresponding to mesh vertices.
Generic GeometricField class.