38#ifndef CentredFitScheme_H
39#define CentredFitScheme_H
53template<
class Type,
class Polynomial,
class Stencil>
63 const scalar linearLimitFactor_;
66 const scalar centralWeight_;
90 linearLimitFactor_(readScalar(is)),
104 linearLimitFactor_(readScalar(is)),
155#define makeCentredFitSurfaceInterpolationTypeScheme\
163typedef CentredFitScheme<TYPE, POLYNOMIAL, STENCIL> \
164 CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
165defineTemplateTypeNameAndDebugWithName \
166 (CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
168surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
169<CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
170 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
172surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
173<CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
174 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
176#define makeCentredFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
178makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
179makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
180makeCentredFitSurfaceInterpolationTypeScheme \
187makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor)\
188makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
Data for the quadratic fit correction interpolation scheme.
const List< scalarList > & coeffs() const
Return reference to fit coefficients.
Centred fit surface interpolation scheme which applies an explicit correction to linear.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
CentredFitScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
CentredFitScheme(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
TypeName("CentredFitScheme")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar > > &stencilWeights) const
Sum vol field contributions to create face values.
Mesh data needed to do the Finite Volume discretisation.
Central-differencing interpolation scheme class.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.