40#ifndef limitedSurfaceInterpolationScheme_H
41#define limitedSurfaceInterpolationScheme_H
55class limitedSurfaceInterpolationScheme
57 public surfaceInterpolationScheme<Type>
62 limitedSurfaceInterpolationScheme
64 const limitedSurfaceInterpolationScheme&
68 void operator=(
const limitedSurfaceInterpolationScheme&) =
delete;
81 TypeName(
"limitedSurfaceInterpolationScheme");
89 limitedSurfaceInterpolationScheme,
101 limitedSurfaceInterpolationScheme,
108 (
mesh, faceFlux, schemeData)
115 limitedSurfaceInterpolationScheme
121 surfaceInterpolationScheme<Type>(
mesh),
135 surfaceInterpolationScheme<Type>(
mesh),
208#define makelimitedSurfaceInterpolationTypeScheme(SS, Type) \
210defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
212surfaceInterpolationScheme<Type>::addMeshConstructorToTable<SS<Type>> \
213 add##SS##Type##MeshConstructorToTable_; \
215surfaceInterpolationScheme<Type>::addMeshFluxConstructorToTable<SS<Type>> \
216 add##SS##Type##MeshFluxConstructorToTable_; \
218limitedSurfaceInterpolationScheme<Type>::addMeshConstructorToTable<SS<Type>> \
219 add##SS##Type##MeshConstructorToLimitedTable_; \
221limitedSurfaceInterpolationScheme<Type>:: \
222 addMeshFluxConstructorToTable<SS<Type>> \
223 add##SS##Type##MeshFluxConstructorToLimitedTable_;
225#define makelimitedSurfaceInterpolationScheme(SS) \
227makelimitedSurfaceInterpolationTypeScheme(SS, scalar) \
228makelimitedSurfaceInterpolationTypeScheme(SS, vector) \
229makelimitedSurfaceInterpolationTypeScheme(SS, sphericalTensor) \
230makelimitedSurfaceInterpolationTypeScheme(SS, symmTensor) \
231makelimitedSurfaceInterpolationTypeScheme(SS, tensor)
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for limited surface interpolation schemes.
const surfaceScalarField & faceFlux_
limitedSurfaceInterpolationScheme(const fvMesh &mesh, const surfaceScalarField &faceFlux)
Construct from mesh and faceFlux.
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
limitedSurfaceInterpolationScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &, const surfaceScalarField &CDweights, tmp< surfaceScalarField > tLimiter) const
Return the interpolation weighting factors for the given field,.
declareRunTimeSelectionTable(tmp, limitedSurfaceInterpolationScheme, MeshFlux,(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData),(mesh, faceFlux, schemeData))
virtual ~limitedSurfaceInterpolationScheme()
Destructor.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
declareRunTimeSelectionTable(tmp, limitedSurfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
TypeName("limitedSurfaceInterpolationScheme")
Runtime type information.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.