Interpolate from cell centres to points (vertices) using inverse distance weighting. More...
Public Member Functions | |
ClassName ("volPointInterpolation") | |
volPointInterpolation (const fvMesh &) | |
Constructor given fvMesh and pointMesh. More... | |
~volPointInterpolation () | |
Destructor. More... | |
void | updateMesh (const mapPolyMesh &) |
Update mesh topology using the morph engine. More... | |
bool | movePoints () |
Correct weighting factors for moving mesh. More... | |
template<class Type > | |
tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &) const |
Interpolate volField using inverse distance weighting. More... | |
template<class Type > | |
tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh > > &) const |
Interpolate tmp<volField> using inverse distance weighting. More... | |
template<class Type > | |
tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const wordList &patchFieldTypes) const |
Interpolate volField using inverse distance weighting. More... | |
template<class Type > | |
tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh > > &, const wordList &patchFieldTypes) const |
Interpolate tmp<volField> using inverse distance weighting. More... | |
template<class Type > | |
tmp< DimensionedField< Type, pointMesh > > | interpolate (const DimensionedField< Type, volMesh > &) const |
Interpolate dimensionedField using inverse distance weighting. More... | |
template<class Type > | |
tmp< DimensionedField< Type, pointMesh > > | interpolate (const tmp< DimensionedField< Type, volMesh > > &) const |
Interpolate tmp<dimensionedField> using inverse distance. More... | |
template<class Type > | |
void | interpolateInternalField (const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const |
Interpolate internal field from volField to pointField. More... | |
template<class Type > | |
void | interpolateBoundaryField (const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const |
Interpolate boundary field without applying constraints/boundary. More... | |
template<class Type > | |
void | interpolateBoundaryField (const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideFixedValue) const |
Interpolate boundary with constraints/boundary conditions. More... | |
template<class Type > | |
void | interpolate (const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const |
Interpolate from volField to pointField. More... | |
template<class Type > | |
tmp< GeometricField< Type, pointPatchField, pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const word &name, const bool cache) const |
Interpolate volField using inverse distance weighting. More... | |
template<class Type > | |
void | interpolateDimensionedInternalField (const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const |
Interpolate dimensioned internal field from cells to points. More... | |
template<class Type > | |
tmp< DimensionedField< Type, pointMesh > > | interpolate (const DimensionedField< Type, volMesh > &, const word &name, const bool cache) const |
Interpolate dimensionedField using inverse distance weighting. More... | |
void | interpolateDisplacement (const volVectorField &, pointVectorField &) const |
Interpolate from volField to pointField. More... | |
template<class Type > | |
Foam::tmp< Foam::Field< Type > > | flatBoundaryField (const GeometricField< Type, fvPatchField, volMesh > &vf) const |
template<class Type > | |
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf, const wordList &patchFieldTypes) const |
template<class Type > | |
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh > > &tvf, const wordList &patchFieldTypes) const |
template<class Type > | |
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name, const bool cache) const |
template<class Type > | |
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf) const |
template<class Type > | |
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > | interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh > > &tvf) const |
template<class Type > | |
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > | interpolate (const DimensionedField< Type, volMesh > &vf, const word &name, const bool cache) const |
template<class Type > | |
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > | interpolate (const DimensionedField< Type, volMesh > &vf) const |
template<class Type > | |
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > | interpolate (const tmp< DimensionedField< Type, volMesh > > &tvf) const |
Public Member Functions inherited from MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation > | |
MeshObject (const fvMesh &mesh) | |
Construct on Mesh type. More... | |
virtual | ~MeshObject ()=default |
Destructor. More... | |
const fvMesh & | mesh () const |
virtual bool | writeData (Ostream &os) const |
Additional Inherited Members | |
Static Public Member Functions inherited from MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation > | |
static const volPointInterpolation & | New (const fvMesh &mesh, Args &&... args) |
Get existing or create a new MeshObject. More... | |
static bool | Delete (const fvMesh &mesh) |
Static destructor. More... | |
Protected Attributes inherited from MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation > | |
const fvMesh & | mesh_ |
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Definition at line 59 of file volPointInterpolation.H.
|
explicit |
Constructor given fvMesh and pointMesh.
Definition at line 475 of file volPointInterpolation.C.
Destructor.
Definition at line 486 of file volPointInterpolation.C.
ClassName | ( | "volPointInterpolation" | ) |
void updateMesh | ( | const mapPolyMesh & | ) |
Update mesh topology using the morph engine.
Definition at line 492 of file volPointInterpolation.C.
References mesh, and Time::New().
bool movePoints | ( | ) |
Correct weighting factors for moving mesh.
Definition at line 501 of file volPointInterpolation.C.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | ) | const |
Interpolate volField using inverse distance weighting.
Referenced by parseDriver::cellToPoint(), shapeSensitivitiesBase::getWallPointSensNormal(), shapeSensitivitiesBase::getWallPointSensNormalVec(), and shapeSensitivitiesBase::getWallPointSensVec().
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | ) | const |
Interpolate tmp<volField> using inverse distance weighting.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | , |
const wordList & | patchFieldTypes | ||
) | const |
Interpolate volField using inverse distance weighting.
returning pointField with the same patchField types. Assigns to any fixedValue boundary conditions to make them consistent with internal field
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | , |
const wordList & | patchFieldTypes | ||
) | const |
Interpolate tmp<volField> using inverse distance weighting.
returning pointField with the same patchField types. Assigns to any fixedValue boundary conditions to make them consistent with internal field
tmp< DimensionedField< Type, pointMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | ) | const |
Interpolate dimensionedField using inverse distance weighting.
tmp< DimensionedField< Type, pointMesh > > interpolate | ( | const tmp< DimensionedField< Type, volMesh > > & | ) | const |
Interpolate tmp<dimensionedField> using inverse distance.
weighting returning pointField
void interpolateInternalField | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
GeometricField< Type, pointPatchField, pointMesh > & | pf | ||
) | const |
Interpolate internal field from volField to pointField.
using inverse distance weighting
Definition at line 129 of file volPointInterpolate.C.
References Foam::endl(), forAll, DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Foam::Pout, and Foam::Zero.
void interpolateBoundaryField | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
GeometricField< Type, pointPatchField, pointMesh > & | pf | ||
) | const |
Interpolate boundary field without applying constraints/boundary.
conditions
Definition at line 273 of file volPointInterpolate.C.
References boundary, forAll, mesh, pFaces, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), pointConstraints::syncUntransformedData(), and Foam::Zero.
void interpolateBoundaryField | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
GeometricField< Type, pointPatchField, pointMesh > & | pf, | ||
const bool | overrideFixedValue | ||
) | const |
Interpolate boundary with constraints/boundary conditions.
Definition at line 340 of file volPointInterpolate.C.
References pointConstraints::constrain(), DimensionedField< Type, GeoMesh >::mesh(), and Time::New().
void interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
GeometricField< Type, pointPatchField, pointMesh > & | pf | ||
) | const |
Interpolate from volField to pointField.
using inverse distance weighting
Definition at line 357 of file volPointInterpolate.C.
References Foam::endl(), IOobject::name(), and Foam::Pout.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | , |
const word & | name, | ||
const bool | cache | ||
) | const |
Interpolate volField using inverse distance weighting.
returning pointField with name. Optionally caches
void interpolateDimensionedInternalField | ( | const DimensionedField< Type, volMesh > & | vf, |
DimensionedField< Type, pointMesh > & | pf | ||
) | const |
Interpolate dimensioned internal field from cells to points.
using inverse distance weighting
Definition at line 166 of file volPointInterpolate.C.
References primitiveMesh::cellCentres(), Foam::endl(), forAll, Foam::mag(), mesh, DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), primitiveMesh::pointCells(), points, polyMesh::points(), Foam::Pout, s(), UList< T >::size(), pointConstraints::syncUntransformedData(), and Foam::Zero.
tmp< DimensionedField< Type, pointMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | , |
const word & | name, | ||
const bool | cache | ||
) | const |
Interpolate dimensionedField using inverse distance weighting.
returning pointField with name. Optionally caches
void interpolateDisplacement | ( | const volVectorField & | vf, |
pointVectorField & | pf | ||
) | const |
Interpolate from volField to pointField.
using inverse distance weighting
Definition at line 509 of file volPointInterpolation.C.
References pointConstraints::constrainDisplacement(), DimensionedField< Type, GeoMesh >::mesh(), and Time::New().
Foam::tmp< Foam::Field< Type > > flatBoundaryField | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) | const |
Definition at line 226 of file volPointInterpolate.C.
References fvMesh::boundary(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), forAll, mesh, DimensionedField< Type, GeoMesh >::mesh(), primitiveMesh::nBoundaryFaces(), primitiveMesh::nInternalFaces(), tmp< T >::ref(), UPtrList< T >::size(), and Foam::Zero.
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
const wordList & | patchFieldTypes | ||
) | const |
Definition at line 381 of file volPointInterpolate.C.
References DimensionedField< Type, GeoMesh >::dimensions(), IOobject::instance(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), Foam::New(), GeometricField< Type, PatchField, GeoMesh >::ref(), and pointMesh::thisDb().
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | tvf, |
const wordList & | patchFieldTypes | ||
) | const |
Definition at line 414 of file volPointInterpolate.C.
References Foam::interpolate().
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
const word & | name, | ||
const bool | cache | ||
) | const |
Definition at line 430 of file volPointInterpolate.C.
References solution::cachePrintMessage(), DimensionedField< Type, GeoMesh >::dimensions(), IOobject::instance(), Foam::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), Time::New(), Foam::New(), GeometricField< Type, PatchField, GeoMesh >::ref(), regIOobject::store(), and pointMesh::thisDb().
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) | const |
Definition at line 500 of file volPointInterpolate.C.
References Foam::interpolate(), and IOobject::name().
Foam::tmp< Foam::GeometricField< Type, Foam::pointPatchField, Foam::pointMesh > > interpolate | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | tvf | ) | const |
Definition at line 511 of file volPointInterpolate.C.
References Foam::interpolate().
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | vf, |
const word & | name, | ||
const bool | cache | ||
) | const |
Definition at line 526 of file volPointInterpolate.C.
References solution::cachePrintMessage(), DimensionedField< Type, GeoMesh >::dimensions(), IOobject::instance(), Foam::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), Time::New(), Foam::New(), regIOobject::store(), and pointMesh::thisDb().
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > interpolate | ( | const DimensionedField< Type, volMesh > & | vf | ) | const |
Definition at line 597 of file volPointInterpolate.C.
References Foam::interpolate(), and IOobject::name().
Foam::tmp< Foam::DimensionedField< Type, Foam::pointMesh > > interpolate | ( | const tmp< DimensionedField< Type, volMesh > > & | tvf | ) | const |
Definition at line 608 of file volPointInterpolate.C.
References Foam::interpolate().