Go to the documentation of this file.
37 #ifndef CentredFitSnGradScheme_H
38 #define CentredFitSnGradScheme_H
58 template<
class Type,
class Polynomial,
class Stencil>
68 const scalar linearLimitFactor_;
71 const scalar centralWeight_;
95 linearLimitFactor_(readScalar(is)),
164 #define makeCentredFitSnGradTypeScheme(SS, POLYNOMIAL, STENCIL, TYPE) \
165 typedef Foam::fv::CentredFitSnGradScheme \
166 <Foam::TYPE, Foam::POLYNOMIAL, Foam::STENCIL> \
167 CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_; \
169 defineTemplateTypeNameAndDebugWithName \
170 (CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
176 snGradScheme<TYPE>::addMeshConstructorToTable \
177 <CentredFitSnGradScheme<TYPE, POLYNOMIAL, STENCIL>> \
178 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
182 #define makeCentredFitSnGradScheme(SS, POLYNOMIAL, STENCIL) \
184 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
185 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
186 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,sphericalTensor) \
187 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
188 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
CentredFitSnGradScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
TypeName("CentredFitSnGradScheme")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Mesh data needed to do the Finite Volume discretisation.
Data for centred fit snGrad schemes.
const fvMesh & mesh() const
Return mesh reference.
Abstract base class for snGrad schemes.
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.
const List< scalarList > & coeffs() const
Return reference to fit coefficients.