Go to the documentation of this file.
32 #ifndef CentredFitSnGradScheme_H
33 #define CentredFitSnGradScheme_H
53 template<
class Type,
class Polynomial,
class Stencil>
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)
const scalar linearLimitFactor_
Factor the fit is allowed to deviate from linear.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
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.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
static const CentredFitSnGradData< Polynomial > & New(const fvMesh &mesh)
TypeName("CentredFitSnGradScheme")
Runtime type information.
const scalar centralWeight_
Weights for central stencil.
CentredFitSnGradScheme(const CentredFitSnGradScheme &)
Disallow default bitwise copy construct.
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.
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.
void operator=(const CentredFitSnGradScheme &)
Disallow default bitwise assignment.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Abstract base class for snGrad schemes.
Generic GeometricField class.
const List< scalarList > & coeffs() const
Return reference to fit coefficients.