32 #ifndef CentredFitSnGradScheme_H
33 #define CentredFitSnGradScheme_H
53 template<
class Type,
class Polynomial,
class Stencil>
63 const scalar linearLimitFactor_;
66 const scalar centralWeight_;
159 #define makeCentredFitSnGradTypeScheme(SS, POLYNOMIAL, STENCIL, TYPE) \
160 typedef Foam::fv::CentredFitSnGradScheme \
161 <Foam::TYPE, Foam::POLYNOMIAL, Foam::STENCIL> \
162 CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_; \
164 defineTemplateTypeNameAndDebugWithName \
165 (CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
171 snGradScheme<TYPE>::addMeshConstructorToTable \
172 <CentredFitSnGradScheme<TYPE, POLYNOMIAL, STENCIL>> \
173 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
177 #define makeCentredFitSnGradScheme(SS, POLYNOMIAL, STENCIL) \
179 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
180 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
181 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,sphericalTensor) \
182 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
183 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
Data for centred fit snGrad schemes.
const List< scalarList > & coeffs() const
Return reference to fit coefficients.
static Type & New(const word &name, const Mesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
tmp< SurfaceField< Type > > weightedSum(const VolField< Type > &fld, const List< List< scalar >> &stencilWeights) const
Sum vol field contributions to create face values.
Mesh data needed to do the Finite Volume discretisation.
Centred fit snGrad scheme which applies an explicit correction to snGrad.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
virtual tmp< surfaceScalarField > deltaCoeffs(const VolField< Type > &) const
Return the interpolation weighting factors for the given field.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &vf) const
Return the explicit correction to the face-interpolate.
TypeName("CentredFitSnGradScheme")
Runtime type information.
Abstract base class for snGrad schemes.
const fvMesh & mesh() const
Return mesh reference.
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionSet dimLength
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.