Go to the documentation of this file.
37 #ifndef volPointInterpolation_H
38 #define volPointInterpolation_H
59 public MeshObject<fvMesh, UpdateableMeshObject, volPointInterpolation>
233 const bool overrideFixedValue
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
A class for handling words, derived from string.
ClassName("volPointInterpolation")
boolList boundaryIsPatchFace_
Per boundary face whether is on non-coupled patch.
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
A class for managing temporary objects.
volPointInterpolation(const volPointInterpolation &)
Disallow default bitwise copy construct.
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
void makeWeights()
Construct all point weighting factors.
Pre-declare SubField and related Field type.
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
void operator=(const volPointInterpolation &)
Disallow default bitwise assignment.
const word & name() const
Return name.
Mesh data needed to do the Finite Volume discretisation.
~volPointInterpolation()
Destructor.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool movePoints()
Correct weighting factors for moving mesh.
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
void makeInternalWeights(scalarField &sumWeights)
Make weights for internal and coupled-only boundarypoints.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
boolList isPatchPoint_
Per mesh(!) point whether is on non-coupled patch (on any.
scalarListList pointWeights_
Interpolation scheme weighting factor array.
Generic GeometricField class.
void interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
Interpolate from cell centres to points (vertices) using inverse distance weighting.
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...