Interpolation of cell-based displacements to the points with additional correction for interpolation inconsistency on patches. More...
Public Member Functions | |
TypeName ("patchCorrected") | |
Runtime type information. More... | |
patchCorrectedInterpolation (const fvMesh &mesh, Istream &entry) | |
Construct from an fvMesh and an Istream. More... | |
virtual | ~patchCorrectedInterpolation () |
Destructor. More... | |
virtual void | interpolate (const volScalarField &, pointScalarField &) const |
Interpolate the given scalar cell displacement. More... | |
virtual void | interpolate (const volVectorField &, pointVectorField &) const |
Interpolate the given vector cell displacement. More... | |
![]() | |
TypeName ("motionInterpolation") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, motionInterpolation, Istream,(const fvMesh &mesh, Istream &entry),(mesh, entry)) | |
motionInterpolation (const fvMesh &mesh) | |
Construct from an fvMesh. More... | |
motionInterpolation (const fvMesh &mesh, Istream &entry) | |
Construct from an fvMesh and an Istream. More... | |
virtual | ~motionInterpolation () |
Destructor. More... | |
const fvMesh & | mesh () const |
Return const-refernce to the mesh. More... | |
Private Member Functions | |
labelListList | getPatchGroups (Istream &entry) const |
Get patch groups from the input stream. More... | |
template<class Type > | |
void | interpolateType (const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const |
Interpolate the given cell displacement. More... | |
template<class Type > | |
void | interpolateDataFromPatchGroups (GeometricField< Type, pointPatchField, pointMesh > &) const |
Interpolate patch data by inverse distance weighting. More... | |
template<class Type > | |
void | propagateDataFromPatchGroup (const label, pointScalarField &, GeometricField< Type, pointPatchField, pointMesh > &) const |
Propagate data from a number of patches into the field. More... | |
Private Attributes | |
const labelListList | patchGroups_ |
Patch groups from which to propagate the displacement. More... | |
Additional Inherited Members | |
![]() | |
static autoPtr< motionInterpolation > | New (const fvMesh &mesh) |
Select default. More... | |
static autoPtr< motionInterpolation > | New (const fvMesh &mesh, Istream &entry) |
Select from stream. More... | |
Interpolation of cell-based displacements to the points with additional correction for interpolation inconsistency on patches.
The default interpolation method interpolates from the cells to all points except those on boundaries with value boundary conditions. The discrepancy across the boundary cell can cause shearing and inversion if the cells are of very high aspect ratio.
This method applies the default interpolation to *all* points, including those on value boundaries. The difference between the interpolated value on the boundary and the desired boundary condition is then propagated into the mesh with a wave. Contributions from different patches are inverse-distance weighted, and the result is added to the default interpolation. The result of this is that thin boundary cells are maintained much more accurately for non-uniform patch displacements.
The user must specify the patch groups from which to propagate the motion. Ideally, these groups will be opposing; i.e. one group with the aerofoil, and one with the far field:
interpolation patchCorrected ( (aerofoilUpper aerofoilLower) (farField) );
Definition at line 74 of file patchCorrectedInterpolation.H.
patchCorrectedInterpolation | ( | const fvMesh & | mesh, |
Istream & | entry | ||
) |
Construct from an fvMesh and an Istream.
Definition at line 83 of file patchCorrectedInterpolation.C.
|
virtual |
Destructor.
Definition at line 95 of file patchCorrectedInterpolation.C.
|
private |
Get patch groups from the input stream.
Definition at line 47 of file patchCorrectedInterpolation.C.
References Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, List::resize(), and List::size().
|
private |
Interpolate the given cell displacement.
Definition at line 36 of file patchCorrectedInterpolationTemplates.C.
References GeometricField::boundaryField(), GeometricField::correctBoundaryConditions(), GeometricField::internalField(), mesh, Foam::compressible::New(), and timeName.
|
private |
Interpolate patch data by inverse distance weighting.
Definition at line 92 of file patchCorrectedInterpolationTemplates.C.
References Foam::dimless, forAll, GeometricField::internalField(), Foam::max(), mesh, Foam::sqr(), and timeName.
|
private |
Propagate data from a number of patches into the field.
Definition at line 152 of file patchCorrectedInterpolationTemplates.C.
References PointEdgeWave< Type, TrackingData >::data(), Foam::distance(), forAll, mesh, pointPatch::meshPoints(), nPoints, pointPatchField::patch(), pointPatchField::patchInternalField(), points, Foam::sqrt(), and pointPatchField::updateCoeffs().
TypeName | ( | "patchCorrected" | ) |
Runtime type information.
|
virtual |
Interpolate the given scalar cell displacement.
Reimplemented from motionInterpolation.
Definition at line 102 of file patchCorrectedInterpolation.C.
|
virtual |
Interpolate the given vector cell displacement.
Reimplemented from motionInterpolation.
Definition at line 112 of file patchCorrectedInterpolation.C.
|
private |
Patch groups from which to propagate the displacement.
Definition at line 81 of file patchCorrectedInterpolation.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.