Go to the documentation of this file.
35 #ifndef multivariateSurfaceInterpolationScheme_H
36 #define multivariateSurfaceInterpolationScheme_H
61 public HashTable<const GeometricField<Type, fvPatchField, volMesh>*>
101 virtual const word&
type()
const = 0;
210 #define makeMultivariateSurfaceInterpolationTypeScheme(SS, Type) \
212 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
214 multivariateSurfaceInterpolationScheme<Type>:: \
215 addIstreamConstructorToTable<SS<Type> > \
216 add##SS##Type##ConstructorToTable_;
219 #define makeMultivariateSurfaceInterpolationScheme(SS) \
221 makeMultivariateSurfaceInterpolationTypeScheme(SS, scalar) \
222 makeMultivariateSurfaceInterpolationTypeScheme(SS, vector) \
223 makeMultivariateSurfaceInterpolationTypeScheme(SS, sphericalTensor) \
224 makeMultivariateSurfaceInterpolationTypeScheme(SS, symmTensor) \
225 makeMultivariateSurfaceInterpolationTypeScheme(SS, tensor)
const fvMesh & mesh_
Hold reference to mesh.
multivariateSurfaceInterpolationScheme(const multivariateSurfaceInterpolationScheme &)
Disallow default bitwise copy construct.
A class for handling words, derived from string.
A class for managing temporary objects.
Reference counter for various OpenFOAM components.
virtual const word & type() const =0
Runtime type information.
fieldScheme(const GeometricField< Type, fvPatchField, volMesh > &field)
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &field) const =0
Return the interpolation weighting factors.
void add(const GeometricField< Type, fvPatchField, volMesh > &f)
surfaceInterpolationScheme sub-class returned by operator(field)
Abstract base class for surface interpolation schemes.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
declareRunTimeSelectionTable(tmp, multivariateSurfaceInterpolationScheme, Istream,(const fvMesh &mesh, const fieldTable &fields, const surfaceScalarField &faceFlux, Istream &is),(mesh, fields, faceFlux, is))
Abstract base class for multi-variate surface interpolation schemes.
virtual ~multivariateSurfaceInterpolationScheme()
Destructor.
const fvMesh & mesh() const
Return mesh reference.
Mesh data needed to do the Finite Volume discretisation.
static tmp< multivariateSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, const fieldTable &fields, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new gradScheme created on freestore.
An STL-conforming hash table.
const fieldTable & fields() const
Return fields to be interpolated.
void operator=(const multivariateSurfaceInterpolationScheme &)
Disallow default bitwise assignment.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
Generic GeometricField class.
const fieldTable & fields_
HashTable of pointers to the field set.