#include <pointVolInterpolation.H>
Public Member Functions | |
| ClassName ("pointVolInterpolation") | |
| pointVolInterpolation (const pointMesh &, const fvMesh &) | |
| Construct given pointMesh and fvMesh. | |
| ~pointVolInterpolation () | |
| Destructor. | |
| const FieldField< Field, scalar > & | volWeights () const |
| Return reference to weights arrays. | |
| void | updateTopology () |
| Update mesh topology using the morph engine. | |
| bool | movePoints () |
| Correct weighting factors for moving mesh. | |
| template<class Type> | |
| void | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const |
| Interpolate from pointField to volField. | |
| template<class Type> | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &) const |
| Interpolate pointField returning volField. | |
| template<class Type> | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | interpolate (const tmp< GeometricField< Type, pointPatchField, pointMesh > > &) const |
| Interpolate tmp<pointField> returning volField. | |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolate (const GeometricField< Type, pointPatchField, pointMesh > &pf) const |
| template<class Type> | |
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > | interpolate (const tmp< GeometricField< Type, pointPatchField, pointMesh > > &tpf) const |
Protected Member Functions | |
| const pointMesh & | pMesh () const noexcept |
| const fvMesh & | vMesh () const noexcept |
Definition at line 59 of file pointVolInterpolation.H.
Construct given pointMesh and fvMesh.
Definition at line 163 of file pointVolInterpolation.C.
Destructor.
Definition at line 176 of file pointVolInterpolation.C.
|
inlineprotectednoexcept |
Definition at line 103 of file pointVolInterpolation.H.
References Foam::noexcept.
|
inlineprotectednoexcept |
Definition at line 105 of file pointVolInterpolation.H.
References Foam::noexcept.
Referenced by interpolate(), and interpolate().

| ClassName | ( | "pointVolInterpolation" | ) |
| const Foam::FieldField< Foam::Field, Foam::scalar > & volWeights | ( | ) | const |
Return reference to weights arrays.
This also constructs the weighting factors if necessary.
Definition at line 187 of file pointVolInterpolation.C.
Referenced by interpolate().

| void updateTopology | ( | ) |
Update mesh topology using the morph engine.
Definition at line 200 of file pointVolInterpolation.C.
| bool movePoints | ( | ) |
Correct weighting factors for moving mesh.
Definition at line 208 of file pointVolInterpolation.C.
| void interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | pf, |
| GeometricField< Type, fvPatchField, volMesh > & | vf ) const |
Interpolate from pointField to volField.
using inverse distance weighting
Definition at line 32 of file pointVolInterpolate.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), DebugInFunction, Foam::endl(), forAll, DimensionedField< Type, GeoMesh >::mesh(), pi(), Foam::type(), vMesh(), and volWeights().
Referenced by volumetricBSplinesDesignVariables::dCdb().


| tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | ) | const |
Interpolate pointField returning volField.
using inverse distance weighting
| tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate | ( | const tmp< GeometricField< Type, pointPatchField, pointMesh > > & | ) | const |
Interpolate tmp<pointField> returning volField.
using inverse distance weighting
| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolate | ( | const GeometricField< Type, pointPatchField, pointMesh > & | pf | ) | const |
Definition at line 88 of file pointVolInterpolate.C.
References IOobject::db(), DimensionedField< Type, GeoMesh >::dimensions(), IOobject::instance(), Foam::interpolate(), IOobject::name(), tmp< T >::ref(), and vMesh().

| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > interpolate | ( | const tmp< GeometricField< Type, pointPatchField, pointMesh > > & | tpf | ) | const |
Definition at line 118 of file pointVolInterpolate.C.
References Foam::interpolate().
