33 #ifndef CentredFitScheme_H 34 #define CentredFitScheme_H 48 template<
class Type,
class Polynomial,
class Stencil>
58 const scalar linearLimitFactor_;
61 const scalar centralWeight_;
150 #define makeCentredFitSurfaceInterpolationTypeScheme\ 158 typedef CentredFitScheme<TYPE, POLYNOMIAL, STENCIL> \ 159 CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \ 160 defineTemplateTypeNameAndDebugWithName \ 161 (CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \ 163 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \ 164 <CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \ 165 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \ 167 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \ 168 <CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \ 169 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_; 171 #define makeCentredFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \ 173 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \ 174 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \ 175 makeCentredFitSurfaceInterpolationTypeScheme \ 182 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor)\ 183 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor) virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Centred interpolation interpolation scheme class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
const List< scalarList > & coeffs() const
Return reference to fit coefficients.
const fvMesh & mesh() const
Return mesh reference.
TypeName("CentredFitScheme")
Runtime type information.
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.
Data for the quadratic fit correction interpolation scheme.
void operator=(const CentredFitScheme &)=delete
Disallow default bitwise assignment.
Centred fit surface interpolation scheme which applies an explicit correction to linear.
CentredFitScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
static CentredFitData< Polynomial > & New(fvMesh &mesh)
Mesh data needed to do the Finite Volume discretisation.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
A class for managing temporary objects.