Go to the documentation of this file.
37 #ifndef surfaceInterpolationScheme_H
38 #define surfaceInterpolationScheme_H
80 TypeName(
"surfaceInterpolationScheme");
107 (
mesh, faceFlux, schemeData)
163 template<
class SFType>
275 const GeometricField<scalar, fvPatchField, volMesh>&
287 #define makeSurfaceInterpolationTypeScheme(SS, Type) \
289 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
291 surfaceInterpolationScheme<Type>::addMeshConstructorToTable<SS<Type>> \
292 add##SS##Type##MeshConstructorToTable_; \
294 surfaceInterpolationScheme<Type>::addMeshFluxConstructorToTable<SS<Type>> \
295 add##SS##Type##MeshFluxConstructorToTable_;
297 #define makeSurfaceInterpolationScheme(SS) \
299 makeSurfaceInterpolationTypeScheme(SS, scalar) \
300 makeSurfaceInterpolationTypeScheme(SS, vector) \
301 makeSurfaceInterpolationTypeScheme(SS, sphericalTensor) \
302 makeSurfaceInterpolationTypeScheme(SS, symmTensor) \
303 makeSurfaceInterpolationTypeScheme(SS, tensor)
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
Return the face-interpolate of the given cell field.
A class for managing temporary objects.
Reference counter for various OpenFOAM components.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
Mesh data needed to do the Finite Volume discretisation.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
TypeName("surfaceInterpolationScheme")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
declareRunTimeSelectionTable(tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
Mesh data needed to do the Finite Volume discretisation.
surfaceInterpolationScheme(const fvMesh &mesh)
Construct from mesh.
Macros to ease declaration of run-time selection tables.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
virtual ~surfaceInterpolationScheme()=default
Destructor.