Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
pointMVCWeight Class Reference

Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coordinates. More...

Collaboration diagram for pointMVCWeight:
Collaboration graph
[legend]

Public Member Functions

 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 scalarFieldweights () 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 int debug
 Debug switch. More...
 
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...
 

Detailed Description

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

Source files

Definition at line 67 of file pointMVCWeight.H.

Constructor & Destructor Documentation

◆ pointMVCWeight()

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().

Here is the call graph for this function:
Here is the caller graph for this function:

Member Function Documentation

◆ calcWeights() [1/2]

void calcWeights ( const Map< label > &  toLocal,
const face f,
const DynamicList< point > &  u,
const scalarField dist,
scalarField weights 
) const
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().

Here is the call graph for this function:

◆ calcWeights() [2/2]

void calcWeights ( const polyMesh mesh,
const labelList toGlobal,
const Map< label > &  toLocal,
const vector position,
const vectorField uVec,
const scalarField dist,
scalarField weights 
) const
protected

◆ cell()

label cell ( ) const
inline

Cell index.

Definition at line 129 of file pointMVCWeight.H.

References pointMVCWeight::cellIndex_.

◆ weights()

const scalarField& weights ( ) const
inline

Interpolation weights (in order of cellPoints)

Definition at line 135 of file pointMVCWeight.H.

References pointMVCWeight::interpolate(), and pointMVCWeight::weights_.

Here is the call graph for this function:

◆ interpolate()

Type interpolate ( const GeometricField< Type, pointPatchField, pointMesh > &  psip) const
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().

Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ cellIndex_

const label cellIndex_
protected

Cell index.

Definition at line 74 of file pointMVCWeight.H.

Referenced by pointMVCWeight::cell().

◆ weights_

scalarField weights_
protected

Weights applied to cell vertices.

Definition at line 77 of file pointMVCWeight.H.

Referenced by pointMVCWeight::weights().

◆ debug

int debug
static

Debug switch.

Definition at line 107 of file pointMVCWeight.H.

◆ tol

Foam::scalar tol
static

Tolerance used in calculating barycentric co-ordinates.

(applied to normalised values)

Definition at line 111 of file pointMVCWeight.H.


The documentation for this class was generated from the following files: